/* Copyright 2020 Remi Lehe * * This file is part of WarpX. * * License: BSD-3-Clause-LBNL */ #include "WarpX.H" #include "Utils/WarpXAlgorithmSelection.H" #include "FiniteDifferenceSolver.H" #ifdef WARPX_DIM_RZ # include "FiniteDifferenceAlgorithms/CylindricalYeeAlgorithm.H" #else # include "FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H" # include "FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H" # include "FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H" #endif #include using namespace amrex; /** * \brief Update the B field, over one timestep */ void FiniteDifferenceSolver::EvolveB ( std::array< std::unique_ptr, 3 >& Bfield, std::array< std::unique_ptr, 3 > const& Efield, int lev, amrex::Real const dt ) { // Select algorithm (The choice of algorithm is a runtime option, // but we compile code for each algorithm, using templates) #ifdef WARPX_DIM_RZ if (m_fdtd_algo == MaxwellSolverAlgo::Yee){ EvolveBCylindrical ( Bfield, Efield, lev, dt ); #else if (m_do_nodal) { EvolveBCartesian ( Bfield, Efield, lev, dt ); } else if (m_fdtd_algo == MaxwellSolverAlgo::Yee) { EvolveBCartesian ( Bfield, Efield, lev, dt ); } else if (m_fdtd_algo == MaxwellSolverAlgo::CKC) { EvolveBCartesian ( Bfield, Efield, lev, dt ); #endif } else { amrex::Abort("EvolveB: Unknown algorithm"); } } #ifndef WARPX_DIM_RZ template void FiniteDifferenceSolver::EvolveBCartesian ( std::array< std::unique_ptr, 3 >& Bfield, std::array< std::unique_ptr, 3 > const& Efield, int lev, amrex::Real const dt ) { amrex::LayoutData* cost = WarpX::getCosts(lev); // Loop through the grids, and over the tiles within each grid #ifdef AMREX_USE_OMP #pragma omp parallel if (amrex::Gpu::notInLaunchRegion()) #endif for ( MFIter mfi(*Bfield[0], TilingIfNotGPU()); mfi.isValid(); ++mfi ) { if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers) { amrex::Gpu::synchronize(); } Real wt = amrex::second(); // Extract field data for this grid/tile Array4 const& Bx = Bfield[0]->array(mfi); Array4 const& By = Bfield[1]->array(mfi); Array4 const& Bz = Bfield[2]->array(mfi); Array4 const& Ex = Efield[0]->array(mfi); Array4 const& Ey = Efield[1]->array(mfi); Array4 const& Ez = Efield[2]->array(mfi); // Extract stencil coefficients Real const * const AMREX_RESTRICT coefs_x = m_stencil_coefs_x.dataPtr(); int const n_coefs_x = m_stencil_coefs_x.size(); Real const * const AMREX_RESTRICT coefs_y = m_stencil_coefs_y.dataPtr(); int const n_coefs_y = m_stencil_coefs_y.size(); Real const * const AMREX_RESTRICT coefs_z = m_stencil_coefs_z.dataPtr(); int const n_coefs_z = m_stencil_coefs_z.size(); // Extract tileboxes for which to loop Box const& tbx = mfi.tilebox(Bfield[0]->ixType().toIntVect()); Box const& tby = mfi.tilebox(Bfield[1]->ixType().toIntVect()); Box const& tbz = mfi.tilebox(Bfield[2]->ixType().toIntVect()); // Loop over the cells and update the fields amrex::ParallelFor(tbx, tby, tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k){ Bx(i, j, k) += dt * T_Algo::UpwardDz(Ey, coefs_z, n_coefs_z, i, j, k) - dt * T_Algo::UpwardDy(Ez, coefs_y, n_coefs_y, i, j, k); }, [=] AMREX_GPU_DEVICE (int i, int j, int k){ By(i, j, k) += dt * T_Algo::UpwardDx(Ez, coefs_x, n_coefs_x, i, j, k) - dt * T_Algo::UpwardDz(Ex, coefs_z, n_coefs_z, i, j, k); }, [=] AMREX_GPU_DEVICE (int i, int j, int k){ Bz(i, j, k) += dt * T_Algo::UpwardDy(Ex, coefs_y, n_coefs_y, i, j, k) - dt * T_Algo::UpwardDx(Ey, coefs_x, n_coefs_x, i, j, k); } ); if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers) { amrex::Gpu::synchronize(); wt = amrex::second() - wt; amrex::HostDevice::Atomic::Add( &(*cost)[mfi.index()], wt); } } } #else // corresponds to ifndef WARPX_DIM_RZ template void FiniteDifferenceSolver::EvolveBCylindrical ( std::array< std::unique_ptr, 3 >& Bfield, std::array< std::unique_ptr, 3 > const& Efield, int lev, amrex::Real const dt ) { amrex::LayoutData* cost = WarpX::getCosts(lev); // Loop through the grids, and over the tiles within each grid #ifdef AMREX_USE_OMP #pragma omp parallel if (amrex::Gpu::notInLaunchRegion()) #endif for ( MFIter mfi(*Bfield[0], TilingIfNotGPU()); mfi.isValid(); ++mfi ) { if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers) { amrex::Gpu::synchronize(); } Real wt = amrex::second(); // Extract field data for this grid/tile Array4 const& Br = Bfield[0]->array(mfi); Array4 const& Bt = Bfield[1]->array(mfi); Array4 const& Bz = Bfield[2]->array(mfi); Array4 const& Er = Efield[0]->array(mfi); Array4 const& Et = Efield[1]->array(mfi); Array4 const& Ez = Efield[2]->array(mfi); // Extract stencil coefficients Real const * const AMREX_RESTRICT coefs_r = m_stencil_coefs_r.dataPtr(); int const n_coefs_r = m_stencil_coefs_r.size(); Real const * const AMREX_RESTRICT coefs_z = m_stencil_coefs_z.dataPtr(); int const n_coefs_z = m_stencil_coefs_z.size(); // Extract cylindrical specific parameters Real const dr = m_dr; int const nmodes = m_nmodes; Real const rmin = m_rmin; // Extract tileboxes for which to loop Box const& tbr = mfi.tilebox(Bfield[0]->ixType().toIntVect()); Box const& tbt = mfi.tilebox(Bfield[1]->ixType().toIntVect()); Box const& tbz = mfi.tilebox(Bfield[2]->ixType().toIntVect()); // Loop over the cells and update the fields amrex::ParallelFor(tbr, tbt, tbz, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/){ Real const r = rmin + i*dr; // r on nodal point (Br is nodal in r) if (r != 0) { // Off-axis, regular Maxwell equations Br(i, j, 0, 0) += dt * T_Algo::UpwardDz(Et, coefs_z, n_coefs_z, i, j, 0, 0); // Mode m=0 for (int m=1; m