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 
15 #include "Evolve/WarpXDtType.H"
18 
19 #include <AMReX_GpuContainers.H>
20 #include <AMReX_REAL.H>
21 
22 #include <AMReX_BaseFwd.H>
23 
24 #include <array>
25 #include <memory>
26 
34 {
35  public:
36 
37  // Constructor
47  int fdtd_algo,
48  std::array<amrex::Real,3> cell_size,
49  short grid_type );
50 
51  void EvolveB ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
52  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
53  std::unique_ptr<amrex::MultiFab> const& Gfield,
54  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
55  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& area_mod,
56  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
57  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Venl,
58  std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
59  std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
60  int lev, amrex::Real dt );
61 
62  void EvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
63  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
64  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
65  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
66  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
67  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
68  std::unique_ptr<amrex::MultiFab> const& Ffield,
69  int lev, amrex::Real dt );
70 
71  void EvolveF ( std::unique_ptr<amrex::MultiFab>& Ffield,
72  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
73  std::unique_ptr<amrex::MultiFab> const& rhofield,
74  int rhocomp,
75  amrex::Real dt );
76 
77  void EvolveG (std::unique_ptr<amrex::MultiFab>& Gfield,
78  std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
79  amrex::Real dt);
80 
81  void EvolveECTRho ( std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
82  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
83  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
84  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
85  int lev );
86 
88  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
89  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
90  amrex::Box domain_box,
91  amrex::Real dt,
92  amrex::Vector<int> field_boundary_lo,
93  amrex::Vector<int> field_boundary_hi);
94 
95  void ComputeDivE ( const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
96  amrex::MultiFab& divE );
97 
110  void MacroscopicEvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3>& Efield,
111  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
112  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
113  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
114  amrex::Real dt,
115  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
116 
117  void EvolveBPML ( std::array< amrex::MultiFab*, 3 > Bfield,
118  std::array< amrex::MultiFab*, 3 > Efield,
119  amrex::Real dt,
120  bool dive_cleaning);
121 
122  void EvolveEPML ( std::array< amrex::MultiFab*, 3 > Efield,
123  std::array< amrex::MultiFab*, 3 > Bfield,
124  std::array< amrex::MultiFab*, 3 > Jfield,
125  std::array< amrex::MultiFab*, 3 > edge_lengths,
126  amrex::MultiFab* Ffield,
127  MultiSigmaBox const& sigba,
128  amrex::Real dt, bool pml_has_particles );
129 
130  void EvolveFPML ( amrex::MultiFab* Ffield,
131  std::array< amrex::MultiFab*, 3 > Efield,
132  amrex::Real dt );
133 
150  void HybridPICSolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3>& Efield,
151  std::array< std::unique_ptr<amrex::MultiFab>, 3>& Jfield,
152  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jifield,
153  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
154  std::unique_ptr<amrex::MultiFab> const& rhofield,
155  std::unique_ptr<amrex::MultiFab> const& Pefield,
156  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
157  int lev, HybridPICModel const* hybrid_pic_model,
158  bool include_resistivity_term );
159 
170  std::array< std::unique_ptr<amrex::MultiFab>, 3>& Jfield,
171  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
172  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
173  int lev );
174 
175  private:
176 
178  short m_grid_type;
179 
180 #ifdef WARPX_DIM_RZ
181  amrex::Real m_dr, m_rmin;
182  int m_nmodes;
183  // host-only
185  // device copy after init
188 #else
189  // host-only
190  amrex::Vector<amrex::Real> m_h_stencil_coefs_x, m_h_stencil_coefs_y, m_h_stencil_coefs_z;
191  // device copy after init
192  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_x;
193  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_y;
195 #endif
196 
197  public:
198  // The member functions below contain extended __device__ lambda.
199  // In order to compile with nvcc, they need to be public.
200 
201 #ifdef WARPX_DIM_RZ
202  template< typename T_Algo >
203  void EvolveBCylindrical (
204  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
205  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
206  int lev,
207  amrex::Real dt );
208 
209  template< typename T_Algo >
210  void EvolveECylindrical (
211  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
212  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
213  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
214  std::unique_ptr<amrex::MultiFab> const& Ffield,
215  int lev,
216  amrex::Real dt );
217 
218  template< typename T_Algo >
219  void EvolveFCylindrical (
220  std::unique_ptr<amrex::MultiFab>& Ffield,
221  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
222  std::unique_ptr<amrex::MultiFab> const& rhofield,
223  int rhocomp,
224  amrex::Real dt );
225 
226  template< typename T_Algo >
228  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
229  amrex::MultiFab& divE );
230 
231  template<typename T_Algo>
233  std::array< std::unique_ptr<amrex::MultiFab>, 3>& Efield,
234  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Jfield,
235  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Jifield,
236  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
237  std::unique_ptr<amrex::MultiFab> const& rhofield,
238  std::unique_ptr<amrex::MultiFab> const& Pefield,
239  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
240  int lev, HybridPICModel const* hybrid_pic_model,
241  bool include_resistivity_term );
242 
243  template<typename T_Algo>
245  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Jfield,
246  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
247  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
248  int lev
249  );
250 
251 #else
252  template< typename T_Algo >
253  void EvolveBCartesian (
254  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
255  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
256  std::unique_ptr<amrex::MultiFab> const& Gfield,
257  int lev, amrex::Real dt );
258 
259  template< typename T_Algo >
260  void EvolveECartesian (
261  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
262  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
263  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
264  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
265  std::unique_ptr<amrex::MultiFab> const& Ffield,
266  int lev, amrex::Real dt );
267 
268  template< typename T_Algo >
269  void EvolveFCartesian (
270  std::unique_ptr<amrex::MultiFab>& Ffield,
271  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
272  std::unique_ptr<amrex::MultiFab> const& rhofield,
273  int rhocomp,
274  amrex::Real dt );
275 
276  template< typename T_Algo >
277  void EvolveGCartesian (
278  std::unique_ptr<amrex::MultiFab>& Gfield,
279  std::array<std::unique_ptr<amrex::MultiFab>,3> const& Bfield,
280  amrex::Real dt);
281 
282  void EvolveRhoCartesianECT (
283  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
284  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
285  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
286  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield, int lev);
287 
288  void EvolveBCartesianECT (
289  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
290  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& face_areas,
291  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& area_mod,
292  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& ECTRhofield,
293  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Venl,
294  std::array< std::unique_ptr<amrex::iMultiFab>, 3 >& flag_info_cell,
295  std::array< std::unique_ptr<amrex::LayoutData<FaceInfoBox> >, 3 >& borrowing,
296  int lev, amrex::Real dt
297  );
298 
299  template< typename T_Algo >
300  void ComputeDivECartesian (
301  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
302  amrex::MultiFab& divE );
303 
304  template< typename T_Algo, typename T_MacroAlgo >
305  void MacroscopicEvolveECartesian (
306  std::array< std::unique_ptr< amrex::MultiFab>, 3>& Efield,
307  std::array< std::unique_ptr< amrex::MultiFab>, 3> const& Bfield,
308  std::array< std::unique_ptr< amrex::MultiFab>, 3> const& Jfield,
309  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
310  amrex::Real dt,
311  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
312 
313  template< typename T_Algo >
314  void EvolveBPMLCartesian (
315  std::array< amrex::MultiFab*, 3 > Bfield,
316  std::array< amrex::MultiFab*, 3 > Efield,
317  amrex::Real dt,
318  bool dive_cleaning);
319 
320  template< typename T_Algo >
321  void EvolveEPMLCartesian (
322  std::array< amrex::MultiFab*, 3 > Efield,
323  std::array< amrex::MultiFab*, 3 > Bfield,
324  std::array< amrex::MultiFab*, 3 > Jfield,
325  std::array< amrex::MultiFab*, 3 > edge_lengths,
326  amrex::MultiFab* Ffield,
327  MultiSigmaBox const& sigba,
328  amrex::Real dt, bool pml_has_particles );
329 
330  template< typename T_Algo >
331  void EvolveFPMLCartesian ( amrex::MultiFab* Ffield,
332  std::array< amrex::MultiFab*, 3 > Efield,
333  amrex::Real dt );
334 
335  template<typename T_Algo>
336  void HybridPICSolveECartesian (
337  std::array< std::unique_ptr<amrex::MultiFab>, 3>& Efield,
338  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Jfield,
339  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Jifield,
340  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
341  std::unique_ptr<amrex::MultiFab> const& rhofield,
342  std::unique_ptr<amrex::MultiFab> const& Pefield,
343  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
344  int lev, HybridPICModel const* hybrid_pic_model,
345  bool include_resistivity_term );
346 
347  template<typename T_Algo>
348  void CalculateCurrentAmpereCartesian (
349  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Jfield,
350  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
351  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& edge_lengths,
352  int lev
353  );
354 #endif
355 
356 };
357 
358 #endif // WARPX_FINITE_DIFFERENCE_SOLVER_H_
Top-level class for the electromagnetic finite-difference solver.
Definition: FiniteDifferenceSolver.H:34
int m_nmodes
Definition: FiniteDifferenceSolver.H:182
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 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 CalculateCurrentAmpere(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Jfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, int lev)
Calculation of total current using Ampere's law (without displacement current): J = (curl x B) / mu0.
Definition: HybridPICSolveE.cpp:25
void EvolveG(std::unique_ptr< amrex::MultiFab > &Gfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, amrex::Real dt)
Definition: EvolveG.cpp:40
void EvolveBPML(std::array< amrex::MultiFab *, 3 > Bfield, std::array< amrex::MultiFab *, 3 > Efield, amrex::Real dt, bool dive_cleaning)
Update the B field, over one timestep.
Definition: EvolveBPML.cpp:43
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 rhocomp, amrex::Real dt)
Definition: EvolveF.cpp:137
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 rhocomp, amrex::Real dt)
Update the F field, over one timestep.
Definition: EvolveF.cpp:46
void ComputeDivECylindrical(const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, amrex::MultiFab &divE)
Definition: ComputeDivE.cpp:125
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_z
Definition: FiniteDifferenceSolver.H:187
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 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
int m_fdtd_algo
Definition: FiniteDifferenceSolver.H:177
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 dt)
Update the E field, over one timestep.
Definition: EvolveE.cpp:48
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 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, int lev)
Update the B field, over one timestep.
Definition: EvolveECTRho.cpp:49
void EvolveFPML(amrex::MultiFab *Ffield, std::array< amrex::MultiFab *, 3 > Efield, amrex::Real dt)
Update the E field, over one timestep.
Definition: EvolveFPML.cpp:42
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, int lev, amrex::Real dt)
Definition: EvolveE.cpp:234
void EvolveEPML(std::array< amrex::MultiFab *, 3 > Efield, std::array< amrex::MultiFab *, 3 > Bfield, std::array< amrex::MultiFab *, 3 > Jfield, std::array< amrex::MultiFab *, 3 > edge_lengths, amrex::MultiFab *Ffield, MultiSigmaBox const &sigba, amrex::Real dt, bool pml_has_particles)
Update the E field, over one timestep.
Definition: EvolveEPML.cpp:46
amrex::Vector< amrex::Real > m_h_stencil_coefs_z
Definition: FiniteDifferenceSolver.H:184
short m_grid_type
Definition: FiniteDifferenceSolver.H:178
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_r
Definition: FiniteDifferenceSolver.H:186
amrex::Real m_dr
Definition: FiniteDifferenceSolver.H:181
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 dt)
Update the B field, over one timestep.
Definition: EvolveB.cpp:50
amrex::Vector< amrex::Real > m_h_stencil_coefs_r
Definition: FiniteDifferenceSolver.H:184
void CalculateCurrentAmpereCylindrical(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Jfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, int lev)
Definition: HybridPICSolveE.cpp:60
FiniteDifferenceSolver(int fdtd_algo, std::array< amrex::Real, 3 > cell_size, short grid_type)
Initialize the finite-difference Maxwell solver (for a given refinement level)
Definition: FiniteDifferenceSolver.cpp:30
amrex::Real m_rmin
Definition: FiniteDifferenceSolver.H:181
void EvolveBCylindrical(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Efield, int lev, amrex::Real dt)
Definition: EvolveB.cpp:372
void HybridPICSolveECylindrical(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Jfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Jifield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::unique_ptr< amrex::MultiFab > const &rhofield, std::unique_ptr< amrex::MultiFab > const &Pefield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, int lev, HybridPICModel const *hybrid_pic_model, bool include_resistivity_term)
Definition: HybridPICSolveE.cpp:405
void HybridPICSolveE(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Jfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Jifield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::unique_ptr< amrex::MultiFab > const &rhofield, std::unique_ptr< amrex::MultiFab > const &Pefield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, int lev, HybridPICModel const *hybrid_pic_model, bool include_resistivity_term)
E-update in the hybrid PIC algorithm as described in Winske et al. (2003) Eq. 10. https://link....
Definition: HybridPICSolveE.cpp:368
This class contains the parameters needed to evaluate hybrid field solutions (kinetic ions with fluid...
Definition: HybridPICModel.H:30
Definition: PML.H:113
cell_size
Definition: compute_domain.py:37
float dt
Definition: stencil.py:442