ARTEMIS
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 EvolveBLondon ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
50  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& current,
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, amrex::Real const penetration_depth );
59 
60  void EvolveB ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
61  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
62  std::unique_ptr<amrex::MultiFab> const& Gfield,
63  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
64  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& area_mod,
65  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
66  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Venl,
67  std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
68  std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
69  int lev, amrex::Real const dt );
70 
71  void EvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
72  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
73  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
74  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
75  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
76  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
77  std::unique_ptr<amrex::MultiFab> const& Ffield,
78  int lev, amrex::Real const dt );
79 
80  void EvolveF ( std::unique_ptr<amrex::MultiFab>& Ffield,
81  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
82  std::unique_ptr<amrex::MultiFab> const& rhofield,
83  int const rhocomp,
84  amrex::Real const dt );
85 
86  void EvolveG (std::unique_ptr<amrex::MultiFab>& Gfield,
87  std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
88  amrex::Real const dt);
89 
90  void EvolveECTRho ( std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
91  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
92  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
93  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
94  const int lev );
95 
97  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
98  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
99  amrex::Box domain_box,
100  amrex::Real const dt,
101  amrex::Vector<int> field_boundary_lo,
102  amrex::Vector<int> field_boundary_hi);
103 
104  void ComputeDivE ( const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
105  amrex::MultiFab& divE );
106 
119  void MacroscopicEvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3>& Efield,
120 #ifndef WARPX_MAG_LLG
121  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
122 #else
123  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Hfield,
124 #endif
125  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
126  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
127  amrex::Real const dt,
128  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
129 #ifndef WARPX_DIM_RZ
130 #ifdef WARPX_MAG_LLG
147  void MacroscopicEvolveHM (
148  std::array<std::unique_ptr<amrex::MultiFab>, 3> &Mfield,
149  std::array<std::unique_ptr<amrex::MultiFab>, 3> &Hfield, // H Maxwell
150  std::array< std::unique_ptr<amrex::MultiFab>, 3>& Bfield,
151  std::array<std::unique_ptr<amrex::MultiFab>, 3> const &H_biasfield, // H bias
152  std::array<std::unique_ptr<amrex::MultiFab>, 3> const &Efield,
153  amrex::Real const dt,
154  std::unique_ptr<MacroscopicProperties> const &macroscopic_properties);
155 
156  void MacroscopicEvolveHM_2nd (
157  int lev,
158  std::array<std::unique_ptr<amrex::MultiFab>, 3> &Mfield,
159  std::array<std::unique_ptr<amrex::MultiFab>, 3> &Hfield, // H Maxwell
160  std::array< std::unique_ptr<amrex::MultiFab>, 3>& Bfield,
161  std::array<std::unique_ptr<amrex::MultiFab>, 3> const &H_biasfield, // H bias
162  std::array<std::unique_ptr<amrex::MultiFab>, 3> const &Efield,
163  amrex::Real const dt,
164  std::unique_ptr<MacroscopicProperties> const &macroscopic_properties);
165 
166 #endif
167 #endif // ifndef WARPX_DIM_RZ
168 
169  void EvolveBPML ( std::array< amrex::MultiFab*, 3 > Bfield,
170  std::array< amrex::MultiFab*, 3 > const Efield,
171  amrex::Real const dt,
172  const bool dive_cleaning);
173 
174  void EvolveEPML ( std::array< amrex::MultiFab*, 3 > Efield,
175  std::array< amrex::MultiFab*, 3 > const Bfield,
176  std::array< amrex::MultiFab*, 3 > const Jfield,
177  std::array< amrex::MultiFab*, 3 > const edge_lengths,
178  amrex::MultiFab* const Ffield,
179  MultiSigmaBox const& sigba,
180  amrex::Real const dt, bool pml_has_particles );
181 
182  void EvolveFPML ( amrex::MultiFab* Ffield,
183  std::array< amrex::MultiFab*, 3 > const Efield,
184  amrex::Real const dt );
185 
186  void MacroscopicEvolveEPML ( std::array< amrex::MultiFab*, 3 > Efield,
187 #ifndef WARPX_MAG_LLG
188  std::array< amrex::MultiFab*, 3 > const Bfield,
189 #else
190  std::array< amrex::MultiFab*, 3 > const Hfield,
191 #endif
192  std::array< amrex::MultiFab*, 3 > const Jfield,
193  amrex::MultiFab* const Ffield,
194  MultiSigmaBox const& sigba,
195  amrex::Real const dt, bool pml_has_particles,
196  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties,
197  amrex::MultiFab* const eps_mf,
198  amrex::MultiFab* const mu_mf,
199  amrex::MultiFab* const sigma_mf);
200 
201 #ifndef WARPX_DIM_RZ
202 #ifdef WARPX_MAG_LLG
203  void EvolveHPML ( std::array< amrex::MultiFab*, 3 > Hfield,
204  std::array< amrex::MultiFab*, 3 > const Efield,
205  amrex::Real const dt,
206  const bool dive_cleaning);
207 #endif
208 #endif // ifndef WARPX_DIM_RZ
209 
210  private:
211 
213  short m_grid_type;
214 
215 #ifdef WARPX_DIM_RZ
216  amrex::Real m_dr, m_rmin;
217  int m_nmodes;
218  // host-only
220  // device copy after init
223 #else
224  // host-only
225  amrex::Vector<amrex::Real> m_h_stencil_coefs_x, m_h_stencil_coefs_y, m_h_stencil_coefs_z;
226  // device copy after init
227  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_x;
228  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_y;
230 #endif
231 
232  public:
233  // The member functions below contain extended __device__ lambda.
234  // In order to compile with nvcc, they need to be public.
235 
236 #ifdef WARPX_DIM_RZ
237  template< typename T_Algo >
238  void EvolveBCylindrical (
239  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
240  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
241  const int lev,
242  amrex::Real const dt );
243 
244  template< typename T_Algo >
245  void EvolveECylindrical (
246  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
247  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
248  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
249  std::unique_ptr<amrex::MultiFab> const& Ffield,
250  const int lev,
251  amrex::Real const dt );
252 
253  template< typename T_Algo >
254  void EvolveFCylindrical (
255  std::unique_ptr<amrex::MultiFab>& Ffield,
256  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
257  std::unique_ptr<amrex::MultiFab> const& rhofield,
258  int const rhocomp,
259  amrex::Real const dt );
260 
261  template< typename T_Algo >
263  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
264  amrex::MultiFab& divE );
265 #else
266  template< typename T_Algo >
267  void EvolveBLondonCartesian (
268  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
269  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& current,
270  std::unique_ptr<amrex::MultiFab> const& Gfield,
271  int lev, amrex::Real const dt, amrex::Real const penetration_depth );
272 
273  template< typename T_Algo >
274  void EvolveBCartesian (
275  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
276  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
277  std::unique_ptr<amrex::MultiFab> const& Gfield,
278  int lev, amrex::Real const dt );
279 
280  template< typename T_Algo >
281  void EvolveECartesian (
282  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
283  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
284  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
285  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
286  std::unique_ptr<amrex::MultiFab> const& Ffield,
287  int lev, amrex::Real const dt );
288 
289  template< typename T_Algo >
290  void EvolveFCartesian (
291  std::unique_ptr<amrex::MultiFab>& Ffield,
292  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
293  std::unique_ptr<amrex::MultiFab> const& rhofield,
294  int const rhocomp,
295  amrex::Real const dt );
296 
297  template< typename T_Algo >
298  void EvolveGCartesian (
299  std::unique_ptr<amrex::MultiFab>& Gfield,
300  std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
301  amrex::Real const dt);
302 
303  void EvolveRhoCartesianECT (
304  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
305  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
306  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
307  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield, int lev);
308 
309  void EvolveBCartesianECT (
310  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
311  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
312  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& area_mod,
313  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
314  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Venl,
315  std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
316  std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
317  int lev, amrex::Real const dt
318  );
319 
320  template< typename T_Algo >
321  void ComputeDivECartesian (
322  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
323  amrex::MultiFab& divE );
324 
325  template< typename T_Algo, typename T_MacroAlgo >
326  void MacroscopicEvolveECartesian (
327  std::array< std::unique_ptr< amrex::MultiFab>, 3>& Efield,
328 #ifndef WARPX_MAG_LLG
329  std::array< std::unique_ptr< amrex::MultiFab>, 3> const &Bfield,
330 #else
331  std::array< std::unique_ptr< amrex::MultiFab>, 3> const &Hfield,
332 #endif
333  std::array< std::unique_ptr< amrex::MultiFab>, 3> const& Jfield,
334  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
335  amrex::Real const dt,
336  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
337 
338 #ifdef WARPX_MAG_LLG
339  template< typename T_Algo >
340  void MacroscopicEvolveHMCartesian(
341  std::array<std::unique_ptr<amrex::MultiFab>, 3> &Mfield,
342  std::array<std::unique_ptr<amrex::MultiFab>, 3> &Hfield, // H Maxwell
343  std::array< std::unique_ptr<amrex::MultiFab>, 3>& Bfield,
344  std::array<std::unique_ptr<amrex::MultiFab>, 3> const &H_biasfield, // H bias
345  std::array<std::unique_ptr<amrex::MultiFab>, 3> const &Efield,
346  amrex::Real const dt,
347  std::unique_ptr<MacroscopicProperties> const &macroscopic_properties);
348 
349  template< typename T_Algo >
350  void MacroscopicEvolveHMCartesian_2nd(
351  int lev,
352  std::array<std::unique_ptr<amrex::MultiFab>, 3> &Mfield,
353  std::array<std::unique_ptr<amrex::MultiFab>, 3> &Hfield, // H Maxwell
354  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
355  std::array<std::unique_ptr<amrex::MultiFab>, 3> const &H_biasfield, // H bias
356  std::array<std::unique_ptr<amrex::MultiFab>, 3> const &Efield,
357  amrex::Real const dt,
358  std::unique_ptr<MacroscopicProperties> const &macroscopic_properties);
359 #endif
360 
361  template< typename T_Algo >
362  void EvolveBPMLCartesian (
363  std::array< amrex::MultiFab*, 3 > Bfield,
364  std::array< amrex::MultiFab*, 3 > const Efield,
365  amrex::Real const dt,
366  const bool dive_cleaning);
367 
368  template< typename T_Algo >
369  void EvolveEPMLCartesian (
370  std::array< amrex::MultiFab*, 3 > Efield,
371  std::array< amrex::MultiFab*, 3 > const Bfield,
372  std::array< amrex::MultiFab*, 3 > const Jfield,
373  std::array< amrex::MultiFab*, 3 > const edge_lengths,
374  amrex::MultiFab* const Ffield,
375  MultiSigmaBox const& sigba,
376  amrex::Real const dt, bool pml_has_particles );
377 
378  template< typename T_Algo >
379  void EvolveFPMLCartesian ( amrex::MultiFab* Ffield,
380  std::array< amrex::MultiFab*, 3 > const Efield,
381  amrex::Real const dt );
382 
383  template< typename T_Algo, typename T_MacroAlgo >
384  void MacroscopicEvolveEPMLCartesian (
385  std::array< amrex::MultiFab*, 3 > Efield,
386 #ifndef WARPX_MAG_LLG
387  std::array< amrex::MultiFab*, 3 > const Bfield,
388 #else
389  std::array< amrex::MultiFab*, 3 > const Hfield,
390 #endif
391  std::array< amrex::MultiFab*, 3 > const Jfield,
392  amrex::MultiFab* const Ffield,
393  MultiSigmaBox const& sigba,
394  amrex::Real const dt, bool pml_has_particles,
395  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties,
396  amrex::MultiFab* const eps_mf,
397  amrex::MultiFab* const mu_mf,
398  amrex::MultiFab* const sigma_mf);
399 
400 
401 #ifdef WARPX_MAG_LLG
402  template< typename T_Algo >
403  void EvolveHPMLCartesian (
404  std::array< amrex::MultiFab*, 3 > Bfield,
405  std::array< amrex::MultiFab*, 3 > const Efield,
406  amrex::Real const dt,
407  const bool dive_cleaning);
408 #endif
409 
410 #endif
411 
412 };
413 
414 #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:217
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 EvolveBLondon(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &current, 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, amrex::Real const penetration_depth)
Update the B field, over one timestep.
Definition: EvolveBLondon.cpp:48
void ComputeDivECylindrical(const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, amrex::MultiFab &divE)
Definition: ComputeDivE.cpp:123
void MacroscopicEvolveEPML(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, std::unique_ptr< MacroscopicProperties > const &macroscopic_properties, amrex::MultiFab *const eps_mf, amrex::MultiFab *const mu_mf, amrex::MultiFab *const sigma_mf)
Update the E field, over one timestep.
Definition: MacroscopicEvolveEPML.cpp:24
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_z
Definition: FiniteDifferenceSolver.H:222
int m_fdtd_algo
Definition: FiniteDifferenceSolver.H:212
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:219
short m_grid_type
Definition: FiniteDifferenceSolver.H:213
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_r
Definition: FiniteDifferenceSolver.H:221
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:216
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:219
amrex::Real m_rmin
Definition: FiniteDifferenceSolver.H:216
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