ARTEMIS
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 #ifdef WARPX_MAG_LLG
155  std::array<amrex::MultiFab*,3> GetH_fp ();
156  std::array<amrex::MultiFab*,3> GetH_cp ();
157 #endif
162  amrex::MultiFab* GetE_fp (int comp);
163 
164  // Used when WarpX::do_pml_dive_cleaning = true
167 
168  // return macroscopic pml multifabs
175  // Used when WarpX::do_pml_divb_cleaning = true
178 
180  { return *sigba_fp; }
181 
183  { return *sigba_cp; }
184 
185 #ifdef WARPX_USE_PSATD
186  void PushPSATD (const int lev);
187 #endif
188 
189  void CopyJtoPMLs (const std::array<amrex::MultiFab*,3>& j_fp,
190  const std::array<amrex::MultiFab*,3>& j_cp);
191 
192  void Exchange (const std::array<amrex::MultiFab*,3>& mf_pml,
193  const std::array<amrex::MultiFab*,3>& mf,
194  const PatchType& patch_type,
195  const int do_pml_in_domain);
196 
197  void CopyJtoPMLs (PatchType patch_type,
198  const std::array<amrex::MultiFab*,3>& jp);
199 
200  void ExchangeF (amrex::MultiFab* F_fp, amrex::MultiFab* F_cp, int do_pml_in_domain);
201  void ExchangeF (PatchType patch_type, amrex::MultiFab* Fp, int do_pml_in_domain);
202 
203 #ifdef WARPX_MAG_LLG
204  void ExchangeH (const std::array<amrex::MultiFab*,3>& H_fp,
205  const std::array<amrex::MultiFab*,3>& H_cp, int do_pml_in_domain);
206  void ExchangeH (PatchType patch_type,
207  const std::array<amrex::MultiFab*,3>& Hp, int do_pml_in_domain);
208 #endif
209  void ExchangeG (amrex::MultiFab* G_fp, amrex::MultiFab* G_cp, int do_pml_in_domain);
210  void ExchangeG (PatchType patch_type, amrex::MultiFab* Gp, int do_pml_in_domain);
211 
212  void FillBoundary ();
213  void FillBoundaryE ();
214  void FillBoundaryB ();
215  void FillBoundaryF ();
216  void FillBoundaryG ();
217  void FillBoundaryE (PatchType patch_type);
218  void FillBoundaryB (PatchType patch_type);
219  void FillBoundaryF (PatchType patch_type);
220 #ifdef WARPX_MAG_LLG
221  void FillBoundaryH ();
222  void FillBoundaryH (PatchType patch_type);
223 #endif
224  void FillBoundaryG (PatchType patch_type);
225 
226  bool ok () const { return m_ok; }
227 
228  void CheckPoint (const std::string& dir) const;
229  void Restart (const std::string& dir);
230 
231  static void Exchange (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom, int do_pml_in_domain);
232 
233  ~PML () = default;
234 
235 private:
236  bool m_ok;
237 
240 
243 
246 
247  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_fp;
248  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_fp;
249  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_fp;
250 
251  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_edge_lengths;
252 
253  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_cp;
254  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_cp;
255  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_cp;
256 
257  // Used when WarpX::do_pml_dive_cleaning = true
258  std::unique_ptr<amrex::MultiFab> pml_F_fp;
259  std::unique_ptr<amrex::MultiFab> pml_F_cp;
260 
261 #ifdef WARPX_MAG_LLG
262  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_H_fp;
263  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_H_cp;
264 #endif
265 
266  // macroproperties in PML, eps-permittivity, mu-permeability, sigma-condictivity
267  std::unique_ptr<amrex::MultiFab> pml_eps_fp;
268  std::unique_ptr<amrex::MultiFab> pml_mu_fp;
269  std::unique_ptr<amrex::MultiFab> pml_sigma_fp;
270  std::unique_ptr<amrex::MultiFab> pml_eps_cp;
271  std::unique_ptr<amrex::MultiFab> pml_mu_cp;
272  std::unique_ptr<amrex::MultiFab> pml_sigma_cp;
273  // Used when WarpX::do_pml_divb_cleaning = true
274  std::unique_ptr<amrex::MultiFab> pml_G_fp;
275  std::unique_ptr<amrex::MultiFab> pml_G_cp;
276 
277  std::unique_ptr<MultiSigmaBox> sigba_fp;
278  std::unique_ptr<MultiSigmaBox> sigba_cp;
279 
280 #ifdef WARPX_USE_PSATD
281  std::unique_ptr<SpectralSolver> spectral_solver_fp;
282  std::unique_ptr<SpectralSolver> spectral_solver_cp;
283 #endif
284 
285  // Factory for field data
286  std::unique_ptr<amrex::FabFactory<amrex::FArrayBox> > pml_field_factory;
287 
289  return *pml_field_factory;
290  }
291 #ifdef AMREX_USE_EB
292  amrex::EBFArrayBoxFactory const& fieldEBFactory () const noexcept {
293  return static_cast<amrex::EBFArrayBoxFactory const&>(*pml_field_factory);
294  }
295 #endif
296 
297  static amrex::BoxArray MakeBoxArray (bool single_box_domain,
298  const amrex::Box& regular_domain,
299  const amrex::Geometry& geom,
300  const amrex::BoxArray& grid_ba,
301  const amrex::IntVect& ncell,
302  int do_pml_in_domain,
303  const amrex::IntVect& do_pml_Lo,
304  const amrex::IntVect& do_pml_Hi);
305 
306  static amrex::BoxArray MakeBoxArray_single (const amrex::Box& regular_domain,
307  const amrex::BoxArray& grid_ba,
308  const amrex::IntVect& ncell,
309  const amrex::IntVect& do_pml_Lo,
310  const amrex::IntVect& do_pml_Hi);
311 
312  static amrex::BoxArray MakeBoxArray_multiple (const amrex::Geometry& geom,
313  const amrex::BoxArray& grid_ba,
314  const amrex::IntVect& ncell,
315  int do_pml_in_domain,
316  const amrex::IntVect& do_pml_Lo,
317  const amrex::IntVect& do_pml_Hi);
318 
319  static void CopyToPML (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom);
320 };
321 
322 #ifdef WARPX_USE_PSATD
323 void PushPMLPSATDSinglePatch( const int lev,
324  SpectralSolver& solver,
325  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_E,
326  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_B,
327  std::unique_ptr<amrex::MultiFab>& pml_F,
328  std::unique_ptr<amrex::MultiFab>& pml_G,
329  const amrex::IntVect& fill_guards);
330 #endif
331 
332 #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:1622
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:533
void ComputePMLFactorsB(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:517
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:508
Definition: PML.H:128
amrex::MultiFab * Getmu_fp()
Definition: PML.cpp:1181
const MultiSigmaBox & GetMultiSigmaBox_fp() const
Definition: PML.H:179
amrex::MultiFab * Getsigma_fp()
Definition: PML.cpp:1187
const amrex::Geometry * m_geom
Definition: PML.H:244
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_edge_lengths
Definition: PML.H:251
std::unique_ptr< amrex::MultiFab > pml_eps_cp
Definition: PML.H:270
void ComputePMLFactors(amrex::Real dt)
Definition: PML.cpp:1084
amrex::MultiFab * GetG_cp()
Definition: PML.cpp:1217
std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > pml_field_factory
Definition: PML.H:286
std::unique_ptr< amrex::MultiFab > pml_F_fp
Definition: PML.H:258
const amrex::IntVect m_fill_guards_fields
Definition: PML.H:241
std::array< amrex::MultiFab *, 3 > GetE_fp()
Definition: PML.cpp:1097
std::array< amrex::MultiFab *, 3 > Getj_cp()
Definition: PML.cpp:1149
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_fp
Definition: PML.H:248
amrex::MultiFab * GetF_fp()
Definition: PML.cpp:1162
std::unique_ptr< MultiSigmaBox > sigba_fp
Definition: PML.H:277
const amrex::Geometry * m_cgeom
Definition: PML.H:245
bool ok() const
Definition: PML.H:226
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_fp
Definition: PML.H:249
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory() const noexcept
Definition: PML.H:288
amrex::MultiFab * Geteps_cp()
Definition: PML.cpp:1193
std::unique_ptr< amrex::MultiFab > pml_mu_fp
Definition: PML.H:268
bool m_divb_cleaning
Definition: PML.H:239
amrex::MultiFab * Geteps_fp()
Definition: PML.cpp:1175
bool m_ok
Definition: PML.H:236
std::array< amrex::MultiFab *, 3 > Getj_fp()
Definition: PML.cpp:1123
std::unique_ptr< amrex::MultiFab > pml_eps_fp
Definition: PML.H:267
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:548
std::array< amrex::MultiFab *, 3 > Get_edge_lengths()
Definition: PML.cpp:1155
bool m_dive_cleaning
Definition: PML.H:238
void PushPSATD(const int lev)
Definition: PML.cpp:1612
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_cp
Definition: PML.H:254
amrex::MultiFab * GetF_cp()
Definition: PML.cpp:1168
amrex::MultiFab * Getsigma_cp()
Definition: PML.cpp:1205
amrex::MultiFab * GetG_fp()
Definition: PML.cpp:1211
const MultiSigmaBox & GetMultiSigmaBox_cp() const
Definition: PML.H:182
std::unique_ptr< amrex::MultiFab > pml_sigma_cp
Definition: PML.H:272
const amrex::IntVect m_fill_guards_current
Definition: PML.H:242
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_cp
Definition: PML.H:253
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_cp
Definition: PML.H:255
std::unique_ptr< MultiSigmaBox > sigba_cp
Definition: PML.H:278
std::array< amrex::MultiFab *, 3 > GetE_cp()
Definition: PML.cpp:1129
std::unique_ptr< amrex::MultiFab > pml_mu_cp
Definition: PML.H:271
~PML()=default
amrex::MultiFab * Getmu_cp()
Definition: PML.cpp:1199
std::unique_ptr< amrex::MultiFab > pml_F_cp
Definition: PML.H:259
std::array< amrex::MultiFab *, 3 > Get_face_areas()
std::unique_ptr< SpectralSolver > spectral_solver_cp
Definition: PML.H:282
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_fp
Definition: PML.H:247
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:1222
std::unique_ptr< amrex::MultiFab > pml_sigma_fp
Definition: PML.H:269
std::unique_ptr< SpectralSolver > spectral_solver_fp
Definition: PML.H:281
std::array< amrex::MultiFab *, 3 > GetB_fp()
Definition: PML.cpp:1109
std::array< amrex::MultiFab *, 3 > GetB_cp()
Definition: PML.cpp:1135
void CopyJtoPMLs(const std::array< amrex::MultiFab *, 3 > &j_fp, const std::array< amrex::MultiFab *, 3 > &j_cp)
Definition: PML.cpp:1283
std::unique_ptr< amrex::MultiFab > pml_G_fp
Definition: PML.H:274
std::unique_ptr< amrex::MultiFab > pml_G_cp
Definition: PML.H:275
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:65
const int[]
std::array< T, N > Array
default
Definition: run_alltests.py:113
int dx
Definition: stencil.py:436
int dt
Definition: stencil.py:440
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:441
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:240
void ComputePMLFactorsE(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:475
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:147
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:197
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