WarpX
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 
19 #include <array>
20 #include <cmath>
21 
22 
28 
30  std::array<amrex::Real,3>& cell_size,
31  amrex::Vector<amrex::Real>& stencil_coefs_x,
32  amrex::Vector<amrex::Real>& stencil_coefs_y,
33  amrex::Vector<amrex::Real>& stencil_coefs_z ) {
34 
35  using namespace amrex;
36  // Store the inverse cell size along each direction in the coefficients
37  stencil_coefs_x.resize(1);
38  stencil_coefs_x[0] = 1._rt/cell_size[0];
39  stencil_coefs_y.resize(1);
40  stencil_coefs_y[0] = 1._rt/cell_size[1];
41  stencil_coefs_z.resize(1);
42  stencil_coefs_z[0] = 1._rt/cell_size[2];
43  }
44 
48  static amrex::Real ComputeMaxDt ( amrex::Real const * const dx ) {
49  using namespace amrex::literals;
50  amrex::Real const delta_t = 1._rt / ( std::sqrt( AMREX_D_TERM(
51  1._rt / (dx[0]*dx[0]),
52  + 1._rt / (dx[1]*dx[1]),
53  + 1._rt / (dx[2]*dx[2])
54  ) ) * PhysConst::c );
55  return delta_t;
56  }
57 
61  static amrex::IntVect GetMaxGuardCell () {
62  // The yee solver requires one guard cell in each dimension
63  return (amrex::IntVect(AMREX_D_DECL(1,1,1)));
64  }
65 
68  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
69  static amrex::Real UpwardDx (
70  amrex::Array4<amrex::Real> const& F,
71  amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
72  int const i, int const j, int const k, int const ncomp=0 ) {
73 
74  amrex::Real const inv_dx = coefs_x[0];
75  return inv_dx*( F(i+1,j,k,ncomp) - F(i,j,k,ncomp) );
76  }
77 
80  template< typename T_Field>
81  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
82  static amrex::Real DownwardDx (
83  T_Field const& F,
84  amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
85  int const i, int const j, int const k, int const ncomp=0 ) {
86 
87  amrex::Real const inv_dx = coefs_x[0];
88  return inv_dx*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
89  }
90 
93  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
94  static amrex::Real UpwardDy (
95  amrex::Array4<amrex::Real> const& F,
96  amrex::Real const * const coefs_y, int const n_coefs_y,
97  int const i, int const j, int const k, int const ncomp=0 ) {
98 
99  using namespace amrex;
100 #if defined WARPX_DIM_3D
101  Real const inv_dy = coefs_y[0];
102  amrex::ignore_unused(n_coefs_y);
103  return inv_dy*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
104 #elif (defined WARPX_DIM_XZ)
105  amrex::ignore_unused(F, coefs_y, n_coefs_y,
106  i, j, k, ncomp);
107  return 0._rt; // 2D Cartesian: derivative along y is 0
108 #else
109  amrex::ignore_unused(F, coefs_y, n_coefs_y,
110  i, j, k, ncomp);
111 #endif
112  }
113 
116  template< typename T_Field>
117  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
118  static amrex::Real DownwardDy (
119  T_Field const& F,
120  amrex::Real const * const coefs_y, int const n_coefs_y,
121  int const i, int const j, int const k, int const ncomp=0 ) {
122 
123  using namespace amrex;
124 #if defined WARPX_DIM_3D
125  Real const inv_dy = coefs_y[0];
126  amrex::ignore_unused(n_coefs_y);
127  return inv_dy*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
128 #elif (defined WARPX_DIM_XZ)
129  amrex::ignore_unused(F, coefs_y, n_coefs_y,
130  i, j, k, ncomp);
131  return 0._rt; // 2D Cartesian: derivative along y is 0
132 #else
133  amrex::ignore_unused(F, coefs_y, n_coefs_y,
134  i, j, k, ncomp);
135 #endif
136  }
137 
140  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
141  static amrex::Real UpwardDz (
142  amrex::Array4<amrex::Real> const& F,
143  amrex::Real const * const coefs_z, int const /*n_coefs_z*/,
144  int const i, int const j, int const k, int const ncomp=0 ) {
145 
146  using namespace amrex;
147  Real const inv_dz = coefs_z[0];
148 #if defined WARPX_DIM_3D
149  return inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k,ncomp) );
150 #elif (defined WARPX_DIM_XZ)
151  return inv_dz*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
152 #endif
153  }
154 
157  template< typename T_Field>
158  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
159  static amrex::Real DownwardDz (
160  T_Field const& F,
161  amrex::Real const * const coefs_z, int const /*n_coefs_z*/,
162  int const i, int const j, int const k, int const ncomp=0 ) {
163 
164  using namespace amrex;
165  Real const inv_dz = coefs_z[0];
166 #if defined WARPX_DIM_3D
167  return inv_dz*( F(i,j,k,ncomp) - F(i,j,k-1,ncomp) );
168 #elif (defined WARPX_DIM_XZ)
169  return inv_dz*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
170 #endif
171  }
172 
173 };
174 
175 #endif // WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
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:29
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:94
int dx
Definition: compute_domain.py:35
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:69
cell_size
Definition: compute_domain.py:37
Definition: CartesianYeeAlgorithm.H:27
i
Definition: check_interp_points_and_weights.py:171
static amrex::IntVect GetMaxGuardCell()
Returns maximum number of guard cells required by the field-solve.
Definition: CartesianYeeAlgorithm.H:61
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:159
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:118
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:82
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition: CartesianYeeAlgorithm.H:48
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:141
Definition: BreitWheelerEngineWrapper.H:35