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

SPADES API: /home/runner/work/spades/spades/Source/MessageParticleContainer.H Source File
SPADES API
MessageParticleContainer.H
Go to the documentation of this file.
1#ifndef MESSAGEPARTICLECONTAINER_H
2#define MESSAGEPARTICLECONTAINER_H
3#include <AMReX.H>
4#include <AMReX_AmrCore.H>
5#include <AMReX_AmrParGDB.H>
6#include <AMReX_Random.H>
8#include "MessageData.H"
9#include "MessageOps.H"
10#include "Utilities.H"
11
15namespace spades::particles {
16
20 MessageTypes::NTYPES,
21 0,
22 0,
23 MessageRealData::ncomps,
24 MessageIntData::ncomps>
25{
26public:
31 static std::string identifier() { return "messages"; }
32
38 explicit MessageParticleContainer(amrex::AmrParGDB* par_gdb, int ngrow = 0);
39
48 const amrex::Vector<amrex::Geometry>& geom,
49 const amrex::Vector<amrex::DistributionMapping>& dmap,
50 const amrex::Vector<amrex::BoxArray>& ba,
51 int ngrow = 0);
52
57 template <typename Model>
58 void initialize_messages(const Model& model);
59
61 void sort() override;
62
64 void update_undefined();
65
67 void resolve_pairs();
68
70 amrex::Real compute_gvt();
71
76 void garbage_collect(const amrex::Real gvt);
77
78 void write_plot_file(const std::string& plt_filename) override
79 {
80 write_plot_file_impl(plt_filename, identifier());
81 };
82
84 void read_parameters() override;
85
87 void initialize_variable_names() override;
88};
89
90template <typename Model>
92{
93 BL_PROFILE("spades::MessageParticleContainer::initialize_messages()");
94
95 const auto& plo = Geom(LEV).ProbLoArray();
96 const auto& dx = Geom(LEV).CellSizeArray();
97 const auto& dom = Geom(LEV).Domain();
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
103 for (amrex::MFIter mfi = MakeMFIter(LEV); mfi.isValid(); ++mfi) {
104 DefineAndReturnParticleTile(LEV, mfi);
105 }
106
107 amrex::iMultiFab num_particles(
108 ParticleBoxArray(LEV), ParticleDistributionMap(LEV), 1, 0,
109 amrex::MFInfo());
110 amrex::iMultiFab init_offsets(
111 ParticleBoxArray(LEV), ParticleDistributionMap(LEV), 1, 0,
112 amrex::MFInfo());
113 num_particles.setVal(total_messages_per_lp);
114 init_offsets.setVal(0);
115
116 for (amrex::MFIter mfi = MakeMFIter(LEV); mfi.isValid(); ++mfi) {
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());
137 auto& pti = GetParticles(LEV)[index];
138 pti.resize(np);
139 const auto parrs = particle_arrays(pti);
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
186 for (MyParIter pti(*this, LEV); pti.isValid(); ++pti) {
187 const size_t np = pti.numParticles();
188 const auto parrs = particle_arrays(pti.GetParticleTile());
189
190 amrex::ParallelFor(np, [=] AMREX_GPU_DEVICE(long pidx) noexcept {
191 bool valid_type = false;
192 for (int typ = 0; typ < MessageTypes::NTYPES; typ++) {
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}
203
204} // namespace spades::particles
205#endif
Main SPADES message container.
Definition MessageParticleContainer.H:25
void initialize_variable_names() override
Initialize variable names.
Definition MessageParticleContainer.cpp:34
void read_parameters() override
Read user parameters.
Definition MessageParticleContainer.cpp:29
void resolve_pairs()
Resolve message pairs (remove message/anti-message pairs)
Definition MessageParticleContainer.cpp:217
amrex::Real compute_gvt()
Compute the minimum time stamp of the messages.
Definition MessageParticleContainer.cpp:350
void write_plot_file(const std::string &plt_filename) override
Write the particles to file.
Definition MessageParticleContainer.H:78
void initialize_messages(const Model &model)
Initialize the messages.
Definition MessageParticleContainer.H:91
MessageParticleContainer(amrex::AmrParGDB *par_gdb, int ngrow=0)
Constructor.
Definition MessageParticleContainer.cpp:6
static std::string identifier()
Class identifier name.
Definition MessageParticleContainer.H:31
void garbage_collect(const amrex::Real gvt)
Perform garbage collection.
Definition MessageParticleContainer.cpp:323
void sort() override
Sort the messages.
Definition MessageParticleContainer.cpp:73
void update_undefined()
Update the undefined messages.
Definition MessageParticleContainer.cpp:80
Main SPADES particle container.
Definition SpadesParticleContainer.H:35
ParticleArrays< NArrayReal, NArrayInt, ParticleType, RealVector, IntVector > particle_arrays(ParticleTileType &pti) const
Definition SpadesParticleContainer.H:218
void write_plot_file_impl(const std::string &plt_filename, const std::string &name)
Definition SpadesParticleContainer.H:760
static constexpr amrex::Real HALF
Definition Constants.H:37
SPADES particles.
Definition EntityData.H:7
@ type_id
Definition ParticleData.H:18
Functor for making a message undefined.
Definition MessageOps.H:60
@ receiver_lp
Definition MessageData.H:23
@ pair_id
Definition MessageData.H:27
@ sender_entity
Definition MessageData.H:24
@ pair_cpu
Definition MessageData.H:26
@ sender_lp
Definition MessageData.H:22
@ receiver_entity
Definition MessageData.H:25
static constexpr int NTYPES
Number of different message types.
Definition MessageData.H:46