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 
11 #include <AMReX_MultiFab.H>
13 #include "BoundaryConditions/PML.H"
14 
22 {
23  public:
24 
25  // Constructor
35  int const fdtd_algo,
36  std::array<amrex::Real,3> cell_size,
37  bool const do_nodal );
38 
39  void EvolveB ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
40  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
41  amrex::Real const dt );
42 
43  void EvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
44  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
45  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
46  std::unique_ptr<amrex::MultiFab> const& Ffield,
47  amrex::Real const dt );
48 
49  void EvolveF ( std::unique_ptr<amrex::MultiFab>& Ffield,
50  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
51  std::unique_ptr<amrex::MultiFab> const& rhofield,
52  int const rhocomp,
53  amrex::Real const dt );
54 
55  void ComputeDivE ( const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
56  amrex::MultiFab& divE );
57 
70  void MacroscopicEvolveE ( std::array< std::unique_ptr<amrex::MultiFab>, 3>& Efield,
71  std::array< std::unique_ptr<amrex::MultiFab>, 3> const& Bfield,
72  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
73  amrex::Real const dt,
74  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
75 
76  void EvolveBPML ( std::array< amrex::MultiFab*, 3 > Bfield,
77  std::array< amrex::MultiFab*, 3 > const Efield,
78  amrex::Real const dt );
79 
80  void EvolveEPML ( std::array< amrex::MultiFab*, 3 > Efield,
81  std::array< amrex::MultiFab*, 3 > const Bfield,
82  std::array< amrex::MultiFab*, 3 > const Jfield,
83  amrex::MultiFab* const Ffield,
84  MultiSigmaBox const& sigba,
85  amrex::Real const dt, bool pml_has_particles );
86 
87  void EvolveFPML ( amrex::MultiFab* Ffield,
88  std::array< amrex::MultiFab*, 3 > const Efield,
89  amrex::Real const dt );
90 
91  private:
92 
94  bool m_do_nodal;
95 
96 #ifdef WARPX_DIM_RZ
97  amrex::Real m_dr, m_rmin;
98  amrex::Real m_nmodes;
99  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_r;
100  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_z;
101 #else
102  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_x;
103  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_y;
104  amrex::Gpu::DeviceVector<amrex::Real> m_stencil_coefs_z;
105 #endif
106 
107  public:
108  // The member functions below contain extended __device__ lambda.
109  // In order to compile with nvcc, they need to be public.
110 
111 #ifdef WARPX_DIM_RZ
112  template< typename T_Algo >
113  void EvolveBCylindrical (
114  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
115  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
116  amrex::Real const dt );
117 
118  template< typename T_Algo >
119  void EvolveECylindrical (
120  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
121  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
122  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
123  std::unique_ptr<amrex::MultiFab> const& Ffield,
124  amrex::Real const dt );
125 
126  template< typename T_Algo >
127  void EvolveFCylindrical (
128  std::unique_ptr<amrex::MultiFab>& Ffield,
129  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
130  std::unique_ptr<amrex::MultiFab> const& rhofield,
131  int const rhocomp,
132  amrex::Real const dt );
133 
134  template< typename T_Algo >
136  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
137  amrex::MultiFab& divE );
138 
139 #else
140  template< typename T_Algo >
141  void EvolveBCartesian (
142  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Bfield,
143  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
144  amrex::Real const dt );
145 
146  template< typename T_Algo >
147  void EvolveECartesian (
148  std::array< std::unique_ptr<amrex::MultiFab>, 3 >& Efield,
149  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Bfield,
150  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Jfield,
151  std::unique_ptr<amrex::MultiFab> const& Ffield,
152  amrex::Real const dt );
153 
154  template< typename T_Algo >
155  void EvolveFCartesian (
156  std::unique_ptr<amrex::MultiFab>& Ffield,
157  std::array< std::unique_ptr<amrex::MultiFab>, 3 > const& Efield,
158  std::unique_ptr<amrex::MultiFab> const& rhofield,
159  int const rhocomp,
160  amrex::Real const dt );
161 
162  template< typename T_Algo >
163  void ComputeDivECartesian (
164  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
165  amrex::MultiFab& divE );
166 
167  template< typename T_Algo, typename T_MacroAlgo >
168  void MacroscopicEvolveECartesian (
169  std::array< std::unique_ptr< amrex::MultiFab>, 3>& Efield,
170  std::array< std::unique_ptr< amrex::MultiFab>, 3> const &Bfield,
171  std::array< std::unique_ptr< amrex::MultiFab>, 3> const& Jfield,
172  amrex::Real const dt,
173  std::unique_ptr<MacroscopicProperties> const& macroscopic_properties);
174 
175  template< typename T_Algo >
176  void EvolveBPMLCartesian (
177  std::array< amrex::MultiFab*, 3 > Bfield,
178  std::array< amrex::MultiFab*, 3 > const Efield,
179  amrex::Real const dt );
180 
181  template< typename T_Algo >
182  void EvolveEPMLCartesian (
183  std::array< amrex::MultiFab*, 3 > Efield,
184  std::array< amrex::MultiFab*, 3 > const Bfield,
185  std::array< amrex::MultiFab*, 3 > const Jfield,
186  amrex::MultiFab* const Ffield,
187  MultiSigmaBox const& sigba,
188  amrex::Real const dt, bool pml_has_particles );
189 
190  template< typename T_Algo >
191  void EvolveFPMLCartesian ( amrex::MultiFab* Ffield,
192  std::array< amrex::MultiFab*, 3 > const Efield,
193  amrex::Real const dt );
194 
195 #endif
196 
197 };
198 
199 #endif // WARPX_FINITE_DIFFERENCE_SOLVER_H_
void EvolveB(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Efield, amrex::Real const dt)
Update the B field, over one timestep.
Definition: EvolveB.cpp:24
void EvolveBPML(std::array< amrex::MultiFab *, 3 > Bfield, std::array< amrex::MultiFab *, 3 > const Efield, amrex::Real const dt)
Update the B field, over one timestep.
Definition: EvolveBPML.cpp:26
void ComputeDivECylindrical(const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, amrex::MultiFab &divE)
Definition: ComputeDivE.cpp:107
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_r
Definition: FiniteDifferenceSolver.H:99
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_z
Definition: FiniteDifferenceSolver.H:100
amrex::Real m_rmin
Definition: FiniteDifferenceSolver.H:97
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:26
int m_fdtd_algo
Definition: FiniteDifferenceSolver.H:93
cell_size
Definition: compute_domain.py:37
bool m_do_nodal
Definition: FiniteDifferenceSolver.H:94
Top-level class for the electromagnetic finite-difference solver.
Definition: FiniteDifferenceSolver.H:21
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:19
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:20
void EvolveBCylindrical(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Efield, amrex::Real const dt)
Definition: EvolveB.cpp:119
Definition: PML.H:76
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, amrex::Real const dt)
Definition: EvolveE.cpp:160
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::unique_ptr< amrex::MultiFab > const &Ffield, amrex::Real const dt)
Update the E field, over one timestep.
Definition: EvolveE.cpp:26
amrex::Real m_nmodes
Definition: FiniteDifferenceSolver.H:98
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:26
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:29
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:26
amrex::Real m_dr
Definition: FiniteDifferenceSolver.H:97
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:117