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, int ncell, int delta);
46 
47  void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
48  void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
49 
50  using SigmaVect = std::array<Sigma,AMREX_SPACEDIM>;
51 
60 
61 };
62 
64  : public amrex::FabFactory<SigmaBox>
65 {
66 public:
67  SigmaBoxFactory (const amrex::BoxArray& grid_ba, const amrex::Real* dx, int ncell, int delta)
68  : m_grids(grid_ba), m_dx(dx), m_ncell(ncell), m_delta(delta) {}
69  virtual ~SigmaBoxFactory () = default;
70 
71  SigmaBoxFactory (const SigmaBoxFactory&) = default;
72  SigmaBoxFactory (SigmaBoxFactory&&) noexcept = default;
73 
74  SigmaBoxFactory () = delete;
75  SigmaBoxFactory& operator= (const SigmaBoxFactory&) = delete;
76  SigmaBoxFactory& operator= (SigmaBoxFactory&&) = delete;
77 
78  virtual SigmaBox* create (const amrex::Box& box, int /*ncomps*/,
79  const amrex::FabInfo& /*info*/, int /*box_index*/) const final
80  { return new SigmaBox(box, m_grids, m_dx, m_ncell, m_delta); }
81  virtual void destroy (SigmaBox* fab) const final {
82  delete fab;
83  }
84  virtual SigmaBoxFactory* clone () const final {
85  return new SigmaBoxFactory(*this);
86  }
87 private:
88  const amrex::BoxArray& m_grids;
89  const amrex::Real* m_dx;
90  int m_ncell;
91  int m_delta;
92 };
93 
95  : public amrex::FabArray<SigmaBox>
96 {
97 public:
98  MultiSigmaBox(const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
99  const amrex::BoxArray& grid_ba, const amrex::Real* dx, int ncell, int delta);
100  void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
101  void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
102 private:
103  amrex::Real dt_B = -1.e10;
104  amrex::Real dt_E = -1.e10;
105 };
106 
107 enum struct PatchType : int;
108 
109 class PML
110 {
111 public:
112  PML (const int lev,
113  const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
114  const amrex::Geometry* geom, const amrex::Geometry* cgeom,
115  int ncell, int delta, amrex::IntVect ref_ratio,
116  amrex::Real dt, int nox_fft, int noy_fft, int noz_fft, bool do_nodal,
117  int do_moving_window, int pml_has_particles, int do_pml_in_domain,
118  const bool J_linear_in_time,
119  const bool do_pml_dive_cleaning, const bool do_pml_divb_cleaning,
120  const amrex::IntVect do_pml_Lo = amrex::IntVect::TheUnitVector(),
121  const amrex::IntVect do_pml_Hi = amrex::IntVect::TheUnitVector());
122 
123  void ComputePMLFactors (amrex::Real dt);
124 
125  std::array<amrex::MultiFab*,3> GetE_fp ();
126  std::array<amrex::MultiFab*,3> GetB_fp ();
127  std::array<amrex::MultiFab*,3> Getj_fp ();
128  std::array<amrex::MultiFab*,3> GetE_cp ();
129  std::array<amrex::MultiFab*,3> GetB_cp ();
130  std::array<amrex::MultiFab*,3> Getj_cp ();
131 
132  // Used when WarpX::do_pml_dive_cleaning = true
133  amrex::MultiFab* GetF_fp ();
134  amrex::MultiFab* GetF_cp ();
135 
136  // Used when WarpX::do_pml_divb_cleaning = true
137  amrex::MultiFab* GetG_fp ();
138  amrex::MultiFab* GetG_cp ();
139 
141  { return *sigba_fp; }
142 
144  { return *sigba_cp; }
145 
146 #ifdef WARPX_USE_PSATD
147  void PushPSATD (const int lev);
148 #endif
149 
150  void ExchangeB (const std::array<amrex::MultiFab*,3>& B_fp,
151  const std::array<amrex::MultiFab*,3>& B_cp, int do_pml_in_domain);
152  void ExchangeE (const std::array<amrex::MultiFab*,3>& E_fp,
153  const std::array<amrex::MultiFab*,3>& E_cp, int do_pml_in_domain);
154  void CopyJtoPMLs (const std::array<amrex::MultiFab*,3>& j_fp,
155  const std::array<amrex::MultiFab*,3>& j_cp);
156 
157  void ExchangeB (PatchType patch_type,
158  const std::array<amrex::MultiFab*,3>& Bp, int do_pml_in_domain);
159  void ExchangeE (PatchType patch_type,
160  const std::array<amrex::MultiFab*,3>& Ep, int do_pml_in_domain);
161 
162  void CopyJtoPMLs (PatchType patch_type,
163  const std::array<amrex::MultiFab*,3>& jp);
164 
165  void ExchangeF (amrex::MultiFab* F_fp, amrex::MultiFab* F_cp, int do_pml_in_domain);
166  void ExchangeF (PatchType patch_type, amrex::MultiFab* Fp, int do_pml_in_domain);
167 
168  void FillBoundary ();
169  void FillBoundaryE ();
170  void FillBoundaryB ();
171  void FillBoundaryF ();
172  void FillBoundaryG ();
173  void FillBoundaryE (PatchType patch_type);
174  void FillBoundaryB (PatchType patch_type);
175  void FillBoundaryF (PatchType patch_type);
176  void FillBoundaryG (PatchType patch_type);
177 
178  bool ok () const { return m_ok; }
179 
180  void CheckPoint (const std::string& dir) const;
181  void Restart (const std::string& dir);
182 
183  static void Exchange (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom, int do_pml_in_domain);
184 
185  ~PML () = default;
186 
187 private:
188  bool m_ok;
189 
192 
193  const amrex::Geometry* m_geom;
194  const amrex::Geometry* m_cgeom;
195 
196  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_fp;
197  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_fp;
198  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_fp;
199 
200  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_E_cp;
201  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_B_cp;
202  std::array<std::unique_ptr<amrex::MultiFab>,3> pml_j_cp;
203 
204  // Used when WarpX::do_pml_dive_cleaning = true
205  std::unique_ptr<amrex::MultiFab> pml_F_fp;
206  std::unique_ptr<amrex::MultiFab> pml_F_cp;
207 
208  // Used when WarpX::do_pml_divb_cleaning = true
209  std::unique_ptr<amrex::MultiFab> pml_G_fp;
210  std::unique_ptr<amrex::MultiFab> pml_G_cp;
211 
212  std::unique_ptr<MultiSigmaBox> sigba_fp;
213  std::unique_ptr<MultiSigmaBox> sigba_cp;
214 
215 #ifdef WARPX_USE_PSATD
216  std::unique_ptr<SpectralSolver> spectral_solver_fp;
217  std::unique_ptr<SpectralSolver> spectral_solver_cp;
218 #endif
219 
220  static amrex::BoxArray MakeBoxArray (const amrex::Geometry& geom,
221  const amrex::BoxArray& grid_ba,
222  int ncell, int do_pml_in_domain,
223  const amrex::IntVect do_pml_Lo = amrex::IntVect::TheUnitVector(),
224  const amrex::IntVect do_pml_Hi = amrex::IntVect::TheUnitVector());
225 
226  static void CopyToPML (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom);
227 };
228 
229 #ifdef WARPX_USE_PSATD
230 void PushPMLPSATDSinglePatch( const int lev,
231  SpectralSolver& solver,
232  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_E,
233  std::array<std::unique_ptr<amrex::MultiFab>,3>& pml_B,
234  std::unique_ptr<amrex::MultiFab>& pml_F,
235  std::unique_ptr<amrex::MultiFab>& pml_G);
236 #endif
237 
238 #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)
Definition: PML.cpp:1211
const amrex::Geometry * m_geom
Definition: PML.H:193
int lo() const
Definition: PML.H:37
bool m_dive_cleaning
Definition: PML.H:190
SigmaVect sigma_star_cumsum
Definition: PML.H:55
SigmaVect sigma_fac
Definition: PML.H:56
std::unique_ptr< SpectralSolver > spectral_solver_cp
Definition: PML.H:217
virtual void destroy(SigmaBox *fab) const final
Definition: PML.H:81
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_fp
Definition: PML.H:198
int dx
Definition: compute_domain.py:35
SigmaVect sigma_cumsum_fac
Definition: PML.H:57
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_fp
Definition: PML.H:197
bool ok() const
Definition: PML.H:178
SigmaVect sigma_star
Definition: PML.H:54
const amrex::BoxArray & m_grids
Definition: PML.H:88
int m_hi
Definition: PML.H:39
std::unique_ptr< amrex::MultiFab > pml_G_fp
Definition: PML.H:209
int m_ncell
Definition: PML.H:90
int m_delta
Definition: PML.H:91
const amrex::Real * m_dx
Definition: PML.H:89
bool m_divb_cleaning
Definition: PML.H:191
const MultiSigmaBox & GetMultiSigmaBox_fp() const
Definition: PML.H:140
Definition: PML.H:35
Definition: PML.H:42
virtual SigmaBoxFactory * clone() const final
Definition: PML.H:84
std::unique_ptr< amrex::MultiFab > pml_G_cp
Definition: PML.H:210
Definition: PML.H:63
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_B_cp
Definition: PML.H:201
SigmaVect sigma_star_cumsum_fac
Definition: PML.H:59
virtual SigmaBox * create(const amrex::Box &box, int, const amrex::FabInfo &, int) const final
Definition: PML.H:78
SigmaVect sigma_cumsum
Definition: PML.H:53
Definition: PML.H:109
SigmaVect sigma_star_fac
Definition: PML.H:58
const amrex::Geometry * m_cgeom
Definition: PML.H:194
SigmaVect sigma
Definition: PML.H:52
bool m_ok
Definition: PML.H:188
SigmaBoxFactory(const amrex::BoxArray &grid_ba, const amrex::Real *dx, int ncell, int delta)
Definition: PML.H:67
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_cp
Definition: PML.H:200
Definition: PML.H:94
std::unique_ptr< amrex::MultiFab > pml_F_fp
Definition: PML.H:205
int m_lo
Definition: PML.H:39
std::unique_ptr< amrex::MultiFab > pml_F_cp
Definition: PML.H:206
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_j_cp
Definition: PML.H:202
std::unique_ptr< SpectralSolver > spectral_solver_fp
Definition: PML.H:216
std::unique_ptr< MultiSigmaBox > sigba_cp
Definition: PML.H:213
int hi() const
Definition: PML.H:38
PatchType
Definition: WarpX.H:71
std::array< std::unique_ptr< amrex::MultiFab >, 3 > pml_E_fp
Definition: PML.H:196
std::array< Sigma, AMREX_SPACEDIM > SigmaVect
Definition: PML.H:50
std::unique_ptr< MultiSigmaBox > sigba_fp
Definition: PML.H:212
const MultiSigmaBox & GetMultiSigmaBox_cp() const
Definition: PML.H:143
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:32