/home/runner/work/spades/spades/Source/ParticleOps.H Source File

SPADES API: /home/runner/work/spades/spades/Source/ParticleOps.H Source File
SPADES API
ParticleOps.H
Go to the documentation of this file.
1#ifndef PARTICLEOPS_H
2#define PARTICLEOPS_H
3
4#include <AMReX.H>
5#include <AMReX_Array4.H>
6#include <AMReX_Print.H>
7#include "ParticleData.H"
8
9namespace spades::particles {
10
12template <size_t NReal, size_t NInt, class P, class RT, class IT>
14{
22 P* aos, std::array<RT, NReal>& rdata, std::array<IT, NInt>& idata)
23 : m_aos(aos)
24 {
25 for (size_t i = 0; i < NReal; i++) {
26 m_rdata[i] = rdata[i].data();
27 }
28 for (size_t i = 0; i < NInt; i++) {
29 m_idata[i] = idata[i].data();
30 }
31 }
32
34 P* m_aos = nullptr;
35
37 amrex::GpuArray<amrex::Real*, NReal> m_rdata;
38
40 amrex::GpuArray<int*, NInt> m_idata;
41};
42
44struct Print
45{
51 template <class PArrs>
52 void operator()(const amrex::Long n, PArrs& parrs) const
53 {
54 auto& p = parrs.m_aos[n];
55 amrex::Print() << "Particle data: \n"
56 << " id = " << p.id() << "\n"
57 << " cpu = " << p.cpu() << "\n"
58 << " type = "
59 << parrs.m_idata[CommonIntData::type_id][n] << "\n";
60 for (int i = 0; i < p.NInt; i++) {
61 amrex::Print() << " int comp(" << i
62 << ") = " << parrs.m_idata[i][n] << "\n";
63 }
64 for (int i = 0; i < p.NReal; i++) {
65 amrex::Print() << " real comp(" << i
66 << ") = " << parrs.m_rdata[i][n] << "\n";
67 }
68 }
69};
70
73{
79 template <class PArrs>
80 AMREX_GPU_DEVICE void operator()(const amrex::Long n, PArrs& parrs) const
81 {
82 auto& p = parrs.m_aos[n];
83 printf(
84 "Particle data: id = %ld, cpu = %ld, type = %d, timestamp = %.8e\n",
85 static_cast<amrex::Long>(p.id()), static_cast<amrex::Long>(p.cpu()),
86 parrs.m_idata[CommonIntData::type_id][n],
87 parrs.m_rdata[CommonRealData::timestamp][n]);
88 }
89};
90
92struct Copy
93{
100 template <class PArrs>
101 AMREX_GPU_DEVICE AMREX_FORCE_INLINE void
102 operator()(const amrex::Long p1, const amrex::Long p2, PArrs& parrs) const
103 {
104 for (int i = 0; i < PArrs::NReal; i++) {
105 parrs.m_rdata[i][p2] = parrs.m_rdata[i][p1];
106 }
107 for (int i = 0; i < PArrs::NInt; i++) {
108 parrs.m_idata[i][p2] = parrs.m_idata[i][p1];
109 }
110 }
111};
112
114template <class PArrs>
115struct Get
116{
124 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
125 Get(const amrex::IntVect& iv,
126 const amrex::Array4<const int>& counts,
127 const amrex::Array4<const int>& offsets,
128 PArrs& parrs)
129 : m_iv(iv), m_counts(counts), m_offsets(offsets), m_parrs(parrs)
130 {}
131
138 AMREX_GPU_DEVICE AMREX_FORCE_INLINE int
139 operator()(const int n, const int typ) const
140 {
141
142 AMREX_ASSERT(m_counts(m_iv, typ) > n);
143 const int idx = m_offsets(m_iv, typ) + n;
144
145 AMREX_ASSERT(m_parrs.m_aos[idx].id() >= 0);
146 AMREX_ASSERT(m_parrs.m_idata[CommonIntData::type_id][idx] == typ);
147
148 return idx;
149 }
150
157 AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool
158 check(const int n, const int typ) const
159 {
160 AMREX_ASSERT(m_counts(m_iv, typ) > n);
161 const int idx = m_offsets(m_iv, typ) + n;
162 return m_parrs.m_idata[CommonIntData::type_id][idx] == typ;
163 }
164
178 AMREX_GPU_DEVICE AMREX_FORCE_INLINE void
179 assert_different(const int n, const int typ, const int expected_type) const
180 {
181#ifdef AMREX_DEBUG
182 const int idx = m_offsets(m_iv, typ) + n;
183 AMREX_ASSERT(
184 m_parrs.m_idata[CommonIntData::type_id][idx] == expected_type);
185#else
186 amrex::ignore_unused(n, typ, expected_type);
187#endif
188 }
189
191 const amrex::IntVect& m_iv;
192
194 const amrex::Array4<const int>& m_counts;
195
197 const amrex::Array4<const int>& m_offsets;
198
200 PArrs& m_parrs;
201};
202
205{
210 template <class P>
211 AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool
212 operator()(const P& /*p1*/, const P& /*p2*/) const
213 {
214 return false;
215 }
216};
217
225{
233 template <class PArrs>
234 AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool
235 operator()(const amrex::Long p1, const amrex::Long p2, PArrs& parrs) const
236 {
237 // sort by iv, then particle type then by timestamp
238 const amrex::IntVect piv1(AMREX_D_DECL(
239 parrs.m_idata[CommonIntData::i][p1],
240 parrs.m_idata[CommonIntData::j][p1],
241 parrs.m_idata[CommonIntData::k][p1]));
242 const amrex::IntVect piv2(AMREX_D_DECL(
243 parrs.m_idata[CommonIntData::i][p2],
244 parrs.m_idata[CommonIntData::j][p2],
245 parrs.m_idata[CommonIntData::k][p2]));
246
247 const auto m1 = parrs.m_idata[CommonIntData::type_id][p1];
248 const auto m2 = parrs.m_idata[CommonIntData::type_id][p2];
249 const auto t1 = parrs.m_rdata[CommonRealData::timestamp][p1];
250 const auto t2 = parrs.m_rdata[CommonRealData::timestamp][p2];
251 return (piv1 < piv2) ||
252 (piv1 == piv2 && ((m1 < m2) || (m1 == m2 && t1 < t2)));
253 }
254};
255
256} // namespace spades::particles
257#endif
SPADES particles.
Definition EntityData.H:7
Functor for pairing particles.
Definition ParticleOps.H:205
AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool operator()(const P &, const P &) const
Compare particles for pairing.
Definition ParticleOps.H:212
@ type_id
Definition ParticleData.H:18
@ timestamp
Definition ParticleData.H:12
Functor for comparing particles.
Definition ParticleOps.H:225
AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool operator()(const amrex::Long p1, const amrex::Long p2, PArrs &parrs) const
Compare particles.
Definition ParticleOps.H:235
Functor for copying a particle's data to another.
Definition ParticleOps.H:93
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void operator()(const amrex::Long p1, const amrex::Long p2, PArrs &parrs) const
Copy particle data.
Definition ParticleOps.H:102
Functor for printing particle data on device.
Definition ParticleOps.H:73
AMREX_GPU_DEVICE void operator()(const amrex::Long n, PArrs &parrs) const
Print particle data to screen.
Definition ParticleOps.H:80
Functor for accessing a particle in a cell.
Definition ParticleOps.H:116
const amrex::Array4< const int > & m_offsets
Cell offsets of particle types.
Definition ParticleOps.H:197
PArrs & m_parrs
Particle array.
Definition ParticleOps.H:200
const amrex::IntVect & m_iv
Cell index.
Definition ParticleOps.H:191
AMREX_GPU_DEVICE AMREX_FORCE_INLINE Get(const amrex::IntVect &iv, const amrex::Array4< const int > &counts, const amrex::Array4< const int > &offsets, PArrs &parrs)
Constructor.
Definition ParticleOps.H:125
AMREX_GPU_DEVICE AMREX_FORCE_INLINE int operator()(const int n, const int typ) const
Get a particle in a cell.
Definition ParticleOps.H:139
AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool check(const int n, const int typ) const
Check validity of the particle type.
Definition ParticleOps.H:158
const amrex::Array4< const int > & m_counts
Cell counts of particle types.
Definition ParticleOps.H:194
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void assert_different(const int n, const int typ, const int expected_type) const
Assert that the particle is of a different but expected type.
Definition ParticleOps.H:179
Particle operations.
Definition ParticleOps.H:14
P * m_aos
AOS data.
Definition ParticleOps.H:34
ParticleArrays(P *aos, std::array< RT, NReal > &rdata, std::array< IT, NInt > &idata)
Constructor.
Definition ParticleOps.H:21
amrex::GpuArray< amrex::Real *, NReal > m_rdata
Pointers to real data.
Definition ParticleOps.H:37
amrex::GpuArray< int *, NInt > m_idata
Pointers to int data.
Definition ParticleOps.H:40
Functor for printing particle data on host.
Definition ParticleOps.H:45
void operator()(const amrex::Long n, PArrs &parrs) const
Print particle data to screen.
Definition ParticleOps.H:52