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

SPADES API: /home/runner/work/spades/spades/Source/SpadesParticleContainer.H Source File
SPADES API
SpadesParticleContainer.H
Go to the documentation of this file.
1#ifndef SPADESPARTICLECONTAINER_H
2#define SPADESPARTICLECONTAINER_H
3#include <AMReX.H>
4#include <AMReX_AmrCore.H>
5#include <AMReX_AmrParGDB.H>
6#include <AMReX_NeighborParticles.H>
7#include <AMReX_Random.H>
8#include "ParticleData.H"
9#include "ParticleOps.H"
10#include "Utilities.H"
11
12#if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
13#include <thrust/sort.h>
14#include <thrust/execution_policy.h>
15#endif
16
20namespace spades::particles {
21
23template <
24 int NType,
25 int NStructReal,
26 int NStructInt,
27 int NArrayReal,
28 int NArrayInt>
30 : public amrex::NeighborParticleContainer<
31 NStructReal,
32 NStructInt,
33 NArrayReal,
34 NArrayInt>
35{
36public:
37 using ParticleType = typename amrex::NeighborParticleContainer<
38 NStructReal,
39 NStructInt,
40 NArrayReal,
41 NArrayInt>::ParticleType;
42 using ParticleTileType = typename amrex::NeighborParticleContainer<
43 NStructReal,
44 NStructInt,
45 NArrayReal,
46 NArrayInt>::ParticleTileType;
47 using IntVector = typename amrex::NeighborParticleContainer<
48 NStructReal,
49 NStructInt,
50 NArrayReal,
51 NArrayInt>::IntVector;
52 using RealVector = typename amrex::NeighborParticleContainer<
53 NStructReal,
54 NStructInt,
55 NArrayReal,
56 NArrayInt>::RealVector;
57
62 static std::string identifier() { return "spades_particles"; }
63
69 explicit SpadesParticleContainer(amrex::AmrParGDB* par_gdb, int ngrow = 0);
70
79 const amrex::Vector<amrex::Geometry>& geom,
80 const amrex::Vector<amrex::DistributionMapping>& dmap,
81 const amrex::Vector<amrex::BoxArray>& ba,
82 int ngrow = 0);
83
86
89
92
95
98
103 const amrex::iMultiFab& counts() const { return m_counts; };
104
109 const amrex::iMultiFab& offsets() const { return m_offsets; };
110
116 amrex::Long total_count(const int typ) const
117 {
118 BL_PROFILE("spades::SpadesParticleContainer::total_count()");
119 return m_counts.sum(typ);
120 }
121
127 amrex::Long min_count(const int typ) const
128 {
129 BL_PROFILE("spades::SpadesParticleContainer::min_count()");
130 return m_counts.min(typ);
131 }
132
138 amrex::Long max_count(const int typ) const
139 {
140 BL_PROFILE("spades::SpadesParticleContainer::max_count()");
141 return m_counts.max(typ);
142 }
143
145 void check_sort(const amrex::MFIter& mfi);
146
148 virtual void sort() = 0;
149
154 template <typename CompareFunctor>
155 void sort_impl(const CompareFunctor& compare);
156
161 template <typename CompareFunctor>
162 void nonencoded_sort_impl(const CompareFunctor& compare);
163
166
168 void print_messages(const std::string& header);
169
172
177 virtual void write_plot_file(const std::string& plt_filename) = 0;
178
180 int ngrow() const { return m_ngrow; }
181
183 static constexpr int LEV{0};
184
186 virtual void initialize_variable_names() = 0;
187
194 const std::string& plt_filename, const std::string& name);
195
197 virtual void read_parameters()
198 {
199 {
200 amrex::ParmParse pp("spades");
201 pp.query("sort_type", m_sort_type);
203 }
204 }
205
207 void check_sort_type(const std::string& sort_type)
208 {
209 const amrex::Vector<std::string> valid_types = {
210 "nonencoded", "encoded"};
211 if (std::find(valid_types.cbegin(), valid_types.cend(), sort_type) ==
212 valid_types.cend()) {
213 amrex::Abort("Invalid sort type. Must be nonencoded or encoded");
214 }
215 }
216
219 {
220 auto& aos = pti.GetArrayOfStructs();
221 auto& soa = pti.GetStructOfArrays();
222 return ParticleArrays<
223 NArrayReal, NArrayInt, ParticleType, RealVector, IntVector>(
224 aos.dataPtr(), soa.GetRealData(), soa.GetIntData());
225 }
226
227protected:
230
232 amrex::Vector<int> m_writeflags_real;
233
235 amrex::Vector<int> m_writeflags_int;
238 amrex::Vector<std::string> m_real_data_names;
239
241 amrex::Vector<std::string> m_int_data_names;
242
244 std::string m_sort_type{"nonencoded"};
245
247 amrex::iMultiFab m_counts;
248
250 amrex::iMultiFab m_offsets;
251
253 amrex::MultiFab m_min_timestamp;
254
256 amrex::MultiFab m_max_timestamp;
257};
258
260
261} // namespace spades::particles
262#endif
Main SPADES particle container.
Definition SpadesParticleContainer.H:35
virtual void read_parameters()
Read user parameters.
Definition SpadesParticleContainer.H:197
static std::string identifier()
Class identifier name.
Definition SpadesParticleContainer.H:62
virtual void write_plot_file(const std::string &plt_filename)=0
Write the particles to file.
amrex::iMultiFab m_offsets
Offsets of particles in each cell.
Definition SpadesParticleContainer.H:250
void sort_impl(const CompareFunctor &compare)
Sort the particles implementation.
amrex::Vector< std::string > m_int_data_names
Names for int data to write to file.
Definition SpadesParticleContainer.H:241
ParticleArrays< NArrayReal, NArrayInt, ParticleType, RealVector, IntVector > particle_arrays(ParticleTileType &pti) const
Definition SpadesParticleContainer.H:218
void clear_state()
Delete particle states (counts and offsets)
Definition SpadesParticleContainer.H:105
amrex::MultiFab m_max_timestamp
Maximum timestamp in each cell for each type.
Definition SpadesParticleContainer.H:256
amrex::Vector< std::string > m_real_data_names
Names for real data to write to file.
Definition SpadesParticleContainer.H:238
amrex::MultiFab m_min_timestamp
Minimum timestamp in each cell for each type.
Definition SpadesParticleContainer.H:253
SpadesParticleContainer(const amrex::Vector< amrex::Geometry > &geom, const amrex::Vector< amrex::DistributionMapping > &dmap, const amrex::Vector< amrex::BoxArray > &ba, int ngrow=0)
Constructor.
Definition SpadesParticleContainer.H:34
amrex::Vector< int > m_writeflags_real
Flags for real data to write to file.
Definition SpadesParticleContainer.H:232
void print_messages(const std::string &header)
Print all the particles to screen.
Definition SpadesParticleContainer.H:658
void count_offsets()
Update the particle offsets.
Definition SpadesParticleContainer.H:167
amrex::Vector< int > m_writeflags_int
Flags for int data to write to file.
Definition SpadesParticleContainer.H:235
void initialize_state()
Initialize particle states (counts and offsets)
Definition SpadesParticleContainer.H:65
typename amrex::NeighborParticleContainer< NStructReal, NStructInt, NArrayReal, NArrayInt >::IntVector IntVector
Definition SpadesParticleContainer.H:47
void count_particles()
Update the particle counts.
Definition SpadesParticleContainer.H:123
amrex::Long total_count(const int typ) const
Get the total number of particles of typ.
Definition SpadesParticleContainer.H:116
std::string m_sort_type
Sort type.
Definition SpadesParticleContainer.H:244
typename amrex::NeighborParticleContainer< NStructReal, NStructInt, NArrayReal, NArrayInt >::RealVector RealVector
Definition SpadesParticleContainer.H:52
typename amrex::NeighborParticleContainer< NStructReal, NStructInt, NArrayReal, NArrayInt >::ParticleType ParticleType
Definition SpadesParticleContainer.H:37
amrex::iMultiFab m_counts
Count of particles in each cell.
Definition SpadesParticleContainer.H:247
static constexpr int LEV
Level index.
Definition SpadesParticleContainer.H:183
amrex::Long max_count(const int typ) const
Get the maximum number of particles of typ.
Definition SpadesParticleContainer.H:138
void write_plot_file_impl(const std::string &plt_filename, const std::string &name)
Write the particles to file (implementation)
Definition SpadesParticleContainer.H:760
void check_sort(const amrex::MFIter &mfi)
Check the result of the sort operation.
Definition SpadesParticleContainer.H:236
typename amrex::NeighborParticleContainer< NStructReal, NStructInt, NArrayReal, NArrayInt >::ParticleTileType ParticleTileType
Definition SpadesParticleContainer.H:42
virtual void sort()=0
Sort the particles.
void check_sort_type(const std::string &sort_type)
Check valid sort type.
Definition SpadesParticleContainer.H:207
SpadesParticleContainer(amrex::AmrParGDB *par_gdb, int ngrow=0)
Constructor.
Definition SpadesParticleContainer.H:9
int ngrow() const
Number of grow cells.
Definition SpadesParticleContainer.H:180
const amrex::iMultiFab & offsets() const
Get the particle offsets.
Definition SpadesParticleContainer.H:109
virtual void initialize_variable_names()=0
Initialize variable names.
void encoded_sort_impl()
Encoded sort the particles implementation.
Definition SpadesParticleContainer.H:423
void nonencoded_sort_impl(const CompareFunctor &compare)
Non-encoded sort the particles implementation.
int m_ngrow
Number of grow cells.
Definition SpadesParticleContainer.H:229
const amrex::iMultiFab & counts() const
Get the particle counts.
Definition SpadesParticleContainer.H:103
void update_counts()
Update the particle counts and offsets.
Definition SpadesParticleContainer.H:218
void reposition_particles()
Reposition the particles inside a cell for visualization.
Definition SpadesParticleContainer.H:691
amrex::Long min_count(const int typ) const
Get the minimum number of particles of typ.
Definition SpadesParticleContainer.H:127
SPADES particles.
Definition EntityData.H:7
Particle operations.
Definition ParticleOps.H:14