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 
15 
16 #ifdef WARPX_USE_FFT
18 #endif
19 
20 #include <AMReX_MultiFab.H>
21 #include <AMReX_BoxArray.H>
22 #include <AMReX_Config.H>
23 #include <AMReX_FabArray.H>
24 #include <AMReX_FabFactory.H>
25 #include <AMReX_GpuContainers.H>
26 #include <AMReX_IntVect.H>
27 #include <AMReX_REAL.H>
28 #include <AMReX_Vector.H>
29 
30 #include <AMReX_BaseFwd.H>
31 
32 #include <array>
33 #include <memory>
34 #include <optional>
35 #include <string>
36 #include <vector>
37 
38 struct Sigma : amrex::Gpu::DeviceVector<amrex::Real>
39 {
40  [[nodiscard]] int lo() const { return m_lo; }
41  [[nodiscard]] int hi() const { return m_hi; }
42  int m_lo, m_hi;
43 };
44 
45 struct SigmaBox
46 {
47  SigmaBox (const amrex::Box& box, const amrex::BoxArray& grids,
48  const amrex::Real* dx, const amrex::IntVect& ncell, const amrex::IntVect& delta,
49  const amrex::Box& regdomain, amrex::Real v_sigma);
50 
51  void define_single (const amrex::Box& regdomain, const amrex::IntVect& ncell,
53  amrex::Real v_sigma);
54  void define_multiple (const amrex::Box& box, const amrex::BoxArray& grids,
55  const amrex::IntVect& ncell,
57  amrex::Real v_sigma);
58 
59  void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
60  void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
61 
62  using SigmaVect = std::array<Sigma,AMREX_SPACEDIM>;
63 
64  using value_type = void; // needed by amrex::FabArray
65 
74  amrex::Real v_sigma;
75 
76 };
77 
79  : public amrex::FabFactory<SigmaBox>
80 {
81 public:
82  SigmaBoxFactory (const amrex::BoxArray& grid_ba, const amrex::Real* dx,
83  const amrex::IntVect& ncell, const amrex::IntVect& delta,
84  const amrex::Box& regular_domain, const amrex::Real v_sigma_sb)
85  : m_grids(grid_ba), m_dx(dx), m_ncell(ncell), m_delta(delta), m_regdomain(regular_domain), m_v_sigma_sb(v_sigma_sb) {}
86  ~SigmaBoxFactory () override = default;
87 
88  SigmaBoxFactory (const SigmaBoxFactory&) = default;
90 
91  SigmaBoxFactory () = delete;
92  SigmaBoxFactory& operator= (const SigmaBoxFactory&) = delete;
93  SigmaBoxFactory& operator= (SigmaBoxFactory&&) = delete;
94 
95  [[nodiscard]] SigmaBox* create (const amrex::Box& box, int /*ncomps*/,
96  const amrex::FabInfo& /*info*/, int /*box_index*/) const final
97  {
99  }
100 
101  void destroy (SigmaBox* fab) const final
102  {
103  delete fab;
104  }
105 
106  [[nodiscard]] SigmaBoxFactory*
107  clone () const final
108  {
109  return new SigmaBoxFactory(*this);
110  }
111 
112 private:
114  const amrex::Real* m_dx;
118  amrex::Real m_v_sigma_sb;
119 };
120 
122  : public amrex::FabArray<SigmaBox>
123 {
124 public:
126  const amrex::BoxArray& grid_ba, const amrex::Real* dx,
127  const amrex::IntVect& ncell, const amrex::IntVect& delta,
128  const amrex::Box& regular_domain, amrex::Real v_sigma_sb);
129  void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
130  void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
131 private:
132  amrex::Real dt_B = -1.e10;
133  amrex::Real dt_E = -1.e10;
134 };
135 
136 class PML
137 {
138 public:
139  PML (int lev, const amrex::BoxArray& ba,
140  const amrex::DistributionMapping& dm, bool do_similar_dm_pml,
141  const amrex::Geometry* geom, const amrex::Geometry* cgeom,
142  int ncell, int delta, amrex::IntVect ref_ratio,
143  amrex::Real dt, int nox_fft, int noy_fft, int noz_fft,
145  int do_moving_window, int pml_has_particles, int do_pml_in_domain,
146  int psatd_solution_type, int J_in_time, int rho_in_time,
147  bool do_pml_dive_cleaning, bool do_pml_divb_cleaning,
148  const amrex::IntVect& fill_guards_fields,
149  const amrex::IntVect& fill_guards_current,
150  int max_guard_EB, amrex::Real v_sigma_sb,
153 
154  void ComputePMLFactors (amrex::Real dt);
155 
156  std::array<amrex::MultiFab*,3> GetE_fp ();
157  std::array<amrex::MultiFab*,3> GetB_fp ();
158  std::array<amrex::MultiFab*,3> Getj_fp ();
159  std::array<amrex::MultiFab*,3> GetE_cp ();
160  std::array<amrex::MultiFab*,3> GetB_cp ();
161  std::array<amrex::MultiFab*,3> Getj_cp ();
162  std::array<amrex::MultiFab*,3> Get_edge_lengths ();
163  std::array<amrex::MultiFab*,3> Get_face_areas ();
164 
165  // Used when WarpX::do_pml_dive_cleaning = true
168 
169  // Used when WarpX::do_pml_divb_cleaning = true
172 
173  [[nodiscard]] const MultiSigmaBox& GetMultiSigmaBox_fp () const
174  {
175  return *sigba_fp;
176  }
177 
178  [[nodiscard]] const MultiSigmaBox& GetMultiSigmaBox_cp () const
179  {
180  return *sigba_cp;
181  }
182 
183 #ifdef WARPX_USE_FFT
184  void PushPSATD (int lev);
185 #endif
186 
187  void CopyJtoPMLs (const std::array<amrex::MultiFab*,3>& j_fp,
188  const std::array<amrex::MultiFab*,3>& j_cp);
189 
190  void Exchange (const std::array<amrex::MultiFab*,3>& mf_pml,
191  const std::array<amrex::MultiFab*,3>& mf,
192  const PatchType& patch_type,
193  int do_pml_in_domain);
194 
195  void CopyJtoPMLs (PatchType patch_type,
196  const std::array<amrex::MultiFab*,3>& jp);
197 
198  void ExchangeF (amrex::MultiFab* F_fp, amrex::MultiFab* F_cp, int do_pml_in_domain);
199  void ExchangeF (PatchType patch_type, amrex::MultiFab* Fp, int do_pml_in_domain);
200 
201  void ExchangeG (amrex::MultiFab* G_fp, amrex::MultiFab* G_cp, int do_pml_in_domain);
202  void ExchangeG (PatchType patch_type, amrex::MultiFab* Gp, int do_pml_in_domain);
203 
204  void FillBoundaryE (PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
205  void FillBoundaryB (PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
206  void FillBoundaryF (PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
207  void FillBoundaryG (PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
208 
209  [[nodiscard]] bool ok () const { return m_ok; }
210 
211  void CheckPoint (const std::string& dir) const;
212  void Restart (const std::string& dir);
213 
214  static void Exchange (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom, int do_pml_in_domain);
215 
216 private:
217  bool m_ok;
218 
221 
224 
227 
228  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_fp;
229  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_fp;
230  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_fp;
231 
232  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_edge_lengths;
233 
234  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_cp;
235  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_cp;
236  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_cp;
237 
238  // Used when WarpX::do_pml_dive_cleaning = true
239  std::unique_ptr<amrex::MultiFab> pml_F_fp;
240  std::unique_ptr<amrex::MultiFab> pml_F_cp;
241 
242  // Used when WarpX::do_pml_divb_cleaning = true
243  std::unique_ptr<amrex::MultiFab> pml_G_fp;
244  std::unique_ptr<amrex::MultiFab> pml_G_cp;
245 
246  std::unique_ptr<MultiSigmaBox> sigba_fp;
247  std::unique_ptr<MultiSigmaBox> sigba_cp;
248 
249 #ifdef WARPX_USE_FFT
250  std::unique_ptr<SpectralSolver> spectral_solver_fp;
251  std::unique_ptr<SpectralSolver> spectral_solver_cp;
252 #endif
253 
254  // Factory for field data
255  std::unique_ptr<amrex::FabFactory<amrex::FArrayBox> > pml_field_factory;
256 
257  [[nodiscard]] amrex::FabFactory<amrex::FArrayBox> const& fieldFactory () const noexcept
258  {
259  return *pml_field_factory;
260  }
261 
262 #ifdef AMREX_USE_EB
263  amrex::EBFArrayBoxFactory const& fieldEBFactory () const noexcept {
264  return static_cast<amrex::EBFArrayBoxFactory const&>(*pml_field_factory);
265  }
266 #endif
267 
268  static amrex::BoxArray MakeBoxArray (bool single_box_domain,
269  const amrex::Box& regular_domain,
270  const amrex::Geometry& geom,
271  const amrex::BoxArray& grid_ba,
272  const amrex::IntVect& ncell,
273  int do_pml_in_domain,
274  const amrex::IntVect& do_pml_Lo,
275  const amrex::IntVect& do_pml_Hi);
276 
277  static amrex::BoxArray MakeBoxArray_single (const amrex::Box& regular_domain,
278  const amrex::BoxArray& grid_ba,
279  const amrex::IntVect& ncell,
280  const amrex::IntVect& do_pml_Lo,
281  const amrex::IntVect& do_pml_Hi);
282 
284  const amrex::BoxArray& grid_ba,
285  const amrex::IntVect& ncell,
286  int do_pml_in_domain,
287  const amrex::IntVect& do_pml_Lo,
288  const amrex::IntVect& do_pml_Hi);
289 
290  static void CopyToPML (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom);
291 };
292 
293 #ifdef WARPX_USE_FFT
294 void PushPMLPSATDSinglePatch( int lev,
295  SpectralSolver& solver,
296  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_E,
297  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_B,
298  std::unique_ptr<amrex::MultiFab>& pml_F,
299  std::unique_ptr<amrex::MultiFab>& pml_G,
300  const amrex::IntVect& fill_guards);
301 #endif
302 
303 #endif
void PushPMLPSATDSinglePatch(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:1386
Definition: PML.H:123
amrex::Real dt_E
Definition: PML.H:133
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:132
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, amrex::Real v_sigma_sb)
Definition: PML.cpp:508
Definition: PML.H:137
static amrex::BoxArray MakeBoxArray(bool single_box_domain, const amrex::Box &regular_domain, const amrex::Geometry &geom, const amrex::BoxArray &grid_ba, const amrex::IntVect &ncell, int do_pml_in_domain, const amrex::IntVect &do_pml_Lo, const amrex::IntVect &do_pml_Hi)
Definition: PML.cpp:891
const MultiSigmaBox & GetMultiSigmaBox_fp() const
Definition: PML.H:173
const amrex::Geometry * m_geom
Definition: PML.H:225
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_edge_lengths
Definition: PML.H:232
void ComputePMLFactors(amrex::Real dt)
Definition: PML.cpp:1021
amrex::MultiFab * GetG_cp()
Definition: PML.cpp:1095
void ExchangeF(amrex::MultiFab *F_fp, amrex::MultiFab *F_cp, int do_pml_in_domain)
Definition: PML.cpp:1138
std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > pml_field_factory
Definition: PML.H:255
std::unique_ptr< amrex::MultiFab > pml_F_fp
Definition: PML.H:239
const amrex::IntVect m_fill_guards_fields
Definition: PML.H:222
void PushPSATD(int lev)
Definition: PML.cpp:1376
std::array< amrex::MultiFab *, 3 > GetE_fp()
Definition: PML.cpp:1034
std::array< amrex::MultiFab *, 3 > Getj_cp()
Definition: PML.cpp:1064
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_fp
Definition: PML.H:229
amrex::MultiFab * GetF_fp()
Definition: PML.cpp:1077
void CheckPoint(const std::string &dir) const
Definition: PML.cpp:1327
static void CopyToPML(amrex::MultiFab &pml, amrex::MultiFab &reg, const amrex::Geometry &geom)
Definition: PML.cpp:1253
std::unique_ptr< MultiSigmaBox > sigba_fp
Definition: PML.H:246
void FillBoundaryG(PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition: PML.cpp:1312
const amrex::Geometry * m_cgeom
Definition: PML.H:226
bool ok() const
Definition: PML.H:209
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_fp
Definition: PML.H:230
void FillBoundaryB(PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition: PML.cpp:1280
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory() const noexcept
Definition: PML.H:257
void Exchange(const std::array< amrex::MultiFab *, 3 > &mf_pml, const std::array< amrex::MultiFab *, 3 > &mf, const PatchType &patch_type, int do_pml_in_domain)
Definition: PML.cpp:1100
bool m_divb_cleaning
Definition: PML.H:220
void FillBoundaryE(PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition: PML.cpp:1263
bool m_ok
Definition: PML.H:217
std::array< amrex::MultiFab *, 3 > Getj_fp()
Definition: PML.cpp:1046
PML(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, bool do_similar_dm_pml, 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, ablastr::utils::enums::GridType grid_type, int do_moving_window, int pml_has_particles, int do_pml_in_domain, int psatd_solution_type, int J_in_time, int rho_in_time, bool do_pml_dive_cleaning, 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, amrex::IntVect do_pml_Lo=amrex::IntVect::TheUnitVector(), amrex::IntVect do_pml_Hi=amrex::IntVect::TheUnitVector())
Definition: PML.cpp:548
std::array< amrex::MultiFab *, 3 > Get_edge_lengths()
Definition: PML.cpp:1070
bool m_dive_cleaning
Definition: PML.H:219
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_cp
Definition: PML.H:235
amrex::MultiFab * GetF_cp()
Definition: PML.cpp:1083
amrex::MultiFab * GetG_fp()
Definition: PML.cpp:1089
const MultiSigmaBox & GetMultiSigmaBox_cp() const
Definition: PML.H:178
static amrex::BoxArray MakeBoxArray_single(const amrex::Box &regular_domain, const amrex::BoxArray &grid_ba, const amrex::IntVect &ncell, const amrex::IntVect &do_pml_Lo, const amrex::IntVect &do_pml_Hi)
Definition: PML.cpp:904
const amrex::IntVect m_fill_guards_current
Definition: PML.H:223
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_cp
Definition: PML.H:234
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_cp
Definition: PML.H:236
std::unique_ptr< MultiSigmaBox > sigba_cp
Definition: PML.H:247
std::array< amrex::MultiFab *, 3 > GetE_cp()
Definition: PML.cpp:1052
static amrex::BoxArray MakeBoxArray_multiple(const amrex::Geometry &geom, const amrex::BoxArray &grid_ba, const amrex::IntVect &ncell, int do_pml_in_domain, const amrex::IntVect &do_pml_Lo, const amrex::IntVect &do_pml_Hi)
Definition: PML.cpp:945
std::unique_ptr< amrex::MultiFab > pml_F_cp
Definition: PML.H:240
std::array< amrex::MultiFab *, 3 > Get_face_areas()
std::unique_ptr< SpectralSolver > spectral_solver_cp
Definition: PML.H:251
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_fp
Definition: PML.H:228
std::unique_ptr< SpectralSolver > spectral_solver_fp
Definition: PML.H:250
void Restart(const std::string &dir)
Definition: PML.cpp:1351
std::array< amrex::MultiFab *, 3 > GetB_fp()
Definition: PML.cpp:1040
std::array< amrex::MultiFab *, 3 > GetB_cp()
Definition: PML.cpp:1058
void ExchangeG(amrex::MultiFab *G_fp, amrex::MultiFab *G_cp, int do_pml_in_domain)
Definition: PML.cpp:1154
void CopyJtoPMLs(const std::array< amrex::MultiFab *, 3 > &j_fp, const std::array< amrex::MultiFab *, 3 > &j_cp)
Definition: PML.cpp:1130
std::unique_ptr< amrex::MultiFab > pml_G_fp
Definition: PML.H:243
std::unique_ptr< amrex::MultiFab > pml_G_cp
Definition: PML.H:244
void FillBoundaryF(PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition: PML.cpp:1297
Definition: PML.H:80
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:82
amrex::IntVect m_delta
Definition: PML.H:116
SigmaBoxFactory(SigmaBoxFactory &&) noexcept=default
amrex::IntVect m_ncell
Definition: PML.H:115
amrex::Box m_regdomain
Definition: PML.H:117
SigmaBoxFactory(const SigmaBoxFactory &)=default
SigmaBoxFactory * clone() const final
Definition: PML.H:107
amrex::Real m_v_sigma_sb
Definition: PML.H:118
SigmaBoxFactory()=delete
void destroy(SigmaBox *fab) const final
Definition: PML.H:101
SigmaBox * create(const amrex::Box &box, int, const amrex::FabInfo &, int) const final
Definition: PML.H:95
const amrex::BoxArray & m_grids
Definition: PML.H:113
~SigmaBoxFactory() override=default
const amrex::Real * m_dx
Definition: PML.H:114
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:35
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVectND< dim > TheUnitVector() noexcept
GridType
Definition: Enums.H:17
PatchType
Definition: Enums.H:28
std::array< T, N > Array
default
Definition: run_alltests.py:113
float dt
Definition: stencil.py:442
tuple dx
lab frame
Definition: stencil.py:429
Definition: PML.H:46
SigmaVect sigma_star_cumsum_fac
Definition: PML.H:73
void define_multiple(const amrex::Box &box, const amrex::BoxArray &grids, const amrex::IntVect &ncell, const amrex::Array< amrex::Real, AMREX_SPACEDIM > &fac, amrex::Real v_sigma)
Definition: PML.cpp:241
SigmaVect sigma_star
Definition: PML.H:68
std::array< Sigma, AMREX_SPACEDIM > SigmaVect
Definition: PML.H:62
SigmaVect sigma_star_fac
Definition: PML.H:72
SigmaVect sigma_cumsum
Definition: PML.H:67
SigmaVect sigma_star_cumsum
Definition: PML.H:69
void ComputePMLFactorsB(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:441
void value_type
Definition: PML.H:64
SigmaVect sigma_cumsum_fac
Definition: PML.H:71
SigmaVect sigma_fac
Definition: PML.H:70
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, amrex::Real v_sigma)
Definition: PML.cpp:148
void ComputePMLFactorsE(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:475
void define_single(const amrex::Box &regdomain, const amrex::IntVect &ncell, const amrex::Array< amrex::Real, AMREX_SPACEDIM > &fac, amrex::Real v_sigma)
Definition: PML.cpp:198
SigmaVect sigma
Definition: PML.H:66
amrex::Real v_sigma
Definition: PML.H:74
Definition: PML.H:39
int m_lo
Definition: PML.H:42
int lo() const
Definition: PML.H:40
int hi() const
Definition: PML.H:41
int m_hi
Definition: PML.H:42