aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
diff options
context:
space:
mode:
Diffstat (limited to 'Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp')
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp26
1 files changed, 23 insertions, 3 deletions
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp b/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
index 2f66796f3..c24b869e0 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
+++ b/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
@@ -25,6 +25,7 @@ using namespace amrex;
void FiniteDifferenceSolver::EvolveB (
std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
+ std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
int lev, amrex::Real const dt ) {
// Select algorithm (The choice of algorithm is a runtime option,
@@ -37,15 +38,15 @@ void FiniteDifferenceSolver::EvolveB (
#else
if (m_do_nodal) {
- EvolveBCartesian <CartesianNodalAlgorithm> ( Bfield, Efield, lev, dt );
+ EvolveBCartesian <CartesianNodalAlgorithm> ( Bfield, Efield, face_areas, lev, dt );
} else if (m_fdtd_algo == MaxwellSolverAlgo::Yee) {
- EvolveBCartesian <CartesianYeeAlgorithm> ( Bfield, Efield, lev, dt );
+ EvolveBCartesian <CartesianYeeAlgorithm> ( Bfield, Efield, face_areas, lev, dt );
} else if (m_fdtd_algo == MaxwellSolverAlgo::CKC) {
- EvolveBCartesian <CartesianCKCAlgorithm> ( Bfield, Efield, lev, dt );
+ EvolveBCartesian <CartesianCKCAlgorithm> ( Bfield, Efield, face_areas, lev, dt );
#endif
} else {
@@ -61,6 +62,7 @@ template<typename T_Algo>
void FiniteDifferenceSolver::EvolveBCartesian (
std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
+ std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
int lev, amrex::Real const dt ) {
amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
@@ -84,6 +86,12 @@ void FiniteDifferenceSolver::EvolveBCartesian (
Array4<Real> const& Ey = Efield[1]->array(mfi);
Array4<Real> const& Ez = Efield[2]->array(mfi);
+#ifdef AMREX_USE_EB
+ amrex::Array4<amrex::Real> const& Sx = face_areas[0]->array(mfi);
+ amrex::Array4<amrex::Real> const& Sy = face_areas[1]->array(mfi);
+ amrex::Array4<amrex::Real> const& Sz = face_areas[2]->array(mfi);
+#endif
+
// 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();
@@ -101,16 +109,28 @@ void FiniteDifferenceSolver::EvolveBCartesian (
amrex::ParallelFor(tbx, tby, tbz,
[=] AMREX_GPU_DEVICE (int i, int j, int k){
+#ifdef AMREX_USE_EB
+ // Skip field push if this cell is fully covered by embedded boundaries
+ if (Sx(i, j, k) <= 0) return;
+#endif
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){
+#ifdef AMREX_USE_EB
+ // Skip field push if this cell is fully covered by embedded boundaries
+ if (Sy(i, j, k) <= 0) return;
+#endif
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){
+#ifdef AMREX_USE_EB
+ // Skip field push if this cell is fully covered by embedded boundaries
+ if (Sz(i, j, k) <= 0) return;
+#endif
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);
}