aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/FiniteDifferenceSolver
diff options
context:
space:
mode:
authorGravatar Remi Lehe <remi.lehe@normalesup.org> 2020-01-27 16:08:21 -0800
committerGravatar Remi Lehe <remi.lehe@normalesup.org> 2020-01-27 16:16:16 -0800
commit0e14e871a29027ded4258b079abd63ada94bc874 (patch)
tree3bd320bc55578f90928af1defa5f7f71f9e79e09 /Source/FieldSolver/FiniteDifferenceSolver
parent6d77161d6e80b943230c2969c15674fe044cdb30 (diff)
downloadWarpX-0e14e871a29027ded4258b079abd63ada94bc874.tar.gz
WarpX-0e14e871a29027ded4258b079abd63ada94bc874.tar.zst
WarpX-0e14e871a29027ded4258b079abd63ada94bc874.zip
Implement nodal solver
Diffstat (limited to 'Source/FieldSolver/FiniteDifferenceSolver')
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp5
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/NodalAlgorithm.H103
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H14
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.cpp21
4 files changed, 129 insertions, 14 deletions
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp b/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
index 72d3a1135..7a3f5aa9b 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
+++ b/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
@@ -5,6 +5,7 @@
#else
#include "FiniteDifferenceAlgorithms/YeeAlgorithm.H"
#include "FiniteDifferenceAlgorithms/CKCAlgorithm.H"
+ #include "FiniteDifferenceAlgorithms/NodalAlgorithm.H"
#endif
#include <AMReX_Gpu.H>
@@ -20,7 +21,9 @@ void FiniteDifferenceSolver::EvolveB ( VectorField& Bfield,
if (m_fdtd_algo == MaxwellSolverAlgo::Yee){
EvolveBCylindrical <CylindricalYeeAlgorithm> ( Bfield, Efield, dt );
#else
- if (m_fdtd_algo == MaxwellSolverAlgo::Yee){
+ if (m_do_nodal) {
+ EvolveBCartesian <NodalAlgorithm> ( Bfield, Efield, dt );
+ } else if (m_fdtd_algo == MaxwellSolverAlgo::Yee) {
EvolveBCartesian <YeeAlgorithm> ( Bfield, Efield, dt );
} else if (m_fdtd_algo == MaxwellSolverAlgo::CKC) {
EvolveBCartesian <CKCAlgorithm> ( Bfield, Efield, dt );
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/NodalAlgorithm.H b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/NodalAlgorithm.H
new file mode 100644
index 000000000..8ffd23b79
--- /dev/null
+++ b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/NodalAlgorithm.H
@@ -0,0 +1,103 @@
+#ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_NODAL_H_
+#define WARPX_FINITE_DIFFERENCE_ALGORITHM_NODAL_H_
+
+#include <AMReX_REAL.H>
+#include <AMReX_Array4.H>
+#include <AMReX_Gpu.H>
+
+struct NodalAlgorithm {
+
+ static void InitializeStencilCoefficients(
+ std::array<amrex::Real,3>& cell_size,
+ amrex::Gpu::ManagedVector<amrex::Real>& stencil_coefs_x,
+ amrex::Gpu::ManagedVector<amrex::Real>& stencil_coefs_y,
+ amrex::Gpu::ManagedVector<amrex::Real>& stencil_coefs_z ) {
+
+ // Store the inverse cell size along each direction in the coefficients
+ stencil_coefs_x.resize(1);
+ stencil_coefs_x[0] = 1./cell_size[0];
+ stencil_coefs_y.resize(1);
+ stencil_coefs_y[0] = 1./cell_size[1];
+ stencil_coefs_z.resize(1);
+ stencil_coefs_z[0] = 1./cell_size[2];
+ }
+
+ AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
+ static amrex::Real UpwardDx(
+ amrex::Array4<amrex::Real> const& F,
+ amrex::Real const* coefs_x, int const n_coefs_x,
+ int const i, int const j, int const k ) {
+
+ amrex::Real inv_dx = coefs_x[0];
+ return 0.5*inv_dx*( F(i+1,j,k) - F(i-1,j,k) );
+ };
+
+ AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
+ static amrex::Real DownwardDx(
+ amrex::Array4<amrex::Real> const& F,
+ amrex::Real const* coefs_x, int const n_coefs_x,
+ int const i, int const j, int const k ) {
+
+ amrex::Real inv_dx = coefs_x[0];
+ return 0.5*inv_dx*( F(i+1,j,k) - F(i-1,j,k) );
+ };
+
+ AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
+ static amrex::Real UpwardDy(
+ amrex::Array4<amrex::Real> const& F,
+ amrex::Real const* coefs_y, int const n_coefs_y,
+ int const i, int const j, int const k ) {
+
+ #if defined WARPX_DIM_3D
+ amrex::Real inv_dy = coefs_y[0];
+ return 0.5*inv_dy*( F(i,j+1,k) - F(i,j-1,k) );
+ #elif (defined WARPX_DIM_XZ)
+ return 0; // 2D Cartesian: derivative along y is 0
+ #endif
+ };
+
+ AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
+ static amrex::Real DownwardDy(
+ amrex::Array4<amrex::Real> const& F,
+ amrex::Real const* coefs_y, int const n_coefs_y,
+ int const i, int const j, int const k ) {
+
+ #if defined WARPX_DIM_3D
+ amrex::Real inv_dy = coefs_y[0];
+ return 0.5*inv_dy*( F(i,j+1,k) - F(i,j-1,k) );
+ #elif (defined WARPX_DIM_XZ)
+ return 0; // 2D Cartesian: derivative along y is 0
+ #endif
+ };
+
+ AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
+ static amrex::Real UpwardDz(
+ amrex::Array4<amrex::Real> const& F,
+ amrex::Real const* coefs_z, int const n_coefs_z,
+ int const i, int const j, int const k ) {
+
+ amrex::Real inv_dz = coefs_z[0];
+ #if defined WARPX_DIM_3D
+ return 0.5*inv_dz*( F(i,j,k+1) - F(i,j,k-1) );
+ #elif (defined WARPX_DIM_XZ)
+ return 0.5*inv_dz*( F(i,j+1,k) - F(i,j-1,k) );
+ #endif
+ };
+
+ AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
+ static amrex::Real DownwardDz(
+ amrex::Array4<amrex::Real> const& F,
+ amrex::Real const* coefs_z, int const n_coefs_z,
+ int const i, int const j, int const k ) {
+
+ amrex::Real inv_dz = coefs_z[0];
+ #if defined WARPX_DIM_3D
+ return 0.5*inv_dz*( F(i,j,k+1) - F(i,j,k-1) );
+ #elif (defined WARPX_DIM_XZ)
+ return 0.5*inv_dz*( F(i,j+1,k) - F(i,j-1,k) );
+ #endif
+ };
+
+};
+
+#endif // WARPX_FINITE_DIFFERENCE_ALGORITHM_NODAL_H_
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H
index 9ce910e3d..0673befca 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H
+++ b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H
@@ -1,13 +1,6 @@
#ifndef WARPX_FINITE_DIFFERENCE_SOLVER_H_
#define WARPX_FINITE_DIFFERENCE_SOLVER_H_
-#include "WarpXAlgorithmSelection.H"
-#ifdef WARPX_DIM_RZ
- #include "FiniteDifferenceAlgorithms/CylindricalYeeAlgorithm.H"
-#else
- #include "FiniteDifferenceAlgorithms/YeeAlgorithm.H"
- #include "FiniteDifferenceAlgorithms/CKCAlgorithm.H"
-#endif
#include <AMReX_MultiFab.H>
/**
@@ -22,8 +15,10 @@ class FiniteDifferenceSolver
using VectorField = std::array< std::unique_ptr<amrex::MultiFab>, 3 >;
// Constructor
- FiniteDifferenceSolver ( int const fdtd_algo,
- std::array<amrex::Real,3> cell_size );
+ FiniteDifferenceSolver (
+ int const fdtd_algo,
+ std::array<amrex::Real,3> cell_size,
+ int const do_nodal );
void EvolveB ( VectorField& Bfield,
VectorField const& Efield,
@@ -31,6 +26,7 @@ class FiniteDifferenceSolver
private:
int m_fdtd_algo;
+ int m_do_nodal;
#ifdef WARPX_DIM_RZ
amrex::Real m_dr, m_rmin;
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.cpp b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.cpp
index d89e6e9d3..94499200a 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.cpp
+++ b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.cpp
@@ -1,9 +1,19 @@
+#include "WarpXAlgorithmSelection.H"
+#ifdef WARPX_DIM_RZ
+ #include "FiniteDifferenceAlgorithms/CylindricalYeeAlgorithm.H"
+#else
+ #include "FiniteDifferenceAlgorithms/YeeAlgorithm.H"
+ #include "FiniteDifferenceAlgorithms/CKCAlgorithm.H"
+ #include "FiniteDifferenceAlgorithms/NodalAlgorithm.H"
+#endif
#include "FiniteDifferenceSolver.H"
#include "WarpX.H"
// Constructor
-FiniteDifferenceSolver::FiniteDifferenceSolver ( int const fdtd_algo,
- std::array<amrex::Real,3> cell_size ) {
+FiniteDifferenceSolver::FiniteDifferenceSolver (
+ int const fdtd_algo,
+ std::array<amrex::Real,3> cell_size,
+ int do_nodal ) {
// Register the type of finite-difference algorithm
m_fdtd_algo = fdtd_algo;
@@ -13,11 +23,14 @@ FiniteDifferenceSolver::FiniteDifferenceSolver ( int const fdtd_algo,
m_dr = cell_size[0];
m_nmodes = WarpX::GetInstance().n_rz_azimuthal_modes;
m_rmin = WarpX::GetInstance().Geom(0).ProbLo(0);
- if (fdtd_algo == MaxwellSolverAlgo::Yee){
+ if (fdtd_algo == MaxwellSolverAlgo::Yee) {
CylindricalYeeAlgorithm::InitializeStencilCoefficients( cell_size,
stencil_coefs_r, stencil_coefs_z );
#else
- if (fdtd_algo == MaxwellSolverAlgo::Yee){
+ if (do_nodal) {
+ NodalAlgorithm::InitializeStencilCoefficients( cell_size,
+ stencil_coefs_x, stencil_coefs_y, stencil_coefs_z );
+ } else if (fdtd_algo == MaxwellSolverAlgo::Yee) {
YeeAlgorithm::InitializeStencilCoefficients( cell_size,
stencil_coefs_x, stencil_coefs_y, stencil_coefs_z );
} else if (fdtd_algo == MaxwellSolverAlgo::CKC) {