ARTEMIS
CartesianYeeAlgorithm.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_ALGORITHM_CARTESIAN_YEE_H_
9 #define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
10 
11 #include "FieldAccessorFunctors.H"
12 #include "Utils/WarpXConst.H"
13 
14 #include <AMReX.H>
15 #include <AMReX_Array4.H>
16 #include <AMReX_Gpu.H>
17 #include <AMReX_REAL.H>
18 #include "Utils/WarpXUtil.H"
19 #include "WarpX.H"
21 
22 #include <array>
23 #include <cmath>
24 
25 
31 
33  std::array<amrex::Real,3>& cell_size,
34  amrex::Vector<amrex::Real>& stencil_coefs_x,
35  amrex::Vector<amrex::Real>& stencil_coefs_y,
36  amrex::Vector<amrex::Real>& stencil_coefs_z ) {
37 
38  using namespace amrex;
39  // Store the inverse cell size along each direction in the coefficients
40  stencil_coefs_x.resize(1);
41  stencil_coefs_x[0] = 1._rt/cell_size[0];
42  stencil_coefs_y.resize(1);
43  stencil_coefs_y[0] = 1._rt/cell_size[1];
44  stencil_coefs_z.resize(1);
45  stencil_coefs_z[0] = 1._rt/cell_size[2];
46  }
47 
51  static amrex::Real ComputeMaxDt ( amrex::Real const * const dx ) {
52  using namespace amrex::literals;
53  amrex::Real const delta_t = 1._rt / ( std::sqrt( AMREX_D_TERM(
54  1._rt / (dx[0]*dx[0]),
55  + 1._rt / (dx[1]*dx[1]),
56  + 1._rt / (dx[2]*dx[2])
57  ) ) * PhysConst::c );
58  return delta_t;
59  }
60 
65  // The yee solver requires one guard cell in each dimension
66  return (amrex::IntVect(AMREX_D_DECL(1,1,1)));
67  }
68 
72  static amrex::Real UpwardDx (
74  amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
75  int const i, int const j, int const k, int const ncomp=0 ) {
76 
77  using namespace amrex;
78 #if (defined WARPX_DIM_1D_Z)
79  amrex::ignore_unused(F, coefs_x, i, j, k, ncomp);
80  return 0._rt; // 1D Cartesian: derivative along x is 0
81 #else
82  amrex::Real const inv_dx = coefs_x[0];
83  return inv_dx*( F(i+1,j,k,ncomp) - F(i,j,k,ncomp) );
84 #endif
85  }
86 
89  template< typename T_Field>
91  static amrex::Real DownwardDx (
92  T_Field const& F,
93  amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
94  int const i, int const j, int const k, int const ncomp=0 ) {
95 
96  using namespace amrex;
97 #if (defined WARPX_DIM_1D_Z)
98  amrex::ignore_unused(F, coefs_x, i, j, k, ncomp);
99  return 0._rt; // 1D Cartesian: derivative along x is 0
100 #else
101  amrex::Real const inv_dx = coefs_x[0];
102  return inv_dx*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
103 #endif
104  }
105 
109  static amrex::Real UpwardDy (
111  amrex::Real const * const coefs_y, int const n_coefs_y,
112  int const i, int const j, int const k, int const ncomp=0 ) {
113 
114  using namespace amrex;
115 #if defined WARPX_DIM_3D
116  Real const inv_dy = coefs_y[0];
117  amrex::ignore_unused(n_coefs_y);
118  return inv_dy*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
119 #elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
120  amrex::ignore_unused(F, coefs_y, n_coefs_y,
121  i, j, k, ncomp);
122  return 0._rt; // 1D and 2D Cartesian: derivative along y is 0
123 #else
124  amrex::ignore_unused(F, coefs_y, n_coefs_y,
125  i, j, k, ncomp);
126  return 0._rt;
127 #endif
128  }
129 
132  template< typename T_Field>
134  static amrex::Real DownwardDy (
135  T_Field const& F,
136  amrex::Real const * const coefs_y, int const n_coefs_y,
137  int const i, int const j, int const k, int const ncomp=0 ) {
138 
139  using namespace amrex;
140 #if defined WARPX_DIM_3D
141  Real const inv_dy = coefs_y[0];
142  amrex::ignore_unused(n_coefs_y);
143  return inv_dy*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
144 #elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
145  amrex::ignore_unused(F, coefs_y, n_coefs_y,
146  i, j, k, ncomp);
147  return 0._rt; // 1D and 2D Cartesian: derivative along y is 0
148 #else
149  amrex::ignore_unused(F, coefs_y, n_coefs_y,
150  i, j, k, ncomp);
151  return 0._rt;
152 #endif
153  }
154 
158  static amrex::Real UpwardDz (
160  amrex::Real const * const coefs_z, int const /*n_coefs_z*/,
161  int const i, int const j, int const k, int const ncomp=0 ) {
162 
163  using namespace amrex;
164  Real const inv_dz = coefs_z[0];
165 #if defined WARPX_DIM_3D
166  return inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k,ncomp) );
167 #elif (defined WARPX_DIM_XZ)
168  return inv_dz*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
169 #elif (defined WARPX_DIM_1D_Z)
170  return inv_dz*( F(i+1,j,k,ncomp) - F(i,j,k,ncomp) );
171 #endif
172  }
173 
176  template< typename T_Field>
178  static amrex::Real DownwardDz (
179  T_Field const& F,
180  amrex::Real const * const coefs_z, int const /*n_coefs_z*/,
181  int const i, int const j, int const k, int const ncomp=0 ) {
182 
183  using namespace amrex;
184  Real const inv_dz = coefs_z[0];
185 #if defined WARPX_DIM_3D
186  return inv_dz*( F(i,j,k,ncomp) - F(i,j,k-1,ncomp) );
187 #elif (defined WARPX_DIM_XZ)
188  return inv_dz*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
189 #elif (defined WARPX_DIM_1D_Z)
190  return inv_dz*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
191 #endif
192  }
193 
194 #ifdef WARPX_MAG_LLG
195 
199  static amrex::Real LaplacianDx (
201  amrex::Real const * const coefs_x,
202  int const i, int const j, int const k, int const ncomp=0 ) {
203 
204  amrex::Real const inv_dx = coefs_x[0];
205  return inv_dx*inv_dx*( F(i+1,j,k,ncomp) - 2.0*F(i,j,k,ncomp) + F(i-1,j,k,ncomp));
206  }
207 
211  static amrex::Real LaplacianDy (
213  amrex::Real const * const coefs_y,
214  int const i, int const j, int const k, int const ncomp=0 ) {
215 
216  amrex::Real const inv_dy = coefs_y[0];
217  return inv_dy*inv_dy*( F(i,j+1,k,ncomp) - 2.0*F(i,j,k,ncomp) + F(i,j-1,k,ncomp));
218  }
219 
223  static amrex::Real LaplacianDz (
225  amrex::Real const * const coefs_z,
226  int const i, int const j, int const k, int const ncomp=0 ) {
227 
228  amrex::Real const inv_dz = coefs_z[0];
229  return inv_dz*inv_dz*( F(i,j,k+1,ncomp) - 2.0*F(i,j,k,ncomp) + F(i,j,k-1,ncomp));
230  }
231 
235  static amrex::Real Laplacian (
237  amrex::Real const * const coefs_x, amrex::Real const * const coefs_y, amrex::Real const * const coefs_z,
238  int const i, int const j, int const k, int const ncomp=0 ) {
239 
240  return LaplacianDx(F, coefs_x, i, j, k, ncomp) + LaplacianDy(F, coefs_y, i, j, k, ncomp) + LaplacianDz(F, coefs_z, i, j, k, ncomp);
241  }
242 
243 #endif
244 
245 #ifdef WARPX_MAG_LLG
246 
250  static amrex::Real LaplacianDx_Mag (
252  amrex::Real const * const coefs_x, int const n_coefs_x, amrex::Real const Ms_lo_x, amrex::Real const Ms_hi_x,
253  int const i, int const j, int const k, int const ncomp=0, int const nodality=0) {
254 
255  amrex::Real const inv_dx = coefs_x[0];
256  if (nodality == 0){ // at x face (normal face). dM/dx = 0
257  if (Ms_hi_x == 0.){
258  return 0.5 * inv_dx * inv_dx * (8. * F(i-1, j, k, ncomp) - F(i-2, j, k, ncomp) - 7. * F(i, j, k, ncomp));
259  } else if (Ms_lo_x == 0){
260  return 0.5 * inv_dx * inv_dx * (8. * F(i+1, j, k, ncomp) - F(i+2, j, k, ncomp) - 7. * F(i, j, k, ncomp));
261  } else {
262  return inv_dx*(UpwardDx(F, coefs_x, n_coefs_x, i, j, k, ncomp) - DownwardDx(F, coefs_x, n_coefs_x, i, j, k, ncomp));
263  }
264  } else { // at y or z faces
265  if (Ms_hi_x == 0.){
266  return inv_dx*(0. - DownwardDx(F, coefs_x, n_coefs_x, i, j, k, ncomp));
267  } else if (Ms_lo_x == 0.){
268  return inv_dx*(UpwardDx(F, coefs_x, n_coefs_x, i, j, k, ncomp) - 0.);
269  } else {
270  return inv_dx*(UpwardDx(F, coefs_x, n_coefs_x, i, j, k, ncomp) - DownwardDx(F, coefs_x, n_coefs_x, i, j, k, ncomp));
271  }
272  }
273  }
274 
278  static amrex::Real LaplacianDy_Mag (
280  amrex::Real const * const coefs_y, int const n_coefs_y, amrex::Real const Ms_lo_y, amrex::Real const Ms_hi_y,
281  int const i, int const j, int const k, int const ncomp=0, int const nodality=0) {
282 
283  amrex::Real const inv_dy = coefs_y[0];
284  if (nodality == 1){ // y face (normal face), dM/dy = 0
285  if (Ms_hi_y == 0.){
286  return 0.5 * inv_dy * inv_dy * (8. * F(i, j-1, k, ncomp) - F(i, j-2, k, ncomp) - 7. * F(i, j, k, ncomp));
287  } else if (Ms_lo_y == 0.) {
288  return 0.5 * inv_dy * inv_dy * (8. * F(i, j+1, k, ncomp) - F(i, j+2, k, ncomp) - 7. * F(i, j, k, ncomp));
289  } else {
290  return inv_dy*(UpwardDy(F, coefs_y, n_coefs_y, i, j, k, ncomp) - DownwardDy(F, coefs_y, n_coefs_y, i, j, k, ncomp));
291  }
292 
293  } else { // x or z faces
294  if (Ms_hi_y == 0.){
295  return inv_dy*(0. - DownwardDy(F, coefs_y, n_coefs_y, i, j, k, ncomp));
296  } else if (Ms_lo_y == 0.) {
297  return inv_dy*(UpwardDy(F, coefs_y, n_coefs_y, i, j, k, ncomp) - 0.);
298  } else {
299  return inv_dy*(UpwardDy(F, coefs_y, n_coefs_y, i, j, k, ncomp) - DownwardDy(F, coefs_y, n_coefs_y, i, j, k, ncomp));
300  }
301  }
302  }
303 
307  static amrex::Real LaplacianDz_Mag (
309  amrex::Real const * const coefs_z, int const n_coefs_z, amrex::Real const Ms_lo_z, amrex::Real const Ms_hi_z,
310  int const i, int const j, int const k, int const ncomp=0, int const nodality=0) {
311 
312  amrex::Real const inv_dz = coefs_z[0];
313  if (nodality == 2){ // z face (normal face) dM/dz = 0
314  if ( Ms_hi_z == 0.) {
315  return 0.5 * inv_dz * inv_dz * (8. * F(i, j, k-1, ncomp) - F(i, j, k-2, ncomp) - 7. * F(i, j, k, ncomp));
316  } else if ( Ms_lo_z == 0.){
317  return 0.5 * inv_dz * inv_dz * (8. * F(i, j, k+1, ncomp) - F(i, j, k+2, ncomp) - 7. * F(i, j, k, ncomp));
318  } else {
319  return inv_dz*(UpwardDz(F, coefs_z, n_coefs_z, i, j, k, ncomp) - DownwardDz(F, coefs_z, n_coefs_z, i, j, k, ncomp));
320  }
321 
322  } else { // x or y face
323  if ( Ms_hi_z == 0.) {
324  return inv_dz*(0. - DownwardDz(F, coefs_z, n_coefs_z, i, j, k, ncomp));
325  } else if ( Ms_lo_z == 0.){
326  return inv_dz*(UpwardDz(F, coefs_z, n_coefs_z, i, j, k, ncomp) - 0.);
327  } else {
328  return inv_dz*(UpwardDz(F, coefs_z, n_coefs_z, i, j, k, ncomp) - DownwardDz(F, coefs_z, n_coefs_z, i, j, k, ncomp));
329  }
330  }
331  }
332 
336  static amrex::Real Laplacian_Mag (
338  amrex::Real const * const coefs_x, amrex::Real const * const coefs_y, amrex::Real const * const coefs_z,
339  int const n_coefs_x, int const n_coefs_y, int const n_coefs_z,
340  amrex::Real const Ms_lo_x, amrex::Real const Ms_hi_x, amrex::Real const Ms_lo_y, amrex::Real const Ms_hi_y, amrex::Real const Ms_lo_z, amrex::Real const Ms_hi_z,
341  int const i, int const j, int const k, int const ncomp=0, int const nodality=0) {
342 
343  return LaplacianDx_Mag(F, coefs_x, n_coefs_x, Ms_lo_x, Ms_hi_x, i, j, k, ncomp, nodality) + LaplacianDy_Mag(F, coefs_y, n_coefs_y, Ms_lo_y, Ms_hi_y, i, j, k, ncomp, nodality) + LaplacianDz_Mag(F, coefs_z, n_coefs_z, Ms_lo_z, Ms_hi_z, i, j, k, ncomp, nodality);
344  }
345 
346 #endif
347 
348 };
349 
350 #endif // WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
#define AMREX_D_TERM(a, b, c)
#define AMREX_D_DECL(a, b, c)
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
i
Definition: check_interp_points_and_weights.py:174
cell_size
Definition: compute_domain.py:37
int dx
Definition: stencil.py:436
Definition: CartesianYeeAlgorithm.H:30
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDx(T_Field const &F, amrex::Real const *const coefs_x, int const, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:91
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDz(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_z, int const, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:158
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDy(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_y, int const n_coefs_y, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:109
static amrex::IntVect GetMaxGuardCell()
Returns maximum number of guard cells required by the field-solve.
Definition: CartesianYeeAlgorithm.H:64
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDx(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_x, int const, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:72
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDy(T_Field const &F, amrex::Real const *const coefs_y, int const n_coefs_y, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:134
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDz(T_Field const &F, amrex::Real const *const coefs_z, int const, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:178
static void InitializeStencilCoefficients(std::array< amrex::Real, 3 > &cell_size, amrex::Vector< amrex::Real > &stencil_coefs_x, amrex::Vector< amrex::Real > &stencil_coefs_y, amrex::Vector< amrex::Real > &stencil_coefs_z)
Definition: CartesianYeeAlgorithm.H:32
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition: CartesianYeeAlgorithm.H:51