WarpX
CartesianNodalAlgorithm.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_NODAL_H_
9 #define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_
10 
11 #include "Utils/WarpXConst.H"
12 
13 #include <AMReX.H>
14 #include <AMReX_Array4.H>
15 #include <AMReX_Gpu.H>
16 #include <AMReX_REAL.H>
17 
18 #include <array>
19 #include <cmath>
20 
21 
27 
29  std::array<amrex::Real,3>& cell_size,
30  amrex::Vector<amrex::Real>& stencil_coefs_x,
31  amrex::Vector<amrex::Real>& stencil_coefs_y,
32  amrex::Vector<amrex::Real>& stencil_coefs_z ) {
33 
34  using namespace amrex;
35 
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 nodal solver requires one guard cell in each dimension
63  return (amrex::IntVect(AMREX_D_DECL(1,1,1)));
64  }
65 
70  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
71  static amrex::Real UpwardDx (
72  amrex::Array4<amrex::Real> const& F,
73  amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
74  int const i, int const j, int const k, int const ncomp=0 ) {
75 
76  using namespace amrex;
77  Real const inv_dx = coefs_x[0];
78  return 0.5_rt*inv_dx*( F(i+1,j,k,ncomp) - F(i-1,j,k,ncomp) );
79  }
80 
85  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
86  static amrex::Real DownwardDx (
87  amrex::Array4<amrex::Real> const& F,
88  amrex::Real const * const coefs_x, int const n_coefs_x,
89  int const i, int const j, int const k, int const ncomp=0 ) {
90 
91  return UpwardDx( F, coefs_x, n_coefs_x, i, j, k ,ncomp);
92  // For CartesianNodalAlgorithm, UpwardDx and DownwardDx are equivalent
93  }
94 
99  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
100  static amrex::Real UpwardDy (
101  amrex::Array4<amrex::Real> const& F,
102  amrex::Real const * const coefs_y, int const /*n_coefs_y*/,
103  int const i, int const j, int const k, int const ncomp=0 ) {
104 
105  using namespace amrex;
106 #if defined WARPX_DIM_3D
107  Real const inv_dy = coefs_y[0];
108  return 0.5_rt*inv_dy*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
109 #elif (defined WARPX_DIM_XZ)
110  ignore_unused(i, j, k, coefs_y, ncomp, F);
111  return 0._rt; // 2D Cartesian: derivative along y is 0
112 #endif
113  }
114 
119  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
120  static amrex::Real DownwardDy (
121  amrex::Array4<amrex::Real> const& F,
122  amrex::Real const * const coefs_y, int const n_coefs_y,
123  int const i, int const j, int const k, int const ncomp=0 ) {
124 
125  return UpwardDy( F, coefs_y, n_coefs_y, i, j, k ,ncomp);
126  // For CartesianNodalAlgorithm, UpwardDy and DownwardDy are equivalent
127  }
128 
133  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
134  static amrex::Real UpwardDz (
135  amrex::Array4<amrex::Real> const& F,
136  amrex::Real const * const coefs_z, int const /*n_coefs_z*/,
137  int const i, int const j, int const k, int const ncomp=0 ) {
138 
139  using namespace amrex;
140  Real const inv_dz = coefs_z[0];
141 #if defined WARPX_DIM_3D
142  return 0.5_rt*inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k-1,ncomp) );
143 #elif (defined WARPX_DIM_XZ)
144  return 0.5_rt*inv_dz*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
145 #endif
146  }
147 
152  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
153  static amrex::Real DownwardDz (
154  amrex::Array4<amrex::Real> const& F,
155  amrex::Real const * const coefs_z, int const n_coefs_z,
156  int const i, int const j, int const k, int const ncomp=0 ) {
157 
158  return UpwardDz( F, coefs_z, n_coefs_z, i, j, k ,ncomp);
159  // For CartesianNodalAlgorithm, UpwardDz and DownwardDz are equivalent
160  }
161 
162 };
163 
164 #endif // WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDx(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_x, int const n_coefs_x, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianNodalAlgorithm.H:86
Definition: CartesianNodalAlgorithm.H:26
static amrex::IntVect GetMaxGuardCell()
Returns maximum number of guard cells required by the field-solve.
Definition: CartesianNodalAlgorithm.H:61
int dx
Definition: compute_domain.py:35
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, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianNodalAlgorithm.H:100
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: CartesianNodalAlgorithm.H:134
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDy(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: CartesianNodalAlgorithm.H:120
cell_size
Definition: compute_domain.py:37
i
Definition: check_interp_points_and_weights.py:171
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDz(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_z, int const n_coefs_z, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianNodalAlgorithm.H:153
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: CartesianNodalAlgorithm.H:28
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: CartesianNodalAlgorithm.H:71
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition: CartesianNodalAlgorithm.H:48
Definition: BreitWheelerEngineWrapper.H:35