aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/ElectrostaticSolver.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/FieldSolver/ElectrostaticSolver.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/FieldSolver/ElectrostaticSolver.cpp')
-rw-r--r--Source/FieldSolver/ElectrostaticSolver.cpp54
1 files changed, 46 insertions, 8 deletions
diff --git a/Source/FieldSolver/ElectrostaticSolver.cpp b/Source/FieldSolver/ElectrostaticSolver.cpp
index 5f33bbfb4..582890aa8 100644
--- a/Source/FieldSolver/ElectrostaticSolver.cpp
+++ b/Source/FieldSolver/ElectrostaticSolver.cpp
@@ -400,7 +400,9 @@ WarpX::computePhiCartesian (const amrex::Vector<std::unique_ptr<amrex::MultiFab>
// Set the value of beta
amrex::Array<amrex::Real,AMREX_SPACEDIM> beta_solver =
-# if (AMREX_SPACEDIM==2)
+# if (AMREX_SPACEDIM==1)
+ {{ beta[2] }}; // beta_x and beta_z
+# elif (AMREX_SPACEDIM==2)
{{ beta[0], beta[2] }}; // beta_x and beta_z
# else
{{ beta[0], beta[1], beta[2] }};
@@ -463,7 +465,13 @@ WarpX::computePhiCartesian (const amrex::Vector<std::unique_ptr<amrex::MultiFab>
if (do_electrostatic == ElectrostaticSolverAlgo::LabFrame)
{
for (int lev = 0; lev <= max_level; ++lev) {
-#if (AMREX_SPACEDIM==2)
+#if (AMREX_SPACEDIM==1)
+ mlmg.getGradSolution(
+ {amrex::Array<amrex::MultiFab*,1>{
+ get_pointer_Efield_fp(lev, 2)
+ }}
+ );
+#elif (AMREX_SPACEDIM==2)
mlmg.getGradSolution(
{amrex::Array<amrex::MultiFab*,2>{
get_pointer_Efield_fp(lev, 0),get_pointer_Efield_fp(lev, 2)
@@ -578,21 +586,28 @@ WarpX::computeE (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >
#endif
for ( MFIter mfi(*phi[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
{
- const Real inv_dx = 1._rt/dx[0];
#if (AMREX_SPACEDIM == 3)
+ const Real inv_dx = 1._rt/dx[0];
const Real inv_dy = 1._rt/dx[1];
const Real inv_dz = 1._rt/dx[2];
-#else
+#elif (AMREX_SPACEDIM == 2)
+ const Real inv_dx = 1._rt/dx[0];
const Real inv_dz = 1._rt/dx[1];
+#else
+ const Real inv_dz = 1._rt/dx[0];
#endif
+#if (AMREX_SPACEDIM >= 2)
const Box& tbx = mfi.tilebox( E[lev][0]->ixType().toIntVect() );
+#endif
#if (AMREX_SPACEDIM == 3)
const Box& tby = mfi.tilebox( E[lev][1]->ixType().toIntVect() );
#endif
const Box& tbz = mfi.tilebox( E[lev][2]->ixType().toIntVect() );
const auto& phi_arr = phi[lev]->array(mfi);
+#if (AMREX_SPACEDIM >= 2)
const auto& Ex_arr = (*E[lev][0])[mfi].array();
+#endif
#if (AMREX_SPACEDIM == 3)
const auto& Ey_arr = (*E[lev][1])[mfi].array();
#endif
@@ -631,7 +646,7 @@ WarpX::computeE (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >
+(beta_y*beta_z-1)*inv_dz*( phi_arr(i,j,k+1)-phi_arr(i,j,k) );
}
);
-#else
+#elif (AMREX_SPACEDIM == 2)
amrex::ParallelFor( tbx, tbz,
[=] AMREX_GPU_DEVICE (int i, int j, int k) {
Ex_arr(i,j,k) +=
@@ -646,6 +661,14 @@ WarpX::computeE (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >
+(beta_y*beta_z-1)*inv_dz*( phi_arr(i,j+1,k)-phi_arr(i,j,k) );
}
);
+#else
+ amrex::ParallelFor( tbz,
+ [=] AMREX_GPU_DEVICE (int i, int j, int k) {
+ Ez_arr(i,j,k) +=
+ +(beta_y*beta_z-1)*inv_dz*( phi_arr(i+1,j,k)-phi_arr(i,j,k) );
+ }
+ );
+ ignore_unused(beta_x);
#endif
}
}
@@ -680,12 +703,15 @@ WarpX::computeB (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >
#endif
for ( MFIter mfi(*phi[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
{
- const Real inv_dx = 1._rt/dx[0];
#if (AMREX_SPACEDIM == 3)
+ const Real inv_dx = 1._rt/dx[0];
const Real inv_dy = 1._rt/dx[1];
const Real inv_dz = 1._rt/dx[2];
-#else
+#elif (AMREX_SPACEDIM == 2)
+ const Real inv_dx = 1._rt/dx[0];
const Real inv_dz = 1._rt/dx[1];
+#else
+ const Real inv_dz = 1._rt/dx[0];
#endif
const Box& tbx = mfi.tilebox( B[lev][0]->ixType().toIntVect() );
const Box& tby = mfi.tilebox( B[lev][1]->ixType().toIntVect() );
@@ -728,7 +754,7 @@ WarpX::computeB (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >
+ phi_arr(i+1,j+1,k)-phi_arr(i,j+1,k)));
}
);
-#else
+#elif (AMREX_SPACEDIM == 2)
amrex::ParallelFor( tbx, tby, tbz,
[=] AMREX_GPU_DEVICE (int i, int j, int k) {
Bx_arr(i,j,k) += inv_c * (
@@ -746,6 +772,18 @@ WarpX::computeB (amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3> >
+beta_y*inv_dx*( phi_arr(i+1,j,k)-phi_arr(i,j,k) ));
}
);
+#else
+ amrex::ParallelFor( tbx, tby,
+ [=] AMREX_GPU_DEVICE (int i, int j, int k) {
+ Bx_arr(i,j,k) += inv_c * (
+ -beta_y*inv_dz*( phi_arr(i+1,j,k)-phi_arr(i,j,k) ));
+ },
+ [=] AMREX_GPU_DEVICE (int i, int j, int k) {
+ By_arr(i,j,k) += inv_c * (
+ +beta_x*inv_dz*(phi_arr(i+1,j,k)-phi_arr(i,j,k)));
+ }
+ );
+ ignore_unused(beta_z,tbz,Bz_arr);
#endif
}
}