WarpX
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 PARTICLEIO_H_
8 #define PARTICLEIO_H_
9 
11 
12 #include <AMReX_AmrParticles.H>
13 #include <AMReX_Gpu.H>
14 #include <AMReX_REAL.H>
15 
17 
33 template< typename T_ParticleContainer >
34 void
35 particlesConvertUnits (ConvertDirection convert_direction, T_ParticleContainer* pc, amrex::ParticleReal const mass )
36 {
37  using namespace amrex;
38 
39  // Compute conversion factor
40  auto factor = 1_rt;
41 
42  if (convert_direction == ConvertDirection::WarpX_to_SI){
43  factor = mass;
44  } else if (convert_direction == ConvertDirection::SI_to_WarpX){
45  factor = 1._rt/mass;
46  }
47 
48  const int nLevels = pc->finestLevel();
49  for (int lev=0; lev<=nLevels; lev++){
50 #ifdef AMREX_USE_OMP
51 #pragma omp parallel if (Gpu::notInLaunchRegion())
52 #endif
53  for (ParIter pti(*pc, lev); pti.isValid(); ++pti)
54  {
55  // - momenta are stored as a struct of array, in `attribs`
56  // The GetStructOfArrays is called directly since the convenience routine GetAttribs
57  // is only available in WarpXParIter. ParIter is used here since the pc passed in
58  // will sometimes be a PinnedMemoryParticleContainer (not derived from a WarpXParticleContainer).
59  auto& attribs = pti.GetStructOfArrays().GetRealData();
60  ParticleReal* AMREX_RESTRICT ux = attribs[PIdx::ux].dataPtr();
61  ParticleReal* AMREX_RESTRICT uy = attribs[PIdx::uy].dataPtr();
62  ParticleReal* AMREX_RESTRICT uz = attribs[PIdx::uz].dataPtr();
63  // Loop over the particles and convert momentum
64  const long np = pti.numParticles();
65  ParallelFor( np,
66  [=] AMREX_GPU_DEVICE (long i) {
67  ux[i] *= factor;
68  uy[i] *= factor;
69  uz[i] *= factor;
70  }
71  );
72  }
73  }
74 }
75 
76 #endif /* PARTICLEIO_H_ */
void particlesConvertUnits(ConvertDirection convert_direction, T_ParticleContainer *pc, amrex::ParticleReal const mass)
Definition: ParticleIO.H:35
ConvertDirection
Definition: ParticleIO.H:16
Definition: NamedComponentParticleContainer.H:25
Definition: NamedComponentParticleContainer.H:25
std::enable_if_t< std::is_integral< T >::value > ParallelFor(TypeList< CTOs... >, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
#define AMREX_GPU_DEVICE
i
Definition: check_interp_points_and_weights.py:174
Definition: NamedComponentParticleContainer.H:25
bool isValid() const noexcept
#define AMREX_RESTRICT