ImpactX
Loading...
Searching...
No Matches
PushAll.H
Go to the documentation of this file.
1/* Copyright 2022-2023 The Regents of the University of California, through Lawrence
2 * Berkeley National Laboratory (subject to receipt of any required
3 * approvals from the U.S. Dept. of Energy). All rights reserved.
4 *
5 * This file is part of ImpactX.
6 *
7 * Authors: Axel Huebl
8 * License: BSD-3-Clause-LBNL
9 */
10#ifndef IMPACTX_PUSH_ALL_H
11#define IMPACTX_PUSH_ALL_H
12
14
15#include <AMReX_BLProfiler.H>
16
17
18namespace impactx
19{
32 template<typename T_Element>
33 void push_all (
35 T_Element & element,
36 [[maybe_unused]] int step,
37 [[maybe_unused]] int period,
38 [[maybe_unused]] bool omp_parallel = true
39 )
40 {
41 // performance profiling per element
42 std::string const profile_name = "impactx::push::" + std::string(T_Element::type);
43 BL_PROFILE(profile_name);
44
45 // is spin tracking enabled?
46 bool spin = false;
47 amrex::ParmParse pp_algo("algo");
48 pp_algo.query("spin", spin);
49
50 // preparing to access reference particle data: RefPart
51 RefPart & ref_part = pc.GetRefParticle();
52
53 // push reference particle in global coordinates
54 {
55 BL_PROFILE("impactx::push::RefPart");
56 element(ref_part);
57 }
58
59 // loop over refinement levels
60 int const nLevel = pc.finestLevel();
61 for (int lev = 0; lev <= nLevel; ++lev)
62 {
63 // loop over all particle boxes
65#ifdef AMREX_USE_OMP
66#pragma omp parallel if (amrex::Gpu::notInLaunchRegion() && omp_parallel)
67#endif
68 for (ParIt pti(pc, lev); pti.isValid(); ++pti) {
69 // push beam particles relative to reference particle
70 element(pti, ref_part, spin);
71 } // end loop over all particle boxes
72 } // env mesh-refinement level loop
73 }
74
75} // namespace impactx
76
78#define IMPACTX_PUSH_EXTERN_TEMPLATE(ElementType) \
79 extern template void impactx::push_all<ElementType>( \
80 impactx::ImpactXParticleContainer &, \
81 ElementType &, \
82 int, int, bool);
83
85#define IMPACTX_PUSH_INSTANTIATE(ElementType) \
86 template void impactx::push_all<ElementType>( \
87 impactx::ImpactXParticleContainer &, \
88 ElementType &, \
89 int, int, bool);
90
91#endif // IMPACTX_PUSH_ALL_H
#define BL_PROFILE(a)
int query(std::string_view name, bool &ref, int ival=FIRST) const
Definition ImpactXParticleContainer.H:136
impactx::ParIterSoA iterator
amrex iterator for particle boxes
Definition ImpactXParticleContainer.H:139
RefPart & GetRefParticle()
Definition ImpactXParticleContainer.cpp:419
Definition CovarianceMatrixMath.H:25
void push_all(ImpactXParticleContainer &pc, T_Element &element, int step, int period, bool omp_parallel=true)
Definition PushAll.H:33
Definition ReferenceParticle.H:33