aboutsummaryrefslogtreecommitdiff
path: root/Source/WarpX.cpp
diff options
context:
space:
mode:
authorGravatar Prabhat Kumar <89051199+prkkumar@users.noreply.github.com> 2021-11-19 00:31:18 -0800
committerGravatar GitHub <noreply@github.com> 2021-11-19 00:31:18 -0800
commit45ef533c39a55d05afbb9872659963d19b79f463 (patch)
tree65d7756277fd82059e8267b92a65e16d7865590f /Source/WarpX.cpp
parentff8b73f58d71762e8d81be797c8afd1519d79c2a (diff)
downloadWarpX-45ef533c39a55d05afbb9872659963d19b79f463.tar.gz
WarpX-45ef533c39a55d05afbb9872659963d19b79f463.tar.zst
WarpX-45ef533c39a55d05afbb9872659963d19b79f463.zip
1D3V Cartesian Support (#2307)
* Build System: Add 1D Geometry * test PR * test PR * 1D cartesian yee algorithm * fixed typo * Fixes for PML * 1D support related multiple changes * Fix compilation * change 1D to 1D_Z * 1D Field Gather and typo fix * 1D Charge Deposition * Particle Pusher * multiple changes related to 1D * 1D diagnostics and initialization * PlasmaInjector and PEC fixes for 1D * clean-up delete diags file * mobility 1D laser particle and bilinear filter * deleted diags files * update laser particle weight formula * delete diags files * Azure: Add 1D Cartesian Runner * 1D fixes for FieldProbe * Update Docs/source/developers/dimensionality.rst Co-authored-by: Remi Lehe <remi.lehe@normalesup.org> * 1d laser injection and langmuir test input files * 1d tests * clean up : delete print statements * analyse simulation result for laser injection and Langmuir tests * EOL * delete input files for which there are no automated tests * delete input files for which there are no automated tests * add ignore_unused to remove warnings * remove space * Fix compilation issues * fix error : macro name must be an identifier * Small bug fix * cleanup Python script for analysis * bug fix * bug fix * Update ParticleProbe: Check 1D in-domain * Update Source/Make.WarpX * Update .azure-pipelines.yml * Add USE_OPENPMD=FALSE to .azure-pipeline.yml * resolve conflict * resolve conflict * fix typo * Correct out-of-bound access * Fix Particle BC in WarpXParticleContainer and correct path to checksumAPI in python analysis scripts * EOL * Fix bug : accessing out of bound index of cell in 1D * remove 1d test for cartesian3d * Fix CI check * Slight style change * Address review comments * Fix GPU compilation Filter.cpp * Fix CI * Fix Indentation * Address review comments * More consistent ifdef for dimension bigger than 1 * Update Examples/Tests/Langmuir/analysis_langmuir_multi_1d.py Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update GNUmakefile Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Regression/prepare_file_ci.py Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianCKCAlgorithm.H Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceAlgorithms/CartesianNodalAlgorithm.H Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Source/Filter/Filter.cpp Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Source/Filter/Filter.cpp Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Source/Filter/Filter.cpp Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Source/Filter/Filter.cpp Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Source/Initialization/PlasmaInjector.cpp Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Update Source/Initialization/PlasmaInjector.cpp Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * add comment inline to explain twice push_back * Add amrex::Abort for NCIGodfreyFilter Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> Co-authored-by: Prabhat Kumar <prabhatkumar@kraken.dhcp.lbl.gov> Co-authored-by: Remi Lehe <remi.lehe@normalesup.org> Co-authored-by: Remi Lehe <rlehe@lbl.gov>
Diffstat (limited to 'Source/WarpX.cpp')
-rw-r--r--Source/WarpX.cpp40
1 files changed, 36 insertions, 4 deletions
diff --git a/Source/WarpX.cpp b/Source/WarpX.cpp
index f41de8b58..5f9106388 100644
--- a/Source/WarpX.cpp
+++ b/Source/WarpX.cpp
@@ -659,7 +659,9 @@ WarpX::ReadParameters ()
Vector<int> parse_filter_npass_each_dir(AMREX_SPACEDIM,1);
queryArrWithParser(pp_warpx, "filter_npass_each_dir", parse_filter_npass_each_dir, 0, AMREX_SPACEDIM);
filter_npass_each_dir[0] = parse_filter_npass_each_dir[0];
+#if (AMREX_SPACEDIM >= 2)
filter_npass_each_dir[1] = parse_filter_npass_each_dir[1];
+#endif
#if (AMREX_SPACEDIM == 3)
filter_npass_each_dir[2] = parse_filter_npass_each_dir[2];
#endif
@@ -1432,7 +1434,9 @@ WarpX::AllocLevelData (int lev, const BoxArray& ba, const DistributionMapping& d
bool aux_is_nodal = (field_gathering_algo == GatheringAlgo::MomentumConserving);
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+ amrex::RealVect dx({WarpX::CellSize(lev)[2]});
+#elif (AMREX_SPACEDIM == 2)
amrex::RealVect dx = {WarpX::CellSize(lev)[0], WarpX::CellSize(lev)[2]};
#elif (AMREX_SPACEDIM == 3)
amrex::RealVect dx = {WarpX::CellSize(lev)[0], WarpX::CellSize(lev)[1], WarpX::CellSize(lev)[2]};
@@ -1493,7 +1497,18 @@ WarpX::AllocLevelMFs (int lev, const BoxArray& ba, const DistributionMapping& dm
amrex::IntVect F_nodal_flag, G_nodal_flag;
// Set nodal flags
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+ // AMReX convention: x = missing dimension, y = missing dimension, z = only dimension
+ Ex_nodal_flag = IntVect(1);
+ Ey_nodal_flag = IntVect(1);
+ Ez_nodal_flag = IntVect(0);
+ Bx_nodal_flag = IntVect(0);
+ By_nodal_flag = IntVect(0);
+ Bz_nodal_flag = IntVect(1);
+ jx_nodal_flag = IntVect(1);
+ jy_nodal_flag = IntVect(1);
+ jz_nodal_flag = IntVect(0);
+#elif (AMREX_SPACEDIM == 2)
// AMReX convention: x = first dimension, y = missing dimension, z = second dimension
Ex_nodal_flag = IntVect(0,1);
Ey_nodal_flag = IntVect(1,1);
@@ -2041,7 +2056,7 @@ WarpX::CellSize (int lev)
#elif (AMREX_SPACEDIM == 2)
return { dx[0], 1.0, dx[1] };
#else
- static_assert(AMREX_SPACEDIM != 1, "1D is not supported");
+ return { 1.0, 1.0, dx[0] };
#endif
}
@@ -2065,6 +2080,9 @@ WarpX::LowerCorner(const Box& bx, std::array<amrex::Real,3> galilean_shift, int
#elif (AMREX_SPACEDIM == 2)
return { xyzmin[0] + galilean_shift[0], std::numeric_limits<Real>::lowest(), xyzmin[1] + galilean_shift[2] };
+
+#elif (AMREX_SPACEDIM == 1)
+ return { std::numeric_limits<Real>::lowest(), std::numeric_limits<Real>::lowest(), xyzmin[0] + galilean_shift[2] };
#endif
}
@@ -2077,6 +2095,8 @@ WarpX::UpperCorner(const Box& bx, int lev)
return { xyzmax[0], xyzmax[1], xyzmax[2] };
#elif (AMREX_SPACEDIM == 2)
return { xyzmax[0], std::numeric_limits<Real>::max(), xyzmax[1] };
+#elif (AMREX_SPACEDIM == 1)
+ return { std::numeric_limits<Real>::max(), std::numeric_limits<Real>::max(), xyzmax[0] };
#endif
}
@@ -2276,7 +2296,19 @@ WarpX::BuildBufferMasksInBox ( const amrex::Box tbx, amrex::IArrayBox &buffer_ma
const auto hi = amrex::ubound( tbx );
Array4<int> msk = buffer_mask.array();
Array4<int const> gmsk = guard_mask.array();
-#if (AMREX_SPACEDIM == 2)
+#if (AMREX_SPACEDIM == 1)
+ int k = lo.z;
+ int j = lo.y;
+ for (int i = lo.x; i <= hi.x; ++i) {
+ setnull = false;
+ // If gmsk=0 for any neighbor within ng cells, current cell is in the buffer region
+ for (int ii = i-ng; ii <= i+ng; ++ii) {
+ if ( gmsk(ii,j,k) == 0 ) setnull = true;
+ }
+ if ( setnull ) msk(i,j,k) = 0;
+ else msk(i,j,k) = 1;
+ }
+#elif (AMREX_SPACEDIM == 2)
int k = lo.z;
for (int j = lo.y; j <= hi.y; ++j) {
for (int i = lo.x; i <= hi.x; ++i) {