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