WarpX
Loading...
Searching...
No Matches
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
21
22#include <AMReX_BoxArray.H>
23#include <AMReX_Config.H>
24#include <AMReX_FabArray.H>
25#include <AMReX_FabFactory.H>
26#include <AMReX_GpuContainers.H>
27#include <AMReX_IntVect.H>
28#include <AMReX_REAL.H>
29#include <AMReX_Vector.H>
30
31#include <AMReX_BaseFwd.H>
32
33#include <array>
34#include <memory>
35#include <optional>
36#include <string>
37#include <vector>
38
39struct Sigma : amrex::Gpu::DeviceVector<amrex::Real>
40{
41 [[nodiscard]] int lo() const { return m_lo; }
42 [[nodiscard]] int hi() const { return m_hi; }
43 int m_lo, m_hi;
44};
45
47{
48 SigmaBox (const amrex::Box& box, const amrex::BoxArray& grids,
49 const amrex::Real* dx, const amrex::IntVect& ncell, const amrex::IntVect& delta,
50 const amrex::Box& regdomain, amrex::Real v_sigma);
51
52 void define_single (const amrex::Box& regdomain, const amrex::IntVect& ncell,
54 amrex::Real v_sigma);
55 void define_multiple (const amrex::Box& box, const amrex::BoxArray& grids,
56 const amrex::IntVect& ncell,
58 amrex::Real v_sigma);
59
60 void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
61 void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
62
63 using SigmaVect = std::array<Sigma,AMREX_SPACEDIM>;
64
65 using value_type = void; // needed by amrex::FabArray
66
75 amrex::Real v_sigma;
76
77};
78
80 : public amrex::FabFactory<SigmaBox>
81{
82public:
83 SigmaBoxFactory (const amrex::BoxArray* grid_ba, const amrex::Real* dx,
84 const amrex::IntVect& ncell, const amrex::IntVect& delta,
85 const amrex::Box& regular_domain, const amrex::Real v_sigma_sb)
86 : m_grids{grid_ba}, m_dx(dx), m_ncell(ncell), m_delta(delta), m_regdomain(regular_domain), m_v_sigma_sb(v_sigma_sb) {}
87 ~SigmaBoxFactory () override = default;
88
89 SigmaBoxFactory (const SigmaBoxFactory&) = default;
90 SigmaBoxFactory (SigmaBoxFactory&&) noexcept = default;
91
92 SigmaBoxFactory () = delete;
93 SigmaBoxFactory& operator= (const SigmaBoxFactory&) = delete;
94 SigmaBoxFactory& operator= (SigmaBoxFactory&&) = delete;
95
96 [[nodiscard]] SigmaBox* create (const amrex::Box& box, int /*ncomps*/,
97 const amrex::FabInfo& /*info*/, int /*box_index*/) const final
98 {
100 }
101
102 void destroy (SigmaBox* fab) const final
103 {
104 delete fab;
105 }
106
107 [[nodiscard]] SigmaBoxFactory*
108 clone () const final
109 {
110 return new SigmaBoxFactory(*this);
111 }
112
113private:
115 const amrex::Real* m_dx;
119 amrex::Real m_v_sigma_sb;
120};
121
123 : public amrex::FabArray<SigmaBox>
124{
125public:
127 const amrex::BoxArray* grid_ba, const amrex::Real* dx,
128 const amrex::IntVect& ncell, const amrex::IntVect& delta,
129 const amrex::Box& regular_domain, amrex::Real v_sigma_sb);
130 void ComputePMLFactorsB (const amrex::Real* dx, amrex::Real dt);
131 void ComputePMLFactorsE (const amrex::Real* dx, amrex::Real dt);
132private:
133 amrex::Real dt_B = -1.e10;
134 amrex::Real dt_E = -1.e10;
135};
136
137class PML
138{
139public:
140 PML (int lev, const amrex::BoxArray& ba,
141 const amrex::DistributionMapping& dm, bool do_similar_dm_pml,
142 const amrex::Geometry* geom, const amrex::Geometry* cgeom,
143 int ncell, int delta, amrex::IntVect ref_ratio,
144 amrex::Real dt, int nox_fft, int noy_fft, int noz_fft,
146 int do_moving_window, int pml_has_particles, int do_pml_in_domain,
147 PSATDSolutionType psatd_solution_type,
148 TimeDependencyJ time_dependency_J, TimeDependencyRho time_dependency_rho,
149 bool do_pml_dive_cleaning, bool do_pml_divb_cleaning,
150 const amrex::IntVect& fill_guards_fields,
151 const amrex::IntVect& fill_guards_current,
152 bool eb_enabled,
153 int max_guard_EB, amrex::Real v_sigma_sb,
157
158 void ComputePMLFactors (amrex::Real dt);
159
160 [[nodiscard]] const MultiSigmaBox& GetMultiSigmaBox_fp () const
161 {
162 return *sigba_fp;
163 }
164
165 [[nodiscard]] const MultiSigmaBox& GetMultiSigmaBox_cp () const
166 {
167 return *sigba_cp;
168 }
169
170#ifdef WARPX_USE_FFT
171 void PushPSATD (ablastr::fields::MultiFabRegister& fields, int lev);
172#endif
173
174 void CopyJtoPMLs (ablastr::fields::MultiFabRegister& fields, int lev);
175
178 const PatchType& patch_type,
179 int do_pml_in_domain);
180 void Exchange (amrex::MultiFab* mf_pml,
181 amrex::MultiFab* mf,
182 const PatchType& patch_type,
183 int do_pml_in_domain);
184
185 void CopyJtoPMLs (
187 PatchType patch_type,
188 int lev
189 );
190
191 void FillBoundary (ablastr::fields::VectorField mf_pml, PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
192 void FillBoundary (amrex::MultiFab & mf_pml, PatchType patch_type, std::optional<bool> nodal_sync=std::nullopt);
193
194 [[nodiscard]] bool ok () const { return m_ok; }
195
196 void CheckPoint (ablastr::fields::MultiFabRegister& fields, const std::string& dir) const;
197 void Restart (ablastr::fields::MultiFabRegister& fields, const std::string& dir);
198
199 static void Exchange (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom, int do_pml_in_domain);
200
201private:
202 bool m_ok;
203
206
209
212
213 std::unique_ptr<MultiSigmaBox> sigba_fp;
214 std::unique_ptr<MultiSigmaBox> sigba_cp;
215
216#ifdef WARPX_USE_FFT
217 std::unique_ptr<SpectralSolver> spectral_solver_fp;
218 std::unique_ptr<SpectralSolver> spectral_solver_cp;
219#endif
220
221 // Factory for field data
222 std::unique_ptr<amrex::FabFactory<amrex::FArrayBox> > pml_field_factory;
223
224#ifdef AMREX_USE_EB
225 [[nodiscard]] amrex::EBFArrayBoxFactory const& fieldEBFactory () const noexcept {
226 return static_cast<amrex::EBFArrayBoxFactory const&>(*pml_field_factory);
227 }
228#endif
229
230 static amrex::BoxArray MakeBoxArray (bool single_box_domain,
231 const amrex::Box& regular_domain,
232 const amrex::Geometry& geom,
233 const amrex::BoxArray& grid_ba,
234 const amrex::IntVect& ncell,
235 int do_pml_in_domain,
236 const amrex::IntVect& do_pml_Lo,
237 const amrex::IntVect& do_pml_Hi);
238
239 static amrex::BoxArray MakeBoxArray_single (const amrex::Box& regular_domain,
240 const amrex::BoxArray& grid_ba,
241 const amrex::IntVect& ncell,
242 const amrex::IntVect& do_pml_Lo,
243 const amrex::IntVect& do_pml_Hi);
244
246 const amrex::BoxArray& grid_ba,
247 const amrex::IntVect& ncell,
248 int do_pml_in_domain,
249 const amrex::IntVect& do_pml_Lo,
250 const amrex::IntVect& do_pml_Hi);
251
252 static void CopyToPML (amrex::MultiFab& pml, amrex::MultiFab& reg, const amrex::Geometry& geom);
253};
254
255#ifdef WARPX_USE_FFT
257 int lev,
263 const amrex::IntVect& fill_guards
264);
265#endif
266
267#endif
HYPRE_Solver solver
void PushPMLPSATDSinglePatch(int lev, SpectralSolver &solver, ablastr::fields::VectorField &pml_E, ablastr::fields::VectorField &pml_B, ablastr::fields::ScalarField pml_F, ablastr::fields::ScalarField pml_G, const amrex::IntVect &fill_guards)
Definition PML.cpp:1321
TimeDependencyJ
Definition WarpXAlgorithmSelection.H:106
@ PML
Definition WarpXAlgorithmSelection.H:138
PSATDSolutionType
Definition WarpXAlgorithmSelection.H:100
TimeDependencyRho
Definition WarpXAlgorithmSelection.H:112
Definition PML.H:124
amrex::Real dt_E
Definition PML.H:134
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:511
void ComputePMLFactorsE(const amrex::Real *dx, amrex::Real dt)
Definition PML.cpp:536
void ComputePMLFactorsB(const amrex::Real *dx, amrex::Real dt)
Definition PML.cpp:520
amrex::Real dt_B
Definition PML.H:133
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:917
const amrex::Geometry * m_geom
Definition PML.H:210
void ComputePMLFactors(amrex::Real dt)
Definition PML.cpp:1047
std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > pml_field_factory
Definition PML.H:222
const MultiSigmaBox & GetMultiSigmaBox_cp() const
Definition PML.H:165
static void CopyToPML(amrex::MultiFab &pml, amrex::MultiFab &reg, const amrex::Geometry &geom)
Definition PML.cpp:1202
std::unique_ptr< MultiSigmaBox > sigba_fp
Definition PML.H:213
const MultiSigmaBox & GetMultiSigmaBox_fp() const
Definition PML.H:160
void CopyJtoPMLs(ablastr::fields::MultiFabRegister &fields, int lev)
Definition PML.cpp:1092
const amrex::Geometry * m_cgeom
Definition PML.H:211
bool ok() const
Definition PML.H:194
void PushPSATD(ablastr::fields::MultiFabRegister &fields, int lev)
Definition PML.cpp:1302
void Exchange(ablastr::fields::VectorField mf_pml, ablastr::fields::VectorField mf, const PatchType &patch_type, int do_pml_in_domain)
Definition PML.cpp:1101
bool m_divb_cleaning
Definition PML.H:205
bool m_ok
Definition PML.H:202
bool m_dive_cleaning
Definition PML.H:204
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:930
amrex::IntVect m_fill_guards_fields
Definition PML.H:207
void CheckPoint(ablastr::fields::MultiFabRegister &fields, const std::string &dir) const
Definition PML.cpp:1235
std::unique_ptr< MultiSigmaBox > sigba_cp
Definition PML.H:214
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, PSATDSolutionType psatd_solution_type, TimeDependencyJ time_dependency_J, TimeDependencyRho time_dependency_rho, bool do_pml_dive_cleaning, bool do_pml_divb_cleaning, const amrex::IntVect &fill_guards_fields, const amrex::IntVect &fill_guards_current, bool eb_enabled, int max_guard_EB, amrex::Real v_sigma_sb, ablastr::fields::MultiFabRegister &fields, amrex::IntVect do_pml_Lo=amrex::IntVect::TheUnitVector(), amrex::IntVect do_pml_Hi=amrex::IntVect::TheUnitVector())
Definition PML.cpp:551
void Restart(ablastr::fields::MultiFabRegister &fields, const std::string &dir)
Definition PML.cpp:1268
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:971
void FillBoundary(ablastr::fields::VectorField mf_pml, PatchType patch_type, std::optional< bool > nodal_sync=std::nullopt)
Definition PML.cpp:1212
amrex::IntVect m_fill_guards_current
Definition PML.H:208
std::unique_ptr< SpectralSolver > spectral_solver_cp
Definition PML.H:218
std::unique_ptr< SpectralSolver > spectral_solver_fp
Definition PML.H:217
Definition PML.H:81
const amrex::BoxArray * m_grids
Definition PML.H:114
SigmaBox * create(const amrex::Box &box, int, const amrex::FabInfo &, int) const final
Definition PML.H:96
amrex::IntVect m_delta
Definition PML.H:117
SigmaBoxFactory(SigmaBoxFactory &&) noexcept=default
amrex::IntVect m_ncell
Definition PML.H:116
amrex::Box m_regdomain
Definition PML.H:118
SigmaBoxFactory(const SigmaBoxFactory &)=default
SigmaBoxFactory * clone() const final
Definition PML.H:108
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:83
amrex::Real m_v_sigma_sb
Definition PML.H:119
SigmaBoxFactory()=delete
void destroy(SigmaBox *fab) const final
Definition PML.H:102
~SigmaBoxFactory() override=default
const amrex::Real * m_dx
Definition PML.H:115
Top-level class for the electromagnetic spectral solver.
Definition SpectralSolver.H:37
__host__ static __device__ constexpr IntVectND< dim > TheUnitVector() noexcept
Definition EffectivePotentialPoissonSolver.H:63
std::array< amrex::MultiFab *, 3 > VectorField
Definition MultiFabRegister.H:191
amrex::MultiFab * ScalarField
Definition MultiFabRegister.H:180
GridType
Definition Enums.H:23
PatchType
Definition Enums.H:30
PODVector< T, ArenaAllocator< T > > DeviceVector
BoxND< 3 > Box
IntVectND< 3 > IntVect
std::array< T, N > Array
Definition PML.H:47
SigmaVect sigma_star_cumsum_fac
Definition PML.H:74
void define_multiple(const amrex::Box &box, const amrex::BoxArray &grids, const amrex::IntVect &ncell, const amrex::Array< amrex::Real, 3 > &fac, amrex::Real v_sigma)
Definition PML.cpp:244
void define_single(const amrex::Box &regdomain, const amrex::IntVect &ncell, const amrex::Array< amrex::Real, 3 > &fac, amrex::Real v_sigma)
Definition PML.cpp:201
SigmaVect sigma_star
Definition PML.H:69
SigmaVect sigma_star_fac
Definition PML.H:73
std::array< Sigma, 3 > SigmaVect
Definition PML.H:63
SigmaVect sigma_cumsum
Definition PML.H:68
SigmaVect sigma_star_cumsum
Definition PML.H:70
void ComputePMLFactorsB(const amrex::Real *dx, amrex::Real dt)
Definition PML.cpp:444
void value_type
Definition PML.H:65
SigmaVect sigma_cumsum_fac
Definition PML.H:72
SigmaVect sigma_fac
Definition PML.H:71
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:151
void ComputePMLFactorsE(const amrex::Real *dx, amrex::Real dt)
Definition PML.cpp:478
SigmaVect sigma
Definition PML.H:67
amrex::Real v_sigma
Definition PML.H:75
Definition PML.H:40
int m_lo
Definition PML.H:43
int lo() const
Definition PML.H:41
int hi() const
Definition PML.H:42
int m_hi
Definition PML.H:43
Definition MultiFabRegister.H:262