WarpX
PML.H
Go to the documentation of this file.
1 /* Copyright 2019 Andrew Myers, Aurore Blelly, Axel Huebl
2  * Maxence Thevenet, Remi Lehe, Weiqun Zhang
3  *
4  *
5  * This file is part of WarpX.
6  *
7  * License: BSD-3-Clause-LBNL
8  */
9 #ifndef WARPX_PML_H_
10 #define WARPX_PML_H_
11 
12 #include "PML_fwd.H"
13 
14 #ifdef WARPX_USE_PSATD
16 #endif
17 
18 #include <AMReX_MultiFab.H>
19 #include <AMReX_BoxArray.H>
20 #include <AMReX_Config.H>
21 #include <AMReX_FabArray.H>
22 #include <AMReX_FabFactory.H>
23 #include <AMReX_GpuContainers.H>
24 #include <AMReX_IntVect.H>
25 #include <AMReX_REAL.H>
26 #include <AMReX_Vector.H>
27 
28 #include <AMReX_BaseFwd.H>
29 
30 #include <array>
31 #include <memory>
32 #include <string>
33 #include <vector>
34 
35 struct Sigma : amrex::Gpu::DeviceVector<amrex::Real>
36 {
37  int lo() const { return m_lo; }
38  int hi() const { return m_hi; }
39  int m_lo, m_hi;
40 };
41 
42 struct SigmaBox
43 {
44  SigmaBox (const amrex::Box& box, const amrex::BoxArray& grids,
45  const amrex::Real* dx, const amrex::IntVect& ncell, const amrex::IntVect& delta,
46  const amrex::Box& regdomain, const amrex::Real v_sigma);
47 
48  void define_single (const amrex::Box& regdomain, const amrex::IntVect& ncell,
50  const amrex::Real v_sigma);
51  void define_multiple (const amrex::Box& box, const amrex::BoxArray& grids,
52  const amrex::IntVect& ncell,
54  const amrex::Real v_sigma);
55 
56  void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
57  void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
58 
59  using SigmaVect = std::array<Sigma,AMREX_SPACEDIM>;
60 
61  using value_type = void; // needed by amrex::FabArray
62 
71  amrex::Real v_sigma;
72 
73 };
74 
76  : public amrex::FabFactory<SigmaBox>
77 {
78 public:
79  SigmaBoxFactory (const amrex::BoxArray& grid_ba, const amrex::Real* dx,
80  const amrex::IntVect& ncell, const amrex::IntVect& delta,
81  const amrex::Box& regular_domain, const amrex::Real v_sigma_sb)
82  : m_grids(grid_ba), m_dx(dx), m_ncell(ncell), m_delta(delta), m_regdomain(regular_domain), m_v_sigma_sb(v_sigma_sb) {}
83  virtual ~SigmaBoxFactory () = default;
84 
85  SigmaBoxFactory (const SigmaBoxFactory&) = default;
87 
88  SigmaBoxFactory () = delete;
89  SigmaBoxFactory& operator= (const SigmaBoxFactory&) = delete;
90  SigmaBoxFactory& operator= (SigmaBoxFactory&&) = delete;
91 
92  virtual SigmaBox* create (const amrex::Box& box, int /*ncomps*/,
93  const amrex::FabInfo& /*info*/, int /*box_index*/) const final
94  { return new SigmaBox(box, m_grids, m_dx, m_ncell, m_delta, m_regdomain, m_v_sigma_sb); }
95  virtual void destroy (SigmaBox* fab) const final {
96  delete fab;
97  }
98  virtual SigmaBoxFactory* clone () const final {
99  return new SigmaBoxFactory(*this);
100  }
101 private:
103  const amrex::Real* m_dx;
107  amrex::Real m_v_sigma_sb;
108 };
109 
111  : public amrex::FabArray<SigmaBox>
112 {
113 public:
115  const amrex::BoxArray& grid_ba, const amrex::Real* dx,
116  const amrex::IntVect& ncell, const amrex::IntVect& delta,
117  const amrex::Box& regular_domain, const amrex::Real v_sigma_sb);
118  void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
119  void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
120 private:
121  amrex::Real dt_B = -1.e10;
122  amrex::Real dt_E = -1.e10;
123 };
124 
125 enum struct PatchType : int;
126 
127 class PML
128 {
129 public:
130  PML (const int lev,
131  const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
132  const amrex::Geometry* geom, const amrex::Geometry* cgeom,
133  int ncell, int delta, amrex::IntVect ref_ratio,
134  amrex::Real dt, int nox_fft, int noy_fft, int noz_fft, short grid_type,
135  int do_moving_window, int pml_has_particles, int do_pml_in_domain,
136  const int psatd_solution_type, const int J_in_time, const int rho_in_time,
137  const bool do_pml_dive_cleaning, const bool do_pml_divb_cleaning,
138  const amrex::IntVect& fill_guards_fields,
139  const amrex::IntVect& fill_guards_current,
140  int max_guard_EB, amrex::Real v_sigma_sb,
141  const amrex::IntVect do_pml_Lo = amrex::IntVect::TheUnitVector(),
142  const amrex::IntVect do_pml_Hi = amrex::IntVect::TheUnitVector());
143 
144  void ComputePMLFactors (amrex::Real dt);
145 
146  std::array<amrex::MultiFab*,3> GetE_fp ();
147  std::array<amrex::MultiFab*,3> GetB_fp ();
148  std::array<amrex::MultiFab*,3> Getj_fp ();
149  std::array<amrex::MultiFab*,3> GetE_cp ();
150  std::array<amrex::MultiFab*,3> GetB_cp ();
151  std::array<amrex::MultiFab*,3> Getj_cp ();
152  std::array<amrex::MultiFab*,3> Get_edge_lengths ();
153  std::array<amrex::MultiFab*,3> Get_face_areas ();
154 
155  // Used when WarpX::do_pml_dive_cleaning = true
158 
159  // Used when WarpX::do_pml_divb_cleaning = true
162 
164  { return *sigba_fp; }
165 
167  { return *sigba_cp; }
168 
169 #ifdef WARPX_USE_PSATD
170  void PushPSATD (const int lev);
171 #endif
172 
173  void CopyJtoPMLs (const std::array<amrex::MultiFab*,3>& j_fp,
174  const std::array<amrex::MultiFab*,3>& j_cp);
175 
176  void Exchange (const std::array<amrex::MultiFab*,3>& mf_pml,
177  const std::array<amrex::MultiFab*,3>& mf,
178  const PatchType& patch_type,
179  const int do_pml_in_domain);
180 
181  void CopyJtoPMLs (PatchType patch_type,
182  const std::array<amrex::MultiFab*,3>& jp);
183 
184  void ExchangeF (amrex::MultiFab* F_fp, amrex::MultiFab* F_cp, int do_pml_in_domain);
185  void ExchangeF (PatchType patch_type, amrex::MultiFab* Fp, int do_pml_in_domain);
186 
187  void ExchangeG (amrex::MultiFab* G_fp, amrex::MultiFab* G_cp, int do_pml_in_domain);
188  void ExchangeG (PatchType patch_type, amrex::MultiFab* Gp, int do_pml_in_domain);
189 
190  void FillBoundary ();
191  void FillBoundaryE ();
192  void FillBoundaryB ();
193  void FillBoundaryF ();
194  void FillBoundaryG ();
195  void FillBoundaryE (PatchType patch_type);
196  void FillBoundaryB (PatchType patch_type);
197  void FillBoundaryF (PatchType patch_type);
198  void FillBoundaryG (PatchType patch_type);
199 
200  bool ok () const { return m_ok; }
201 
202  void CheckPoint (const std::string& dir) const;
203  void Restart (const std::string& dir);
204 
205  static void Exchange (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom, int do_pml_in_domain);
206 
207  ~PML () = default;
208 
209 private:
210  bool m_ok;
211 
214 
217 
220 
221  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_fp;
222  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_fp;
223  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_fp;
224 
225  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_edge_lengths;
226 
227  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_cp;
228  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_cp;
229  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_cp;
230 
231  // Used when WarpX::do_pml_dive_cleaning = true
232  std::unique_ptr<amrex::MultiFab> pml_F_fp;
233  std::unique_ptr<amrex::MultiFab> pml_F_cp;
234 
235  // Used when WarpX::do_pml_divb_cleaning = true
236  std::unique_ptr<amrex::MultiFab> pml_G_fp;
237  std::unique_ptr<amrex::MultiFab> pml_G_cp;
238 
239  std::unique_ptr<MultiSigmaBox> sigba_fp;
240  std::unique_ptr<MultiSigmaBox> sigba_cp;
241 
242 #ifdef WARPX_USE_PSATD
243  std::unique_ptr<SpectralSolver> spectral_solver_fp;
244  std::unique_ptr<SpectralSolver> spectral_solver_cp;
245 #endif
246 
247  // Factory for field data
248  std::unique_ptr<amrex::FabFactory<amrex::FArrayBox> > pml_field_factory;
249 
251  return *pml_field_factory;
252  }
253 #ifdef AMREX_USE_EB
254  amrex::EBFArrayBoxFactory const& fieldEBFactory () const noexcept {
255  return static_cast<amrex::EBFArrayBoxFactory const&>(*pml_field_factory);
256  }
257 #endif
258 
259  static amrex::BoxArray MakeBoxArray (bool single_box_domain,
260  const amrex::Box& regular_domain,
261  const amrex::Geometry& geom,
262  const amrex::BoxArray& grid_ba,
263  const amrex::IntVect& ncell,
264  int do_pml_in_domain,
265  const amrex::IntVect& do_pml_Lo,
266  const amrex::IntVect& do_pml_Hi);
267 
268  static amrex::BoxArray MakeBoxArray_single (const amrex::Box& regular_domain,
269  const amrex::BoxArray& grid_ba,
270  const amrex::IntVect& ncell,
271  const amrex::IntVect& do_pml_Lo,
272  const amrex::IntVect& do_pml_Hi);
273 
274  static amrex::BoxArray MakeBoxArray_multiple (const amrex::Geometry& geom,
275  const amrex::BoxArray& grid_ba,
276  const amrex::IntVect& ncell,
277  int do_pml_in_domain,
278  const amrex::IntVect& do_pml_Lo,
279  const amrex::IntVect& do_pml_Hi);
280 
281  static void CopyToPML (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom);
282 };
283 
284 #ifdef WARPX_USE_PSATD
285 void PushPMLPSATDSinglePatch( const int lev,
286  SpectralSolver& solver,
287  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_E,
288  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_B,
289  std::unique_ptr<amrex::MultiFab>& pml_F,
290  std::unique_ptr<amrex::MultiFab>& pml_G,
291  const amrex::IntVect& fill_guards);
292 #endif
293 
294 #endif
void PushPMLPSATDSinglePatch(const int lev, SpectralSolver &solver, std::array< std::unique_ptr< amrex::MultiFab >, 3 > &pml_E, std::array< std::unique_ptr< amrex::MultiFab >, 3 > &pml_B, std::unique_ptr< amrex::MultiFab > &pml_F, std::unique_ptr< amrex::MultiFab > &pml_G, const amrex::IntVect &fill_guards)
Definition: PML.cpp:1392
PatchType
Definition: WarpX.H:73
Definition: PML.H:112
amrex::Real dt_E
Definition: PML.H:122
void ComputePMLFactorsE(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:530
void ComputePMLFactorsB(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:514
amrex::Real dt_B
Definition: PML.H:121
MultiSigmaBox(const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const amrex::BoxArray &grid_ba, const amrex::Real *dx, const amrex::IntVect &ncell, const amrex::IntVect &delta, const amrex::Box &regular_domain, const amrex::Real v_sigma_sb)
Definition: PML.cpp:505
Definition: PML.H:128
const MultiSigmaBox & GetMultiSigmaBox_fp() const
Definition: PML.H:163
const amrex::Geometry * m_geom
Definition: PML.H:218
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_edge_lengths
Definition: PML.H:225
void ComputePMLFactors(amrex::Real dt)
Definition: PML.cpp:990
amrex::MultiFab * GetG_cp()
Definition: PML.cpp:1064
std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > pml_field_factory
Definition: PML.H:248
std::unique_ptr< amrex::MultiFab > pml_F_fp
Definition: PML.H:232
const amrex::IntVect m_fill_guards_fields
Definition: PML.H:215
std::array< amrex::MultiFab *, 3 > GetE_fp()
Definition: PML.cpp:1003
std::array< amrex::MultiFab *, 3 > Getj_cp()
Definition: PML.cpp:1033
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_fp
Definition: PML.H:222
amrex::MultiFab * GetF_fp()
Definition: PML.cpp:1046
std::unique_ptr< MultiSigmaBox > sigba_fp
Definition: PML.H:239
const amrex::Geometry * m_cgeom
Definition: PML.H:219
bool ok() const
Definition: PML.H:200
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_fp
Definition: PML.H:223
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory() const noexcept
Definition: PML.H:250
bool m_divb_cleaning
Definition: PML.H:213
bool m_ok
Definition: PML.H:210
std::array< amrex::MultiFab *, 3 > Getj_fp()
Definition: PML.cpp:1015
PML(const int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const amrex::Geometry *geom, const amrex::Geometry *cgeom, int ncell, int delta, amrex::IntVect ref_ratio, amrex::Real dt, int nox_fft, int noy_fft, int noz_fft, short grid_type, int do_moving_window, int pml_has_particles, int do_pml_in_domain, const int psatd_solution_type, const int J_in_time, const int rho_in_time, const bool do_pml_dive_cleaning, const bool do_pml_divb_cleaning, const amrex::IntVect &fill_guards_fields, const amrex::IntVect &fill_guards_current, int max_guard_EB, amrex::Real v_sigma_sb, const amrex::IntVect do_pml_Lo=amrex::IntVect::TheUnitVector(), const amrex::IntVect do_pml_Hi=amrex::IntVect::TheUnitVector())
Definition: PML.cpp:545
std::array< amrex::MultiFab *, 3 > Get_edge_lengths()
Definition: PML.cpp:1039
bool m_dive_cleaning
Definition: PML.H:212
void PushPSATD(const int lev)
Definition: PML.cpp:1382
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_cp
Definition: PML.H:228
amrex::MultiFab * GetF_cp()
Definition: PML.cpp:1052
amrex::MultiFab * GetG_fp()
Definition: PML.cpp:1058
const MultiSigmaBox & GetMultiSigmaBox_cp() const
Definition: PML.H:166
const amrex::IntVect m_fill_guards_current
Definition: PML.H:216
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_cp
Definition: PML.H:227
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_cp
Definition: PML.H:229
std::unique_ptr< MultiSigmaBox > sigba_cp
Definition: PML.H:240
std::array< amrex::MultiFab *, 3 > GetE_cp()
Definition: PML.cpp:1021
~PML()=default
std::unique_ptr< amrex::MultiFab > pml_F_cp
Definition: PML.H:233
std::array< amrex::MultiFab *, 3 > Get_face_areas()
std::unique_ptr< SpectralSolver > spectral_solver_cp
Definition: PML.H:244
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_fp
Definition: PML.H:221
void Exchange(const std::array< amrex::MultiFab *, 3 > &mf_pml, const std::array< amrex::MultiFab *, 3 > &mf, const PatchType &patch_type, const int do_pml_in_domain)
Definition: PML.cpp:1069
std::unique_ptr< SpectralSolver > spectral_solver_fp
Definition: PML.H:243
std::array< amrex::MultiFab *, 3 > GetB_fp()
Definition: PML.cpp:1009
std::array< amrex::MultiFab *, 3 > GetB_cp()
Definition: PML.cpp:1027
void CopyJtoPMLs(const std::array< amrex::MultiFab *, 3 > &j_fp, const std::array< amrex::MultiFab *, 3 > &j_cp)
Definition: PML.cpp:1099
std::unique_ptr< amrex::MultiFab > pml_G_fp
Definition: PML.H:236
std::unique_ptr< amrex::MultiFab > pml_G_cp
Definition: PML.H:237
Definition: PML.H:77
SigmaBoxFactory(const amrex::BoxArray &grid_ba, const amrex::Real *dx, const amrex::IntVect &ncell, const amrex::IntVect &delta, const amrex::Box &regular_domain, const amrex::Real v_sigma_sb)
Definition: PML.H:79
virtual SigmaBoxFactory * clone() const final
Definition: PML.H:98
amrex::IntVect m_delta
Definition: PML.H:105
SigmaBoxFactory(SigmaBoxFactory &&) noexcept=default
amrex::IntVect m_ncell
Definition: PML.H:104
amrex::Box m_regdomain
Definition: PML.H:106
SigmaBoxFactory(const SigmaBoxFactory &)=default
virtual void destroy(SigmaBox *fab) const final
Definition: PML.H:95
amrex::Real m_v_sigma_sb
Definition: PML.H:107
virtual ~SigmaBoxFactory()=default
SigmaBoxFactory()=delete
const amrex::BoxArray & m_grids
Definition: PML.H:102
virtual SigmaBox * create(const amrex::Box &box, int, const amrex::FabInfo &, int) const final
Definition: PML.H:92
const amrex::Real * m_dx
Definition: PML.H:103
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:33
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVect TheUnitVector() noexcept
void FillBoundary(amrex::MultiFab &mf, amrex::IntVect ng, bool do_single_precision_comms, const amrex::Periodicity &period, std::optional< bool > nodal_sync)
Definition: Communication.cpp:71
const int[]
std::array< T, N > Array
default
Definition: run_alltests.py:113
float dt
Definition: stencil.py:440
tuple dx
lab frame
Definition: stencil.py:427
Definition: PML.H:43
SigmaVect sigma_star_cumsum_fac
Definition: PML.H:70
SigmaVect sigma_star
Definition: PML.H:65
std::array< Sigma, AMREX_SPACEDIM > SigmaVect
Definition: PML.H:59
SigmaVect sigma_star_fac
Definition: PML.H:69
SigmaVect sigma_cumsum
Definition: PML.H:64
SigmaVect sigma_star_cumsum
Definition: PML.H:66
void ComputePMLFactorsB(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:438
void value_type
Definition: PML.H:61
SigmaVect sigma_cumsum_fac
Definition: PML.H:68
SigmaVect sigma_fac
Definition: PML.H:67
void define_multiple(const amrex::Box &box, const amrex::BoxArray &grids, const amrex::IntVect &ncell, const amrex::Array< amrex::Real, AMREX_SPACEDIM > &fac, const amrex::Real v_sigma)
Definition: PML.cpp:238
void ComputePMLFactorsE(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:472
SigmaVect sigma
Definition: PML.H:63
SigmaBox(const amrex::Box &box, const amrex::BoxArray &grids, const amrex::Real *dx, const amrex::IntVect &ncell, const amrex::IntVect &delta, const amrex::Box &regdomain, const amrex::Real v_sigma)
Definition: PML.cpp:145
void define_single(const amrex::Box &regdomain, const amrex::IntVect &ncell, const amrex::Array< amrex::Real, AMREX_SPACEDIM > &fac, const amrex::Real v_sigma)
Definition: PML.cpp:195
amrex::Real v_sigma
Definition: PML.H:71
Definition: PML.H:36
int m_lo
Definition: PML.H:39
int lo() const
Definition: PML.H:37
int hi() const
Definition: PML.H:38
int m_hi
Definition: PML.H:39