aboutsummaryrefslogtreecommitdiff
path: root/Source/Parallelization/WarpXComm.cpp
diff options
context:
space:
mode:
authorGravatar Andrew Myers <atmyers2@gmail.com> 2020-01-31 10:13:04 -0800
committerGravatar Andrew Myers <atmyers2@gmail.com> 2020-01-31 10:13:04 -0800
commit58e64aafe0103b6644048d7480a3e62fe01bd5cb (patch)
tree8f0d7ed234d780b929d493d7a217fa504647fa8c /Source/Parallelization/WarpXComm.cpp
parente45710b641c01970a1cc9eb2eae244899091106b (diff)
parent3f5bcb4a798862e0a5aa604e5dce162bb0e291b3 (diff)
downloadWarpX-58e64aafe0103b6644048d7480a3e62fe01bd5cb.tar.gz
WarpX-58e64aafe0103b6644048d7480a3e62fe01bd5cb.tar.zst
WarpX-58e64aafe0103b6644048d7480a3e62fe01bd5cb.zip
Merge branch 'dev' into elementary_process
Diffstat (limited to 'Source/Parallelization/WarpXComm.cpp')
-rw-r--r--Source/Parallelization/WarpXComm.cpp107
1 files changed, 79 insertions, 28 deletions
diff --git a/Source/Parallelization/WarpXComm.cpp b/Source/Parallelization/WarpXComm.cpp
index b61ae4fc7..31bb802f7 100644
--- a/Source/Parallelization/WarpXComm.cpp
+++ b/Source/Parallelization/WarpXComm.cpp
@@ -1,3 +1,11 @@
+/* Copyright 2019 Andrew Myers, Aurore Blelly, Axel Huebl
+ * David Grote, Maxence Thevenet, Remi Lehe
+ * Revathi Jambunathan, Weiqun Zhang
+ *
+ * This file is part of WarpX.
+ *
+ * License: BSD-3-Clause-LBNL
+ */
#include <WarpXComm.H>
#include <WarpXComm_K.H>
#include <WarpX.H>
@@ -321,41 +329,41 @@ WarpX::UpdateAuxilaryDataSameType ()
}
void
-WarpX::FillBoundaryB ()
+WarpX::FillBoundaryB (IntVect ng, IntVect ng_extra_fine)
{
for (int lev = 0; lev <= finest_level; ++lev)
{
- FillBoundaryB(lev);
+ FillBoundaryB(lev, ng, ng_extra_fine);
}
}
void
-WarpX::FillBoundaryE ()
+WarpX::FillBoundaryE (IntVect ng, IntVect ng_extra_fine)
{
for (int lev = 0; lev <= finest_level; ++lev)
{
- FillBoundaryE(lev);
+ FillBoundaryE(lev, ng, ng_extra_fine);
}
}
void
-WarpX::FillBoundaryF ()
+WarpX::FillBoundaryF (IntVect ng)
{
for (int lev = 0; lev <= finest_level; ++lev)
{
- FillBoundaryF(lev);
+ FillBoundaryF(lev, ng);
}
}
void
-WarpX::FillBoundaryE(int lev)
+WarpX::FillBoundaryE(int lev, IntVect ng, IntVect ng_extra_fine)
{
- FillBoundaryE(lev, PatchType::fine);
- if (lev > 0) FillBoundaryE(lev, PatchType::coarse);
+ FillBoundaryE(lev, PatchType::fine, ng+ng_extra_fine);
+ if (lev > 0) FillBoundaryE(lev, PatchType::coarse, ng);
}
void
-WarpX::FillBoundaryE (int lev, PatchType patch_type)
+WarpX::FillBoundaryE (int lev, PatchType patch_type, IntVect ng)
{
if (patch_type == PatchType::fine)
{
@@ -370,8 +378,12 @@ WarpX::FillBoundaryE (int lev, PatchType patch_type)
}
const auto& period = Geom(lev).periodicity();
- Vector<MultiFab*> mf{Efield_fp[lev][0].get(),Efield_fp[lev][1].get(),Efield_fp[lev][2].get()};
- amrex::FillBoundary(mf, period);
+ AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
+ ng <= Efield_fp[lev][0]->nGrowVect(),
+ "Error: in FillBoundaryE, requested more guard cells than allocated");
+ Efield_fp[lev][0]->FillBoundary(ng, period);
+ Efield_fp[lev][1]->FillBoundary(ng, period);
+ Efield_fp[lev][2]->FillBoundary(ng, period);
}
else if (patch_type == PatchType::coarse)
{
@@ -386,20 +398,24 @@ WarpX::FillBoundaryE (int lev, PatchType patch_type)
}
const auto& cperiod = Geom(lev-1).periodicity();
- Vector<MultiFab*> mf{Efield_cp[lev][0].get(),Efield_cp[lev][1].get(),Efield_cp[lev][2].get()};
- amrex::FillBoundary(mf, cperiod);
+ AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
+ ng <= Efield_cp[lev][0]->nGrowVect(),
+ "Error: in FillBoundaryE, requested more guard cells than allocated");
+ Efield_cp[lev][0]->FillBoundary(ng, cperiod);
+ Efield_cp[lev][1]->FillBoundary(ng, cperiod);
+ Efield_cp[lev][2]->FillBoundary(ng, cperiod);
}
}
void
-WarpX::FillBoundaryB (int lev)
+WarpX::FillBoundaryB (int lev, IntVect ng, IntVect ng_extra_fine)
{
- FillBoundaryB(lev, PatchType::fine);
- if (lev > 0) FillBoundaryB(lev, PatchType::coarse);
+ FillBoundaryB(lev, PatchType::fine, ng + ng_extra_fine);
+ if (lev > 0) FillBoundaryB(lev, PatchType::coarse, ng);
}
void
-WarpX::FillBoundaryB (int lev, PatchType patch_type)
+WarpX::FillBoundaryB (int lev, PatchType patch_type, IntVect ng)
{
if (patch_type == PatchType::fine)
{
@@ -413,8 +429,12 @@ WarpX::FillBoundaryB (int lev, PatchType patch_type)
pml[lev]->FillBoundaryB(patch_type);
}
const auto& period = Geom(lev).periodicity();
- Vector<MultiFab*> mf{Bfield_fp[lev][0].get(),Bfield_fp[lev][1].get(),Bfield_fp[lev][2].get()};
- amrex::FillBoundary(mf, period);
+ AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
+ ng <= Bfield_fp[lev][0]->nGrowVect(),
+ "Error: in FillBoundaryB, requested more guard cells than allocated");
+ Bfield_fp[lev][0]->FillBoundary(ng, period);
+ Bfield_fp[lev][1]->FillBoundary(ng, period);
+ Bfield_fp[lev][2]->FillBoundary(ng, period);
}
else if (patch_type == PatchType::coarse)
{
@@ -428,20 +448,24 @@ WarpX::FillBoundaryB (int lev, PatchType patch_type)
pml[lev]->FillBoundaryB(patch_type);
}
const auto& cperiod = Geom(lev-1).periodicity();
- Vector<MultiFab*> mf{Bfield_cp[lev][0].get(),Bfield_cp[lev][1].get(),Bfield_cp[lev][2].get()};
- amrex::FillBoundary(mf, cperiod);
+ AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
+ ng <= Bfield_cp[lev][0]->nGrowVect(),
+ "Error: in FillBoundaryB, requested more guard cells than allocated");
+ Bfield_cp[lev][0]->FillBoundary(ng, cperiod);
+ Bfield_cp[lev][1]->FillBoundary(ng, cperiod);
+ Bfield_cp[lev][2]->FillBoundary(ng, cperiod);
}
}
void
-WarpX::FillBoundaryF (int lev)
+WarpX::FillBoundaryF (int lev, IntVect ng)
{
- FillBoundaryF(lev, PatchType::fine);
- if (lev > 0) FillBoundaryF(lev, PatchType::coarse);
+ FillBoundaryF(lev, PatchType::fine, ng);
+ if (lev > 0) FillBoundaryF(lev, PatchType::coarse, ng);
}
void
-WarpX::FillBoundaryF (int lev, PatchType patch_type)
+WarpX::FillBoundaryF (int lev, PatchType patch_type, IntVect ng)
{
if (patch_type == PatchType::fine && F_fp[lev])
{
@@ -453,7 +477,10 @@ WarpX::FillBoundaryF (int lev, PatchType patch_type)
}
const auto& period = Geom(lev).periodicity();
- F_fp[lev]->FillBoundary(period);
+ AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
+ ng <= F_fp[lev]->nGrowVect(),
+ "Error: in FillBoundaryF, requested more guard cells than allocated");
+ F_fp[lev]->FillBoundary(ng, period);
}
else if (patch_type == PatchType::coarse && F_cp[lev])
{
@@ -465,11 +492,35 @@ WarpX::FillBoundaryF (int lev, PatchType patch_type)
}
const auto& cperiod = Geom(lev-1).periodicity();
- F_cp[lev]->FillBoundary(cperiod);
+ AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
+ ng <= F_cp[lev]->nGrowVect(),
+ "Error: in FillBoundaryF, requested more guard cells than allocated");
+ F_cp[lev]->FillBoundary(ng, cperiod);
}
}
void
+WarpX::FillBoundaryAux (IntVect ng)
+{
+ for (int lev = 0; lev <= finest_level-1; ++lev)
+ {
+ FillBoundaryAux(lev, ng);
+ }
+}
+
+void
+WarpX::FillBoundaryAux (int lev, IntVect ng)
+{
+ const auto& period = Geom(lev).periodicity();
+ Efield_aux[lev][0]->FillBoundary(ng, period);
+ Efield_aux[lev][1]->FillBoundary(ng, period);
+ Efield_aux[lev][2]->FillBoundary(ng, period);
+ Bfield_aux[lev][0]->FillBoundary(ng, period);
+ Bfield_aux[lev][1]->FillBoundary(ng, period);
+ Bfield_aux[lev][2]->FillBoundary(ng, period);
+}
+
+void
WarpX::SyncCurrent ()
{
BL_PROFILE("SyncCurrent()");