Skip to content

Commit a7c1afc

Browse files
authored
Merge pull request #259 from ablelly/choose_pml_direction
Choose pml direction
2 parents 01cff40 + 9aac67c commit a7c1afc

File tree

6 files changed

+62
-14
lines changed

6 files changed

+62
-14
lines changed

Docs/source/running_cpp/parameters.rst

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -661,6 +661,12 @@ Boundary conditions
661661
The characteristic depth, in number of cells, over which
662662
the absorption coefficients of the PML increases.
663663

664+
* ``warpx.do_pml_Lo`` (`2 ints in 2D`, `3 ints in 3D`; default: `1 1 1`)
665+
The directions along which one wants a pml boundary condition for lower boundaries on mother grid.
666+
667+
* ``warpx.do_pml_Hi`` (`2 floats in 2D`, `3 floats in 3D`; default: `1 1 1`)
668+
The directions along which one wants a pml boundary condition for upper boundaries on mother grid.
669+
664670
Diagnostics and output
665671
----------------------
666672

Source/BoundaryConditions/PML.H

Lines changed: 5 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -101,7 +101,8 @@ public:
101101
#ifdef WARPX_USE_PSATD
102102
amrex::Real dt, int nox_fft, int noy_fft, int noz_fft, bool do_nodal,
103103
#endif
104-
int do_dive_cleaning, int do_moving_window);
104+
int do_dive_cleaning, int do_moving_window,
105+
const amrex::IntVect do_pml_Lo, const amrex::IntVect do_pml_Hi);
105106

106107
void ComputePMLFactors (amrex::Real dt);
107108

@@ -172,7 +173,9 @@ private:
172173
#endif
173174

174175
static amrex::BoxArray MakeBoxArray (const amrex::Geometry& geom,
175-
const amrex::BoxArray& grid_ba, int ncell);
176+
const amrex::BoxArray& grid_ba, int ncell,
177+
const amrex::IntVect do_pml_Lo,
178+
const amrex::IntVect do_pml_Hi);
176179

177180
static void Exchange (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom);
178181
};

Source/BoundaryConditions/PML.cpp

Lines changed: 12 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -319,11 +319,12 @@ PML::PML (const BoxArray& grid_ba, const DistributionMapping& grid_dm,
319319
#ifdef WARPX_USE_PSATD
320320
Real dt, int nox_fft, int noy_fft, int noz_fft, bool do_nodal,
321321
#endif
322-
int do_dive_cleaning, int do_moving_window)
322+
int do_dive_cleaning, int do_moving_window,
323+
const amrex::IntVect do_pml_Lo, const amrex::IntVect do_pml_Hi)
323324
: m_geom(geom),
324325
m_cgeom(cgeom)
325326
{
326-
const BoxArray& ba = MakeBoxArray(*geom, grid_ba, ncell);
327+
const BoxArray& ba = MakeBoxArray(*geom, grid_ba, ncell, do_pml_Lo, do_pml_Hi);
327328
if (ba.size() == 0) {
328329
m_ok = false;
329330
return;
@@ -399,7 +400,7 @@ PML::PML (const BoxArray& grid_ba, const DistributionMapping& grid_dm,
399400

400401
BoxArray grid_cba = grid_ba;
401402
grid_cba.coarsen(ref_ratio);
402-
const BoxArray& cba = MakeBoxArray(*cgeom, grid_cba, ncell);
403+
const BoxArray& cba = MakeBoxArray(*cgeom, grid_cba, ncell, do_pml_Lo, do_pml_Hi);
403404

404405
DistributionMapping cdm{cba};
405406

@@ -438,12 +439,18 @@ PML::PML (const BoxArray& grid_ba, const DistributionMapping& grid_dm,
438439
}
439440

440441
BoxArray
441-
PML::MakeBoxArray (const amrex::Geometry& geom, const amrex::BoxArray& grid_ba, int ncell)
442+
PML::MakeBoxArray (const amrex::Geometry& geom, const amrex::BoxArray& grid_ba, int ncell,
443+
const amrex::IntVect do_pml_Lo, const amrex::IntVect do_pml_Hi)
442444
{
443445
Box domain = geom.Domain();
444446
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
445447
if ( ! geom.isPeriodic(idim) ) {
446-
domain.grow(idim, ncell);
448+
if (do_pml_Lo[idim]){
449+
domain.growLo(idim, ncell);
450+
}
451+
if (do_pml_Hi[idim]){
452+
domain.growHi(idim, ncell);
453+
}
447454
}
448455
}
449456

Source/Initialization/WarpXInitData.cpp

Lines changed: 16 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -131,21 +131,35 @@ WarpX::InitPML ()
131131
{
132132
if (do_pml)
133133
{
134+
amrex::IntVect do_pml_Lo_corrected = do_pml_Lo;
135+
136+
#ifdef WARPX_DIM_RZ
137+
do_pml_Lo_corrected[0] = 0; // no PML at r=0, in cylindrical geometry
138+
#endif
134139
pml[0].reset(new PML(boxArray(0), DistributionMap(0), &Geom(0), nullptr,
135140
pml_ncell, pml_delta, 0,
136141
#ifdef WARPX_USE_PSATD
137142
dt[0], nox_fft, noy_fft, noz_fft, do_nodal,
138143
#endif
139-
do_dive_cleaning, do_moving_window));
144+
do_dive_cleaning, do_moving_window,
145+
do_pml_Lo_corrected, do_pml_Hi));
140146
for (int lev = 1; lev <= finest_level; ++lev)
141147
{
148+
amrex::IntVect do_pml_Lo_MR = amrex::IntVect::TheUnitVector();
149+
#ifdef WARPX_DIM_RZ
150+
//In cylindrical geometry, if the edge of the patch is at r=0, do not add PML
151+
if ((max_level > 0) && (fine_tag_lo[0]==0.)) {
152+
do_pml_Lo_MR[0] = 0;
153+
}
154+
#endif
142155
pml[lev].reset(new PML(boxArray(lev), DistributionMap(lev),
143156
&Geom(lev), &Geom(lev-1),
144157
pml_ncell, pml_delta, refRatio(lev-1)[0],
145158
#ifdef WARPX_USE_PSATD
146159
dt[lev], nox_fft, noy_fft, noz_fft, do_nodal,
147160
#endif
148-
do_dive_cleaning, do_moving_window));
161+
do_dive_cleaning, do_moving_window,
162+
do_pml_Lo_MR, amrex::IntVect::TheUnitVector()));
149163
}
150164
}
151165
}

Source/WarpX.H

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -496,6 +496,8 @@ private:
496496
int do_pml = 1;
497497
int pml_ncell = 10;
498498
int pml_delta = 10;
499+
amrex::IntVect do_pml_Lo = amrex::IntVect::TheUnitVector();
500+
amrex::IntVect do_pml_Hi = amrex::IntVect::TheUnitVector();
499501
amrex::Vector<std::unique_ptr<PML> > pml;
500502

501503
amrex::Real moving_window_x = std::numeric_limits<amrex::Real>::max();

Source/WarpX.cpp

Lines changed: 21 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -383,6 +383,22 @@ WarpX::ReadParameters ()
383383
pp.query("pml_ncell", pml_ncell);
384384
pp.query("pml_delta", pml_delta);
385385

386+
Vector<int> parse_do_pml_Lo(AMREX_SPACEDIM,1);
387+
pp.queryarr("do_pml_Lo", parse_do_pml_Lo);
388+
do_pml_Lo[0] = parse_do_pml_Lo[0];
389+
do_pml_Lo[1] = parse_do_pml_Lo[1];
390+
#if (AMREX_SPACEDIM == 3)
391+
do_pml_Lo[2] = parse_do_pml_Lo[2];
392+
#endif
393+
Vector<int> parse_do_pml_Hi(AMREX_SPACEDIM,1);
394+
pp.queryarr("do_pml_Hi", parse_do_pml_Hi);
395+
do_pml_Hi[0] = parse_do_pml_Hi[0];
396+
do_pml_Hi[1] = parse_do_pml_Hi[1];
397+
#if (AMREX_SPACEDIM == 3)
398+
do_pml_Hi[2] = parse_do_pml_Hi[2];
399+
#endif
400+
401+
386402
pp.query("dump_openpmd", dump_openpmd);
387403
pp.query("dump_plotfiles", dump_plotfiles);
388404
pp.query("plot_raw_fields", plot_raw_fields);
@@ -393,7 +409,7 @@ WarpX::ReadParameters ()
393409
if (not user_fields_to_plot){
394410
// If not specified, set default values
395411
fields_to_plot = {"Ex", "Ey", "Ez", "Bx", "By",
396-
"Bz", "jx", "jy", "jz",
412+
"Bz", "jx", "jy", "jz",
397413
"part_per_cell"};
398414
}
399415
// set plot_rho to true of the users requests it, so that
@@ -411,9 +427,9 @@ WarpX::ReadParameters ()
411427
// If user requests to plot proc_number for a serial run,
412428
// delete proc_number from fields_to_plot
413429
if (ParallelDescriptor::NProcs() == 1){
414-
fields_to_plot.erase(std::remove(fields_to_plot.begin(),
415-
fields_to_plot.end(),
416-
"proc_number"),
430+
fields_to_plot.erase(std::remove(fields_to_plot.begin(),
431+
fields_to_plot.end(),
432+
"proc_number"),
417433
fields_to_plot.end());
418434
}
419435

@@ -497,7 +513,7 @@ WarpX::ReadParameters ()
497513
{
498514
ParmParse pp("algo");
499515
// If not in RZ mode, read use_picsar_deposition
500-
// In RZ mode, use_picsar_deposition is on, as the C++ version
516+
// In RZ mode, use_picsar_deposition is on, as the C++ version
501517
// of the deposition does not support RZ
502518
pp.query("use_picsar_deposition", use_picsar_deposition);
503519
current_deposition_algo = GetAlgorithmInteger(pp, "current_deposition");

0 commit comments

Comments
 (0)