aboutsummaryrefslogtreecommitdiff
path: root/Source
diff options
context:
space:
mode:
Diffstat (limited to 'Source')
-rw-r--r--Source/Evolve/WarpXEvolve.cpp34
-rw-r--r--Source/WarpX.cpp5
2 files changed, 31 insertions, 8 deletions
diff --git a/Source/Evolve/WarpXEvolve.cpp b/Source/Evolve/WarpXEvolve.cpp
index b3b2462b1..1a6932696 100644
--- a/Source/Evolve/WarpXEvolve.cpp
+++ b/Source/Evolve/WarpXEvolve.cpp
@@ -10,6 +10,7 @@
*/
#include "WarpX.H"
+#include "BoundaryConditions/PML.H"
#include "Diagnostics/BackTransformedDiagnostic.H"
#include "Diagnostics/MultiDiagnostics.H"
#include "Diagnostics/ReducedDiags/MultiReducedDiags.H"
@@ -610,23 +611,50 @@ WarpX::OneStep_multiJ (const amrex::Real cur_time)
}
}
- // Transform fields back to real space and exchange guard cells
+ // Transform fields back to real space
if (WarpX::fft_do_time_averaging)
{
// We summed the integral of the field over 2*dt
PSATDScaleAverageFields(1._rt / (2._rt*dt[0]));
PSATDBackwardTransformEBavg();
}
+
+ // Evolve fields in PML
+ for (int lev = 0; lev <= finest_level; ++lev)
+ {
+ if (do_pml && pml[lev]->ok())
+ {
+ pml[lev]->PushPSATD(lev);
+ }
+ ApplyEfieldBoundary(lev, PatchType::fine);
+ if (lev > 0) ApplyEfieldBoundary(lev, PatchType::coarse);
+ ApplyBfieldBoundary(lev, PatchType::fine, DtType::FirstHalf);
+ if (lev > 0) ApplyBfieldBoundary(lev, PatchType::coarse, DtType::FirstHalf);
+ }
+
+ // Damp fields in PML before exchanging guard cells
+ if (do_pml)
+ {
+ DampPML();
+ }
+
+ // Exchange guard cells
FillBoundaryE(guard_cells.ng_alloc_EB);
FillBoundaryB(guard_cells.ng_alloc_EB);
- if (WarpX::do_dive_cleaning) FillBoundaryF(guard_cells.ng_alloc_F);
- if (WarpX::do_divb_cleaning) FillBoundaryG(guard_cells.ng_alloc_G);
+ if (WarpX::do_dive_cleaning || WarpX::do_pml_dive_cleaning) FillBoundaryF(guard_cells.ng_alloc_F);
+ if (WarpX::do_divb_cleaning || WarpX::do_pml_divb_cleaning) FillBoundaryG(guard_cells.ng_alloc_G);
// Synchronize E, B, F, G fields on nodal points
NodalSync(Efield_fp, Efield_cp);
NodalSync(Bfield_fp, Bfield_cp);
if (WarpX::do_dive_cleaning) NodalSync(F_fp, F_cp);
if (WarpX::do_divb_cleaning) NodalSync(G_fp, G_cp);
+
+ // Synchronize fields on nodal points in PML
+ if (do_pml)
+ {
+ NodalSyncPML();
+ }
}
else
{
diff --git a/Source/WarpX.cpp b/Source/WarpX.cpp
index b616ff534..57c053386 100644
--- a/Source/WarpX.cpp
+++ b/Source/WarpX.cpp
@@ -745,11 +745,6 @@ WarpX::ReadParameters ()
pp_warpx.query("do_pml_j_damping", do_pml_j_damping);
pp_warpx.query("do_pml_in_domain", do_pml_in_domain);
- if (do_multi_J && isAnyBoundaryPML())
- {
- amrex::Abort("Multi-J algorithm not implemented with PMLs");
- }
-
// Default values of WarpX::do_pml_dive_cleaning and WarpX::do_pml_divb_cleaning:
// false for FDTD solver, true for PSATD solver.
if (maxwell_solver_id != MaxwellSolverAlgo::PSATD)