MessageParticleContainer Class Reference
SPADES API
|
spades::particles::MessageParticleContainer Class Reference
Main SPADES message container. More...
#include <MessageParticleContainer.H>
Inheritance diagram for spades::particles::MessageParticleContainer:
Collaboration diagram for spades::particles::MessageParticleContainer:
Public Member Functions | |
MessageParticleContainer (amrex::AmrParGDB *par_gdb, int ngrow=0) | |
Constructor. | |
MessageParticleContainer (const amrex::Vector< amrex::Geometry > &geom, const amrex::Vector< amrex::DistributionMapping > &dmap, const amrex::Vector< amrex::BoxArray > &ba, int ngrow=0) | |
Constructor. | |
template<typename Model > | |
void | initialize_messages (const Model &model) |
Initialize the messages. | |
void | sort () override |
Sort the messages. | |
void | update_undefined () |
Update the undefined messages. | |
void | resolve_pairs () |
Resolve message pairs (remove message/anti-message pairs) | |
amrex::Real | compute_gvt () |
Compute the minimum time stamp of the messages. | |
void | garbage_collect (const amrex::Real gvt) |
Perform garbage collection. | |
void | write_plot_file (const std::string &plt_filename) override |
Write the particles to file. | |
void | read_parameters () override |
Read user parameters. | |
void | initialize_variable_names () override |
Initialize variable names. | |
![]() | |
SpadesParticleContainer (amrex::AmrParGDB *par_gdb, int ngrow=0) | |
Constructor. | |
SpadesParticleContainer (const amrex::Vector< amrex::Geometry > &geom, const amrex::Vector< amrex::DistributionMapping > &dmap, const amrex::Vector< amrex::BoxArray > &ba, int ngrow=0) | |
Constructor. | |
void | initialize_state () |
Initialize particle states (counts and offsets) | |
void | clear_state () |
Delete particle states (counts and offsets) | |
void | update_counts () |
Update the particle counts and offsets. | |
void | count_particles () |
Update the particle counts. | |
void | count_offsets () |
Update the particle offsets. | |
const amrex::iMultiFab & | counts () const |
Get the particle counts. | |
const amrex::iMultiFab & | offsets () const |
Get the particle offsets. | |
amrex::Long | total_count (const int typ) const |
Get the total number of particles of typ . | |
amrex::Long | min_count (const int typ) const |
Get the minimum number of particles of typ . | |
amrex::Long | max_count (const int typ) const |
Get the maximum number of particles of typ . | |
void | check_sort (const amrex::MFIter &mfi) |
Check the result of the sort operation. | |
void | sort_impl (const CompareFunctor &compare) |
Sort the particles implementation. | |
void | sort_impl (const CompareFunc &compare) |
void | nonencoded_sort_impl (const CompareFunctor &compare) |
Non-encoded sort the particles implementation. | |
void | nonencoded_sort_impl (const CompareFunc &compare) |
void | encoded_sort_impl () |
Encoded sort the particles implementation. | |
void | print_messages (const std::string &header) |
Print all the particles to screen. | |
void | reposition_particles () |
Reposition the particles inside a cell for visualization. | |
int | ngrow () const |
Number of grow cells. | |
void | write_plot_file_impl (const std::string &plt_filename, const std::string &name) |
Write the particles to file (implementation) | |
void | check_sort_type (const std::string &sort_type) |
Check valid sort type. | |
ParticleArrays< NArrayReal, NArrayInt, ParticleType, RealVector, IntVector > | particle_arrays (ParticleTileType &pti) const |
Static Public Member Functions | |
static std::string | identifier () |
Class identifier name. | |
![]() | |
static std::string | identifier () |
Class identifier name. | |
Additional Inherited Members | |
![]() | |
using | ParticleType |
using | ParticleTileType |
using | IntVector |
using | RealVector |
![]() | |
static constexpr int | LEV |
Level index. | |
![]() | |
int | m_ngrow |
Number of grow cells. | |
amrex::Vector< int > | m_writeflags_real |
Flags for real data to write to file. | |
amrex::Vector< int > | m_writeflags_int |
Flags for int data to write to file. | |
amrex::Vector< std::string > | m_real_data_names |
Names for real data to write to file. | |
amrex::Vector< std::string > | m_int_data_names |
Names for int data to write to file. | |
std::string | m_sort_type |
Sort type. | |
amrex::iMultiFab | m_counts |
Count of particles in each cell. | |
amrex::iMultiFab | m_offsets |
Offsets of particles in each cell. | |
amrex::MultiFab | m_min_timestamp |
Minimum timestamp in each cell for each type. | |
amrex::MultiFab | m_max_timestamp |
Maximum timestamp in each cell for each type. | |
Detailed Description
Main SPADES message container.
Constructor & Destructor Documentation
◆ MessageParticleContainer() [1/2]
|
explicit |
Constructor.
- Parameters
-
par_gdb [in] particle database ngrow [in] number of grow cells
10 0,
11 0,
14{}
SpadesParticleContainer(amrex::AmrParGDB *par_gdb, int ngrow=0)
Definition SpadesParticleContainer.H:9
int ngrow() const
Definition SpadesParticleContainer.H:180
static constexpr int NTYPES
Number of different message types.
Definition MessageData.H:46
◆ MessageParticleContainer() [2/2]
|
explicit |
Constructor.
- Parameters
-
geom [in] geometry dmap [in] distribution map ba [in] box array ngrow [in] number of grow cells
Member Function Documentation
◆ compute_gvt()
amrex::Real spades::particles::MessageParticleContainer::compute_gvt | ( | ) |
Compute the minimum time stamp of the messages.
351{
352 BL_PROFILE("spades::MessageParticleContainer::compute_gvt()");
353 // If this becomes a performance bottleneck it could be sped up by
354 // making a vector of just the message time stamps before the min op
355
356 amrex::Real gvt = constants::LARGE_NUM;
357
358#ifdef AMREX_USE_OMP
359#pragma omp parallel if (amrex::Gpu::notInLaunchRegion()) reduction(min : gvt)
360#endif
362 const size_t np = pti.numParticles();
364
365 amrex::Gpu::DeviceVector<amrex::Real> ts(np, constants::LARGE_NUM);
366 auto* p_ts = ts.data();
367 amrex::ParallelFor(np, [=] AMREX_GPU_DEVICE(long pidx) noexcept {
369 MessageTypes::MESSAGE) {
370 p_ts[pidx] = parrs.m_rdata[CommonRealData::timestamp][pidx];
371 }
372 });
373 gvt = std::min(amrex::Reduce::Min(np, ts.data()), gvt);
374 }
375 amrex::ParallelDescriptor::ReduceRealMin(gvt);
376
377 return gvt;
378}
ParticleArrays< NArrayReal, NArrayInt, ParticleType, RealVector, IntVector > particle_arrays(ParticleTileType &pti) const
Definition SpadesParticleContainer.H:218
static constexpr int LEV
Definition SpadesParticleContainer.H:183
static constexpr amrex::Real LARGE_NUM
A large number.
Definition Constants.H:34
static constexpr int MESSAGE
Message.
Definition MessageData.H:36
Here is the call graph for this function:
◆ garbage_collect()
void spades::particles::MessageParticleContainer::garbage_collect | ( | const amrex::Real | gvt | ) |
Perform garbage collection.
- Parameters
-
gvt [in] global virtual time
324{
325 BL_PROFILE("spades::MessageParticleContainer::garbage_collect()");
326
327#ifdef AMREX_USE_OMP
328#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
329#endif
331 const size_t np = pti.numParticles();
333
334 amrex::ParallelFor(np, [=] AMREX_GPU_DEVICE(long pidx) noexcept {
336 (parrs.m_idata[CommonIntData::type_id][pidx] !=
337 MessageTypes::UNDEFINED)) ||
338 ((parrs.m_idata[CommonIntData::type_id][pidx] ==
339 MessageTypes::CONJUGATE) &&
340 (parrs.m_rdata[MessageRealData::creation_time][pidx] < gvt))) {
341 AMREX_ASSERT(
342 parrs.m_idata[CommonIntData::type_id][pidx] !=
343 MessageTypes::MESSAGE);
344 MarkMessageUndefined()(pidx, parrs);
345 }
346 });
347 }
348}
static constexpr int UNDEFINED
Undefined message.
Definition MessageData.H:44
static constexpr int CONJUGATE
Conjugate message (unsent anti-message)
Definition MessageData.H:42
Here is the call graph for this function:
◆ identifier()
|
inlinestatic |
Class identifier name.
- Returns
- class identifier
31{ return "messages"; }
Here is the caller graph for this function:
◆ initialize_messages()
template<typename Model >
void spades::particles::MessageParticleContainer::initialize_messages | ( | const Model & | model | ) |
Initialize the messages.
- Parameters
-
model [in] model
92{
93 BL_PROFILE("spades::MessageParticleContainer::initialize_messages()");
94
98 const auto init_message_op = model.init_message_op();
99 const auto messages_per_lp = init_message_op.m_messages_per_lp;
100 const int total_messages_per_lp = 3 * messages_per_lp;
101 AMREX_ALWAYS_ASSERT(total_messages_per_lp > messages_per_lp);
102
104 DefineAndReturnParticleTile(LEV, mfi);
105 }
106
107 amrex::iMultiFab num_particles(
109 amrex::MFInfo());
110 amrex::iMultiFab init_offsets(
112 amrex::MFInfo());
113 num_particles.setVal(total_messages_per_lp);
114 init_offsets.setVal(0);
115
117 const amrex::Box& box = mfi.tilebox();
118
119 const auto ncells = static_cast<int>(box.numPts());
120 const int* in = num_particles[mfi].dataPtr();
121 int* out = init_offsets[mfi].dataPtr();
122 const auto np = amrex::Scan::PrefixSum<int>(
123 ncells, [=] AMREX_GPU_DEVICE(int i) -> int { return in[i]; },
124 [=] AMREX_GPU_DEVICE(int i, int const& xi) { out[i] = xi; },
125 amrex::Scan::Type::exclusive, amrex::Scan::retSum);
126
127 const amrex::Long pid = ParticleType::NextID();
128 ParticleType::NextID(pid + np);
129 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
130 static_cast<amrex::Long>(pid + np) < amrex::LastParticleID,
131 "Error: overflow on particle id numbers!");
132
133 const auto my_proc = amrex::ParallelDescriptor::MyProc();
134 const auto& offset_arr = init_offsets[mfi].const_array();
135 const auto& num_particles_arr = num_particles[mfi].const_array();
136 const auto index = std::make_pair(mfi.index(), mfi.LocalTileIndex());
138 pti.resize(np);
140 amrex::ParallelForRNG(
141 box, [=] AMREX_GPU_DEVICE(
142 int i, int j, int AMREX_D_PICK(, , k),
143 amrex::RandomEngine const& engine) noexcept {
144 const amrex::IntVect iv(AMREX_D_DECL(i, j, k));
145 const int start = offset_arr(iv);
146 for (int n = start; n < start + num_particles_arr(iv); n++) {
147 auto& p = parrs.m_aos[n];
148 p.id() = pid + n;
149 p.cpu() = my_proc;
150
151 MarkMessageUndefined()(n, parrs);
152 parrs.m_idata[MessageIntData::sender_lp][n] =
153 static_cast<int>(dom.index(iv));
154 parrs.m_idata[MessageIntData::sender_entity][n] = 0;
155 parrs.m_idata[MessageIntData::receiver_lp][n] =
156 static_cast<int>(dom.index(iv));
157 parrs.m_idata[MessageIntData::receiver_entity][n] = 0;
158
159 AMREX_D_TERM(
160 p.pos(0) = plo[0] + (iv[0] + constants::HALF) * dx[0];
161 , p.pos(1) = plo[1] + (iv[1] + constants::HALF) * dx[1];
162 ,
163 p.pos(2) = plo[2] + (iv[2] + constants::HALF) * dx[2];)
164
165 AMREX_D_TERM(parrs.m_idata[CommonIntData::i][n] = iv[0];
166 , parrs.m_idata[CommonIntData::j][n] = iv[1];
167 , parrs.m_idata[CommonIntData::k][n] = iv[2];)
168 }
169
170 for (int n = start; n < start + messages_per_lp; n++) {
171 init_message_op(parrs, n, engine);
172 parrs.m_idata[MessageIntData::pair_id][n] = -1;
173 parrs.m_idata[MessageIntData::pair_cpu][n] = -1;
174 }
175 });
176
177 // This is necessary
178 amrex::Gpu::streamSynchronize();
179 }
180 Redistribute();
181
182 // Sanity check all initial particles
183#ifdef AMREX_USE_OMP
184#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
185#endif
187 const size_t np = pti.numParticles();
189
190 amrex::ParallelFor(np, [=] AMREX_GPU_DEVICE(long pidx) noexcept {
191 bool valid_type = false;
193 valid_type = parrs.m_idata[CommonIntData::type_id][pidx] == typ;
194 if (valid_type) {
195 break;
196 }
197 }
198 AMREX_ASSERT(valid_type);
199 AMREX_ASSERT(parrs.m_aos[pidx].id() >= 0);
200 });
201 }
202}
Here is the call graph for this function:
◆ initialize_variable_names()
|
overridevirtual |
Initialize variable names.
35{
36 BL_PROFILE("spades::MessageParticleContainer::initialize_variable_names()");
37
42
49
50 AMREX_D_TERM(
52 m_writeflags_int[CommonIntData::i] = 1;
54 m_writeflags_int[CommonIntData::j] = 1;
56 m_writeflags_int[CommonIntData::k] = 1;)
71}
amrex::Vector< std::string > m_int_data_names
Definition SpadesParticleContainer.H:241
amrex::Vector< std::string > m_real_data_names
Definition SpadesParticleContainer.H:238
amrex::Vector< int > m_writeflags_real
Definition SpadesParticleContainer.H:232
amrex::Vector< int > m_writeflags_int
Definition SpadesParticleContainer.H:235
◆ read_parameters()
|
overridevirtual |
Read user parameters.
Reimplemented from spades::particles::SpadesParticleContainer< MessageTypes::NTYPES, 0, 0, MessageRealData::ncomps, MessageIntData::ncomps >.
30{
32}
virtual void read_parameters()
Read user parameters.
Definition SpadesParticleContainer.H:197
Here is the call graph for this function:
◆ resolve_pairs()
void spades::particles::MessageParticleContainer::resolve_pairs | ( | ) |
Resolve message pairs (remove message/anti-message pairs)
218{
219 BL_PROFILE("spades::MessageParticleContainer::resolve_pairs()");
220
221#ifdef AMREX_DEBUG
223#endif
224
225#ifdef AMREX_USE_OMP
226#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
227#endif
229 const amrex::Box& box = mfi.tilebox();
232 const auto index = std::make_pair(mfi.index(), mfi.LocalTileIndex());
235
236 amrex::ParallelFor(
237 box, [=] AMREX_GPU_DEVICE(
238 int i, int j, int AMREX_D_PICK(, , k)) noexcept {
239 const amrex::IntVect iv(AMREX_D_DECL(i, j, k));
240 const auto getter = Get(iv, cnt_arr, offsets_arr, parrs);
241
244 AMREX_ASSERT(
245 parrs.m_idata[MessageIntData::pair_id][pant_soa] != -1);
246 AMREX_ASSERT(
247 parrs.m_idata[MessageIntData::pair_cpu][pant_soa] !=
248 -1);
249 AMREX_ASSERT(
250 parrs.m_idata[MessageIntData::receiver_lp][pant_soa] ==
251 dom.index(iv));
252
253#ifdef AMREX_DEBUG
254 bool found_pair = false;
255#endif
257 m++) {
258 // This is a message that was already treated,
259 // expect it to be undefined
261 getter.assert_different(
262 m, MessageTypes::MESSAGE,
263 MessageTypes::UNDEFINED);
264 continue;
265 }
268 parrs
269 .m_idata[MessageIntData::pair_id][pant_soa]) &&
270 (parrs
271 .m_idata[MessageIntData::pair_cpu][pmsg_soa] ==
272 parrs.m_idata[MessageIntData::pair_cpu]
273 [pant_soa]) &&
274 (std::abs(
275 parrs.m_rdata[CommonRealData::timestamp]
276 [pmsg_soa] -
277 parrs.m_rdata[CommonRealData::timestamp]
278 [pant_soa]) < constants::EPS)) {
279 AMREX_ASSERT(
280 parrs.m_idata[MessageIntData::sender_lp]
281 [pmsg_soa] ==
282 parrs.m_idata[MessageIntData::sender_lp]
283 [pant_soa]);
284 AMREX_ASSERT(
285 parrs.m_idata[MessageIntData::sender_entity]
286 [pmsg_soa] ==
287 parrs.m_idata[MessageIntData::sender_entity]
288 [pant_soa]);
289 AMREX_ASSERT(
290 parrs.m_idata[MessageIntData::receiver_lp]
291 [pmsg_soa] ==
292 parrs.m_idata[MessageIntData::receiver_lp]
293 [pant_soa]);
294 AMREX_ASSERT(
295 parrs.m_idata[MessageIntData::receiver_entity]
296 [pmsg_soa] ==
297 parrs.m_idata[MessageIntData::receiver_entity]
298 [pant_soa]);
299 AMREX_ASSERT(
300 std::abs(
301 parrs
302 .m_rdata[MessageRealData::creation_time]
303 [pmsg_soa] -
304 parrs
305 .m_rdata[MessageRealData::creation_time]
306 [pant_soa]) < constants::EPS);
307 MarkMessageUndefined()(pant_soa, parrs);
308 MarkMessageUndefined()(pmsg_soa, parrs);
309#ifdef AMREX_DEBUG
310 found_pair = true;
311#endif
312 break;
313 }
314 }
315 AMREX_ASSERT(found_pair);
316 }
317 });
318 }
319
320 sort();
321}
void sort() override
Sort the messages.
Definition MessageParticleContainer.cpp:73
amrex::iMultiFab m_offsets
Definition SpadesParticleContainer.H:250
amrex::iMultiFab m_counts
Definition SpadesParticleContainer.H:247
static constexpr amrex::Real EPS
A number very close to zero.
Definition Constants.H:27
static constexpr int ANTI
Anti-message.
Definition MessageData.H:40
Here is the call graph for this function:
◆ sort()
|
overridevirtual |
Sort the messages.
74{
75 BL_PROFILE("spades::MessageParticleContainer::sort()");
76
77 sort_impl(CompareParticle());
78}
void sort_impl(const CompareFunctor &compare)
Here is the call graph for this function:
Here is the caller graph for this function:
◆ update_undefined()
void spades::particles::MessageParticleContainer::update_undefined | ( | ) |
Update the undefined messages.
81{
82 BL_PROFILE("spades::MessageParticleContainer::update_undefined()");
83
87 int n_removals = 0;
88
89 MessageParticleContainer pc_adds(
90 m_gdb->Geom(), m_gdb->DistributionMap(), m_gdb->boxArray(), ngrow());
92 pc_adds.DefineAndReturnParticleTile(LEV, mfi);
93 }
94
96 const amrex::Box& box = mfi.tilebox();
97 const auto index = std::make_pair(mfi.index(), mfi.LocalTileIndex());
102
103 // remove particles
104 const auto ncells = static_cast<int>(box.numPts());
105 amrex::Gpu::DeviceVector<int> removals(ncells, 0);
106 auto* p_removals = removals.data();
107 amrex::ParallelFor(ncells, [=] AMREX_GPU_DEVICE(long icell) noexcept {
108 const auto iv = box.atOffset(icell);
110 const int target_count = std::max(2 * msg_count, 2);
112 if (current_count > target_count) {
113 p_removals[icell] = current_count - target_count;
114 }
115 });
116
117 n_removals += amrex::Reduce::Sum(removals.size(), p_removals);
118
119 amrex::ParallelFor(
120 box, [=] AMREX_GPU_DEVICE(
121 int i, int j, int AMREX_D_PICK(, , k)) noexcept {
122 const amrex::IntVect iv(AMREX_D_DECL(i, j, k));
123 const auto idx = box.index(iv);
124 const auto np = p_removals[idx];
125 const auto getter = Get(iv, cnt_arr, offsets_arr, parrs);
126 for (int m = 0; m < np; m++) {
128 parrs.m_aos[pidx].id() = -1;
129 }
130 });
131
132 amrex::Gpu::DeviceVector<int> additions(ncells, 0);
133 auto* p_additions = additions.data();
134 amrex::ParallelFor(ncells, [=] AMREX_GPU_DEVICE(long icell) noexcept {
135 const auto iv = box.atOffset(icell);
137 const int target_count = std::max(2 * msg_count, 2);
139 if (target_count > current_count) {
140 p_additions[icell] = target_count - current_count;
141 }
142#ifdef AMREX_DEBUG
143 if (p_removals[icell] > 0) {
144 AMREX_ASSERT(p_additions[icell] == 0);
145 }
146 if (p_additions[icell] > 0) {
147 AMREX_ASSERT(p_removals[icell] == 0);
148 }
149#endif
150 });
151
152 amrex::Gpu::DeviceVector<int> update_offsets(ncells, 0);
153 auto* p_update_offsets = update_offsets.data();
154 const auto np = amrex::Scan::PrefixSum<int>(
155 ncells,
156 [=] AMREX_GPU_DEVICE(int i) -> int { return p_additions[i]; },
157 [=] AMREX_GPU_DEVICE(int i, int const& xi) {
158 p_update_offsets[i] = xi;
159 },
160 amrex::Scan::Type::exclusive, amrex::Scan::retSum);
161
162 const amrex::Long pid = ParticleType::NextID();
163 ParticleType::NextID(pid + np);
164 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
165 static_cast<amrex::Long>(pid + np) < amrex::LastParticleID,
166 "Error: overflow on particle id numbers!");
167
169 pti_adds.resize(np);
170 const auto my_proc = amrex::ParallelDescriptor::MyProc();
171 const auto parrs_adds = pc_adds.particle_arrays(pti_adds);
172 amrex::ParallelFor(ncells, [=] AMREX_GPU_DEVICE(long icell) noexcept {
173 const int start = p_update_offsets[icell];
174 const auto iv = box.atOffset(icell);
175 for (int n = start; n < start + p_additions[icell]; n++) {
176 auto& p = parrs_adds.m_aos[n];
177 p.id() = pid + n;
178 p.cpu() = my_proc;
179
180 AMREX_D_TERM(
181 p.pos(0) = plo[0] + (iv[0] + constants::HALF) * dx[0];
182 , p.pos(1) = plo[1] + (iv[1] + constants::HALF) * dx[1];
183 , p.pos(2) = plo[2] + (iv[2] + constants::HALF) * dx[2];)
184
185 MarkMessageUndefined()(n, parrs_adds);
186 parrs_adds.m_idata[MessageIntData::sender_lp][n] =
187 static_cast<int>(dom.index(iv));
188 parrs_adds.m_idata[MessageIntData::sender_entity][n] = 0;
189 parrs_adds.m_idata[MessageIntData::receiver_lp][n] =
190 static_cast<int>(dom.index(iv));
191 parrs_adds.m_idata[MessageIntData::receiver_entity][n] = 0;
192
193 AMREX_D_TERM(
194 p.pos(0) = plo[0] + (iv[0] + constants::HALF) * dx[0];
195 , p.pos(1) = plo[1] + (iv[1] + constants::HALF) * dx[1];
196 , p.pos(2) = plo[2] + (iv[2] + constants::HALF) * dx[2];)
197
198 AMREX_D_TERM(parrs_adds.m_idata[CommonIntData::i][n] = iv[0];
199 , parrs_adds.m_idata[CommonIntData::j][n] = iv[1];
200 , parrs_adds.m_idata[CommonIntData::k][n] = iv[2];)
201 }
202 });
203
204 // This is necessary
205 amrex::Gpu::streamSynchronize();
206 }
207
208 addParticles(pc_adds, true);
209
210 amrex::ParallelAllReduce::Sum<int>(
211 n_removals, amrex::ParallelContext::CommunicatorSub());
212 if (n_removals > 0) {
213 Redistribute(); // kind of a bummer
214 }
215}
MessageParticleContainer(amrex::AmrParGDB *par_gdb, int ngrow=0)
Constructor.
Definition MessageParticleContainer.cpp:6
Here is the call graph for this function:
◆ write_plot_file()
|
inlineoverridevirtual |
Write the particles to file.
- Parameters
-
plt_filename [in] file name for the plot file
79 {
81 };
static std::string identifier()
Class identifier name.
Definition MessageParticleContainer.H:31
void write_plot_file_impl(const std::string &plt_filename, const std::string &name)
Definition SpadesParticleContainer.H:760
Here is the call graph for this function:
The documentation for this class was generated from the following files:
- /home/runner/work/spades/spades/Source/MessageParticleContainer.H
- /home/runner/work/spades/spades/Source/MessageParticleContainer.cpp
Generated by