aboutsummaryrefslogtreecommitdiff
path: root/Source/Particles/PhysicalParticleContainer.cpp
diff options
context:
space:
mode:
Diffstat (limited to 'Source/Particles/PhysicalParticleContainer.cpp')
-rw-r--r--Source/Particles/PhysicalParticleContainer.cpp107
1 files changed, 83 insertions, 24 deletions
diff --git a/Source/Particles/PhysicalParticleContainer.cpp b/Source/Particles/PhysicalParticleContainer.cpp
index 936a64c9c..21c8ddd31 100644
--- a/Source/Particles/PhysicalParticleContainer.cpp
+++ b/Source/Particles/PhysicalParticleContainer.cpp
@@ -967,6 +967,80 @@ PhysicalParticleContainer::EvolveES (const Vector<std::array<std::unique_ptr<Mul
#endif // WARPX_DO_ELECTROSTATIC
void
+PhysicalParticleContainer::AssignExternalFieldOnParticles(WarpXParIter& pti,
+ RealVector& Exp, RealVector& Eyp, RealVector& Ezp,
+ RealVector& Bxp, RealVector& Byp, RealVector& Bzp,
+ Gpu::ManagedDeviceVector<ParticleReal> xp,
+ Gpu::ManagedDeviceVector<ParticleReal> yp,
+ Gpu::ManagedDeviceVector<ParticleReal> zp, int lev)
+{
+ const long np = pti.numParticles();
+ /// get WarpX class object
+ auto & warpx = WarpX::GetInstance();
+ /// get MultiParticleContainer class object
+ auto & mypc = warpx.GetPartContainer();
+ if (mypc.m_E_ext_particle_s=="constant" ||
+ mypc.m_E_ext_particle_s=="default") {
+ Exp.assign(np,mypc.m_E_external_particle[0]);
+ Eyp.assign(np,mypc.m_E_external_particle[1]);
+ Ezp.assign(np,mypc.m_E_external_particle[2]);
+ }
+ if (mypc.m_B_ext_particle_s=="constant" ||
+ mypc.m_B_ext_particle_s=="default") {
+ Bxp.assign(np,mypc.m_B_external_particle[0]);
+ Byp.assign(np,mypc.m_B_external_particle[1]);
+ Bzp.assign(np,mypc.m_B_external_particle[2]);
+ }
+ if (mypc.m_E_ext_particle_s=="parse_e_ext_particle_function") {
+ Real* const AMREX_RESTRICT xp_data = xp.dataPtr();
+ Real* const AMREX_RESTRICT yp_data = yp.dataPtr();
+ Real* const AMREX_RESTRICT zp_data = zp.dataPtr();
+ Real* const AMREX_RESTRICT Exp_data = Exp.dataPtr();
+ Real* const AMREX_RESTRICT Eyp_data = Eyp.dataPtr();
+ Real* const AMREX_RESTRICT Ezp_data = Ezp.dataPtr();
+ ParserWrapper *xfield_partparser = mypc.m_Ex_particle_parser.get();
+ ParserWrapper *yfield_partparser = mypc.m_Ey_particle_parser.get();
+ ParserWrapper *zfield_partparser = mypc.m_Ez_particle_parser.get();
+ Real time = warpx.gett_new(lev);
+ amrex::ParallelFor(pti.numParticles(),
+ [=] AMREX_GPU_DEVICE (long i) {
+ Exp_data[i] = xfield_partparser->getField(xp_data[i],yp_data[i],zp_data[i],time);
+ Eyp_data[i] = yfield_partparser->getField(xp_data[i],yp_data[i],zp_data[i],time);
+ Ezp_data[i] = zfield_partparser->getField(xp_data[i],yp_data[i],zp_data[i],time);
+ },
+ /* To allocate shared memory for the GPU threads. */
+ /* But, for now only 4 doubles (x,y,z,t) are allocated. */
+ amrex::Gpu::numThreadsPerBlockParallelFor() * sizeof(double) * 4
+ );
+ }
+ if (mypc.m_B_ext_particle_s=="parse_b_ext_particle_function") {
+ Real* const AMREX_RESTRICT xp_data = xp.dataPtr();
+ Real* const AMREX_RESTRICT yp_data = yp.dataPtr();
+ Real* const AMREX_RESTRICT zp_data = zp.dataPtr();
+ Real* const AMREX_RESTRICT Bxp_data = Bxp.dataPtr();
+ Real* const AMREX_RESTRICT Byp_data = Byp.dataPtr();
+ Real* const AMREX_RESTRICT Bzp_data = Bzp.dataPtr();
+ ParserWrapper *xfield_partparser = mypc.m_Bx_particle_parser.get();
+ ParserWrapper *yfield_partparser = mypc.m_By_particle_parser.get();
+ ParserWrapper *zfield_partparser = mypc.m_Bz_particle_parser.get();
+ Real time = warpx.gett_new(lev);
+ amrex::ParallelFor(pti.numParticles(),
+ [=] AMREX_GPU_DEVICE (long i) {
+ Bxp_data[i] = xfield_partparser->getField(xp_data[i],yp_data[i],zp_data[i],time);
+ Byp_data[i] = yfield_partparser->getField(xp_data[i],yp_data[i],zp_data[i],time);
+ Bzp_data[i] = zfield_partparser->getField(xp_data[i],yp_data[i],zp_data[i],time);
+ },
+ /* To allocate shared memory for the GPU threads. */
+ /* But, for now only 4 doubles (x,y,z,t) are allocated. */
+ amrex::Gpu::numThreadsPerBlockParallelFor() * sizeof(double) * 4
+ );
+ }
+
+}
+
+
+
+void
PhysicalParticleContainer::FieldGather (int lev,
const amrex::MultiFab& Ex,
const amrex::MultiFab& Ey,
@@ -1015,13 +1089,6 @@ PhysicalParticleContainer::FieldGather (int lev,
const FArrayBox& byfab = By[pti];
const FArrayBox& bzfab = Bz[pti];
- Exp.assign(np,WarpX::E_external_particle[0]);
- Eyp.assign(np,WarpX::E_external_particle[1]);
- Ezp.assign(np,WarpX::E_external_particle[2]);
- Bxp.assign(np,WarpX::B_external_particle[0]);
- Byp.assign(np,WarpX::B_external_particle[1]);
- Bzp.assign(np,WarpX::B_external_particle[2]);
-
//
// copy data from particle container to temp arrays
//
@@ -1152,14 +1219,6 @@ PhysicalParticleContainer::Evolve (int lev,
exfab, eyfab, ezfab, bxfab, byfab, bzfab);
}
- Exp.assign(np,WarpX::E_external_particle[0]);
- Eyp.assign(np,WarpX::E_external_particle[1]);
- Ezp.assign(np,WarpX::E_external_particle[2]);
-
- Bxp.assign(np,WarpX::B_external_particle[0]);
- Byp.assign(np,WarpX::B_external_particle[1]);
- Bzp.assign(np,WarpX::B_external_particle[2]);
-
// Determine which particles deposit/gather in the buffer, and
// which particles deposit/gather in the fine patch
long nfine_current = np;
@@ -1799,14 +1858,6 @@ PhysicalParticleContainer::PushP (int lev, Real dt,
const FArrayBox& byfab = By[pti];
const FArrayBox& bzfab = Bz[pti];
- Exp.assign(np,WarpX::E_external_particle[0]);
- Eyp.assign(np,WarpX::E_external_particle[1]);
- Ezp.assign(np,WarpX::E_external_particle[2]);
-
- Bxp.assign(np,WarpX::B_external_particle[0]);
- Byp.assign(np,WarpX::B_external_particle[1]);
- Bzp.assign(np,WarpX::B_external_particle[2]);
-
//
// copy data from particle container to temp arrays
//
@@ -2176,9 +2227,16 @@ PhysicalParticleContainer::FieldGather (WarpXParIter& pti,
AMREX_ALWAYS_ASSERT_WITH_MESSAGE((gather_lev==(lev-1)) ||
(gather_lev==(lev )),
"Gather buffers only work for lev-1");
-
// If no particles, do not do anything
if (np_to_gather == 0) return;
+
+ // initializing the field value to the externally applied field before
+ // gathering fields from the grid to the particles.
+ AssignExternalFieldOnParticles(pti, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
+ m_xp[thread_num], m_yp[thread_num],
+ m_zp[thread_num], lev);
+
+
// Get cell size on gather_lev
const std::array<Real,3>& dx = WarpX::CellSize(std::max(gather_lev,0));
@@ -2424,4 +2482,5 @@ set_quantum_sync_engine_ptr(std::shared_ptr<QuantumSynchrotronEngine> ptr)
{
m_shr_p_qs_engine = ptr;
}
+
#endif