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 <optional>
33 #include <string>
34 #include <vector>
35 
36 struct Sigma : amrex::Gpu::DeviceVector<amrex::Real>
37 {
38  [[nodiscard]] int lo() const { return m_lo; }
39  [[nodiscard]] int hi() const { return m_hi; }
40  int m_lo, m_hi;
41 };
42 
43 struct SigmaBox
44 {
45  SigmaBox (const amrex::Box& box, const amrex::BoxArray& grids,
46  const amrex::Real* dx, const amrex::IntVect& ncell, const amrex::IntVect& delta,
47  const amrex::Box& regdomain, amrex::Real v_sigma);
48 
49  void define_single (const amrex::Box& regdomain, const amrex::IntVect& ncell,
51  amrex::Real v_sigma);
52  void define_multiple (const amrex::Box& box, const amrex::BoxArray& grids,
53  const amrex::IntVect& ncell,
55  amrex::Real v_sigma);
56 
57  void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
58  void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
59 
60  using SigmaVect = std::array<Sigma,AMREX_SPACEDIM>;
61 
62  using value_type = void; // needed by amrex::FabArray
63 
72  amrex::Real v_sigma;
73 
74 };
75 
77  : public amrex::FabFactory<SigmaBox>
78 {
79 public:
80  SigmaBoxFactory (const amrex::BoxArray& grid_ba, const amrex::Real* dx,
81  const amrex::IntVect& ncell, const amrex::IntVect& delta,
82  const amrex::Box& regular_domain, const amrex::Real v_sigma_sb)
83  : m_grids(grid_ba), m_dx(dx), m_ncell(ncell), m_delta(delta), m_regdomain(regular_domain), m_v_sigma_sb(v_sigma_sb) {}
84  ~SigmaBoxFactory () override = default;
85 
86  SigmaBoxFactory (const SigmaBoxFactory&) = default;
88 
89  SigmaBoxFactory () = delete;
90  SigmaBoxFactory& operator= (const SigmaBoxFactory&) = delete;
91  SigmaBoxFactory& operator= (SigmaBoxFactory&&) = delete;
92 
93  [[nodiscard]] SigmaBox* create (const amrex::Box& box, int /*ncomps*/,
94  const amrex::FabInfo& /*info*/, int /*box_index*/) const final
95  {
97  }
98 
99  void destroy (SigmaBox* fab) const final
100  {
101  delete fab;
102  }
103 
104  [[nodiscard]] SigmaBoxFactory*
105  clone () const final
106  {
107  return new SigmaBoxFactory(*this);
108  }
109 
110 private:
112  const amrex::Real* m_dx;
116  amrex::Real m_v_sigma_sb;
117 };
118 
120  : public amrex::FabArray<SigmaBox>
121 {
122 public:
124  const amrex::BoxArray& grid_ba, const amrex::Real* dx,
125  const amrex::IntVect& ncell, const amrex::IntVect& delta,
126  const amrex::Box& regular_domain, amrex::Real v_sigma_sb);
127  void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
128  void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
129 private:
130  amrex::Real dt_B = -1.e10;
131  amrex::Real dt_E = -1.e10;
132 };
133 
134 enum struct PatchType : int;
135 
136 class PML
137 {
138 public:
139  PML (int lev,
140  const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
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, short grid_type,
144  int do_moving_window, int pml_has_particles, int do_pml_in_domain,
145  int psatd_solution_type, int J_in_time, int rho_in_time,
146  bool do_pml_dive_cleaning, bool do_pml_divb_cleaning,
147  const amrex::IntVect& fill_guards_fields,
148  const amrex::IntVect& fill_guards_current,
149  int max_guard_EB, amrex::Real v_sigma_sb,
152 
153  void ComputePMLFactors (amrex::Real dt);
154 
155  std::array<amrex::MultiFab*,3> GetE_fp ();
156  std::array<amrex::MultiFab*,3> GetB_fp ();
157  std::array<amrex::MultiFab*,3> Getj_fp ();
158  std::array<amrex::MultiFab*,3> GetE_cp ();
159  std::array<amrex::MultiFab*,3> GetB_cp ();
160  std::array<amrex::MultiFab*,3> Getj_cp ();
161  std::array<amrex::MultiFab*,3> Get_edge_lengths ();
162  std::array<amrex::MultiFab*,3> Get_face_areas ();
163 
164  // Used when WarpX::do_pml_dive_cleaning = true
167 
168  // Used when WarpX::do_pml_divb_cleaning = true
171 
172  [[nodiscard]] const MultiSigmaBox& GetMultiSigmaBox_fp () const
173  {
174  return *sigba_fp;
175  }
176 
177  [[nodiscard]] const MultiSigmaBox& GetMultiSigmaBox_cp () const
178  {
179  return *sigba_cp;
180  }
181 
182 #ifdef WARPX_USE_PSATD
183  void PushPSATD (int lev);
184 #endif
185 
186  void CopyJtoPMLs (const std::array<amrex::MultiFab*,3>& j_fp,
187  const std::array<amrex::MultiFab*,3>& j_cp);
188 
189  void Exchange (const std::array<amrex::MultiFab*,3>& mf_pml,
190  const std::array<amrex::MultiFab*,3>& mf,
191  const PatchType& patch_type,
192  int do_pml_in_domain);
193 
194  void CopyJtoPMLs (PatchType patch_type,
195  const std::array<amrex::MultiFab*,3>& jp);
196 
197  void ExchangeF (amrex::MultiFab* F_fp, amrex::MultiFab* F_cp, int do_pml_in_domain);
198  void ExchangeF (PatchType patch_type, amrex::MultiFab* Fp, int do_pml_in_domain);
199 
200  void ExchangeG (amrex::MultiFab* G_fp, amrex::MultiFab* G_cp, int do_pml_in_domain);
201  void ExchangeG (PatchType patch_type, amrex::MultiFab* Gp, int do_pml_in_domain);
202 
203  void FillBoundaryE (PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
204  void FillBoundaryB (PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
205  void FillBoundaryF (PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
206  void FillBoundaryG (PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
207 
208  [[nodiscard]] bool ok () const { return m_ok; }
209 
210  void CheckPoint (const std::string& dir) const;
211  void Restart (const std::string& dir);
212 
213  static void Exchange (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom, int do_pml_in_domain);
214 
215 private:
216  bool m_ok;
217 
220 
223 
226 
227  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_fp;
228  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_fp;
229  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_fp;
230 
231  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_edge_lengths;
232 
233  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_cp;
234  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_cp;
235  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_cp;
236 
237  // Used when WarpX::do_pml_dive_cleaning = true
238  std::unique_ptr<amrex::MultiFab> pml_F_fp;
239  std::unique_ptr<amrex::MultiFab> pml_F_cp;
240 
241  // Used when WarpX::do_pml_divb_cleaning = true
242  std::unique_ptr<amrex::MultiFab> pml_G_fp;
243  std::unique_ptr<amrex::MultiFab> pml_G_cp;
244 
245  std::unique_ptr<MultiSigmaBox> sigba_fp;
246  std::unique_ptr<MultiSigmaBox> sigba_cp;
247 
248 #ifdef WARPX_USE_PSATD
249  std::unique_ptr<SpectralSolver> spectral_solver_fp;
250  std::unique_ptr<SpectralSolver> spectral_solver_cp;
251 #endif
252 
253  // Factory for field data
254  std::unique_ptr<amrex::FabFactory<amrex::FArrayBox> > pml_field_factory;
255 
256  [[nodiscard]] amrex::FabFactory<amrex::FArrayBox> const& fieldFactory () const noexcept
257  {
258  return *pml_field_factory;
259  }
260 
261 #ifdef AMREX_USE_EB
262  amrex::EBFArrayBoxFactory const& fieldEBFactory () const noexcept {
263  return static_cast<amrex::EBFArrayBoxFactory const&>(*pml_field_factory);
264  }
265 #endif
266 
267  static amrex::BoxArray MakeBoxArray (bool single_box_domain,
268  const amrex::Box& regular_domain,
269  const amrex::Geometry& geom,
270  const amrex::BoxArray& grid_ba,
271  const amrex::IntVect& ncell,
272  int do_pml_in_domain,
273  const amrex::IntVect& do_pml_Lo,
274  const amrex::IntVect& do_pml_Hi);
275 
276  static amrex::BoxArray MakeBoxArray_single (const amrex::Box& regular_domain,
277  const amrex::BoxArray& grid_ba,
278  const amrex::IntVect& ncell,
279  const amrex::IntVect& do_pml_Lo,
280  const amrex::IntVect& do_pml_Hi);
281 
282  static amrex::BoxArray MakeBoxArray_multiple (const amrex::Geometry& geom,
283  const amrex::BoxArray& grid_ba,
284  const amrex::IntVect& ncell,
285  int do_pml_in_domain,
286  const amrex::IntVect& do_pml_Lo,
287  const amrex::IntVect& do_pml_Hi);
288 
289  static void CopyToPML (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom);
290 };
291 
292 #ifdef WARPX_USE_PSATD
293 void PushPMLPSATDSinglePatch( int lev,
294  SpectralSolver& solver,
295  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_E,
296  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_B,
297  std::unique_ptr<amrex::MultiFab>& pml_F,
298  std::unique_ptr<amrex::MultiFab>& pml_G,
299  const amrex::IntVect& fill_guards);
300 #endif
301 
302 #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:1359
PatchType
Definition: WarpX.H:80
Definition: PML.H:121
amrex::Real dt_E
Definition: PML.H:131
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:130
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:505
Definition: PML.H:137
const MultiSigmaBox & GetMultiSigmaBox_fp() const
Definition: PML.H:172
const amrex::Geometry * m_geom
Definition: PML.H:224
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_edge_lengths
Definition: PML.H:231
void ComputePMLFactors(amrex::Real dt)
Definition: PML.cpp:994
amrex::MultiFab * GetG_cp()
Definition: PML.cpp:1068
std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > pml_field_factory
Definition: PML.H:254
std::unique_ptr< amrex::MultiFab > pml_F_fp
Definition: PML.H:238
const amrex::IntVect m_fill_guards_fields
Definition: PML.H:221
void PushPSATD(int lev)
Definition: PML.cpp:1349
std::array< amrex::MultiFab *, 3 > GetE_fp()
Definition: PML.cpp:1007
std::array< amrex::MultiFab *, 3 > Getj_cp()
Definition: PML.cpp:1037
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_fp
Definition: PML.H:228
amrex::MultiFab * GetF_fp()
Definition: PML.cpp:1050
std::unique_ptr< MultiSigmaBox > sigba_fp
Definition: PML.H:245
const amrex::Geometry * m_cgeom
Definition: PML.H:225
bool ok() const
Definition: PML.H:208
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_fp
Definition: PML.H:229
amrex::FabFactory< amrex::FArrayBox > const & fieldFactory() const noexcept
Definition: PML.H:256
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:1073
bool m_divb_cleaning
Definition: PML.H:219
bool m_ok
Definition: PML.H:216
std::array< amrex::MultiFab *, 3 > Getj_fp()
Definition: PML.cpp:1019
std::array< amrex::MultiFab *, 3 > Get_edge_lengths()
Definition: PML.cpp:1043
bool m_dive_cleaning
Definition: PML.H:218
PML(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, 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:545
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_cp
Definition: PML.H:234
amrex::MultiFab * GetF_cp()
Definition: PML.cpp:1056
amrex::MultiFab * GetG_fp()
Definition: PML.cpp:1062
const MultiSigmaBox & GetMultiSigmaBox_cp() const
Definition: PML.H:177
const amrex::IntVect m_fill_guards_current
Definition: PML.H:222
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_cp
Definition: PML.H:233
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_cp
Definition: PML.H:235
std::unique_ptr< MultiSigmaBox > sigba_cp
Definition: PML.H:246
std::array< amrex::MultiFab *, 3 > GetE_cp()
Definition: PML.cpp:1025
std::unique_ptr< amrex::MultiFab > pml_F_cp
Definition: PML.H:239
std::array< amrex::MultiFab *, 3 > Get_face_areas()
std::unique_ptr< SpectralSolver > spectral_solver_cp
Definition: PML.H:250
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_fp
Definition: PML.H:227
std::unique_ptr< SpectralSolver > spectral_solver_fp
Definition: PML.H:249
std::array< amrex::MultiFab *, 3 > GetB_fp()
Definition: PML.cpp:1013
std::array< amrex::MultiFab *, 3 > GetB_cp()
Definition: PML.cpp:1031
void CopyJtoPMLs(const std::array< amrex::MultiFab *, 3 > &j_fp, const std::array< amrex::MultiFab *, 3 > &j_cp)
Definition: PML.cpp:1103
std::unique_ptr< amrex::MultiFab > pml_G_fp
Definition: PML.H:242
std::unique_ptr< amrex::MultiFab > pml_G_cp
Definition: PML.H:243
Definition: PML.H:78
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:80
amrex::IntVect m_delta
Definition: PML.H:114
SigmaBoxFactory(SigmaBoxFactory &&) noexcept=default
amrex::IntVect m_ncell
Definition: PML.H:113
amrex::Box m_regdomain
Definition: PML.H:115
SigmaBoxFactory(const SigmaBoxFactory &)=default
SigmaBoxFactory * clone() const final
Definition: PML.H:105
amrex::Real m_v_sigma_sb
Definition: PML.H:116
SigmaBoxFactory()=delete
void destroy(SigmaBox *fab) const final
Definition: PML.H:99
SigmaBox * create(const amrex::Box &box, int, const amrex::FabInfo &, int) const final
Definition: PML.H:93
const amrex::BoxArray & m_grids
Definition: PML.H:111
~SigmaBoxFactory() override=default
const amrex::Real * m_dx
Definition: PML.H:112
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:33
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVect TheUnitVector() noexcept
const int[]
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:44
SigmaVect sigma_star_cumsum_fac
Definition: PML.H:71
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:238
SigmaVect sigma_star
Definition: PML.H:66
std::array< Sigma, AMREX_SPACEDIM > SigmaVect
Definition: PML.H:60
SigmaVect sigma_star_fac
Definition: PML.H:70
SigmaVect sigma_cumsum
Definition: PML.H:65
SigmaVect sigma_star_cumsum
Definition: PML.H:67
void ComputePMLFactorsB(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:438
void value_type
Definition: PML.H:62
SigmaVect sigma_cumsum_fac
Definition: PML.H:69
SigmaVect sigma_fac
Definition: PML.H:68
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:145
void ComputePMLFactorsE(const amrex::Real *dx, amrex::Real dt)
Definition: PML.cpp:472
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:195
SigmaVect sigma
Definition: PML.H:64
amrex::Real v_sigma
Definition: PML.H:72
Definition: PML.H:37
int m_lo
Definition: PML.H:40
int lo() const
Definition: PML.H:38
int hi() const
Definition: PML.H:39
int m_hi
Definition: PML.H:40