diff options
author | 2021-01-03 14:10:24 -0800 | |
---|---|---|
committer | 2021-01-03 14:10:24 -0800 | |
commit | 2fe64fb007abd3b55b320ba1acf2f285345248ae (patch) | |
tree | c9bac6b32c296889717598445917e8704b8c485e /Source/Particles/MultiParticleContainer.cpp | |
parent | 4a62e7c9c62f3ffee3f61902139fd302fefba521 (diff) | |
download | WarpX-2fe64fb007abd3b55b320ba1acf2f285345248ae.tar.gz WarpX-2fe64fb007abd3b55b320ba1acf2f285345248ae.tar.zst WarpX-2fe64fb007abd3b55b320ba1acf2f285345248ae.zip |
Reconfigured the collision classes to allow for generalization (#1583)
* Reconfigured the collision classes to allow for generalization
* Added literals to PairWiseCoulombCollision
* Minor cleanup of PairWiseCoulombCollision.cpp, removing query of ndt
* Add boilerplate to CollisionBase class
* Fixed white space in CollisionBase.H
Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja>
Diffstat (limited to '')
-rw-r--r-- | Source/Particles/MultiParticleContainer.cpp | 52 |
1 files changed, 5 insertions, 47 deletions
diff --git a/Source/Particles/MultiParticleContainer.cpp b/Source/Particles/MultiParticleContainer.cpp index dacc7ba87..9205e14da 100644 --- a/Source/Particles/MultiParticleContainer.cpp +++ b/Source/Particles/MultiParticleContainer.cpp @@ -11,7 +11,6 @@ */ #include "MultiParticleContainer.H" #include "SpeciesPhysicalProperties.H" -#include "Utils/WarpXUtil.H" #include "WarpX.H" #ifdef WARPX_QED #include "Particles/ElementaryProcess/QEDInternals/SchwingerProcessWrapper.H" @@ -71,13 +70,8 @@ MultiParticleContainer::MultiParticleContainer (AmrCore* amr_core) } } - // collision - auto const ncollisions = collision_names.size(); - allcollisions.resize(ncollisions); - for (int i = 0; i < static_cast<int>(ncollisions); ++i) { - allcollisions[i] = - std::make_unique<CollisionType>(species_names, collision_names[i]); - } + // Setup particle collisions + collisionhandler = std::make_unique<CollisionHandler>(); } @@ -241,10 +235,6 @@ MultiParticleContainer::ReadParameters () } } - // binary collisions - ParmParse pc("collisions"); - pc.queryarr("collision_names", collision_names); - } pp.query("use_fdtd_nci_corr", WarpX::use_fdtd_nci_corr); #ifdef WARPX_DIM_RZ @@ -730,42 +720,10 @@ MultiParticleContainer::doFieldIonization (int lev, } void -MultiParticleContainer::doCoulombCollisions ( Real cur_time ) +MultiParticleContainer::doCollisions ( 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 ]; - - // Enable tiling - MFItInfo info; - if (Gpu::notInLaunchRegion()) info.EnableTiling(species1->tile_size); - - // Loop over refinement levels - for (int lev = 0; lev <= species1->finestLevel(); ++lev){ - - // Loop over all grids/tiles at this level -#ifdef _OPENMP - info.SetDynamic(true); -#pragma omp parallel if (Gpu::notInLaunchRegion()) -#endif - for (MFIter mfi = species1->MakeMFIter(lev, info); mfi.isValid(); ++mfi){ - - CollisionType::doCoulombCollisionsWithinTile - ( lev, mfi, species1, species2, - collision->m_isSameSpecies, - collision->m_CoulombLog, - collision->m_ndt ); - - } - } - } + WARPX_PROFILE("MultiParticleContainer::doCollisions()"); + collisionhandler->doCollisions(cur_time, this); } void MultiParticleContainer::doResampling (const int timestep) |