aboutsummaryrefslogtreecommitdiff
path: root/Source
diff options
context:
space:
mode:
Diffstat (limited to 'Source')
-rw-r--r--Source/Evolve/WarpXEvolve.cpp2
-rw-r--r--Source/Particles/Collision/CollisionType.H4
-rw-r--r--Source/Particles/Collision/CollisionType.cpp13
-rw-r--r--Source/Particles/MultiParticleContainer.H2
-rw-r--r--Source/Particles/MultiParticleContainer.cpp9
5 files changed, 20 insertions, 10 deletions
diff --git a/Source/Evolve/WarpXEvolve.cpp b/Source/Evolve/WarpXEvolve.cpp
index 78f4376eb..3dda534f1 100644
--- a/Source/Evolve/WarpXEvolve.cpp
+++ b/Source/Evolve/WarpXEvolve.cpp
@@ -272,7 +272,7 @@ WarpX::OneStep_nosub (Real cur_time)
// product species.
doFieldIonization();
- mypc->doCoulombCollisions();
+ mypc->doCoulombCollisions( cur_time );
#ifdef WARPX_QED
mypc->doQEDSchwinger();
#endif
diff --git a/Source/Particles/Collision/CollisionType.H b/Source/Particles/Collision/CollisionType.H
index 962530973..f1ffbe439 100644
--- a/Source/Particles/Collision/CollisionType.H
+++ b/Source/Particles/Collision/CollisionType.H
@@ -17,6 +17,7 @@ class CollisionType
public:
int m_species1_index;
int m_species2_index;
+ int m_ndt;
bool m_isSameSpecies;
amrex::Real m_CoulombLog;
@@ -38,7 +39,8 @@ public:
int const lev, amrex::MFIter const& mfi,
std::unique_ptr<WarpXParticleContainer>& species1,
std::unique_ptr<WarpXParticleContainer>& species2,
- bool const isSameSpecies, amrex::Real const CoulombLog );
+ bool const isSameSpecies, amrex::Real const CoulombLog,
+ int const ndt );
};
diff --git a/Source/Particles/Collision/CollisionType.cpp b/Source/Particles/Collision/CollisionType.cpp
index 79962b803..36ce7952d 100644
--- a/Source/Particles/Collision/CollisionType.cpp
+++ b/Source/Particles/Collision/CollisionType.cpp
@@ -26,6 +26,10 @@ CollisionType::CollisionType(
m_CoulombLog = -1.0;
pp.query("CoulombLog", m_CoulombLog);
+ // number of time steps between collisions
+ m_ndt = 1;
+ pp.query("ndt", m_ndt);
+
for (int i=0; i<static_cast<int>(species_names.size()); i++)
{
if (species_names[i] == collision_species[0])
@@ -57,13 +61,14 @@ using namespace ParticleUtils;
* @param species1/2 pointer to species container
* @param isSameSpecies true if collision is between same species
* @param CoulombLog user input Coulomb logrithm
+ * @param ndt user input number of time stpes between collisions
*
*/
void CollisionType::doCoulombCollisionsWithinTile
( int const lev, MFIter const& mfi,
std::unique_ptr<WarpXParticleContainer>& species_1,
std::unique_ptr<WarpXParticleContainer>& species_2,
- bool const isSameSpecies, Real const CoulombLog )
+ bool const isSameSpecies, Real const CoulombLog, int const ndt )
{
if ( isSameSpecies ) // species_1 == species_2
@@ -139,8 +144,7 @@ void CollisionType::doCoulombCollisionsWithinTile
indices_1, indices_1,
ux_1, uy_1, uz_1, ux_1, uy_1, uz_1, w_1, w_1,
q1, q1, m1, m1, Real(-1.0), Real(-1.0),
- dt, CoulombLog, dV,
- engine);
+ dt*ndt, CoulombLog, dV, engine );
}
}
);
@@ -236,8 +240,7 @@ void CollisionType::doCoulombCollisionsWithinTile
indices_1, indices_2,
ux_1, uy_1, uz_1, ux_2, uy_2, uz_2, w_1, w_2,
q1, q2, m1, m2, Real(-1.0), Real(-1.0),
- dt, CoulombLog, dV,
- engine);
+ dt*ndt, CoulombLog, dV, engine );
}
}
);
diff --git a/Source/Particles/MultiParticleContainer.H b/Source/Particles/MultiParticleContainer.H
index 73a8207bd..4a515acf8 100644
--- a/Source/Particles/MultiParticleContainer.H
+++ b/Source/Particles/MultiParticleContainer.H
@@ -122,7 +122,7 @@ public:
const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
- void doCoulombCollisions ();
+ void doCoulombCollisions (amrex::Real cur_time);
/**
* \brief This function loops over all species and performs resampling if appropriate.
diff --git a/Source/Particles/MultiParticleContainer.cpp b/Source/Particles/MultiParticleContainer.cpp
index d1dde1837..f986e7dba 100644
--- a/Source/Particles/MultiParticleContainer.cpp
+++ b/Source/Particles/MultiParticleContainer.cpp
@@ -683,12 +683,16 @@ MultiParticleContainer::doFieldIonization (int lev,
}
void
-MultiParticleContainer::doCoulombCollisions ()
+MultiParticleContainer::doCoulombCollisions ( Real cur_time )
{
WARPX_PROFILE("MultiParticleContainer::doCoulombCollisions()");
for( auto const& collision : allcollisions )
{
+
+ const Real dt = WarpX::GetInstance().getdt(0);
+ if ( int(std::floor(cur_time/dt)) % collision->m_ndt != 0 ) continue;
+
auto& species1 = allcontainers[ collision->m_species1_index ];
auto& species2 = allcontainers[ collision->m_species2_index ];
@@ -709,7 +713,8 @@ MultiParticleContainer::doCoulombCollisions ()
CollisionType::doCoulombCollisionsWithinTile
( lev, mfi, species1, species2,
collision->m_isSameSpecies,
- collision->m_CoulombLog );
+ collision->m_CoulombLog,
+ collision->m_ndt );
}
}