WarpX
Loading...
Searching...
No Matches
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
62 // The nodal solver requires one guard cell in each dimension
63 return amrex::IntVect{AMREX_D_DECL(1,1,1)};
64 }
65
71 static amrex::Real UpwardDx (
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#if (defined WARPX_DIM_1D_Z)
78 ignore_unused(i, j, k, coefs_x, ncomp, F);
79 return 0._rt; // 1D Cartesian: derivative along x is 0
80#else
81 Real const inv_dx = coefs_x[0];
82 return 0.5_rt*inv_dx*( F(i+1,j,k,ncomp) - F(i-1,j,k,ncomp) );
83#endif
84 }
85
91 static amrex::Real DownwardDx (
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 ignore_unused(i, j, k, coefs_x, n_coefs_x, ncomp, F);
99 return 0._rt; // 1D Cartesian: derivative along x is 0
100#else
101 return UpwardDx( F, coefs_x, n_coefs_x, i, j, k ,ncomp);
102 // For CartesianNodalAlgorithm, UpwardDx and DownwardDx are equivalent
103#endif
104 }
105
108 template< typename T_Field>
110 static amrex::Real Dxx (
111 T_Field const& F,
112 amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
113 int const i, int const j, int const k, int const ncomp=0 ) {
114
115 using namespace amrex;
116#if (defined WARPX_DIM_1D_Z)
117 amrex::ignore_unused(F, coefs_x, i, j, k, ncomp);
118 return 0._rt; // 1D Cartesian: derivative along x is 0
119#else
120 amrex::Real const inv_dx2 = coefs_x[0]*coefs_x[0];
121 return inv_dx2*( F(i-1,j,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i+1,j,k,ncomp) );
122#endif
123 }
124
130 static amrex::Real UpwardDy (
132 amrex::Real const * const coefs_y, int const /*n_coefs_y*/,
133 int const i, int const j, int const k, int const ncomp=0 ) {
134
135 using namespace amrex;
136#if defined WARPX_DIM_3D
137 Real const inv_dy = coefs_y[0];
138 return 0.5_rt*inv_dy*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
139#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z)
140 ignore_unused(i, j, k, coefs_y, ncomp, F);
141 return 0._rt; // 1D and 2D Cartesian: derivative along y is 0
142#endif
143 }
144
150 static amrex::Real DownwardDy (
152 amrex::Real const * const coefs_y, int const n_coefs_y,
153 int const i, int const j, int const k, int const ncomp=0 ) {
154
155 return UpwardDy( F, coefs_y, n_coefs_y, i, j, k ,ncomp);
156 // For CartesianNodalAlgorithm, UpwardDy and DownwardDy are equivalent
157 }
158
161 template< typename T_Field>
163 static amrex::Real Dyy (
164 T_Field const& F,
165 amrex::Real const * const coefs_y, int const /*n_coefs_y*/,
166 int const i, int const j, int const k, int const ncomp=0 ) {
167
168 using namespace amrex;
169#if defined WARPX_DIM_3D
170 Real const inv_dy2 = coefs_y[0]*coefs_y[0];
171 return inv_dy2*( F(i,j-1,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j+1,k,ncomp) );
172#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
173 amrex::ignore_unused(F, coefs_y, i, j, k, ncomp);
174 return 0._rt; // 1D and 2D Cartesian: derivative along y is 0
175#endif
176 }
177
183 static amrex::Real UpwardDz (
185 amrex::Real const * const coefs_z, int const /*n_coefs_z*/,
186 int const i, int const j, int const k, int const ncomp=0 ) {
187
188 using namespace amrex;
189 Real const inv_dz = coefs_z[0];
190#if defined WARPX_DIM_3D
191 return 0.5_rt*inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k-1,ncomp) );
192#elif (defined WARPX_DIM_XZ)
193 return 0.5_rt*inv_dz*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
194#elif (defined WARPX_DIM_1D_Z)
195 return 0.5_rt*inv_dz*( F(i+1,j,k,ncomp) - F(i-1,j,k,ncomp) );
196#endif
197 }
198
204 static amrex::Real DownwardDz (
206 amrex::Real const * const coefs_z, int const n_coefs_z,
207 int const i, int const j, int const k, int const ncomp=0 ) {
208
209 return UpwardDz( F, coefs_z, n_coefs_z, i, j, k ,ncomp);
210 // For CartesianNodalAlgorithm, UpwardDz and DownwardDz are equivalent
211 }
212
213
216 template< typename T_Field>
218 static amrex::Real Dzz (
219 T_Field const& F,
220 amrex::Real const * const coefs_z, int const /*n_coefs_z*/,
221 int const i, int const j, int const k, int const ncomp=0 ) {
222
223 using namespace amrex;
224 Real const inv_dz2 = coefs_z[0]*coefs_z[0];
225#if defined WARPX_DIM_3D
226 return inv_dz2*( F(i,j,k-1,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j,k+1,ncomp) );
227#elif (defined WARPX_DIM_XZ)
228 return inv_dz2*( F(i,j-1,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j+1,k,ncomp) );
229#elif (defined WARPX_DIM_1D_Z)
230 return inv_dz2*( F(i-1,j,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i+1,j,k,ncomp) );
231#endif
232 }
233};
234
235#endif // WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_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
__host__ __device__ void ignore_unused(const Ts &...)
IntVectND< 3 > IntVect
Definition CartesianNodalAlgorithm.H:26
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDy(amrex::Array4< amrex::Real const > 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:130
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDz(amrex::Array4< amrex::Real const > 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:204
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDz(amrex::Array4< amrex::Real const > 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:183
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dzz(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 CartesianNodalAlgorithm.H:218
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition CartesianNodalAlgorithm.H:48
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dxx(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 CartesianNodalAlgorithm.H:110
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 > 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
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDx(amrex::Array4< amrex::Real const > 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:91
static amrex::IntVect GetMaxGuardCell()
Returns maximum number of guard cells required by the field-solve.
Definition CartesianNodalAlgorithm.H:61
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDy(amrex::Array4< amrex::Real const > 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:150
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dyy(T_Field 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:163