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