aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/FiniteDifferenceSolver
diff options
context:
space:
mode:
Diffstat (limited to 'Source/FieldSolver/FiniteDifferenceSolver')
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/ApplySilverMuellerBoundary.cpp28
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/EvolveE.cpp1
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H36
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H17
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H24
5 files changed, 93 insertions, 13 deletions
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/ApplySilverMuellerBoundary.cpp b/Source/FieldSolver/FiniteDifferenceSolver/ApplySilverMuellerBoundary.cpp
index 02b712679..46b32a8ff 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/ApplySilverMuellerBoundary.cpp
+++ b/Source/FieldSolver/FiniteDifferenceSolver/ApplySilverMuellerBoundary.cpp
@@ -165,9 +165,11 @@ void FiniteDifferenceSolver::ApplySilverMuellerBoundary (
#else
// Calculate relevant coefficients
+#if (defined WARPX_DIM_3D || WARPX_DIM_XZ)
amrex::Real const cdt_over_dx = PhysConst::c*dt*m_h_stencil_coefs_x[0];
amrex::Real const coef1_x = (1._rt - cdt_over_dx)/(1._rt + cdt_over_dx);
amrex::Real const coef2_x = 2._rt*cdt_over_dx/(1._rt + cdt_over_dx) / PhysConst::c;
+#endif
#ifdef WARPX_DIM_3D
amrex::Real const cdt_over_dy = PhysConst::c*dt*m_h_stencil_coefs_y[0];
amrex::Real const coef1_y = (1._rt - cdt_over_dy)/(1._rt + cdt_over_dy);
@@ -177,8 +179,10 @@ void FiniteDifferenceSolver::ApplySilverMuellerBoundary (
amrex::Real const coef1_z = (1._rt - cdt_over_dz)/(1._rt + cdt_over_dz);
amrex::Real const coef2_z = 2._rt*cdt_over_dz/(1._rt + cdt_over_dz) / PhysConst::c;
+#if (defined WARPX_DIM_3D || WARPX_DIM_XZ)
bool const apply_lo_x = (field_boundary_lo[0] == FieldBoundaryType::Absorbing_SilverMueller);
bool const apply_hi_x = (field_boundary_hi[0] == FieldBoundaryType::Absorbing_SilverMueller);
+#endif
#ifdef WARPX_DIM_3D
bool const apply_lo_y = (field_boundary_lo[1] == FieldBoundaryType::Absorbing_SilverMueller);
bool const apply_hi_y = (field_boundary_hi[1] == FieldBoundaryType::Absorbing_SilverMueller);
@@ -201,10 +205,14 @@ void FiniteDifferenceSolver::ApplySilverMuellerBoundary (
// Extract field data for this grid/tile
Array4<Real> const& Ex = Efield[0]->array(mfi);
Array4<Real> const& Ey = Efield[1]->array(mfi);
+#ifndef WARPX_DIM_1D_Z
Array4<Real> const& Ez = Efield[2]->array(mfi);
+#endif
Array4<Real> const& Bx = Bfield[0]->array(mfi);
Array4<Real> const& By = Bfield[1]->array(mfi);
+#ifndef WARPX_DIM_1D_Z
Array4<Real> const& Bz = Bfield[2]->array(mfi);
+#endif
// Extract the tileboxes for which to loop
Box tbx = mfi.tilebox(Bfield[0]->ixType().toIntVect());
@@ -244,18 +252,27 @@ void FiniteDifferenceSolver::ApplySilverMuellerBoundary (
// At the -z boundary (innermost guard cell)
if ( apply_lo_z && ( j==domain_box.smallEnd(1)-1 ) )
Bx(i,j,k) = coef1_z * Bx(i,j,k) + coef2_z * Ey(i,j+1,k);
+#elif WARPX_DIM_1D_Z
+ // At the +z boundary (innermost guard cell)
+ if ( apply_hi_z && ( i==domain_box.bigEnd(0)+1 ) )
+ Bx(i,j,k) = coef1_z * Bx(i,j,k) - coef2_z * Ey(i,j,k);
+ // At the -z boundary (innermost guard cell)
+ if ( apply_lo_z && ( i==domain_box.smallEnd(0)-1 ) )
+ Bx(i,j,k) = coef1_z * Bx(i,j,k) + coef2_z * Ey(i+1,j,k);
#endif
},
// Apply Boundary condition to By
[=] AMREX_GPU_DEVICE (int i, int j, int k){
+#if (defined WARPX_DIM_3D || WARPX_DIM_XZ)
// At the +x boundary (innermost guard cell)
if ( apply_hi_x && ( i==domain_box.bigEnd(0)+1 ) )
By(i,j,k) = coef1_x * By(i,j,k) - coef2_x * Ez(i,j,k);
// At the -x boundary (innermost guard cell)
if ( apply_lo_x && ( i==domain_box.smallEnd(0)-1 ) )
By(i,j,k) = coef1_x * By(i,j,k) + coef2_x * Ez(i+1,j,k);
+#endif
#ifdef WARPX_DIM_3D
// At the +z boundary (innermost guard cell)
if ( apply_hi_z && ( k==domain_box.bigEnd(2)+1 ) )
@@ -270,18 +287,27 @@ void FiniteDifferenceSolver::ApplySilverMuellerBoundary (
// At the -z boundary (innermost guard cell)
if ( apply_lo_z && ( j==domain_box.smallEnd(1)-1 ) )
By(i,j,k) = coef1_z * By(i,j,k) - coef2_z * Ex(i,j+1,k);
+#elif WARPX_DIM_1D_Z
+ // At the +z boundary (innermost guard cell)
+ if ( apply_hi_z && ( i==domain_box.bigEnd(0)+1 ) )
+ By(i,j,k) = coef1_z * By(i,j,k) + coef2_z * Ex(i,j,k);
+ // At the -z boundary (innermost guard cell)
+ if ( apply_lo_z && ( i==domain_box.smallEnd(0)-1 ) )
+ By(i,j,k) = coef1_z * By(i,j,k) - coef2_z * Ex(i+1,j,k);
#endif
},
// Apply Boundary condition to Bz
[=] AMREX_GPU_DEVICE (int i, int j, int k){
+#if (defined WARPX_DIM_3D || WARPX_DIM_XZ)
// At the +x boundary (innermost guard cell)
if ( apply_hi_x && ( i==domain_box.bigEnd(0)+1 ) )
Bz(i,j,k) = coef1_x * Bz(i,j,k) + coef2_x * Ey(i,j,k);
// At the -x boundary (innermost guard cell)
if ( apply_lo_x && ( i==domain_box.smallEnd(0)-1 ) )
Bz(i,j,k) = coef1_x * Bz(i,j,k) - coef2_x * Ey(i+1,j,k);
+#endif
#ifdef WARPX_DIM_3D
// At the +y boundary (innermost guard cell)
if ( apply_hi_y && ( j==domain_box.bigEnd(1)+1 ) )
@@ -289,6 +315,8 @@ void FiniteDifferenceSolver::ApplySilverMuellerBoundary (
// At the -y boundary (innermost guard cell)
if ( apply_lo_y && ( j==domain_box.smallEnd(1)-1 ) )
Bz(i,j,k) = coef1_y * Bz(i,j,k) + coef2_y * Ex(i,j+1,k);
+#elif WARPX_DIM_1D_Z
+ ignore_unused(i,j,k);
#endif
}
);
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/EvolveE.cpp b/Source/FieldSolver/FiniteDifferenceSolver/EvolveE.cpp
index ac19984cd..318081159 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/EvolveE.cpp
+++ b/Source/FieldSolver/FiniteDifferenceSolver/EvolveE.cpp
@@ -176,7 +176,6 @@ void FiniteDifferenceSolver::EvolveECartesian (
- T_Algo::DownwardDx(Bz, coefs_x, n_coefs_x, i, j, k)
+ T_Algo::DownwardDz(Bx, coefs_z, n_coefs_z, i, j, k)
- PhysConst::mu0 * jy(i, j, k) );
-
},
[=] AMREX_GPU_DEVICE (int i, int j, int k){
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H
index 216723553..d7946c97e 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H
+++ b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H
@@ -71,6 +71,13 @@ struct CartesianCKCAlgorithm {
constexpr Real gammax=0._rt, gammay=0._rt, gammaz=0._rt;
constexpr Real betaxy=0._rt, betazy=0._rt, betayx=0._rt, betayz=0._rt;
constexpr Real alphay=0._rt;
+#elif defined WARPX_DIM_1D_Z
+ Real const alphaz = inv_dz;
+ // Other coefficients are 0 in 1D Cartesian
+ // (and will actually not be used in the stencil)
+ constexpr Real gammax=0._rt, gammay=0._rt, gammaz=0._rt;
+ constexpr Real betaxy=0._rt, betazy=0._rt, betayx=0._rt, betayz=0._rt, betaxz=0._rt, betazx=0._rt;
+ constexpr Real alphay=0._rt, alphax=0._rt;
#endif
// Store the coefficients in array `stencil_coefs`, in prescribed order
@@ -124,14 +131,16 @@ struct CartesianCKCAlgorithm {
amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
int const i, int const j, int const k, int const ncomp=0 ) {
+ using namespace amrex;
+#if (defined WARPX_DIM_3D || WARPX_DIM_XZ)
amrex::Real const alphax = coefs_x[1];
-#if defined WARPX_DIM_3D
- amrex::Real const betaxy = coefs_x[2];
-#endif
amrex::Real const betaxz = coefs_x[3];
+#endif
#if defined WARPX_DIM_3D
+ amrex::Real const betaxy = coefs_x[2];
amrex::Real const gammax = coefs_x[4];
#endif
+
#if defined WARPX_DIM_3D
return alphax * (F(i+1,j ,k ,ncomp) - F(i, j, k ,ncomp))
+ betaxy * (F(i+1,j+1,k ,ncomp) - F(i ,j+1,k ,ncomp)
@@ -146,6 +155,9 @@ struct CartesianCKCAlgorithm {
return alphax * (F(i+1,j ,k ,ncomp) - F(i, j, k ,ncomp))
+ betaxz * (F(i+1,j+1,k ,ncomp) - F(i ,j+1,k ,ncomp)
+ F(i+1,j-1,k ,ncomp) - F(i ,j-1,k ,ncomp));
+#elif (defined WARPX_DIM_1D_Z)
+ amrex::ignore_unused(F, coefs_x, i, j, k, ncomp);
+ return 0._rt; // 1D Cartesian: derivative along x is 0
#endif
}
@@ -158,8 +170,14 @@ struct CartesianCKCAlgorithm {
amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
int const i, int const j, int const k, int const ncomp=0 ) {
+ using namespace amrex;
+#if (defined WARPX_DIM_1D_Z)
+ amrex::ignore_unused(F, coefs_x, i, j, k, ncomp);
+ return 0._rt; // 1D Cartesian: derivative along x is 0
+#else
amrex::Real const inv_dx = coefs_x[0];
return inv_dx*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
+#endif
}
/**
@@ -186,10 +204,10 @@ struct CartesianCKCAlgorithm {
+ F(i-1,j+1,k+1,ncomp) - F(i-1,j ,k+1,ncomp)
+ F(i+1,j+1,k-1,ncomp) - F(i+1,j ,k-1,ncomp)
+ F(i-1,j+1,k-1,ncomp) - F(i-1,j ,k-1,ncomp));
-#elif (defined WARPX_DIM_XZ)
+#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
amrex::ignore_unused(F, coefs_y, n_coefs_y,
i, j, k, ncomp);
- return 0._rt; // 2D Cartesian: derivative along y is 0
+ return 0._rt; // 1D and 2D Cartesian: derivative along y is 0
#else
amrex::ignore_unused(F, coefs_y, n_coefs_y,
i, j, k, ncomp);
@@ -210,7 +228,7 @@ struct CartesianCKCAlgorithm {
amrex::Real const inv_dy = coefs_y[0];
amrex::ignore_unused(n_coefs_y);
return inv_dy*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
-#elif (defined WARPX_DIM_XZ)
+#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
amrex::ignore_unused(F, coefs_y, n_coefs_y,
i, j, k, ncomp);
return 0._rt; // 2D Cartesian: derivative along y is 0
@@ -233,7 +251,9 @@ struct CartesianCKCAlgorithm {
amrex::ignore_unused(n_coefs_z);
Real const alphaz = coefs_z[1];
+#if (defined WARPX_DIM_3D || WARPX_DIM_XZ)
Real const betazx = coefs_z[2];
+#endif
#if defined WARPX_DIM_3D
Real const betazy = coefs_z[3];
Real const gammaz = coefs_z[4];
@@ -252,6 +272,8 @@ struct CartesianCKCAlgorithm {
return alphaz * (F(i ,j+1,k ,ncomp) - F(i ,j ,k ,ncomp))
+ betazx * (F(i+1,j+1,k ,ncomp) - F(i+1,j ,k ,ncomp)
+ F(i-1,j+1,k ,ncomp) - F(i-1,j ,k ,ncomp));
+#elif (defined WARPX_DIM_1D_Z)
+ return alphaz * (F(i+1 ,j,k ,ncomp) - F(i ,j ,k ,ncomp));
#endif
}
@@ -269,6 +291,8 @@ struct CartesianCKCAlgorithm {
return inv_dz*( F(i,j,k,ncomp) - F(i,j,k-1,ncomp) );
#elif (defined WARPX_DIM_XZ)
return inv_dz*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
+#elif (defined WARPX_DIM_1D_Z)
+ return inv_dz*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
#endif
}
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H
index 5c1d68687..01d83804b 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H
+++ b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H
@@ -74,8 +74,13 @@ struct CartesianNodalAlgorithm {
int const i, int const j, int const k, int const ncomp=0 ) {
using namespace amrex;
+#if (defined WARPX_DIM_1D_Z)
+ ignore_unused(i, j, k, coefs_x, ncomp, F);
+ return 0._rt; // 1D Cartesian: derivative along x is 0
+#else
Real const inv_dx = coefs_x[0];
return 0.5_rt*inv_dx*( F(i+1,j,k,ncomp) - F(i-1,j,k,ncomp) );
+#endif
}
/**
@@ -88,8 +93,14 @@ struct CartesianNodalAlgorithm {
amrex::Real const * const coefs_x, int const n_coefs_x,
int const i, int const j, int const k, int const ncomp=0 ) {
+ using namespace amrex;
+#if (defined WARPX_DIM_1D_Z)
+ ignore_unused(i, j, k, coefs_x, n_coefs_x, ncomp, F);
+ return 0._rt; // 1D Cartesian: derivative along x is 0
+#else
return UpwardDx( F, coefs_x, n_coefs_x, i, j, k ,ncomp);
// For CartesianNodalAlgorithm, UpwardDx and DownwardDx are equivalent
+#endif
}
/**
@@ -106,9 +117,9 @@ struct CartesianNodalAlgorithm {
#if defined WARPX_DIM_3D
Real const inv_dy = coefs_y[0];
return 0.5_rt*inv_dy*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
-#elif (defined WARPX_DIM_XZ)
+#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z)
ignore_unused(i, j, k, coefs_y, ncomp, F);
- return 0._rt; // 2D Cartesian: derivative along y is 0
+ return 0._rt; // 1D and 2D Cartesian: derivative along y is 0
#endif
}
@@ -142,6 +153,8 @@ struct CartesianNodalAlgorithm {
return 0.5_rt*inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k-1,ncomp) );
#elif (defined WARPX_DIM_XZ)
return 0.5_rt*inv_dz*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
+#elif (defined WARPX_DIM_1D_Z)
+ return 0.5_rt*inv_dz*( F(i+1,j,k,ncomp) - F(i-1,j,k,ncomp) );
#endif
}
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H
index 4fe813754..f2c8553d5 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H
+++ b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianYeeAlgorithm.H
@@ -71,8 +71,14 @@ struct CartesianYeeAlgorithm {
amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
int const i, int const j, int const k, int const ncomp=0 ) {
+ using namespace amrex;
+#if (defined WARPX_DIM_1D_Z)
+ amrex::ignore_unused(F, coefs_x, i, j, k, ncomp);
+ return 0._rt; // 1D Cartesian: derivative along x is 0
+#else
amrex::Real const inv_dx = coefs_x[0];
return inv_dx*( F(i+1,j,k,ncomp) - F(i,j,k,ncomp) );
+#endif
}
/**
@@ -84,8 +90,14 @@ struct CartesianYeeAlgorithm {
amrex::Real const * const coefs_x, int const /*n_coefs_x*/,
int const i, int const j, int const k, int const ncomp=0 ) {
+ using namespace amrex;
+#if (defined WARPX_DIM_1D_Z)
+ amrex::ignore_unused(F, coefs_x, i, j, k, ncomp);
+ return 0._rt; // 1D Cartesian: derivative along x is 0
+#else
amrex::Real const inv_dx = coefs_x[0];
return inv_dx*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
+#endif
}
/**
@@ -101,10 +113,10 @@ struct CartesianYeeAlgorithm {
Real const inv_dy = coefs_y[0];
amrex::ignore_unused(n_coefs_y);
return inv_dy*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
-#elif (defined WARPX_DIM_XZ)
+#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
amrex::ignore_unused(F, coefs_y, n_coefs_y,
i, j, k, ncomp);
- return 0._rt; // 2D Cartesian: derivative along y is 0
+ return 0._rt; // 1D and 2D Cartesian: derivative along y is 0
#else
amrex::ignore_unused(F, coefs_y, n_coefs_y,
i, j, k, ncomp);
@@ -125,10 +137,10 @@ struct CartesianYeeAlgorithm {
Real const inv_dy = coefs_y[0];
amrex::ignore_unused(n_coefs_y);
return inv_dy*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
-#elif (defined WARPX_DIM_XZ)
+#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
amrex::ignore_unused(F, coefs_y, n_coefs_y,
i, j, k, ncomp);
- return 0._rt; // 2D Cartesian: derivative along y is 0
+ return 0._rt; // 1D and 2D Cartesian: derivative along y is 0
#else
amrex::ignore_unused(F, coefs_y, n_coefs_y,
i, j, k, ncomp);
@@ -149,6 +161,8 @@ struct CartesianYeeAlgorithm {
return inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k,ncomp) );
#elif (defined WARPX_DIM_XZ)
return inv_dz*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
+#elif (defined WARPX_DIM_1D_Z)
+ return inv_dz*( F(i+1,j,k,ncomp) - F(i,j,k,ncomp) );
#endif
}
@@ -167,6 +181,8 @@ struct CartesianYeeAlgorithm {
return inv_dz*( F(i,j,k,ncomp) - F(i,j,k-1,ncomp) );
#elif (defined WARPX_DIM_XZ)
return inv_dz*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
+#elif (defined WARPX_DIM_1D_Z)
+ return inv_dz*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
#endif
}