aboutsummaryrefslogtreecommitdiff
path: root/Source/Parallelization
diff options
context:
space:
mode:
Diffstat (limited to 'Source/Parallelization')
-rw-r--r--Source/Parallelization/GuardCellManager.H2
-rw-r--r--Source/Parallelization/GuardCellManager.cpp5
2 files changed, 7 insertions, 0 deletions
diff --git a/Source/Parallelization/GuardCellManager.H b/Source/Parallelization/GuardCellManager.H
index d07c5c99a..4568db45a 100644
--- a/Source/Parallelization/GuardCellManager.H
+++ b/Source/Parallelization/GuardCellManager.H
@@ -80,6 +80,8 @@ public:
amrex::IntVect ng_UpdateAux = amrex::IntVect::TheZeroVector();
// Number of guard cells of all MultiFabs that must exchanged before moving window
amrex::IntVect ng_MovingWindow = amrex::IntVect::TheZeroVector();
+ // Number of guard cells of E and B that are exchanged immediatly after the main PSATD push
+ amrex::IntVect ng_afterPushPSATD = amrex::IntVect::TheZeroVector();
// Number of guard cells for local deposition of J and rho
amrex::IntVect ng_depos_J = amrex::IntVect::TheZeroVector();
diff --git a/Source/Parallelization/GuardCellManager.cpp b/Source/Parallelization/GuardCellManager.cpp
index 79b613cf6..8dd44f706 100644
--- a/Source/Parallelization/GuardCellManager.cpp
+++ b/Source/Parallelization/GuardCellManager.cpp
@@ -232,6 +232,10 @@ guardCellManager::Init (
ng_alloc_F.max( ng_FieldSolverF );
ng_alloc_G.max( ng_FieldSolverG );
+ if (do_moving_window && maxwell_solver_id == MaxwellSolverAlgo::PSATD) {
+ ng_afterPushPSATD = ng_alloc_EB;
+ }
+
if (safe_guard_cells){
// Run in safe mode: exchange all allocated guard cells at each
// call of FillBoundary
@@ -240,6 +244,7 @@ guardCellManager::Init (
ng_FieldSolverG = ng_alloc_G;
ng_FieldGather = ng_alloc_EB;
ng_UpdateAux = ng_alloc_EB;
+ ng_afterPushPSATD = ng_alloc_EB;
if (do_moving_window){
ng_MovingWindow = ng_alloc_EB;
}