WarpX
Loading...
Searching...
No Matches
ParticleIO.H
Go to the documentation of this file.
1/* Copyright 2021 Axel Huebl
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef WARPX_PARTICLEIO_H_
8#define WARPX_PARTICLEIO_H_
9
12
14
15#include <AMReX_AmrParticles.H>
16#include <AMReX_ParIter.H>
17#include <AMReX_Gpu.H>
18#include <AMReX_REAL.H>
19
20
22
38template< typename T_ParticleContainer >
39void
40particlesConvertUnits (ConvertDirection convert_direction, T_ParticleContainer* pc, amrex::ParticleReal const mass )
41{
42 using namespace amrex;
43
44 // Compute conversion factor
45 auto factor = 1_rt;
46
47 if (convert_direction == ConvertDirection::WarpX_to_SI){
48 factor = mass;
49 } else if (convert_direction == ConvertDirection::SI_to_WarpX){
50 factor = 1._rt/mass;
51 }
52
53 using PinnedParIter = typename T_ParticleContainer::ParIterType;
54
55 const int nLevels = pc->finestLevel();
56 for (int lev=0; lev<=nLevels; lev++){
57#ifdef AMREX_USE_OMP
58#pragma omp parallel if (Gpu::notInLaunchRegion())
59#endif
60 for (PinnedParIter pti(*pc, lev); pti.isValid(); ++pti)
61 {
62 // - momenta are stored as a struct of array, in `attribs`
63 // The GetStructOfArrays is called directly since the convenience routine GetAttribs
64 // is only available in WarpXParIter. ParIter is used here since the pc passed in
65 // will sometimes be a WarpXParticleContainer (not derived from a WarpXParticleContainer).
66 auto& attribs = pti.GetStructOfArrays().GetRealData();
67 ParticleReal* AMREX_RESTRICT ux = attribs[PIdx::ux].dataPtr();
68 ParticleReal* AMREX_RESTRICT uy = attribs[PIdx::uy].dataPtr();
69 ParticleReal* AMREX_RESTRICT uz = attribs[PIdx::uz].dataPtr();
70 // Loop over the particles and convert momentum
71 const long np = pti.numParticles();
72 ParallelFor( np,
73 [=] AMREX_GPU_DEVICE (long i) {
74 ux[i] *= factor;
75 uy[i] *= factor;
76 uz[i] *= factor;
77 }
78 );
79 }
80 }
81}
82
90void
92 ElectrostaticSolverAlgo electrostatic_solver_id, bool is_full_diagnostic );
93
106void
107storeFieldOnParticles (WarpXParticleContainer::Base& tmp, bool is_full_diagnostic,
108 bool save_Ex, bool save_Ey, bool save_Ez,
109 bool save_Bx, bool save_By, bool save_Bz);
110
111#endif /* WARPX_PARTICLEIO_H_ */
#define AMREX_RESTRICT
#define AMREX_GPU_DEVICE
void storeFieldOnParticles(WarpXParticleContainer::Base &tmp, bool is_full_diagnostic, bool save_Ex, bool save_Ey, bool save_Ez, bool save_Bx, bool save_By, bool save_Bz)
Gathers fields from a MultiFab to the macroparticles. Adds a runtime component of the particle contai...
Definition ParticleIO.cpp:308
void particlesConvertUnits(ConvertDirection convert_direction, T_ParticleContainer *pc, amrex::ParticleReal const mass)
Definition ParticleIO.H:40
void storePhiOnParticles(WarpXParticleContainer::Base &tmp, ElectrostaticSolverAlgo electrostatic_solver_id, bool is_full_diagnostic)
Definition ParticleIO.cpp:253
ConvertDirection
Definition ParticleIO.H:21
@ WarpX_to_SI
Definition ParticleIO.H:21
@ SI_to_WarpX
Definition ParticleIO.H:21
ElectrostaticSolverAlgo
Definition WarpXAlgorithmSelection.H:66
amrex::ParticleContainerPureSoA< PIdx::nattribs, 0, amrex::PolymorphicArenaAllocator > Base
Definition WarpXParticleContainer.H:199
amrex_particle_real ParticleReal
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
@ uz
Definition WarpXParticleContainer.H:70
@ uy
Definition WarpXParticleContainer.H:70
@ ux
Definition WarpXParticleContainer.H:70