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  short const grid_type );
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 
79  void EvolveECTRho ( std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
80  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
81  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
82  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
83  const int lev );
84 
86  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
87  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
88  amrex::Box domain_box,
89  amrex::Real const dt,
90  amrex::Vector<int> field_boundary_lo,
91  amrex::Vector<int> field_boundary_hi);
92 
93  void ComputeDivE ( const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
94  amrex::MultiFab& divE );
95 
108  void MacroscopicEvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3>& Efield,
109  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
110  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
111  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
112  amrex::Real const dt,
113  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
114 
115  void EvolveBPML ( std::array< amrex::MultiFab*, 3 > Bfield,
116  std::array< amrex::MultiFab*, 3 > const Efield,
117  amrex::Real const dt,
118  const bool dive_cleaning);
119 
120  void EvolveEPML ( std::array< amrex::MultiFab*, 3 > Efield,
121  std::array< amrex::MultiFab*, 3 > const Bfield,
122  std::array< amrex::MultiFab*, 3 > const Jfield,
123  std::array< amrex::MultiFab*, 3 > const edge_lengths,
124  amrex::MultiFab* const Ffield,
125  MultiSigmaBox const& sigba,
126  amrex::Real const dt, bool pml_has_particles );
127 
128  void EvolveFPML ( amrex::MultiFab* Ffield,
129  std::array< amrex::MultiFab*, 3 > const Efield,
130  amrex::Real const dt );
131 
132  private:
133 
135  short m_grid_type;
136 
137 #ifdef WARPX_DIM_RZ
138  amrex::Real m_dr, m_rmin;
139  int m_nmodes;
140  // host-only
142  // device copy after init
145 #else
146  // host-only
147  amrex::Vector<amrex::Real> m_h_stencil_coefs_x, m_h_stencil_coefs_y, m_h_stencil_coefs_z;
148  // device copy after init
149  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_x;
150  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_y;
152 #endif
153 
154  public:
155  // The member functions below contain extended __device__ lambda.
156  // In order to compile with nvcc, they need to be public.
157 
158 #ifdef WARPX_DIM_RZ
159  template< typename T_Algo >
160  void EvolveBCylindrical (
161  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
162  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
163  const int lev,
164  amrex::Real const dt );
165 
166  template< typename T_Algo >
167  void EvolveECylindrical (
168  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
169  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
170  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
171  std::unique_ptr<amrex::MultiFab> const& Ffield,
172  const int lev,
173  amrex::Real const dt );
174 
175  template< typename T_Algo >
176  void EvolveFCylindrical (
177  std::unique_ptr<amrex::MultiFab>& Ffield,
178  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
179  std::unique_ptr<amrex::MultiFab> const& rhofield,
180  int const rhocomp,
181  amrex::Real const dt );
182 
183  template< typename T_Algo >
185  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
186  amrex::MultiFab& divE );
187 #else
188  template< typename T_Algo >
189  void EvolveBCartesian (
190  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
191  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
192  std::unique_ptr<amrex::MultiFab> const& Gfield,
193  int lev, amrex::Real const dt );
194 
195  template< typename T_Algo >
196  void EvolveECartesian (
197  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
198  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
199  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
200  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
201  std::unique_ptr<amrex::MultiFab> const& Ffield,
202  int lev, amrex::Real const dt );
203 
204  template< typename T_Algo >
205  void EvolveFCartesian (
206  std::unique_ptr<amrex::MultiFab>& Ffield,
207  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
208  std::unique_ptr<amrex::MultiFab> const& rhofield,
209  int const rhocomp,
210  amrex::Real const dt );
211 
212  template< typename T_Algo >
213  void EvolveGCartesian (
214  std::unique_ptr<amrex::MultiFab>& Gfield,
215  std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
216  amrex::Real const dt);
217 
218  void EvolveRhoCartesianECT (
219  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
220  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
221  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
222  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield, int lev);
223 
224  void EvolveBCartesianECT (
225  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
226  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
227  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& area_mod,
228  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
229  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Venl,
230  std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
231  std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
232  int lev, amrex::Real const dt
233  );
234 
235  template< typename T_Algo >
236  void ComputeDivECartesian (
237  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
238  amrex::MultiFab& divE );
239 
240  template< typename T_Algo, typename T_MacroAlgo >
241  void MacroscopicEvolveECartesian (
242  std::array< std::unique_ptr< amrex::MultiFab>, 3>& Efield,
243  std::array< std::unique_ptr< amrex::MultiFab>, 3> const &Bfield,
244  std::array< std::unique_ptr< amrex::MultiFab>, 3> const& Jfield,
245  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
246  amrex::Real const dt,
247  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
248 
249  template< typename T_Algo >
250  void EvolveBPMLCartesian (
251  std::array< amrex::MultiFab*, 3 > Bfield,
252  std::array< amrex::MultiFab*, 3 > const Efield,
253  amrex::Real const dt,
254  const bool dive_cleaning);
255 
256  template< typename T_Algo >
257  void EvolveEPMLCartesian (
258  std::array< amrex::MultiFab*, 3 > Efield,
259  std::array< amrex::MultiFab*, 3 > const Bfield,
260  std::array< amrex::MultiFab*, 3 > const Jfield,
261  std::array< amrex::MultiFab*, 3 > const edge_lengths,
262  amrex::MultiFab* const Ffield,
263  MultiSigmaBox const& sigba,
264  amrex::Real const dt, bool pml_has_particles );
265 
266  template< typename T_Algo >
267  void EvolveFPMLCartesian ( amrex::MultiFab* Ffield,
268  std::array< amrex::MultiFab*, 3 > const Efield,
269  amrex::Real const dt );
270 #endif
271 
272 };
273 
274 #endif // WARPX_FINITE_DIFFERENCE_SOLVER_H_
Top-level class for the electromagnetic finite-difference solver.
Definition: FiniteDifferenceSolver.H:32
void EvolveECTRho(std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Efield, 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, const int lev)
Update the B field, over one timestep.
Definition: EvolveECTRho.cpp:49
int m_nmodes
Definition: FiniteDifferenceSolver.H:139
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:48
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:234
FiniteDifferenceSolver(int const fdtd_algo, std::array< amrex::Real, 3 > cell_size, short const grid_type)
Initialize the finite-difference Maxwell solver (for a given refinement level)
Definition: FiniteDifferenceSolver.cpp:30
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:37
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:50
void EvolveEPML(std::array< amrex::MultiFab *, 3 > Efield, std::array< amrex::MultiFab *, 3 > const Bfield, std::array< amrex::MultiFab *, 3 > const Jfield, std::array< amrex::MultiFab *, 3 > const edge_lengths, 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:46
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, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, 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:37
void ComputeDivECylindrical(const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, amrex::MultiFab &divE)
Definition: ComputeDivE.cpp:123
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_z
Definition: FiniteDifferenceSolver.H:144
int m_fdtd_algo
Definition: FiniteDifferenceSolver.H:134
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:42
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:137
amrex::Vector< amrex::Real > m_h_stencil_coefs_z
Definition: FiniteDifferenceSolver.H:141
short m_grid_type
Definition: FiniteDifferenceSolver.H:135
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_r
Definition: FiniteDifferenceSolver.H:143
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:46
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:40
amrex::Real m_dr
Definition: FiniteDifferenceSolver.H:138
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:43
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:42
amrex::Vector< amrex::Real > m_h_stencil_coefs_r
Definition: FiniteDifferenceSolver.H:141
amrex::Real m_rmin
Definition: FiniteDifferenceSolver.H:138
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:370
Definition: PML.H:112
cell_size
Definition: compute_domain.py:37
int dt
Definition: stencil.py:440