aboutsummaryrefslogtreecommitdiff
path: root/Source/Particles/PhysicalParticleContainer.cpp
diff options
context:
space:
mode:
authorGravatar Luca Fedeli <luca.fedeli@cea.fr> 2019-10-31 11:19:20 +0100
committerGravatar Luca Fedeli <luca.fedeli@cea.fr> 2019-10-31 11:19:20 +0100
commitecaedd97cb868a266971c81545c68930a2802126 (patch)
tree437edd87aa107ddb97e6f20d4dde563f47e11929 /Source/Particles/PhysicalParticleContainer.cpp
parenta2fd7d66ad2f25ef5777d94655da9aeb56f410b6 (diff)
downloadWarpX-ecaedd97cb868a266971c81545c68930a2802126.tar.gz
WarpX-ecaedd97cb868a266971c81545c68930a2802126.tar.zst
WarpX-ecaedd97cb868a266971c81545c68930a2802126.zip
removed QED from PushP
Diffstat (limited to 'Source/Particles/PhysicalParticleContainer.cpp')
-rw-r--r--Source/Particles/PhysicalParticleContainer.cpp250
1 files changed, 48 insertions, 202 deletions
diff --git a/Source/Particles/PhysicalParticleContainer.cpp b/Source/Particles/PhysicalParticleContainer.cpp
index ee0711372..2fa80a7fc 100644
--- a/Source/Particles/PhysicalParticleContainer.cpp
+++ b/Source/Particles/PhysicalParticleContainer.cpp
@@ -1911,214 +1911,60 @@ PhysicalParticleContainer::PushP (int lev, Real dt,
ion_lev = pti.GetiAttribs(particle_icomps["ionization_level"]).dataPtr();
}
-#ifdef WARPX_QED
- if(has_quantum_sync()){
- Real* AMREX_RESTRICT p_tau =
- pti.GetAttribs(particle_comps["tau"]).dataPtr();
- PushP_QedQuantumSynchrotron(
- ux, uy, uz, Expp, Eypp, Ezpp, Bxpp, Bypp, Bzpp,
- q, m, p_tau, dt, pti.numParticles());
- }
- else{
- PushP_classical(
- ux, uy, uz, Expp, Eypp, Ezpp, Bxpp, Bypp, Bzpp,
- q, m, ion_lev, dt, pti.numParticles());
- }
-#else
- PushP_classical(
- ux, uy, uz, Expp, Eypp, Ezpp, Bxpp, Bypp, Bzpp,
- q, m, ion_lev, dt, pti.numParticles());
-#endif
-
-
- }
- }
-}
-
-void PhysicalParticleContainer::PushP_classical(
- ParticleReal* const AMREX_RESTRICT ux,
- ParticleReal* const AMREX_RESTRICT uy,
- ParticleReal* const AMREX_RESTRICT uz,
- const ParticleReal* const AMREX_RESTRICT Expp,
- const ParticleReal* const AMREX_RESTRICT Eypp,
- const ParticleReal* const AMREX_RESTRICT Ezpp,
- const ParticleReal* const AMREX_RESTRICT Bxpp,
- const ParticleReal* const AMREX_RESTRICT Bypp,
- const ParticleReal* const AMREX_RESTRICT Bzpp,
- Real q, Real m, int* AMREX_RESTRICT ion_lev,
- Real dt, long num_particles
-)
-{
- //Assumes that all consistency checks have been done at initialization
- if(do_classical_radiation_reaction){
- amrex::ParallelFor(num_particles,
- [=] AMREX_GPU_DEVICE (long i){
- Real qp = q;
- if (ion_lev){ qp *= ion_lev[i]; }
- UpdateMomentumBorisWithRadiationReaction(
- ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- qp, m, dt);
- }
- );
- } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::Boris){
- amrex::ParallelFor(num_particles,
- [=] AMREX_GPU_DEVICE (long i) {
- Real qp = q;
- if (ion_lev){ qp *= ion_lev[i]; }
- UpdateMomentumBoris(
- ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- qp, m, dt);
- }
- );
- } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::Vay){
- amrex::ParallelFor(num_particles,
- [=] AMREX_GPU_DEVICE (long i) {
- Real qp = q;
- if (ion_lev){ qp *= ion_lev[i]; }
- UpdateMomentumVay(
- ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- qp, m, dt);
- }
- );
- } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::HigueraCary){
- amrex::ParallelFor(num_particles,
- [=] AMREX_GPU_DEVICE (long i) {
- UpdateMomentumHigueraCary( ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i], Bxpp[i], Bypp[i], Bzpp[i], q, m, dt);
- }
- );
- } else {
+ //Assumes that all consistency checks have been done at initialization
+ if(do_classical_radiation_reaction){
+ amrex::ParallelFor(pti.numParticles(),
+ [=] AMREX_GPU_DEVICE (long i){
+ Real qp = q;
+ if (ion_lev){ qp *= ion_lev[i]; }
+ UpdateMomentumBorisWithRadiationReaction(
+ ux[i], uy[i], uz[i],
+ Expp[i], Eypp[i], Ezpp[i],
+ Bxpp[i], Bypp[i], Bzpp[i],
+ qp, m, dt);
+ }
+ );
+ } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::Boris){
+ amrex::ParallelFor(pti.numParticles(),
+ [=] AMREX_GPU_DEVICE (long i) {
+ Real qp = q;
+ if (ion_lev){ qp *= ion_lev[i]; }
+ UpdateMomentumBoris(
+ ux[i], uy[i], uz[i],
+ Expp[i], Eypp[i], Ezpp[i],
+ Bxpp[i], Bypp[i], Bzpp[i],
+ qp, m, dt);
+ }
+ );
+ } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::Vay){
+ amrex::ParallelFor(pti.numParticles(),
+ [=] AMREX_GPU_DEVICE (long i) {
+ Real qp = q;
+ if (ion_lev){ qp *= ion_lev[i]; }
+ UpdateMomentumVay(
+ ux[i], uy[i], uz[i],
+ Expp[i], Eypp[i], Ezpp[i],
+ Bxpp[i], Bypp[i], Bzpp[i],
+ qp, m, dt);
+ }
+ );
+ } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::HigueraCary){
+ amrex::ParallelFor(pti.numParticles(),
+ [=] AMREX_GPU_DEVICE (long i) {
+ UpdateMomentumHigueraCary( ux[i], uy[i], uz[i],
+ Expp[i], Eypp[i], Ezpp[i],
+ Bxpp[i], Bypp[i], Bzpp[i],
+ q, m, dt);
+ }
+ );
+ } else {
amrex::Abort("Unknown particle pusher");
- }
-}
-
-#ifdef WARPX_QED
-
-void PhysicalParticleContainer::PushP_QedQuantumSynchrotron(
- ParticleReal* const AMREX_RESTRICT ux,
- ParticleReal* const AMREX_RESTRICT uy,
- ParticleReal* const AMREX_RESTRICT uz,
- const ParticleReal* const AMREX_RESTRICT Expp,
- const ParticleReal* const AMREX_RESTRICT Eypp,
- const ParticleReal* const AMREX_RESTRICT Ezpp,
- const ParticleReal* const AMREX_RESTRICT Bxpp,
- const ParticleReal* const AMREX_RESTRICT Bypp,
- const ParticleReal* const AMREX_RESTRICT Bzpp,
- Real q, Real m,
- ParticleReal* const AMREX_RESTRICT tau,
- Real dt, long num_particles
-)
-{
- QuantumSynchrotronEvolveOpticalDepth evolve_opt =
- m_shr_p_qs_engine->build_evolve_functor();
-
- const auto chi_min = m_shr_p_qs_engine->get_ref_ctrl().chi_part_min;
-
- //Assumes that all consistency checks have been done at initialization
- if(do_classical_radiation_reaction){
- amrex::ParallelFor(num_particles,
- [=] AMREX_GPU_DEVICE (long i) {
- const ParticleReal px = m * ux[i];
- const ParticleReal py = m * uy[i];
- const ParticleReal pz = m * uz[i];
-
- const ParticleReal chi = QedUtils::chi_lepton(px, py, pz,
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i]);
-
- if(chi < chi_min){
- UpdateMomentumBorisWithRadiationReaction(
- ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- q, m, dt);
- }
- else
- {
- bool has_event_happened = evolve_opt(
- px, py, pz,
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- dt, tau[i]);
-
- UpdateMomentumBoris(
- ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- q, m, dt);
- }
}
- );
- } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::Boris){
- amrex::ParallelFor(num_particles,
- [=] AMREX_GPU_DEVICE (long i) {
- const ParticleReal px = m * ux[i];
- const ParticleReal py = m * uy[i];
- const ParticleReal pz = m * uz[i];
-
- bool has_event_happened = evolve_opt(
- px, py, pz,
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- dt, tau[i]);
- UpdateMomentumBoris(
- ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- q, m, dt);
- }
- );
- } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::Vay){
- amrex::ParallelFor(num_particles,
- [=] AMREX_GPU_DEVICE (long i) {
- const ParticleReal px = m * ux[i];
- const ParticleReal py = m * uy[i];
- const ParticleReal pz = m * uz[i];
-
- bool has_event_happened = evolve_opt(
- px, py, pz,
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- dt, tau[i]);
-
- UpdateMomentumVay(
- ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- q, m, dt);
- }
- );
- } else if (WarpX::particle_pusher_algo == ParticlePusherAlgo::HigueraCary){
- amrex::ParallelFor(num_particles,
- [=] AMREX_GPU_DEVICE (long i) {
- const ParticleReal px = m * ux[i];
- const ParticleReal py = m * uy[i];
- const ParticleReal pz = m * uz[i];
-
- bool has_event_happened = evolve_opt(
- px, py, pz,
- Expp[i], Eypp[i], Ezpp[i],
- Bxpp[i], Bypp[i], Bzpp[i],
- dt, tau[i]);
-
- UpdateMomentumHigueraCary( ux[i], uy[i], uz[i],
- Expp[i], Eypp[i], Ezpp[i], Bxpp[i], Bypp[i], Bzpp[i], q, m, dt);
- }
- );
- } else {
- amrex::Abort("Unknown particle pusher");
+ }
}
}
-#endif
-
void PhysicalParticleContainer::copy_attribs(WarpXParIter& pti,const ParticleReal* xp,
const ParticleReal* yp, const ParticleReal* zp)
{