aboutsummaryrefslogtreecommitdiff
path: root/Source/Particles
diff options
context:
space:
mode:
Diffstat (limited to 'Source/Particles')
-rw-r--r--Source/Particles/LaserParticleContainer.H3
-rw-r--r--Source/Particles/LaserParticleContainer.cpp8
-rw-r--r--Source/Particles/MultiParticleContainer.H2
-rw-r--r--Source/Particles/MultiParticleContainer.cpp20
-rw-r--r--Source/Particles/PhotonParticleContainer.H3
-rw-r--r--Source/Particles/PhotonParticleContainer.cpp4
-rw-r--r--Source/Particles/PhysicalParticleContainer.H3
-rw-r--r--Source/Particles/PhysicalParticleContainer.cpp10
-rw-r--r--Source/Particles/RigidInjectedParticleContainer.H3
-rw-r--r--Source/Particles/RigidInjectedParticleContainer.cpp4
-rw-r--r--Source/Particles/WarpXParticleContainer.H2
11 files changed, 34 insertions, 28 deletions
diff --git a/Source/Particles/LaserParticleContainer.H b/Source/Particles/LaserParticleContainer.H
index 98273b58a..eabd64bbd 100644
--- a/Source/Particles/LaserParticleContainer.H
+++ b/Source/Particles/LaserParticleContainer.H
@@ -46,7 +46,8 @@ public:
amrex::MultiFab* rho, amrex::MultiFab* crho,
const amrex::MultiFab*, const amrex::MultiFab*, const amrex::MultiFab*,
const amrex::MultiFab*, const amrex::MultiFab*, const amrex::MultiFab*,
- amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full) final;
+ amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full,
+ bool skip_deposition=false) final;
virtual void PushP (int lev, amrex::Real dt,
const amrex::MultiFab& ,
diff --git a/Source/Particles/LaserParticleContainer.cpp b/Source/Particles/LaserParticleContainer.cpp
index 03b623a6d..c5fe3e4fb 100644
--- a/Source/Particles/LaserParticleContainer.cpp
+++ b/Source/Particles/LaserParticleContainer.cpp
@@ -418,7 +418,7 @@ LaserParticleContainer::Evolve (int lev,
MultiFab* rho, MultiFab* crho,
const MultiFab*, const MultiFab*, const MultiFab*,
const MultiFab*, const MultiFab*, const MultiFab*,
- Real t, Real dt, DtType /*a_dt_type*/)
+ Real t, Real dt, DtType /*a_dt_type*/, bool skip_deposition)
{
WARPX_PROFILE("LaserParticleContainer::Evolve()");
WARPX_PROFILE_VAR_NS("LaserParticleContainer::Evolve::ParticlePush", blp_pp);
@@ -475,7 +475,7 @@ LaserParticleContainer::Evolve (int lev,
plane_Yp.resize(np);
amplitude_E.resize(np);
- if (rho) {
+ if (rho && ! skip_deposition) {
int* AMREX_RESTRICT ion_lev = nullptr;
DepositCharge(pti, wp, ion_lev, rho, 0, 0,
np_current, thread_num, lev, lev);
@@ -510,7 +510,7 @@ LaserParticleContainer::Evolve (int lev,
// Current Deposition
//
// Deposit inside domains
- {
+ if (! skip_deposition ) {
int* ion_lev = nullptr;
DepositCurrent(pti, wp, uxp, uyp, uzp, ion_lev, &jx, &jy, &jz,
0, np_current, thread_num,
@@ -525,7 +525,7 @@ LaserParticleContainer::Evolve (int lev,
}
}
- if (rho) {
+ if (rho && ! skip_deposition) {
int* AMREX_RESTRICT ion_lev = nullptr;
DepositCharge(pti, wp, ion_lev, rho, 1, 0,
np_current, thread_num, lev, lev);
diff --git a/Source/Particles/MultiParticleContainer.H b/Source/Particles/MultiParticleContainer.H
index 58bfd0e79..d6a7967b5 100644
--- a/Source/Particles/MultiParticleContainer.H
+++ b/Source/Particles/MultiParticleContainer.H
@@ -100,7 +100,7 @@ public:
amrex::MultiFab* rho, amrex::MultiFab* crho,
const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
- amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full);
+ amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false);
///
/// This pushes the particle positions by one half time step for all the species in the
diff --git a/Source/Particles/MultiParticleContainer.cpp b/Source/Particles/MultiParticleContainer.cpp
index 4218e9fdf..34d12a2b3 100644
--- a/Source/Particles/MultiParticleContainer.cpp
+++ b/Source/Particles/MultiParticleContainer.cpp
@@ -328,16 +328,18 @@ MultiParticleContainer::Evolve (int lev,
MultiFab* rho, MultiFab* crho,
const MultiFab* cEx, const MultiFab* cEy, const MultiFab* cEz,
const MultiFab* cBx, const MultiFab* cBy, const MultiFab* cBz,
- Real t, Real dt, DtType a_dt_type)
+ Real t, Real dt, DtType a_dt_type, bool skip_deposition)
{
- jx.setVal(0.0);
- jy.setVal(0.0);
- jz.setVal(0.0);
- if (cjx) cjx->setVal(0.0);
- if (cjy) cjy->setVal(0.0);
- if (cjz) cjz->setVal(0.0);
- if (rho) rho->setVal(0.0);
- if (crho) crho->setVal(0.0);
+ if (! skip_deposition) {
+ jx.setVal(0.0);
+ jy.setVal(0.0);
+ jz.setVal(0.0);
+ if (cjx) cjx->setVal(0.0);
+ if (cjy) cjy->setVal(0.0);
+ if (cjz) cjz->setVal(0.0);
+ if (rho) rho->setVal(0.0);
+ if (crho) crho->setVal(0.0);
+ }
for (auto& pc : allcontainers) {
pc->Evolve(lev, Ex, Ey, Ez, Bx, By, Bz, jx, jy, jz, cjx, cjy, cjz,
rho, crho, cEx, cEy, cEz, cBx, cBy, cBz, t, dt, a_dt_type);
diff --git a/Source/Particles/PhotonParticleContainer.H b/Source/Particles/PhotonParticleContainer.H
index 4d5ec590a..0cfb37729 100644
--- a/Source/Particles/PhotonParticleContainer.H
+++ b/Source/Particles/PhotonParticleContainer.H
@@ -52,7 +52,8 @@ public:
const amrex::MultiFab* cBz,
amrex::Real t,
amrex::Real dt,
- DtType a_dt_type=DtType::Full) override;
+ DtType a_dt_type=DtType::Full,
+ bool skip_deposition=false) override;
virtual void PushPX(WarpXParIter& pti,
amrex::FArrayBox const * exfab,
diff --git a/Source/Particles/PhotonParticleContainer.cpp b/Source/Particles/PhotonParticleContainer.cpp
index e8ba363a3..bda0f984e 100644
--- a/Source/Particles/PhotonParticleContainer.cpp
+++ b/Source/Particles/PhotonParticleContainer.cpp
@@ -198,7 +198,7 @@ PhotonParticleContainer::Evolve (int lev,
MultiFab* rho, MultiFab* crho,
const MultiFab* cEx, const MultiFab* cEy, const MultiFab* cEz,
const MultiFab* cBx, const MultiFab* cBy, const MultiFab* cBz,
- Real t, Real dt, DtType /*a_dt_type*/)
+ Real t, Real dt, DtType a_dt_type, bool skip_deposition)
{
// This does gather, push and depose.
// Push and depose have been re-written for photons
@@ -210,6 +210,6 @@ PhotonParticleContainer::Evolve (int lev,
rho, crho,
cEx, cEy, cEz,
cBx, cBy, cBz,
- t, dt);
+ t, dt, a_dt_type, skip_deposition);
}
diff --git a/Source/Particles/PhysicalParticleContainer.H b/Source/Particles/PhysicalParticleContainer.H
index 250f3fb02..000caf790 100644
--- a/Source/Particles/PhysicalParticleContainer.H
+++ b/Source/Particles/PhysicalParticleContainer.H
@@ -120,7 +120,8 @@ public:
const amrex::MultiFab* cBz,
amrex::Real t,
amrex::Real dt,
- DtType a_dt_type=DtType::Full) override;
+ DtType a_dt_type=DtType::Full,
+ bool skip_deposition=false ) override;
virtual void PushPX (WarpXParIter& pti,
amrex::FArrayBox const * exfab,
diff --git a/Source/Particles/PhysicalParticleContainer.cpp b/Source/Particles/PhysicalParticleContainer.cpp
index 6fbe442e3..681559fbc 100644
--- a/Source/Particles/PhysicalParticleContainer.cpp
+++ b/Source/Particles/PhysicalParticleContainer.cpp
@@ -953,7 +953,7 @@ PhysicalParticleContainer::Evolve (int lev,
MultiFab* rho, MultiFab* crho,
const MultiFab* cEx, const MultiFab* cEy, const MultiFab* cEz,
const MultiFab* cBx, const MultiFab* cBy, const MultiFab* cBz,
- Real /*t*/, Real dt, DtType a_dt_type)
+ Real /*t*/, Real dt, DtType a_dt_type, bool skip_deposition)
{
WARPX_PROFILE("PhysicalParticleContainer::Evolve()");
@@ -1055,7 +1055,7 @@ PhysicalParticleContainer::Evolve (int lev,
const long np_current = (cjx) ? nfine_current : np;
- if (rho) {
+ if (rho && ! skip_deposition) {
// Deposit charge before particle push, in component 0 of MultiFab rho.
int* AMREX_RESTRICT ion_lev;
if (do_field_ionization){
@@ -1125,9 +1125,9 @@ PhysicalParticleContainer::Evolve (int lev,
WARPX_PROFILE_VAR_STOP(blp_fg);
//
- // Current Deposition (only needed for electromagnetic solver)
+ // Current Deposition
//
- if (WarpX::do_electrostatic == ElectrostaticSolverAlgo::None) {
+ if (! skip_deposition) {
int* AMREX_RESTRICT ion_lev;
if (do_field_ionization){
ion_lev = pti.GetiAttribs(particle_icomps["ionization_level"]).dataPtr();
@@ -1147,7 +1147,7 @@ PhysicalParticleContainer::Evolve (int lev,
} // end of "if do_electrostatic == ElectrostaticSolverAlgo::None"
} // end of "if do_not_push"
- if (rho) {
+ if (rho && ! skip_deposition) {
// Deposit charge after particle push, in component 1 of MultiFab rho.
// (Skipped for electrostatic solver, as this may lead to out-of-bounds)
if (WarpX::do_electrostatic == ElectrostaticSolverAlgo::None) {
diff --git a/Source/Particles/RigidInjectedParticleContainer.H b/Source/Particles/RigidInjectedParticleContainer.H
index 70ac7677e..1f52c0861 100644
--- a/Source/Particles/RigidInjectedParticleContainer.H
+++ b/Source/Particles/RigidInjectedParticleContainer.H
@@ -66,7 +66,8 @@ public:
const amrex::MultiFab* cBz,
amrex::Real t,
amrex::Real dt,
- DtType a_dt_type=DtType::Full) override;
+ DtType a_dt_type=DtType::Full,
+ bool skip_deposition=false ) override;
virtual void PushPX (WarpXParIter& pti,
amrex::FArrayBox const * exfab,
diff --git a/Source/Particles/RigidInjectedParticleContainer.cpp b/Source/Particles/RigidInjectedParticleContainer.cpp
index dc2f88239..5cb49f69c 100644
--- a/Source/Particles/RigidInjectedParticleContainer.cpp
+++ b/Source/Particles/RigidInjectedParticleContainer.cpp
@@ -347,7 +347,7 @@ RigidInjectedParticleContainer::Evolve (int lev,
MultiFab* rho, MultiFab* crho,
const MultiFab* cEx, const MultiFab* cEy, const MultiFab* cEz,
const MultiFab* cBx, const MultiFab* cBy, const MultiFab* cBz,
- Real t, Real dt, DtType a_dt_type)
+ Real t, Real dt, DtType a_dt_type, bool skip_deposition)
{
// Update location of injection plane in the boosted frame
@@ -374,7 +374,7 @@ RigidInjectedParticleContainer::Evolve (int lev,
rho, crho,
cEx, cEy, cEz,
cBx, cBy, cBz,
- t, dt, a_dt_type);
+ t, dt, a_dt_type, skip_deposition);
}
void
diff --git a/Source/Particles/WarpXParticleContainer.H b/Source/Particles/WarpXParticleContainer.H
index 91d5e64ee..07a9ddcd6 100644
--- a/Source/Particles/WarpXParticleContainer.H
+++ b/Source/Particles/WarpXParticleContainer.H
@@ -164,7 +164,7 @@ public:
amrex::MultiFab* rho, amrex::MultiFab* crho,
const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
- amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full) = 0;
+ amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false) = 0;
virtual void PostRestart () = 0;