aboutsummaryrefslogtreecommitdiff
path: root/Source/Parallelization
diff options
context:
space:
mode:
Diffstat (limited to 'Source/Parallelization')
-rw-r--r--Source/Parallelization/GuardCellManager.cpp5
-rw-r--r--Source/Parallelization/WarpXComm.cpp4
-rw-r--r--Source/Parallelization/WarpXComm_K.H140
3 files changed, 134 insertions, 15 deletions
diff --git a/Source/Parallelization/GuardCellManager.cpp b/Source/Parallelization/GuardCellManager.cpp
index 44848c915..b72e13cdc 100644
--- a/Source/Parallelization/GuardCellManager.cpp
+++ b/Source/Parallelization/GuardCellManager.cpp
@@ -115,6 +115,9 @@ guardCellManager::Init (
#elif (AMREX_SPACEDIM == 2)
ng_alloc_EB = IntVect(ngx,ngz);
ng_alloc_J = IntVect(ngJx,ngJz);
+#elif (AMREX_SPACEDIM == 1)
+ ng_alloc_EB = IntVect(ngz);
+ ng_alloc_J = IntVect(ngJz);
#endif
// TODO Adding one cell for rho should not be necessary, given that the number of guard cells
@@ -195,6 +198,8 @@ guardCellManager::Init (
IntVect ngFFT = IntVect(ngFFt_x, ngFFt_y, ngFFt_z);
#elif (AMREX_SPACEDIM == 2)
IntVect ngFFT = IntVect(ngFFt_x, ngFFt_z);
+#elif (AMREX_SPACEDIM == 1)
+ IntVect ngFFT = IntVect(ngFFt_z);
#endif
// All boxes should have the same number of guard cells, to avoid temporary parallel copies:
diff --git a/Source/Parallelization/WarpXComm.cpp b/Source/Parallelization/WarpXComm.cpp
index 2d1e1f04a..ecaddb570 100644
--- a/Source/Parallelization/WarpXComm.cpp
+++ b/Source/Parallelization/WarpXComm.cpp
@@ -979,7 +979,9 @@ WarpX::ApplyFilterandSumBoundaryJ (int lev, PatchType patch_type)
IntVect ng_depos_J = get_ng_depos_J();
if (WarpX::do_current_centering)
{
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+ ng_depos_J[0] += WarpX::current_centering_noz / 2;
+#elif (AMREX_SPACEDIM == 2)
ng_depos_J[0] += WarpX::current_centering_nox / 2;
ng_depos_J[1] += WarpX::current_centering_noz / 2;
#elif (AMREX_SPACEDIM == 3)
diff --git a/Source/Parallelization/WarpXComm_K.H b/Source/Parallelization/WarpXComm_K.H
index 04226e971..5995243f5 100644
--- a/Source/Parallelization/WarpXComm_K.H
+++ b/Source/Parallelization/WarpXComm_K.H
@@ -20,17 +20,17 @@ void warpx_interp (int j, int k, int l,
{
using namespace amrex;
- // NOTE Indices (j,k,l) in the following refer to (x,z,-) in 2D and (x,y,z) in 3D
+ // NOTE Indices (j,k,l) in the following refer to (z,-,-) in 1D, (x,z,-) in 2D, and (x,y,z) in 3D
// Refinement ratio
const int rj = rr[0];
- const int rk = rr[1];
- const int rl = (AMREX_SPACEDIM == 2) ? 1 : rr[2];
+ const int rk = (AMREX_SPACEDIM == 1) ? 1 : rr[1];
+ const int rl = (AMREX_SPACEDIM <= 2) ? 1 : rr[2];
// Staggering (0: cell-centered; 1: nodal)
const int sj = arr_stag[0];
- const int sk = arr_stag[1];
- const int sl = (AMREX_SPACEDIM == 2) ? 0 : arr_stag[2];
+ const int sk = (AMREX_SPACEDIM == 1) ? 0 : arr_stag[1];
+ const int sl = (AMREX_SPACEDIM <= 2) ? 0 : arr_stag[2];
// Number of points used for interpolation from coarse grid to fine grid
const int nj = (sj == 0) ? 1 : 2;
@@ -87,7 +87,21 @@ void warpx_interp_nd_bfield_x (int j, int k, int l,
Real wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
Real owy = 1.0_rt-wy;
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+
+ // interp from coarse nodal to fine nodal
+ Real bg = owx * Bxg(jg ,0,0)
+ + wx * Bxg(jg+1,0,0);
+
+ // interp from coarse staggered to fine nodal
+ Real bc = owx * Bxc(jg ,0,0)
+ + wx * Bxc(jg+1,0,0);
+
+ // interp from fine staggered to fine nodal
+ Real bf = 0.5_rt*(Bxf_zeropad(j-1,0,0) + Bxf_zeropad(j,0,0));
+ amrex::ignore_unused(owy);
+
+#elif (AMREX_SPACEDIM == 2)
// interp from coarse nodal to fine nodal
Real bg = owx * owy * Bxg(jg ,kg ,0)
@@ -163,7 +177,21 @@ void warpx_interp_nd_bfield_y (int j, int k, int l,
Real wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
Real owy = 1.0_rt-wy;
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+
+ // interp from coarse nodal to fine nodal
+ Real bg = owx * Byg(jg ,0,0)
+ + wx * Byg(jg+1,0,0);
+
+ // interp from coarse staggered to fine nodal
+ Real bc = owx * Byc(jg ,0,0)
+ + wx * Byc(jg+1,0,0);
+
+ // interp from fine staggered to fine nodal
+ Real bf = 0.5_rt*(Byf_zeropad(j-1,0,0) + Byf_zeropad(j,0,0));
+ amrex::ignore_unused(owy);
+
+#elif (AMREX_SPACEDIM == 2)
// interp from coarse nodal to fine nodal
Real bg = owx * owy * Byg(jg ,kg ,0)
@@ -241,7 +269,21 @@ void warpx_interp_nd_bfield_z (int j, int k, int l,
Real wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
Real owy = 1.0_rt-wy;
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+
+ // interp from coarse nodal to fine nodal
+ Real bg = owx * Bzg(jg ,0,0)
+ + wx * Bzg(jg+1,0,0);
+
+ // interp from coarse staggered to fine nodal
+ Real bc = owx * Bzc(jg ,0,0)
+ + wx * Bzc(jg+1,0,0);
+
+ // interp from fine staggered to fine nodal
+ Real bf = Bzf_zeropad(j,0,0);
+ amrex::ignore_unused(owy);
+
+#elif (AMREX_SPACEDIM == 2)
// interp from coarse nodal to fine nodal
Real bg = owx * owy * Bzg(jg ,kg ,0)
@@ -318,7 +360,21 @@ void warpx_interp_nd_efield_x (int j, int k, int l,
Real wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
Real owy = 1.0_rt-wy;
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+
+ // interp from coarse nodal to fine nodal
+ Real eg = owx * Exg(jg ,0,0)
+ + wx * Exg(jg+1,0,0);
+
+ // interp from coarse staggered to fine nodal
+ Real ec = owx * Exc(jg ,0,0)
+ + wx * Exc(jg+1,0,0);
+
+ // interp from fine staggered to fine nodal
+ Real ef = Exf_zeropad(j,0,0);
+ amrex::ignore_unused(owy);
+
+#elif (AMREX_SPACEDIM == 2)
// interp from coarse nodal to fine nodal
Real eg = owx * owy * Exg(jg ,kg ,0)
@@ -394,7 +450,23 @@ void warpx_interp_nd_efield_y (int j, int k, int l,
Real wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
Real owy = 1.0_rt-wy;
-#if (AMREX_SPACEDIM == 2)
+
+
+#if (AMREX_SPACEDIM == 1)
+
+ // interp from coarse nodal to fine nodal
+ Real eg = owx * Eyg(jg ,0,0)
+ + wx * Eyg(jg+1,0,0);
+
+ // interp from coarse staggered to fine nodal
+ Real ec = owx * Eyc(jg ,0,0)
+ + wx * Eyc(jg+1,0,0);
+
+ // interp from fine staggered to fine nodal
+ Real ef = Eyf_zeropad(j,0,0);
+ amrex::ignore_unused(owy);
+
+#elif (AMREX_SPACEDIM == 2)
// interp from coarse nodal to fine nodal
Real eg = owx * owy * Eyg(jg ,kg ,0)
@@ -469,7 +541,21 @@ void warpx_interp_nd_efield_z (int j, int k, int l,
Real wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
Real owy = 1.0_rt-wy;
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+
+ // interp from coarse nodal to fine nodal
+ Real eg = owx * Ezg(jg ,0,0)
+ + wx * Ezg(jg+1,0,0);
+
+ // interp from coarse staggered to fine nodal
+ Real ec = owx * Ezc(jg ,0,0)
+ + wx * Ezc(jg+1,0,0);
+
+ // interp from fine staggered to fine nodal
+ Real ef = 0.5_rt*(Ezf_zeropad(j-1,0,0) + Ezf_zeropad(j,0,0));
+ amrex::ignore_unused(owy);
+
+#elif (AMREX_SPACEDIM == 2)
// interp from coarse nodal to fine nodal
Real eg = owx * owy * Ezg(jg ,kg ,0)
@@ -568,6 +654,9 @@ void warpx_interp (const int j,
};
// Avoid compiler warnings
+#if (AMREX_SPACEDIM == 1)
+ amrex::ignore_unused(nox, noy, stencil_coeffs_x, stencil_coeffs_y);
+#endif
#if (AMREX_SPACEDIM == 2)
amrex::ignore_unused(noy, stencil_coeffs_y);
#endif
@@ -586,32 +675,43 @@ void warpx_interp (const int j,
// Staggering (s = 0 if cell-centered, s = 1 if nodal)
const int sj = (dst_nodal) ? src_stag[0] : dst_stag[0];
+#if (AMREX_SPACEDIM >= 2)
const int sk = (dst_nodal) ? src_stag[1] : dst_stag[1];
+#endif
#if (AMREX_SPACEDIM == 3)
const int sl = (dst_nodal) ? src_stag[2] : dst_stag[2];
#endif
// Interpolate along j,k,l only if source MultiFab is staggered along j,k,l
const bool interp_j = (sj == 0);
+#if (AMREX_SPACEDIM >= 2)
const bool interp_k = (sk == 0);
+#endif
#if (AMREX_SPACEDIM == 3)
const bool interp_l = (sl == 0);
#endif
+#if (AMREX_SPACEDIM == 1)
+ const int noj = noz;
+#elif (AMREX_SPACEDIM == 2)
const int noj = nox;
-#if (AMREX_SPACEDIM == 2)
const int nok = noz;
#elif (AMREX_SPACEDIM == 3)
+ const int noj = nox;
const int nok = noy;
const int nol = noz;
#endif
// Additional normalization factor
const amrex::Real wj = (interp_j) ? 0.5_rt : 1.0_rt;
+#if (AMREX_SPACEDIM == 1)
+ constexpr amrex::Real wk = 1.0_rt;
+ constexpr amrex::Real wl = 1.0_rt;
+#elif (AMREX_SPACEDIM == 2)
const amrex::Real wk = (interp_k) ? 0.5_rt : 1.0_rt;
-#if (AMREX_SPACEDIM == 2)
constexpr amrex::Real wl = 1.0_rt;
#elif (AMREX_SPACEDIM == 3)
+ const amrex::Real wk = (interp_k) ? 0.5_rt : 1.0_rt;
const amrex::Real wl = (interp_l) ? 0.5_rt : 1.0_rt;
#endif
@@ -620,11 +720,17 @@ void warpx_interp (const int j,
const int jmax = (interp_j) ? j + noj/2 + shift - 1 : j;
// Min and max for interpolation loop along k
+#if (AMREX_SPACEDIM == 1)
+ // k = 0 always
+ const int kmin = k;
+ const int kmax = k;
+#else
const int kmin = (interp_k) ? k - nok/2 + shift : k;
const int kmax = (interp_k) ? k + nok/2 + shift - 1 : k;
+#endif
// Min and max for interpolation loop along l
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM <= 2)
// l = 0 always
const int lmin = l;
const int lmax = l;
@@ -696,6 +802,11 @@ void warpx_interp (const int j,
else // PSATD (finite-order interpolation)
{
+#if (AMREX_SPACEDIM == 1)
+ amrex::Abort("PSATD not implemented in 1D");
+#endif
+
+#if (AMREX_SPACEDIM >= 2) // 1D not implemented for PSATD
amrex::Real cj = 1.0_rt;
amrex::Real ck = 1.0_rt;
amrex::Real cl = 1.0_rt;
@@ -725,6 +836,7 @@ void warpx_interp (const int j,
}
}
}
+#endif //1D
}
dst_arr(j,k,l) = wj * wk * wl * res;