WarpX
FiniteDifferenceSolver.H
Go to the documentation of this file.
1 /* Copyright 2020 Remi Lehe
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 
8 #ifndef WARPX_FINITE_DIFFERENCE_SOLVER_H_
9 #define WARPX_FINITE_DIFFERENCE_SOLVER_H_
10 
13 
16 
17 #include <AMReX_GpuContainers.H>
18 #include <AMReX_REAL.H>
19 
20 #include <AMReX_BaseFwd.H>
21 
22 #include <array>
23 #include <memory>
24 
32 {
33  public:
34 
35  // Constructor
45  int const fdtd_algo,
46  std::array<amrex::Real,3> cell_size,
47  bool const do_nodal );
48 
49  void EvolveB ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
50  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
51  std::unique_ptr<amrex::MultiFab> const& Gfield,
52  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
53  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& area_mod,
54  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
55  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Venl,
56  std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
57  std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
58  int lev, amrex::Real const dt );
59 
60  void EvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
61  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
62  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
63  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
64  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
65  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
66  std::unique_ptr<amrex::MultiFab> const& Ffield,
67  int lev, amrex::Real const dt );
68 
69  void EvolveF ( std::unique_ptr<amrex::MultiFab>& Ffield,
70  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
71  std::unique_ptr<amrex::MultiFab> const& rhofield,
72  int const rhocomp,
73  amrex::Real const dt );
74 
75  void EvolveG (std::unique_ptr<amrex::MultiFab>& Gfield,
76  std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
77  amrex::Real const dt);
78 
80  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
81  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
82  amrex::Box domain_box,
83  amrex::Real const dt,
84  amrex::Vector<int> field_boundary_lo,
85  amrex::Vector<int> field_boundary_hi);
86 
87  void ComputeDivE ( const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
88  amrex::MultiFab& divE );
89 
102  void MacroscopicEvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3>& Efield,
103  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
104  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
105  amrex::Real const dt,
106  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
107 
108  void EvolveBPML ( std::array< amrex::MultiFab*, 3 > Bfield,
109  std::array< amrex::MultiFab*, 3 > const Efield,
110  amrex::Real const dt,
111  const bool dive_cleaning);
112 
113  void EvolveEPML ( std::array< amrex::MultiFab*, 3 > Efield,
114  std::array< amrex::MultiFab*, 3 > const Bfield,
115  std::array< amrex::MultiFab*, 3 > const Jfield,
116  amrex::MultiFab* const Ffield,
117  MultiSigmaBox const& sigba,
118  amrex::Real const dt, bool pml_has_particles );
119 
120  void EvolveFPML ( amrex::MultiFab* Ffield,
121  std::array< amrex::MultiFab*, 3 > const Efield,
122  amrex::Real const dt );
123 
124  private:
125 
128 
129 #ifdef WARPX_DIM_RZ
130  amrex::Real m_dr, m_rmin;
131  int m_nmodes;
132  // host-only
133  amrex::Vector<amrex::Real> m_h_stencil_coefs_r, m_h_stencil_coefs_z;
134  // device copy after init
135  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_r;
136  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_z;
137 #else
138  // host-only
139  amrex::Vector<amrex::Real> m_h_stencil_coefs_x, m_h_stencil_coefs_y, m_h_stencil_coefs_z;
140  // device copy after init
141  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_x;
142  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_y;
143  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_z;
144 #endif
145 
146  public:
147  // The member functions below contain extended __device__ lambda.
148  // In order to compile with nvcc, they need to be public.
149 
150 #ifdef WARPX_DIM_RZ
151  template< typename T_Algo >
152  void EvolveBCylindrical (
153  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
154  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
155  const int lev,
156  amrex::Real const dt );
157 
158  template< typename T_Algo >
159  void EvolveECylindrical (
160  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
161  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
162  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
163  std::unique_ptr<amrex::MultiFab> const& Ffield,
164  const int lev,
165  amrex::Real const dt );
166 
167  template< typename T_Algo >
168  void EvolveFCylindrical (
169  std::unique_ptr<amrex::MultiFab>& Ffield,
170  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
171  std::unique_ptr<amrex::MultiFab> const& rhofield,
172  int const rhocomp,
173  amrex::Real const dt );
174 
175  template< typename T_Algo >
177  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
178  amrex::MultiFab& divE );
179 #else
180  template< typename T_Algo >
181  void EvolveBCartesian (
182  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
183  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
184  std::unique_ptr<amrex::MultiFab> const& Gfield,
185  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
186  int lev, amrex::Real const dt );
187 
188  template< typename T_Algo >
189  void EvolveECartesian (
190  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
191  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
192  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
193  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
194  std::unique_ptr<amrex::MultiFab> const& Ffield,
195  int lev, amrex::Real const dt );
196 
197  template< typename T_Algo >
198  void EvolveFCartesian (
199  std::unique_ptr<amrex::MultiFab>& Ffield,
200  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
201  std::unique_ptr<amrex::MultiFab> const& rhofield,
202  int const rhocomp,
203  amrex::Real const dt );
204 
205  template< typename T_Algo >
206  void EvolveGCartesian (
207  std::unique_ptr<amrex::MultiFab>& Gfield,
208  std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
209  amrex::Real const dt);
210 
211  void EvolveRhoCartesianECT (
212  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
213  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
214  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
215  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield, int lev);
216 
217  void EvolveBCartesianECT (
218  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
219  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
220  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& area_mod,
221  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
222  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Venl,
223  std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
224  std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
225  int lev, amrex::Real const dt
226  );
227 
228  template< typename T_Algo >
229  void ComputeDivECartesian (
230  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
231  amrex::MultiFab& divE );
232 
233  template< typename T_Algo, typename T_MacroAlgo >
234  void MacroscopicEvolveECartesian (
235  std::array< std::unique_ptr< amrex::MultiFab>, 3>& Efield,
236  std::array< std::unique_ptr< amrex::MultiFab>, 3> const &Bfield,
237  std::array< std::unique_ptr< amrex::MultiFab>, 3> const& Jfield,
238  amrex::Real const dt,
239  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
240 
241  template< typename T_Algo >
242  void EvolveBPMLCartesian (
243  std::array< amrex::MultiFab*, 3 > Bfield,
244  std::array< amrex::MultiFab*, 3 > const Efield,
245  amrex::Real const dt,
246  const bool dive_cleaning);
247 
248  template< typename T_Algo >
249  void EvolveEPMLCartesian (
250  std::array< amrex::MultiFab*, 3 > Efield,
251  std::array< amrex::MultiFab*, 3 > const Bfield,
252  std::array< amrex::MultiFab*, 3 > const Jfield,
253  amrex::MultiFab* const Ffield,
254  MultiSigmaBox const& sigba,
255  amrex::Real const dt, bool pml_has_particles );
256 
257  template< typename T_Algo >
258  void EvolveFPMLCartesian ( amrex::MultiFab* Ffield,
259  std::array< amrex::MultiFab*, 3 > const Efield,
260  amrex::Real const dt );
261 #endif
262 
263 };
264 
265 #endif // WARPX_FINITE_DIFFERENCE_SOLVER_H_
void EvolveECylindrical(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Jfield, std::unique_ptr< amrex::MultiFab > const &Ffield, const int lev, amrex::Real const dt)
Definition: EvolveE.cpp:233
void EvolveBCylindrical(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Efield, const int lev, amrex::Real const dt)
Definition: EvolveB.cpp:436
void ComputeDivECylindrical(const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, amrex::MultiFab &divE)
Definition: ComputeDivE.cpp:122
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_r
Definition: FiniteDifferenceSolver.H:135
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_z
Definition: FiniteDifferenceSolver.H:136
amrex::Real m_rmin
Definition: FiniteDifferenceSolver.H:130
void EvolveE(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Jfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &face_areas, std::array< std::unique_ptr< amrex::MultiFab >, 3 > &ECTRhofield, std::unique_ptr< amrex::MultiFab > const &Ffield, int lev, amrex::Real const dt)
Update the E field, over one timestep.
Definition: EvolveE.cpp:47
void EvolveF(std::unique_ptr< amrex::MultiFab > &Ffield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Efield, std::unique_ptr< amrex::MultiFab > const &rhofield, int const rhocomp, amrex::Real const dt)
Update the F field, over one timestep.
Definition: EvolveF.cpp:45
void EvolveB(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Efield, std::unique_ptr< amrex::MultiFab > const &Gfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &face_areas, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &area_mod, std::array< std::unique_ptr< amrex::MultiFab >, 3 > &ECTRhofield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Venl, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > &flag_info_cell, std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > &borrowing, int lev, amrex::Real const dt)
Update the B field, over one timestep.
Definition: EvolveB.cpp:48
int m_fdtd_algo
Definition: FiniteDifferenceSolver.H:126
cell_size
Definition: compute_domain.py:37
void EvolveG(std::unique_ptr< amrex::MultiFab > &Gfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, amrex::Real const dt)
Definition: EvolveG.cpp:39
bool m_do_nodal
Definition: FiniteDifferenceSolver.H:127
Top-level class for the electromagnetic finite-difference solver.
Definition: FiniteDifferenceSolver.H:31
amrex::Vector< amrex::Real > m_h_stencil_coefs_z
Definition: FiniteDifferenceSolver.H:133
void MacroscopicEvolveE(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Jfield, amrex::Real const dt, std::unique_ptr< MacroscopicProperties > const &macroscopic_properties)
Macroscopic E-update for non-vacuum medium using the user-selected finite-difference algorithm and ma...
Definition: MacroscopicEvolveE.cpp:35
FiniteDifferenceSolver(int const fdtd_algo, std::array< amrex::Real, 3 > cell_size, bool const do_nodal)
Initialize the finite-difference Maxwell solver (for a given refinement level)
Definition: FiniteDifferenceSolver.cpp:29
void ApplySilverMuellerBoundary(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Bfield, amrex::Box domain_box, amrex::Real const dt, amrex::Vector< int > field_boundary_lo, amrex::Vector< int > field_boundary_hi)
Update the B field at the boundary, using the Silver-Mueller condition.
Definition: ApplySilverMuellerBoundary.cpp:36
Definition: PML.H:94
void ComputeDivE(const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, amrex::MultiFab &divE)
Update the F field, over one timestep.
Definition: ComputeDivE.cpp:41
void EvolveEPML(std::array< amrex::MultiFab *, 3 > Efield, std::array< amrex::MultiFab *, 3 > const Bfield, std::array< amrex::MultiFab *, 3 > const Jfield, amrex::MultiFab *const Ffield, MultiSigmaBox const &sigba, amrex::Real const dt, bool pml_has_particles)
Update the E field, over one timestep.
Definition: EvolveEPML.cpp:45
void EvolveBPML(std::array< amrex::MultiFab *, 3 > Bfield, std::array< amrex::MultiFab *, 3 > const Efield, amrex::Real const dt, const bool dive_cleaning)
Update the B field, over one timestep.
Definition: EvolveBPML.cpp:42
amrex::Vector< amrex::Real > m_h_stencil_coefs_r
Definition: FiniteDifferenceSolver.H:133
void EvolveFPML(amrex::MultiFab *Ffield, std::array< amrex::MultiFab *, 3 > const Efield, amrex::Real const dt)
Update the E field, over one timestep.
Definition: EvolveFPML.cpp:41
amrex::Real m_dr
Definition: FiniteDifferenceSolver.H:130
int m_nmodes
Definition: FiniteDifferenceSolver.H:131
void EvolveFCylindrical(std::unique_ptr< amrex::MultiFab > &Ffield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Efield, std::unique_ptr< amrex::MultiFab > const &rhofield, int const rhocomp, amrex::Real const dt)
Definition: EvolveF.cpp:136