/* Copyright 2020 * * This file is part of WarpX. * * License: BSD-3-Clause-LBNL */ #include "FiniteDifferenceSolver.H" #ifndef WARPX_DIM_RZ # include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H" # include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H" # include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H" #else # include "FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CylindricalYeeAlgorithm.H" #endif #include "Utils/TextMsg.H" #include "Utils/WarpXAlgorithmSelection.H" #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include using namespace amrex; void FiniteDifferenceSolver::EvolveG ( std::unique_ptr& Gfield, std::array,3> const& Bfield, amrex::Real const dt) { #ifdef WARPX_DIM_RZ // TODO Implement G update equation in RZ geometry amrex::ignore_unused(Gfield, Bfield, dt); #else // Select algorithm if (m_do_nodal) { EvolveGCartesian(Gfield, Bfield, dt); } else if (m_fdtd_algo == ElectromagneticSolverAlgo::Yee) { EvolveGCartesian(Gfield, Bfield, dt); } else if (m_fdtd_algo == ElectromagneticSolverAlgo::CKC) { EvolveGCartesian(Gfield, Bfield, dt); } else { amrex::Abort(Utils::TextMsg::Err("EvolveG: unknown FDTD algorithm")); } #endif } #ifndef WARPX_DIM_RZ template void FiniteDifferenceSolver::EvolveGCartesian ( std::unique_ptr& Gfield, std::array,3> const& Bfield, amrex::Real const dt) { amrex::Real constexpr c2 = PhysConst::c * PhysConst::c; #ifdef AMREX_USE_OMP #pragma omp parallel if (amrex::Gpu::notInLaunchRegion()) #endif // Loop over grids and over tiles within each grid for (amrex::MFIter mfi(*Gfield, TilingIfNotGPU()); mfi.isValid(); ++mfi) { // Extract field data for this grid/tile amrex::Array4 const& G = Gfield->array(mfi); amrex::Array4 const& Bx = Bfield[0]->array(mfi); amrex::Array4 const& By = Bfield[1]->array(mfi); amrex::Array4 const& Bz = Bfield[2]->array(mfi); // Extract stencil coefficients amrex::Real const* const AMREX_RESTRICT coefs_x = m_stencil_coefs_x.dataPtr(); amrex::Real const* const AMREX_RESTRICT coefs_y = m_stencil_coefs_y.dataPtr(); amrex::Real const* const AMREX_RESTRICT coefs_z = m_stencil_coefs_z.dataPtr(); const int n_coefs_x = m_stencil_coefs_x.size(); const int n_coefs_y = m_stencil_coefs_y.size(); const int n_coefs_z = m_stencil_coefs_z.size(); // Extract tilebox to loop over amrex::Box const& tf = mfi.tilebox(Gfield->ixType().toIntVect()); // Loop over cells and update G amrex::ParallelFor(tf, [=] AMREX_GPU_DEVICE (int i, int j, int k) { G(i,j,k) += c2 * dt * (T_Algo::UpwardDx(Bx, coefs_x, n_coefs_x, i, j, k) + T_Algo::UpwardDy(By, coefs_y, n_coefs_y, i, j, k) + T_Algo::UpwardDz(Bz, coefs_z, n_coefs_z, i, j, k)); }); } } #endif