aboutsummaryrefslogtreecommitdiff
path: root/Source/Particles/LaserParticleContainer.cpp
diff options
context:
space:
mode:
authorGravatar Kevin Z. Zhu <86268612+KZhu-ME@users.noreply.github.com> 2022-05-24 17:57:51 -0700
committerGravatar GitHub <noreply@github.com> 2022-05-24 17:57:51 -0700
commit0344c42e27c61b286e885889028c7f77725e7f74 (patch)
treeec116dbe570b342c7df3da3e6d620cbdd4032db5 /Source/Particles/LaserParticleContainer.cpp
parenta681d1994aa5748d4faceb658199c1a61d418a79 (diff)
downloadWarpX-0344c42e27c61b286e885889028c7f77725e7f74.tar.gz
WarpX-0344c42e27c61b286e885889028c7f77725e7f74.tar.zst
WarpX-0344c42e27c61b286e885889028c7f77725e7f74.zip
Specify particle precision (#3065)
* Add precision to printed PIC parameters * Added WarpX_PARTICLE_PRECISION as a build option * Update types to ParticleReal * Updated libwarpx to inject particles with correct ParticleReal type * Fix syntax error * Add logic to avoid duplicate definitions * Use correct ParticleReal type in add_particles * Cleaned up code, addressed comments * Update Python/pywarpx/_libwarpx.py Co-authored-by: Peter Scherpelz <31747262+peterscherpelz@users.noreply.github.com> * Removed redundant functions, fixed some typing * Modified template functions * Cast d_w to Real * Fixed failing tests * Cast types to be consistent * removed in-tree-build from pip command * Added GPU device macros to PDim3 methods * rerun tests * Removed unecessary casting, update calls to use PDim3 instead of XDim3 * Refactored comments * Added mcc fields double precision, particles single precision test * Updated casting and formatting * Removed cast, updated declaration Co-authored-by: Peter Scherpelz <peter.scherpelz@modernelectron.com> Co-authored-by: Peter Scherpelz <31747262+peterscherpelz@users.noreply.github.com>
Diffstat (limited to 'Source/Particles/LaserParticleContainer.cpp')
-rw-r--r--Source/Particles/LaserParticleContainer.cpp8
1 files changed, 4 insertions, 4 deletions
diff --git a/Source/Particles/LaserParticleContainer.cpp b/Source/Particles/LaserParticleContainer.cpp
index 12d7e78cd..fbe2b9710 100644
--- a/Source/Particles/LaserParticleContainer.cpp
+++ b/Source/Particles/LaserParticleContainer.cpp
@@ -453,7 +453,7 @@ LaserParticleContainer::InitData (int lev)
BoxArray plane_ba { Box {IntVect(0), IntVect(0)} };
#endif
- amrex::Vector<amrex::Real> particle_x, particle_y, particle_z, particle_w;
+ amrex::Vector<amrex::ParticleReal> particle_x, particle_y, particle_z, particle_w;
const DistributionMapping plane_dm {plane_ba, nprocs};
const Vector<int>& procmap = plane_dm.ProcessorMap();
@@ -506,9 +506,9 @@ LaserParticleContainer::InitData (int lev)
}
}
const int np = particle_z.size();
- amrex::Vector<amrex::Real> particle_ux(np, 0.0);
- amrex::Vector<amrex::Real> particle_uy(np, 0.0);
- amrex::Vector<amrex::Real> particle_uz(np, 0.0);
+ amrex::Vector<amrex::ParticleReal> particle_ux(np, 0.0);
+ amrex::Vector<amrex::ParticleReal> particle_uy(np, 0.0);
+ amrex::Vector<amrex::ParticleReal> particle_uz(np, 0.0);
if (Verbose()) amrex::Print() << Utils::TextMsg::Info("Adding laser particles");
// Add particles on level 0. They will be redistributed afterwards