SPADES Class Reference

SPADES API: spades::SPADES Class Reference
SPADES API

Main SPADES class. More...

#include <SPADES.H>

Inheritance diagram for spades::SPADES:
[legend]
Collaboration diagram for spades::SPADES:
[legend]

Public Member Functions

 SPADES ()
 
 ~SPADES () override
 
void init_data ()
 Initializes data.
 
void evolve ()
 Advance solution to final time.
 
void MakeNewLevelFromCoarse (int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
 Make a new level.
 
void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
 Remake an existing level.
 
void ClearLevel (int lev) override
 Delete level data.
 
void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
 Make a level from scratch.
 
void ErrorEst (int, amrex::TagBoxArray &, amrex::Real, int) override
 Tag cells for refinement.
 
std::unique_ptr< amrex::MultiFab > get_field (const std::string &name, const int ngrow)
 Get a field based on a variable name.
 
void plot_file_mf ()
 Put together the MultiFab for output.
 
void initialize_state ()
 Initialize state.
 
void init_particle_containers ()
 Initialize the particle containers.
 
void process_messages ()
 Process messages.
 
void rollback ()
 Perform rollback.
 
void rollback_statistics ()
 Print rollback statistics.
 
void update_gvt ()
 Update the global virtual time.
 
void update_lbts ()
 Update the Lower Bound on Incoming Time Stamp.
 

Static Public Member Functions

static amrex::Real est_time_step ()
 Compute the time step.
 

Static Public Attributes

static constexpr int LEV {0}
 Level index.
 

Private Member Functions

void read_parameters ()
 Read parameters.
 
void compute_dt ()
 Wrapper for EstTimeStep.
 
void time_step (const amrex::Real time)
 Advance by the time step.
 
void advance (const amrex::Real time, const amrex::Real dt)
 Advance for a single time step.
 
bool check_field_existence (const std::string &name)
 Check if a field exists.
 
void set_ics ()
 Set the user defined IC functions.
 
std::string plot_file_name (const int step) const
 Get plotfile name.
 
std::string chk_file_name (const int step) const
 Get checkpoint file name.
 
amrex::Vector< std::string > plot_file_var_names () const
 Set plotfile variables names.
 
void write_plot_file ()
 Write plotfile to disk.
 
void write_checkpoint_file () const
 Write checkpoint file to disk.
 
void read_checkpoint_file ()
 Read checkpoint file from disk.
 
void write_info_file (const std::string &path) const
 Write job info to disk.
 
void write_rng_file (const std::string &path) const
 Write random number generator seed info.
 
void read_rng_file (const std::string &path) const
 Read random number generator seed info.
 
void init_rng () const
 Initialize the random number generator.
 
void write_data_file (const bool is_init) const
 Write simulation information.
 
void summary ()
 Compute summary data.
 

Static Private Member Functions

static void post_time_step ()
 Perform work after a time step.
 
static int get_field_component (const std::string &name, const amrex::Vector< std::string > &varnames)
 Get field component.
 

Private Attributes

int m_istep {0}
 Current step.
 
amrex::Real m_gvt {constants::LOW_NUM}
 Global virtual time.
 
amrex::Real m_t_new {0.0}
 New time.
 
amrex::Real m_t_old {constants::LOW_NUM}
 Old time.
 
amrex::Real m_dt {constants::LARGE_NUM}
 Time step.
 
amrex::Real m_lbts {constants::LOW_NUM}
 Lower bound on incoming time stamp.
 
amrex::Long m_ncells {0}
 Cell count.
 
amrex::Long m_ntotal_messages {0}
 Total message count.
 
amrex::Vector< amrex::Long > m_nmessages
 Message counts.
 
amrex::Vector< amrex::Long > m_min_messages
 Min of message counts.
 
amrex::Vector< amrex::Long > m_max_messages
 Max of message counts.
 
amrex::Long m_nprocessed_messages {0}
 Count of processed messages.
 
amrex::Long m_ntotal_entities {0}
 Total entity count.
 
amrex::Long m_nentities {0}
 Entity count.
 
amrex::Vector< int > m_nrollbacks
 Number of rollbacks.
 
amrex::Vector< amrex::Real > m_min_timings
 Min timings for each step.
 
amrex::Vector< amrex::Real > m_max_timings
 Max timings for each step.
 
amrex::Vector< amrex::Real > m_avg_timings
 Average timings for each step.
 
int m_max_step = std::numeric_limits<int>::max()
 Maximum number of steps.
 
amrex::Real m_stop_time = std::numeric_limits<amrex::Real>::max()
 Stop time.
 
amrex::Vector< std::string > m_spades_varnames
 All variable names for output.
 
amrex::Vector< std::string > m_state_varnames
 State variable names for output.
 
amrex::Vector< std::string > m_message_counts_varnames
 Message count names for output.
 
amrex::Vector< std::string > m_entity_counts_varnames
 Entity count names for output.
 
amrex::MultiFab m_state
 Multifabs to store the solution.
 
amrex::MultiFab m_plt_mf
 Multifabs for output.
 
const int m_state_ngrow = 0
 Number of ghost cells.
 
std::string m_restart_chkfile
 Restart file name, restart from this checkpoint if it is not empty.
 
std::string m_plot_file {"plt"}
 Plotfile prefix (optional user input)
 
int m_plot_int = -1
 Plotfile frequency (optional user input)
 
std::string m_chk_file {"chk"}
 Checkpoint prefix (optional user input)
 
int m_chk_int = -1
 Checkpoint frequency (optional user input)
 
int m_nfiles {constants::TWO_TO_EIGHT}
 Number of plot and checkpoint data files per write.
 
int m_file_name_digits {constants::FILE_NAME_DIGITS}
 Digits used in the plotfile and checkpoint file names.
 
int m_rng_file_name_digits {constants::RNG_FILE_NAME_DIGITS}
 Digits used in the rng seed file names.
 
bool m_write_messages {false}
 Boolean to output messages (optional user input)
 
bool m_write_entities {false}
 Boolean to output entities (optional user input)
 
std::string m_data_fname {"data.csv"}
 Filename for simulation data.
 
const int m_data_precision {6}
 Data precision for data output.
 
amrex::Real m_lookahead {1.0}
 Lookahead value (optional user input)
 
amrex::Real m_window_size {constants::LARGE_NUM}
 Window size for processing messages (optional user input)
 
int m_messages_per_step {1}
 Number of messages to process in each step (optional user input)
 
int m_entities_per_lp {1}
 Number of entities per logical process (optional user input)
 
int m_messages_per_lp {1}
 Number of messages per logical process (optional user input)
 
int m_seed {0}
 Random number generator seed (optional user input)
 
amrex::Real m_lambda {1.0}
 
std::unique_ptr< ic::InitializerBasem_ic_op
 Initial condition operator.
 
std::string m_ic_type
 Initial condition type (optional user input)
 
std::unique_ptr< particles::MessageParticleContainerm_message_pc
 Message particle container.
 
std::unique_ptr< particles::EntityParticleContainerm_entity_pc
 Entity article container.
 

Detailed Description

Main SPADES class.

Constructor & Destructor Documentation

◆ SPADES()

spades::SPADES::SPADES ( )
7{
8 BL_PROFILE("spades::SPADES::SPADES()");
10
11 if (max_level > 0) {
12 amrex::Abort(
13 "spades::SPADES::SPADES(): not supporting multilevel right now");
14 }
15
16 m_state_varnames.push_back("lvt");
17 m_state_varnames.push_back("rollback");
18 AMREX_ALWAYS_ASSERT(m_state_varnames.size() == constants::N_STATES);
19
23 "processed_message";
26 "conjugate_message";
28 "undefined_message";
29 for (const auto& mm : m_message_counts_varnames) {
30 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
31 !mm.empty(), "Message variable name needs to be specified");
32 }
36
41 "undefined_entity";
42 for (const auto& mm : m_entity_counts_varnames) {
43 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
44 !mm.empty(), "Entity variable name needs to be specified");
45 }
46
47 for (const auto& vname : m_state_varnames) {
48 m_spades_varnames.push_back(vname);
49 }
50 for (const auto& vname : m_message_counts_varnames) {
51 m_spades_varnames.push_back(vname);
52 }
53 for (const auto& vname : m_entity_counts_varnames) {
54 m_spades_varnames.push_back(vname);
55 }
56 AMREX_ALWAYS_ASSERT(
57 m_spades_varnames.size() ==
60
62 m_min_timings.resize(2, 0);
63 m_max_timings.resize(2, 0);
64 m_avg_timings.resize(2, 0);
65}
amrex::Vector< std::string > m_entity_counts_varnames
Entity count names for output.
Definition SPADES.H:342
amrex::Vector< amrex::Long > m_nmessages
Message counts.
Definition SPADES.H:297
amrex::Vector< std::string > m_spades_varnames
All variable names for output.
Definition SPADES.H:333
amrex::Vector< amrex::Real > m_max_timings
Max timings for each step.
Definition SPADES.H:321
amrex::Vector< int > m_nrollbacks
Number of rollbacks.
Definition SPADES.H:315
amrex::Vector< amrex::Real > m_avg_timings
Average timings for each step.
Definition SPADES.H:324
amrex::Vector< amrex::Long > m_min_messages
Min of message counts.
Definition SPADES.H:300
amrex::Vector< amrex::Real > m_min_timings
Min timings for each step.
Definition SPADES.H:318
amrex::Vector< std::string > m_state_varnames
State variable names for output.
Definition SPADES.H:336
void read_parameters()
Read parameters.
Definition SPADES.cpp:123
amrex::Vector< amrex::Long > m_max_messages
Max of message counts.
Definition SPADES.H:303
amrex::Vector< std::string > m_message_counts_varnames
Message count names for output.
Definition SPADES.H:339
static constexpr int MAX_ROLLBACK_ITERS
Maximum number of rollback iterations allowed.
Definition Constants.H:51
static constexpr int N_STATES
Number of states.
Definition Constants.H:10
static constexpr int BACKUP
Processed entity.
Definition EntityData.H:27
static constexpr int ENTITY
Entity.
Definition EntityData.H:25
static constexpr int NTYPES
Number of different entity types.
Definition EntityData.H:31
static constexpr int UNDEFINED
Undefined entity.
Definition EntityData.H:29
static constexpr int ANTI
Anti-message.
Definition MessageData.H:40
static constexpr int UNDEFINED
Undefined message.
Definition MessageData.H:44
static constexpr int MESSAGE
Message.
Definition MessageData.H:36
static constexpr int NTYPES
Number of different message types.
Definition MessageData.H:46
static constexpr int CONJUGATE
Conjugate message (unsent anti-message)
Definition MessageData.H:42
static constexpr int PROCESSED
Processed message.
Definition MessageData.H:38
Here is the call graph for this function:

◆ ~SPADES()

spades::SPADES::~SPADES ( )
overridedefault

Member Function Documentation

◆ advance()

void spades::SPADES::advance ( const amrex::Real time,
const amrex::Real dt )
private

Advance for a single time step.

Parameters
time[in] time
dt[in] time step
281{
282 BL_PROFILE("spades::SPADES::advance()");
283
284 m_message_pc->update_undefined();
285
286 update_gvt();
287 m_message_pc->garbage_collect(m_gvt);
288 m_message_pc->sort();
289
290 update_lbts();
291
292 const auto n_processed_messages_pre =
294
296 m_message_pc->Redistribute();
297 m_message_pc->sort();
298
299 rollback();
300
301 const auto n_processed_messages_post =
303
305 static_cast<int>(n_processed_messages_post - n_processed_messages_pre);
306 AMREX_ALWAYS_ASSERT(m_nprocessed_messages >= 0);
307
308 m_t_old = m_t_new; // old time is now current time (time)
309 m_t_new += dt; // new time is ahead
310}
amrex::Real m_t_old
Old time.
Definition SPADES.H:282
std::unique_ptr< particles::MessageParticleContainer > m_message_pc
Message particle container.
Definition SPADES.H:419
amrex::Real m_gvt
Global virtual time.
Definition SPADES.H:276
amrex::Real m_t_new
New time.
Definition SPADES.H:279
void update_lbts()
Update the Lower Bound on Incoming Time Stamp.
Definition SPADES.cpp:830
void process_messages()
Process messages.
Definition SPADES.cpp:359
void rollback()
Perform rollback.
Definition SPADES.cpp:539
amrex::Long m_nprocessed_messages
Count of processed messages.
Definition SPADES.H:306
void update_gvt()
Update the global virtual time.
Definition SPADES.cpp:821
Here is the call graph for this function:
Here is the caller graph for this function:

◆ check_field_existence()

bool spades::SPADES::check_field_existence ( const std::string & name)
private

Check if a field exists.

Parameters
name[in] field name
Returns
True if field of name exists
979{
980 BL_PROFILE("spades::SPADES::check_field_existence()");
981 const auto vnames = {
983 return std::any_of(vnames.begin(), vnames.end(), [=](const auto& vn) {
984 return get_field_component(name, vn) != -1;
985 });
986}
Here is the caller graph for this function:

◆ chk_file_name()

std::string spades::SPADES::chk_file_name ( const int step) const
nodiscardprivate

Get checkpoint file name.

Parameters
step[in] current time step
Returns
chkfile name
1053{
1054 return amrex::Concatenate(m_chk_file, step, m_file_name_digits);
1055}
std::string m_chk_file
Checkpoint prefix (optional user input)
Definition SPADES.H:363
int m_file_name_digits
Digits used in the plotfile and checkpoint file names.
Definition SPADES.H:372
Here is the caller graph for this function:

◆ ClearLevel()

void spades::SPADES::ClearLevel ( int lev)
override

Delete level data.

Overrides the pure virtual function in AmrCore.

Parameters
lev[in] level
957{
958 BL_PROFILE("spades::SPADES::ClearLevel()");
959 m_message_pc->clear_state();
960 m_entity_pc->clear_state();
961 m_state.clear();
962 m_plt_mf.clear();
963}
amrex::MultiFab m_state
Multifabs to store the solution.
Definition SPADES.H:345
std::unique_ptr< particles::EntityParticleContainer > m_entity_pc
Entity article container.
Definition SPADES.H:422
amrex::MultiFab m_plt_mf
Multifabs for output.
Definition SPADES.H:348

◆ compute_dt()

void spades::SPADES::compute_dt ( )
private

Wrapper for EstTimeStep.

872{
873 BL_PROFILE("spades::SPADES::compute_dt()");
874 amrex::Real dt_tmp = est_time_step();
875
876 amrex::ParallelDescriptor::ReduceRealMin(dt_tmp);
877
878 constexpr amrex::Real change_max = 1.1;
879 amrex::Real dt_0 = std::min(dt_tmp, change_max * m_dt);
880
881 // Limit dt's by the value of stop_time.
882 const amrex::Real eps = 1.e-3 * dt_0;
883 if (m_t_new + dt_0 > m_stop_time - eps) {
884 dt_0 = m_stop_time - m_t_new;
885 }
886
887 m_dt = dt_0;
888}
static amrex::Real est_time_step()
Compute the time step.
Definition SPADES.cpp:890
amrex::Real m_stop_time
Stop time.
Definition SPADES.H:330
amrex::Real m_dt
Time step.
Definition SPADES.H:285
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ErrorEst()

void spades::SPADES::ErrorEst ( int ,
amrex::TagBoxArray & ,
amrex::Real ,
int  )
inlineoverride

Tag cells for refinement.

Overrides the pure virtual function in AmrCore

112 {}

◆ est_time_step()

amrex::Real spades::SPADES::est_time_step ( )
static

Compute the time step.

Returns
the time step
891{
892 BL_PROFILE("spades::SPADES::est_time_step()");
893 return 1.0;
894}
Here is the caller graph for this function:

◆ evolve()

void spades::SPADES::evolve ( )

Advance solution to final time.

175{
176 BL_PROFILE("spades::SPADES::evolve()");
177
178 amrex::Real cur_time = m_t_new;
179 int last_plot_file_step = 0;
180 int last_chk_file_step = 0;
181
182 for (int step = m_istep; step < m_max_step && cur_time < m_stop_time;
183 ++step) {
184 compute_dt();
185 const amrex::Real start_time = amrex::ParallelDescriptor::second();
186
187 amrex::Print() << "\n==============================================="
188 "================================================="
189 << std::endl;
190 amrex::Print() << "Starting from step " << step
191 << ", advancing from time " << cur_time << " to time "
192 << cur_time + m_dt << " (dt = " << m_dt << ")"
193 << std::endl;
194
195 time_step(cur_time);
196
198
199 cur_time += m_dt;
200
201 // sync up time
202 m_t_new = cur_time;
203
204 if (m_plot_int > 0 && (step + 1) % m_plot_int == 0) {
205 last_plot_file_step = step + 1;
207 }
208
209 if (m_chk_int > 0 && (step + 1) % m_chk_int == 0) {
210 last_chk_file_step = step + 1;
212 }
213
214 const amrex::Real delta_time =
215 amrex::ParallelDescriptor::second() - start_time;
216 const auto n_processed_messages = m_nprocessed_messages;
217 const amrex::Vector<amrex::Real> timings{
218 delta_time,
219 static_cast<amrex::Real>(n_processed_messages) / delta_time};
220 for (int i = 0; i < timings.size(); i++) {
221 m_min_timings[i] = timings[i];
222 m_max_timings[i] = timings[i];
223 m_avg_timings[i] = timings[i];
224 }
225
226 amrex::ParallelReduce::Min(
227 m_min_timings.data(), static_cast<int>(m_min_timings.size()),
228 amrex::ParallelDescriptor::IOProcessorNumber(),
229 amrex::ParallelDescriptor::Communicator());
230 amrex::ParallelReduce::Max(
231 m_max_timings.data(), static_cast<int>(m_max_timings.size()),
232 amrex::ParallelDescriptor::IOProcessorNumber(),
233 amrex::ParallelDescriptor::Communicator());
234 amrex::ParallelReduce::Sum(
235 m_avg_timings.data(), static_cast<int>(m_avg_timings.size()),
236 amrex::ParallelDescriptor::IOProcessorNumber(),
237 amrex::ParallelDescriptor::Communicator());
238
239 if (amrex::ParallelDescriptor::IOProcessor()) {
240 for (auto& tm : m_avg_timings) {
241 tm /= amrex::ParallelDescriptor::NProcs();
242 }
243 }
244
245 write_data_file(false);
246 amrex::Print() << "Wallclock time for this step (min, avg, max) [s]: "
247 << m_min_timings[0] << " " << m_avg_timings[0] << " "
248 << m_max_timings[0] << std::endl;
249 amrex::Print() << "Current rate (min, avg, max) [msg/s]: "
250 << m_min_timings[1] << " " << m_avg_timings[1] << " "
251 << m_max_timings[1] << std::endl;
252
253 if (cur_time >= m_stop_time - constants::TOL * m_dt) {
254 break;
255 }
256 }
257 if (m_plot_int > 0 && m_istep > last_plot_file_step) {
259 }
260 if (m_chk_int > 0 && m_istep > last_chk_file_step) {
262 }
263}
int m_max_step
Maximum number of steps.
Definition SPADES.H:327
static void post_time_step()
Perform work after a time step.
Definition SPADES.cpp:354
void compute_dt()
Wrapper for EstTimeStep.
Definition SPADES.cpp:871
void time_step(const amrex::Real time)
Advance by the time step.
Definition SPADES.cpp:265
int m_chk_int
Checkpoint frequency (optional user input)
Definition SPADES.H:366
void write_checkpoint_file() const
Write checkpoint file to disk.
Definition SPADES.cpp:1107
void write_plot_file()
Write plotfile to disk.
Definition SPADES.cpp:1084
int m_plot_int
Plotfile frequency (optional user input)
Definition SPADES.H:360
void write_data_file(const bool is_init) const
Write simulation information.
Definition SPADES.cpp:1381
int m_istep
Current step.
Definition SPADES.H:273
static constexpr amrex::Real TOL
A tolerance.
Definition Constants.H:31
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_field()

std::unique_ptr< amrex::MultiFab > spades::SPADES::get_field ( const std::string & name,
const int ngrow )

Get a field based on a variable name.

Parameters
name[in] field name
ngrow[in] number of grow cells
Returns
the requested field
1001{
1002 BL_PROFILE("spades::SPADES::get_field()");
1003
1004 if (!check_field_existence(name)) {
1005 amrex::Abort(
1006 "spades::SPADES::get_field(): this field was not found: " + name);
1007 }
1008
1009 const int nc = 1;
1010 std::unique_ptr<amrex::MultiFab> mf = std::make_unique<amrex::MultiFab>(
1011 boxArray(LEV), DistributionMap(LEV), nc, ngrow);
1012
1013 const int srccomp_sd = get_field_component(name, m_state_varnames);
1014 if (srccomp_sd != -1) {
1015 amrex::MultiFab::Copy(*mf, m_state, srccomp_sd, 0, nc, ngrow);
1016 }
1017 const int srccomp_mid =
1019 if (srccomp_mid != -1) {
1020 auto const& msg_cnt_arrs = m_message_pc->counts().const_arrays();
1021 auto const& mf_arrs = mf->arrays();
1022 amrex::ParallelFor(
1023 *mf, mf->nGrowVect(), m_message_pc->counts().nComp(),
1024 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int n) noexcept {
1025 mf_arrs[nbx](i, j, k, n) = msg_cnt_arrs[nbx](i, j, k, n);
1026 });
1027 }
1028 const int srccomp_eid = get_field_component(name, m_entity_counts_varnames);
1029 if (srccomp_eid != -1) {
1030 auto const& msg_cnt_arrs = m_entity_pc->counts().const_arrays();
1031 auto const& mf_arrs = mf->arrays();
1032 amrex::ParallelFor(
1033 *mf, mf->nGrowVect(), m_entity_pc->counts().nComp(),
1034 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int n) noexcept {
1035 mf_arrs[nbx](i, j, k, n) = msg_cnt_arrs[nbx](i, j, k, n);
1036 });
1037 }
1038 amrex::Gpu::streamSynchronize();
1039 return mf;
1040}
bool check_field_existence(const std::string &name)
Check if a field exists.
Definition SPADES.cpp:978
static constexpr int LEV
Level index.
Definition SPADES.H:166
static int get_field_component(const std::string &name, const amrex::Vector< std::string > &varnames)
Get field component.
Definition SPADES.cpp:988
Here is the call graph for this function:

◆ get_field_component()

int spades::SPADES::get_field_component ( const std::string & name,
const amrex::Vector< std::string > & varnames )
staticprivate

Get field component.

Parameters
name[in] field name
varnames[in] vector of the field names
Returns
field component
990{
991 BL_PROFILE("spades::SPADES::get_field_component()");
992 const auto itr = std::find(varnames.begin(), varnames.end(), name);
993 if (itr != varnames.cend()) {
994 return static_cast<int>(std::distance(varnames.begin(), itr));
995 }
996 return -1;
997}
Here is the caller graph for this function:

◆ init_data()

void spades::SPADES::init_data ( )

Initializes data.

70{
71 BL_PROFILE("spades::SPADES::init_data()");
72
73 if (m_restart_chkfile.empty()) {
74
75 init_rng();
76
77 // start simulation from the beginning
78 const amrex::Real time = 0.0;
79 set_ics();
80
82
83 InitFromScratch(time);
84
85 update_gvt();
87
88 compute_dt();
89
90 if (m_chk_int > 0) {
92 }
93 } else {
94 // restart from a checkpoint
96 }
97
98 if (m_plot_int > 0) {
100 }
101
102 amrex::Print() << "Grid summary: " << std::endl;
103 if (amrex::ParallelDescriptor::IOProcessor()) {
104 printGridSummary(amrex::OutStream(), 0, finest_level);
105 }
106
107 summary();
108
109 write_data_file(true);
110}
void set_ics()
Set the user defined IC functions.
Definition SPADES.cpp:965
std::string m_restart_chkfile
Restart file name, restart from this checkpoint if it is not empty.
Definition SPADES.H:354
void summary()
Compute summary data.
Definition SPADES.cpp:312
void read_checkpoint_file()
Read checkpoint file from disk.
Definition SPADES.cpp:1194
void init_rng() const
Initialize the random number generator.
Definition SPADES.cpp:1318
void init_particle_containers()
Initialize the particle containers.
Definition SPADES.cpp:112
Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_particle_containers()

void spades::SPADES::init_particle_containers ( )

Initialize the particle containers.

113{
114 BL_PROFILE("spades::SPADES::init_particle_containers()");
116 std::make_unique<particles::MessageParticleContainer>(GetParGDB());
118 std::make_unique<particles::EntityParticleContainer>(GetParGDB());
119 m_message_pc->read_parameters();
120 m_entity_pc->read_parameters();
121}
Here is the caller graph for this function:

◆ init_rng()

void spades::SPADES::init_rng ( ) const
private

Initialize the random number generator.

1319{
1320 BL_PROFILE("spades::SPADES::init_rng()");
1321
1322 if (m_seed > 0) {
1323 amrex::InitRandom(
1324 m_seed + amrex::ParallelDescriptor::MyProc(),
1325 amrex::ParallelDescriptor::NProcs(),
1326 m_seed + amrex::ParallelDescriptor::MyProc());
1327 } else if (m_seed == 0) {
1328 auto now = std::chrono::system_clock::now();
1329 auto now_ns =
1330 std::chrono::time_point_cast<std::chrono::nanoseconds>(now);
1331 auto rand_seed = now_ns.time_since_epoch().count();
1332 amrex::ParallelDescriptor::Bcast(
1333 &rand_seed, 1, amrex::ParallelDescriptor::IOProcessorNumber());
1334 amrex::InitRandom(
1335 rand_seed + amrex::ParallelDescriptor::MyProc(),
1336 amrex::ParallelDescriptor::NProcs(),
1337 rand_seed + amrex::ParallelDescriptor::MyProc());
1338 } else {
1339 amrex::Abort("RNG seed must be non-negative");
1340 }
1341}
int m_seed
Random number generator seed (optional user input)
Definition SPADES.H:405
Here is the caller graph for this function:

◆ initialize_state()

void spades::SPADES::initialize_state ( )

Initialize state.

938{
939 BL_PROFILE("spades::SPADES::initialize_state()");
940
941 m_ic_op->initialize(Geom(LEV).data());
942
943 m_state.FillBoundary(Geom(LEV).periodicity());
944}
std::unique_ptr< ic::InitializerBase > m_ic_op
Initial condition operator.
Definition SPADES.H:413
Here is the caller graph for this function:

◆ MakeNewLevelFromCoarse()

void spades::SPADES::MakeNewLevelFromCoarse ( int lev,
amrex::Real time,
const amrex::BoxArray & ba,
const amrex::DistributionMapping & dm )
override

Make a new level.

Make a new level using provided BoxArray and DistributionMapping and fill with interpolated coarse level data. Overrides the pure virtual function in AmrCore.

Parameters
lev[in] level
time[in] time
ba[in] box array
dm[in] distribution map
901{
902 BL_PROFILE("spades::SPADES::MakeNewLevelFromCoarse()");
903 amrex::Abort("spades::SPADES::MakeNewLevelFromCoarse(): not implemented");
904}

◆ MakeNewLevelFromScratch()

void spades::SPADES::MakeNewLevelFromScratch ( int lev,
amrex::Real time,
const amrex::BoxArray & ba,
const amrex::DistributionMapping & dm )
override

Make a level from scratch.

Make a new level from scratch using provided BoxArray and DistributionMapping. Only used during initialization. Overrides the pure virtual function in AmrCore

Parameters
lev[in] level
time[in] time
ba[in] box array
dm[in] distribution map
911{
912 BL_PROFILE("spades::SPADES::MakeNewLevelFromScratch()");
913
914 m_state.define(ba, dm, constants::N_STATES, m_state_ngrow, amrex::MFInfo());
915
916 m_t_new = time;
918
919 // Initialize the data
921
922 // Update message particle container
923 m_message_pc->Define(Geom(LEV), dm, ba);
924 m_message_pc->initialize_variable_names();
925 m_message_pc->initialize_messages(m_lookahead);
926 m_message_pc->initialize_state();
927 m_message_pc->sort();
928
929 // Update entity particle container
930 m_entity_pc->Define(Geom(LEV), dm, ba);
931 m_entity_pc->initialize_variable_names();
932 m_entity_pc->initialize_entities();
933 m_entity_pc->initialize_state();
934 m_entity_pc->sort();
935}
void initialize_state()
Initialize state.
Definition SPADES.cpp:937
const int m_state_ngrow
Number of ghost cells.
Definition SPADES.H:351
amrex::Real m_lookahead
Lookahead value (optional user input)
Definition SPADES.H:390
static constexpr amrex::Real LOW_NUM
A large negative number.
Definition Constants.H:19
Here is the call graph for this function:

◆ plot_file_mf()

void spades::SPADES::plot_file_mf ( )

Put together the MultiFab for output.

1058{
1059 m_plt_mf.clear();
1060 m_plt_mf.define(
1061 boxArray(LEV), DistributionMap(LEV),
1062 static_cast<int>(plot_file_var_names().size()), 0, amrex::MFInfo());
1063 m_plt_mf.setVal(0.0);
1064 int cnt = 0;
1065 amrex::MultiFab::Copy(m_plt_mf, m_state, 0, cnt, m_state.nComp(), 0);
1066 cnt += m_state.nComp();
1067 auto const& plt_mf_arrs = m_plt_mf.arrays();
1068 auto const& msg_cnt_arrs = m_message_pc->counts().const_arrays();
1069 amrex::ParallelFor(
1070 m_plt_mf, m_plt_mf.nGrowVect(), m_message_pc->counts().nComp(),
1071 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int n) noexcept {
1072 plt_mf_arrs[nbx](i, j, k, n + cnt) = msg_cnt_arrs[nbx](i, j, k, n);
1073 });
1074 cnt += m_message_pc->counts().nComp();
1075 auto const& ent_cnt_arrs = m_entity_pc->counts().const_arrays();
1076 amrex::ParallelFor(
1077 m_plt_mf, m_plt_mf.nGrowVect(), m_entity_pc->counts().nComp(),
1078 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k, int n) noexcept {
1079 plt_mf_arrs[nbx](i, j, k, n + cnt) = ent_cnt_arrs[nbx](i, j, k, n);
1080 });
1081 amrex::Gpu::streamSynchronize();
1082}
amrex::Vector< std::string > plot_file_var_names() const
Set plotfile variables names.
Definition SPADES.cpp:1042
Here is the call graph for this function:
Here is the caller graph for this function:

◆ plot_file_name()

std::string spades::SPADES::plot_file_name ( const int step) const
nodiscardprivate

Get plotfile name.

Parameters
step[in] current time step
Returns
plotfile name
1048{
1049 return amrex::Concatenate(m_plot_file, step, m_file_name_digits);
1050}
std::string m_plot_file
Plotfile prefix (optional user input)
Definition SPADES.H:357
Here is the caller graph for this function:

◆ plot_file_var_names()

amrex::Vector< std::string > spades::SPADES::plot_file_var_names ( ) const
nodiscardprivate

Set plotfile variables names.

Returns
vector of variable names
1043{
1044 return m_spades_varnames;
1045}
Here is the caller graph for this function:

◆ post_time_step()

void spades::SPADES::post_time_step ( )
staticprivate

Perform work after a time step.

355{
356 BL_PROFILE("spades::SPADES::post_time_step()");
357}
Here is the caller graph for this function:

◆ process_messages()

void spades::SPADES::process_messages ( )

Process messages.

360{
361 BL_PROFILE("spades::SPADES::process_messages()");
362 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
363 m_state_ngrow == m_message_pc->ngrow(),
364 "Particle and state cells must be equal for now");
365
366 const auto& plo = Geom(LEV).ProbLoArray();
367 const auto& dx = Geom(LEV).CellSizeArray();
368 const auto& dom = Geom(LEV).Domain();
369 const auto& dlo = dom.smallEnd();
370 const auto& dhi = dom.bigEnd();
371 const auto lbts = m_lbts;
372 const auto lookahead = m_lookahead;
373 const auto window_size = m_window_size;
374 const auto messages_per_step = m_messages_per_step;
375 const auto lambda = m_lambda;
376 const auto entities_per_lp = m_entities_per_lp;
377
378#ifdef AMREX_USE_OMP
379#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
380#endif
381 for (amrex::MFIter mfi = m_message_pc->MakeMFIter(LEV); mfi.isValid();
382 ++mfi) {
383 const amrex::Box& box = mfi.tilebox();
384 const auto index = std::make_pair(mfi.index(), mfi.LocalTileIndex());
385 const auto& sarr = m_state.array(mfi);
386 const auto& msg_cnt_arr = m_message_pc->counts().const_array(mfi);
387 const auto& msg_offsets_arr = m_message_pc->offsets().const_array(mfi);
388 const auto& ent_cnt_arr = m_entity_pc->counts().const_array(mfi);
389 const auto& ent_offsets_arr = m_entity_pc->offsets().const_array(mfi);
390 auto& msg_pti = m_message_pc->GetParticles(LEV)[index];
391 const auto msg_parrs = m_message_pc->particle_arrays(msg_pti);
392 auto& ent_pti = m_entity_pc->GetParticles(LEV)[index];
393 const auto ent_parrs = m_entity_pc->particle_arrays(ent_pti);
394
395 amrex::ParallelForRNG(
396 box, [=] AMREX_GPU_DEVICE(
397 int i, int j, int AMREX_D_PICK(, , k),
398 amrex::RandomEngine const& engine) noexcept {
399 const amrex::IntVect iv(AMREX_D_DECL(i, j, k));
400 const auto msg_getter =
401 particles::Get(iv, msg_cnt_arr, msg_offsets_arr, msg_parrs);
402 const auto ent_getter =
403 particles::Get(iv, ent_cnt_arr, ent_offsets_arr, ent_parrs);
404
405 for (int n = 0;
406 n < amrex::min<int>(
407 msg_cnt_arr(iv, particles::MessageTypes::MESSAGE),
408 messages_per_step);
409 n++) {
410
411 const auto prcv_soa =
413 const auto ts =
414 msg_parrs.m_rdata[particles::CommonRealData::timestamp]
415 [prcv_soa];
416 if (ts >= lbts + lookahead + window_size) {
417 break;
418 }
419
420 AMREX_ASSERT(sarr(iv, constants::LVT_IDX) < ts);
421
422 const auto ent =
423 msg_parrs
425 [prcv_soa];
426 const auto pe_soa =
427 ent_getter(ent, particles::EntityTypes::ENTITY);
428 AMREX_ASSERT(
429 dom.atOffset(
430 ent_parrs.m_idata[particles::EntityIntData::owner]
431 [pe_soa]) == iv);
432 AMREX_ASSERT(
433 ent_parrs.m_rdata[particles::CommonRealData::timestamp]
434 [pe_soa] < ts);
435
436 // process the event
437 AMREX_ASSERT(
438 dom.atOffset(
439 msg_parrs
441 [prcv_soa]) == iv);
443 [prcv_soa] =
444 ent_parrs.m_rdata[particles::CommonRealData::timestamp]
445 [pe_soa];
446 ent_parrs
447 .m_rdata[particles::CommonRealData::timestamp][pe_soa] =
448 ts;
449 msg_parrs
450 .m_idata[particles::CommonIntData::type_id][prcv_soa] =
452
453 // Create a new message to send
454 const auto ent_lvt =
455 ent_parrs.m_rdata[particles::CommonRealData::timestamp]
456 [pe_soa];
457 const amrex::IntVect iv_dest(AMREX_D_DECL(
458 amrex::Random_int(dhi[0] - dlo[0] + 1, engine) + dlo[0],
459 amrex::Random_int(dhi[1] - dlo[1] + 1, engine) + dlo[1],
460 amrex::Random_int(dhi[2] - dlo[2] + 1, engine) +
461 dlo[2]));
462 const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM> pos = {
463 AMREX_D_DECL(
464 plo[0] + (iv_dest[0] + constants::HALF) * dx[0],
465 plo[1] + (iv_dest[1] + constants::HALF) * dx[1],
466 plo[2] + (iv_dest[2] + constants::HALF) * dx[2])};
467 const int rcv_ent = static_cast<int>(
468 amrex::Random_int(entities_per_lp, engine));
469 const amrex::Real next_ts =
470 ent_lvt + random_exponential(lambda, engine) +
471 lookahead;
472
473 const auto psnd_soa =
474 msg_getter(2 * n, particles::MessageTypes::UNDEFINED);
475 particles::CreateMessage()(
476 psnd_soa, msg_parrs, next_ts, pos, iv_dest,
477 static_cast<int>(dom.index(iv)), ent,
478 static_cast<int>(dom.index(iv_dest)), rcv_ent);
479 auto& prcv = msg_parrs.m_aos[prcv_soa];
480 AMREX_ALWAYS_ASSERT(
481 prcv.id() < std::numeric_limits<int>::max());
482 msg_parrs
483 .m_idata[particles::MessageIntData::pair_id][psnd_soa] =
484 static_cast<int>(prcv.id());
485 msg_parrs.m_idata[particles::MessageIntData::pair_cpu]
486 [psnd_soa] = prcv.cpu();
488 [psnd_soa] = ent_lvt;
489
490 // Create the conjugate message
491 const auto pcnj_soa = msg_getter(
493
494 // This is weird. The conjugate
495 // position is iv but the receiver_lp is
496 // still updated (we need to know who to
497 // send this to)
498 const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>
499 conj_pos = {AMREX_D_DECL(
500 plo[0] + (iv[0] + constants::HALF) * dx[0],
501 plo[1] + (iv[1] + constants::HALF) * dx[1],
502 plo[2] + (iv[2] + constants::HALF) * dx[2])};
503
504 particles::CreateMessage()(
505 pcnj_soa, msg_parrs, next_ts, conj_pos, iv,
506 static_cast<int>(dom.index(iv)), ent,
507 static_cast<int>(dom.index(iv_dest)), rcv_ent);
508 msg_parrs
509 .m_idata[particles::MessageIntData::pair_id][pcnj_soa] =
510 static_cast<int>(prcv.id());
511 msg_parrs.m_idata[particles::MessageIntData::pair_cpu]
512 [pcnj_soa] = prcv.cpu();
514 [pcnj_soa] = ent_lvt;
515 msg_parrs
516 .m_idata[particles::CommonIntData::type_id][pcnj_soa] =
518 }
519
520 // Update LVT for the logical process
521 amrex::Real max_ent_lvt = constants::LOW_NUM;
522 for (int n = 0;
523 n < ent_cnt_arr(iv, particles::EntityTypes::ENTITY); n++) {
524 const auto pe =
525 ent_getter(n, particles::EntityTypes::ENTITY);
526 const auto ent_lvt =
527 ent_parrs
529 if (ent_lvt > max_ent_lvt) {
530 max_ent_lvt = ent_lvt;
531 }
532 }
533 AMREX_ASSERT(sarr(iv, constants::LVT_IDX) <= max_ent_lvt);
534 sarr(iv, constants::LVT_IDX) = max_ent_lvt;
535 });
536 }
537}
int m_entities_per_lp
Number of entities per logical process (optional user input)
Definition SPADES.H:399
amrex::Real m_lambda
Definition SPADES.H:410
amrex::Real m_lbts
Lower bound on incoming time stamp.
Definition SPADES.H:288
int m_messages_per_step
Number of messages to process in each step (optional user input)
Definition SPADES.H:396
amrex::Real m_window_size
Window size for processing messages (optional user input)
Definition SPADES.H:393
static constexpr int LVT_IDX
Index of the local virtual time.
Definition Constants.H:13
static constexpr amrex::Real HALF
Definition Constants.H:37
amrex::Real random_exponential(const amrex::Real lambda)
Exponential distribution.
Definition Utilities.cpp:4
@ type_id
Definition ParticleData.H:18
@ timestamp
Definition ParticleData.H:12
@ owner
Definition EntityData.H:18
@ receiver_lp
Definition MessageData.H:23
@ pair_id
Definition MessageData.H:27
@ pair_cpu
Definition MessageData.H:26
@ receiver_entity
Definition MessageData.H:25
@ old_timestamp
Definition MessageData.H:12
@ creation_time
Definition MessageData.H:13
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_checkpoint_file()

void spades::SPADES::read_checkpoint_file ( )
private

Read checkpoint file from disk.

1195{
1196 BL_PROFILE("spades::SPADES::read_checkpoint_file()");
1197 const auto& varnames = m_state_varnames;
1198
1199 amrex::Print() << "Restarting from checkpoint file " << m_restart_chkfile
1200 << std::endl;
1201
1202 // Header
1203 const std::string file(m_restart_chkfile + "/Header");
1204
1205 amrex::VisMF::IO_Buffer io_buffer(amrex::VisMF::GetIOBufferSize());
1206
1207 amrex::Vector<char> file_char_ptr;
1208 amrex::ParallelDescriptor::ReadAndBcastFile(file, file_char_ptr);
1209 std::string file_char_ptr_string(file_char_ptr.dataPtr());
1210 std::istringstream is(file_char_ptr_string, std::istringstream::in);
1211
1212 std::string line;
1213 std::string word;
1214
1215 // read in title line
1216 std::getline(is, line);
1217
1218 // read in finest_level
1219 is >> finest_level;
1220 goto_next_line(is);
1221
1222 // read in array of istep
1223 std::getline(is, line);
1224 {
1225 std::istringstream lis(line);
1226 while (lis >> word) {
1227 m_istep = std::stoi(word);
1228 }
1229 }
1230
1231 // read in array of dt
1232 std::getline(is, line);
1233 {
1234 std::istringstream lis(line);
1235 while (lis >> word) {
1236 m_dt = std::stod(word);
1237 }
1238 }
1239
1240 // read in array of t_new
1241 std::getline(is, line);
1242 {
1243 std::istringstream lis(line);
1244 while (lis >> word) {
1245 m_t_new = std::stod(word);
1246 }
1247 }
1248
1249 // read in level 'lev' BoxArray from Header
1250 amrex::BoxArray ba;
1251 ba.readFrom(is);
1252 goto_next_line(is);
1253
1254 // create a distribution mapping
1255 amrex::DistributionMapping dm{ba, amrex::ParallelDescriptor::NProcs()};
1256
1257 // set BoxArray grids and DistributionMapping dmap in
1258 // AMReX_AmrMesh.H class
1259 SetBoxArray(LEV, ba);
1260 SetDistributionMap(LEV, dm);
1261
1262 // build MultiFabs
1263 const int ncomp = static_cast<int>(varnames.size());
1264 AMREX_ALWAYS_ASSERT(ncomp == constants::N_STATES);
1265 m_state.define(ba, dm, constants::N_STATES, m_state_ngrow, amrex::MFInfo());
1266
1267 // read in the MultiFab data
1268 amrex::VisMF::Read(
1269 m_state, amrex::MultiFabFileFullPrefix(
1270 LEV, m_restart_chkfile, "Level_", varnames[0]));
1271 update_gvt();
1272
1274
1276
1277 m_message_pc->initialize_variable_names();
1278 m_message_pc->Restart(m_restart_chkfile, m_message_pc->identifier());
1279 m_message_pc->initialize_state();
1280 m_message_pc->sort();
1281
1282 m_entity_pc->initialize_variable_names();
1283 m_entity_pc->Restart(m_restart_chkfile, m_entity_pc->identifier());
1284 m_entity_pc->initialize_state();
1285 m_entity_pc->sort();
1286}
void read_rng_file(const std::string &path) const
Read random number generator seed info.
Definition SPADES.cpp:1358
void goto_next_line(std::istream &is)
Skip to the next line in stream.
Definition Utilities.cpp:10
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_parameters()

void spades::SPADES::read_parameters ( )
private

Read parameters.

124{
125 BL_PROFILE("spades::SPADES::read_parameters()");
126
127 {
128 amrex::ParmParse pp;
129 pp.query("max_step", m_max_step);
130 pp.query("stop_time", m_stop_time);
131 }
132
133 {
134 amrex::ParmParse pp("amr");
135
136 pp.query("plot_file", m_plot_file);
137 pp.query("plot_int", m_plot_int);
138 pp.query("chk_file", m_chk_file);
139 pp.query("chk_int", m_chk_int);
140 pp.query("nfiles", m_nfiles);
141 pp.query("restart", m_restart_chkfile);
142 pp.query("file_name_digits", m_file_name_digits);
143 pp.query("rng_file_name_digits", m_rng_file_name_digits);
144 AMREX_ALWAYS_ASSERT(
146 std::ceil(std::log10(amrex::ParallelDescriptor::NProcs())));
147 }
148
149 {
150 amrex::ParmParse pp("spades");
151 pp.get("ic_type", m_ic_type);
152 pp.query("write_messages", m_write_messages);
153 pp.query("write_entities", m_write_entities);
154 pp.query("seed", m_seed);
155 pp.query("lookahead", m_lookahead);
156 pp.query("window_size", m_window_size);
157 pp.query("messages_per_step", m_messages_per_step);
158 pp.query("lambda", m_lambda);
159 pp.query("data_fname", m_data_fname);
160 pp.query("entities_per_lp", m_entities_per_lp);
161 pp.query("messages_per_lp", m_messages_per_lp);
162 }
163
164 // force periodic bcs
165 for (int dir = 0; dir < AMREX_SPACEDIM; dir++) {
166 if (!amrex::DefaultGeometry().isPeriodic(dir)) {
167 amrex::Abort(
168 "spades::SPADES::read_parameters(): geometry needs to be "
169 "periodic");
170 }
171 }
172}
std::string m_ic_type
Initial condition type (optional user input)
Definition SPADES.H:416
int m_messages_per_lp
Number of messages per logical process (optional user input)
Definition SPADES.H:402
int m_rng_file_name_digits
Digits used in the rng seed file names.
Definition SPADES.H:375
std::string m_data_fname
Filename for simulation data.
Definition SPADES.H:384
int m_nfiles
Number of plot and checkpoint data files per write.
Definition SPADES.H:369
bool m_write_entities
Boolean to output entities (optional user input)
Definition SPADES.H:381
bool m_write_messages
Boolean to output messages (optional user input)
Definition SPADES.H:378
Here is the caller graph for this function:

◆ read_rng_file()

void spades::SPADES::read_rng_file ( const std::string & path) const
private

Read random number generator seed info.

Parameters
path[in] path for file
1359{
1360 BL_PROFILE("spades::SPADES::read_rng_file()");
1361
1362 if (m_seed < 0) {
1363 const std::string base(path + "/rng");
1364 const std::string& fname = amrex::Concatenate(
1365 base, amrex::ParallelDescriptor::MyProc(), m_rng_file_name_digits);
1366 amrex::Vector<char> file_char_ptr;
1367 read_file(fname, file_char_ptr, true);
1368 std::string file_str(file_char_ptr.dataPtr());
1369 std::istringstream is(file_str, std::istringstream::in);
1370 amrex::RestoreRandomState(is, 1, 0);
1371 } else {
1372 amrex::Warning(
1373 "Warning: Restarting without using the saved RNG seeds. Run with "
1374 "spades.seed=-1 to use those");
1375 AMREX_ASSERT_WITH_MESSAGE(
1376 m_seed < 0,
1377 "Use spades.seed=-1 when debugging and using restart files");
1378 }
1379}
void read_file(const std::string &filename, amrex::Vector< char > &charBuf, bool bExitOnError)
Read and broadcast file to all ranks.
Definition Utilities.cpp:16
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RemakeLevel()

void spades::SPADES::RemakeLevel ( int lev,
amrex::Real time,
const amrex::BoxArray & ba,
const amrex::DistributionMapping & dm )
override

Remake an existing level.

Remake an existing level using provided BoxArray and DistributionMapping and fill with existing fine and coarse data. Overrides the pure virtual function in AmrCore.

Parameters
lev[in] level
time[in] time
ba[in] box array
dm[in] distribution map
951{
952 BL_PROFILE("spades::SPADES::RemakeLevel()");
953 amrex::Abort("spades::SPADES::RemakeLevel(): not implemented");
954}

◆ rollback()

void spades::SPADES::rollback ( )

Perform rollback.

540{
541 BL_PROFILE("spades::SPADES::rollback()");
542
543 if (m_window_size <= 0) {
544 return;
545 }
546
547 const auto& plo = Geom(LEV).ProbLoArray();
548 const auto& dx = Geom(LEV).CellSizeArray();
549 const auto& dom = Geom(LEV).Domain();
550
551 amrex::iMultiFab rollback(boxArray(LEV), DistributionMap(LEV), 1, 0);
552 rollback.setVal(1.0);
553 bool require_rollback = true;
554
555 int iter = 0;
556 m_state.setVal(0.0, constants::RLB_IDX, 1);
557 while (require_rollback && (iter < constants::MAX_ROLLBACK_ITERS)) {
558
559#ifdef AMREX_USE_OMP
560#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
561#endif
562 for (amrex::MFIter mfi = m_message_pc->MakeMFIter(LEV); mfi.isValid();
563 ++mfi) {
564 const amrex::Box& box = mfi.tilebox();
565 const auto& sarr = m_state.array(mfi);
566 const auto& rarr = rollback.array(mfi);
567 const auto& msg_cnt_arr = m_message_pc->counts().const_array(mfi);
568 const auto& msg_offsets_arr =
569 m_message_pc->offsets().const_array(mfi);
570 const auto& ent_cnt_arr = m_entity_pc->counts().const_array(mfi);
571 const auto& ent_offsets_arr =
572 m_entity_pc->offsets().const_array(mfi);
573 const auto index =
574 std::make_pair(mfi.index(), mfi.LocalTileIndex());
575 auto& msg_pti = m_message_pc->GetParticles(LEV)[index];
576 const auto msg_parrs = m_message_pc->particle_arrays(msg_pti);
577 auto& ent_pti = m_entity_pc->GetParticles(LEV)[index];
578 const auto ent_parrs = m_entity_pc->particle_arrays(ent_pti);
579
580 amrex::ParallelFor(
581 box, [=] AMREX_GPU_DEVICE(
582 int i, int j, int AMREX_D_PICK(, , k)) noexcept {
583 const amrex::IntVect iv(AMREX_D_DECL(i, j, k));
584 const auto msg_getter = particles::Get(
585 iv, msg_cnt_arr, msg_offsets_arr, msg_parrs);
586
587 const amrex::Real msg_lvt =
588 msg_cnt_arr(iv, particles::MessageTypes::MESSAGE) > 0
589 ? msg_parrs
591 [msg_getter(
592 0,
595 const amrex::Real ant_lvt =
596 msg_cnt_arr(iv, particles::MessageTypes::ANTI) > 0
597 ? msg_parrs
599 [msg_getter(
602 const amrex::Real rollback_timestamp =
603 amrex::min<amrex::Real>(msg_lvt, ant_lvt);
604
605 rarr(iv) = 0;
606 if (rollback_timestamp <= sarr(iv, constants::LVT_IDX)) {
607 AMREX_ASSERT(
608 msg_cnt_arr(
610 rarr(iv) = 1;
611 sarr(iv, constants::RLB_IDX) += 1;
612
613 for (int n = 0;
614 n < msg_cnt_arr(
616 n++) {
617 const auto pprd_soa = msg_getter(
619 if (msg_parrs.m_rdata
621 [pprd_soa] >= rollback_timestamp) {
622
623 auto& pprd = msg_parrs.m_aos[pprd_soa];
624 for (int m = 0;
625 m <
626 msg_cnt_arr(
627 iv,
629 m++) {
630
631 // This is a conjugate message that was
632 // already treated, expect it to be an anti
633 // message
634 if (!msg_getter.check(
635 m, particles::MessageTypes::
636 CONJUGATE)) {
637 msg_getter.assert_different(
638 m,
641 continue;
642 }
643
644 const auto pcnj_soa = msg_getter(
646
647 if ((pprd.id() ==
648 msg_parrs.m_idata
649 [particles::MessageIntData::
650 pair_id][pcnj_soa]) &&
651 (pprd.cpu() ==
652 msg_parrs.m_idata
653 [particles::MessageIntData::
654 pair_cpu][pcnj_soa])) {
655 AMREX_ASSERT(
656 std::abs(
657 msg_parrs.m_rdata
658 [particles::CommonRealData::
659 timestamp][pprd_soa] -
660 msg_parrs.m_rdata
661 [particles::
662 MessageRealData::
663 creation_time]
664 [pcnj_soa]) <
666 msg_parrs.m_idata
668 [pcnj_soa] =
670 const auto piv_soa = dom.atOffset(
671 msg_parrs.m_idata
672 [particles::MessageIntData::
673 receiver_lp][pcnj_soa]);
674 auto& pcnj = msg_parrs.m_aos[pcnj_soa];
675 AMREX_D_TERM(
676 pcnj.pos(0) =
677 plo[0] +
678 (piv_soa[0] + constants::HALF) *
679 dx[0];
680 , pcnj.pos(1) =
681 plo[1] + (piv_soa[1] +
683 dx[1];
684 , pcnj.pos(2) =
685 plo[2] + (piv_soa[2] +
687 dx[2];)
688 AMREX_D_TERM(
689 msg_parrs.m_idata
690 [particles::CommonIntData::i]
691 [pcnj_soa] = piv_soa[0];
692 , msg_parrs.m_idata
693 [particles::CommonIntData::j]
694 [pcnj_soa] = piv_soa[1];
695 , msg_parrs.m_idata
696 [particles::CommonIntData::k]
697 [pcnj_soa] = piv_soa[2];)
698 }
699 }
700
701 msg_parrs
703 [pprd_soa] =
705
706 // restore the state
707 const auto ent_getter = particles::Get(
708 iv, ent_cnt_arr, ent_offsets_arr,
709 ent_parrs);
710 auto pe_soa = ent_getter(
711 msg_parrs
712 .m_idata[particles::MessageIntData::
713 receiver_entity][pprd_soa],
715 ent_parrs.m_rdata[particles::CommonRealData::
716 timestamp][pe_soa] =
717 ent_parrs.m_rdata
719 timestamp][pe_soa] >
720 msg_parrs.m_rdata
722 old_timestamp][pprd_soa]
723 ? msg_parrs.m_rdata
725 old_timestamp][pprd_soa]
726 : ent_parrs.m_rdata
728 timestamp][pe_soa];
729
730 amrex::Real max_ent_lvt = constants::LOW_NUM;
731 for (int ne = 0;
732 ne <
733 ent_cnt_arr(
735 ne++) {
736 const auto pe = ent_getter(
738 const auto ent_lvt =
739 ent_parrs
741 timestamp][pe];
742 if (ent_lvt > max_ent_lvt) {
743 max_ent_lvt = ent_lvt;
744 }
745 }
746 sarr(iv, constants::LVT_IDX) =
747 sarr(iv, constants::LVT_IDX) > max_ent_lvt
748 ? max_ent_lvt
749 : sarr(iv, constants::LVT_IDX);
750 }
751 }
752 }
753 });
754 }
755
756 m_message_pc->Redistribute();
757 m_message_pc->sort();
758 require_rollback = rollback.max(0) > 0;
759 iter++;
760 }
761 if (iter == constants::MAX_ROLLBACK_ITERS) {
762 amrex::Abort(
763 "Maximum number of rollbacks reached "
764 "(constants::MAX_ROLLBACK_ITERS = " +
765 std::to_string(constants::MAX_ROLLBACK_ITERS) + ")");
766 } else {
768 }
769
770 m_message_pc->resolve_pairs();
771
772 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
774 "There should be no anti-messages left after rollback");
775}
void rollback_statistics()
Print rollback statistics.
Definition SPADES.cpp:777
static constexpr amrex::Real EPS
A number very close to zero.
Definition Constants.H:27
static constexpr amrex::Real LARGE_NUM
A large number.
Definition Constants.H:34
static constexpr int RLB_IDX
Index of the rollback counts.
Definition Constants.H:16
Here is the call graph for this function:
Here is the caller graph for this function:

◆ rollback_statistics()

void spades::SPADES::rollback_statistics ( )

Print rollback statistics.

778{
779 BL_PROFILE("spades::SPADES::rollback_statistics()");
780
781 auto max_rlb = static_cast<int>(m_state.max(constants::RLB_IDX));
782 amrex::ParallelAllReduce::Max<int>(
783 max_rlb, amrex::ParallelContext::CommunicatorSub());
784
785 amrex::Vector<amrex::Real> nrlbks(max_rlb + 1, 0);
786 AMREX_ALWAYS_ASSERT(nrlbks.size() <= m_nrollbacks.size());
787 for (int n = 0; n < nrlbks.size(); n++) {
788 nrlbks[n] = amrex::ReduceSum(
789 m_state, 0,
790 [=] AMREX_GPU_HOST_DEVICE(
791 const amrex::Box& bx,
792 const amrex::Array4<amrex::Real const>& sarr) -> amrex::Real {
793 amrex::Real nrlbk = 0;
794 amrex::Loop(
795 bx, [=, &nrlbk](
796 int i, int j, int AMREX_D_PICK(, , k)) noexcept {
797 const amrex::IntVect iv(AMREX_D_DECL(i, j, k));
798 nrlbk += (sarr(iv, constants::RLB_IDX) == n) ? 1 : 0;
799 });
800 return nrlbk;
801 });
802 }
803 amrex::ParallelDescriptor::ReduceRealSum(
804 nrlbks.data(), static_cast<int>(nrlbks.size()),
805 amrex::ParallelDescriptor::IOProcessorNumber());
806
807 if (amrex::ParallelDescriptor::IOProcessor()) {
808 m_nrollbacks.assign(m_nrollbacks.size(), 0);
809 for (int n = 0; n < nrlbks.size(); n++) {
810 m_nrollbacks[n] = static_cast<int>(nrlbks[n]);
811 }
812 const auto nta = static_cast<amrex::Long>(
813 std::accumulate(nrlbks.begin(), nrlbks.end(), 0.0));
814 const auto ntb = static_cast<amrex::Long>(
815 std::accumulate(m_nrollbacks.begin(), m_nrollbacks.end(), 0.0));
816 AMREX_ALWAYS_ASSERT(nta == boxArray(LEV).numPts());
817 AMREX_ALWAYS_ASSERT(nta == ntb);
818 }
819}
Here is the caller graph for this function:

◆ set_ics()

void spades::SPADES::set_ics ( )
private

Set the user defined IC functions.

966{
967 BL_PROFILE("spades::SPADES::set_ics()");
968 if (m_ic_type == "constant") {
969 m_ic_op = std::make_unique<ic::Initializer<ic::Constant>>(
970 ic::Constant(ic::Constant()), m_state);
971 } else {
972 amrex::Abort(
973 "spades::SPADES::set_ics(): User must specify a valid initial "
974 "condition");
975 }
976}
Here is the caller graph for this function:

◆ summary()

void spades::SPADES::summary ( )
private

Compute summary data.

313{
314 BL_PROFILE("spades::SPADES::summary()");
315
316 m_ntotal_messages = m_message_pc->TotalNumberOfParticles(LEV != 0);
317 for (int typ = 0; typ < particles::MessageTypes::NTYPES; typ++) {
318 m_nmessages[typ] = m_message_pc->total_count(typ);
319 m_min_messages[typ] = m_message_pc->min_count(typ);
320 m_max_messages[typ] = m_message_pc->max_count(typ);
321 }
322 m_ntotal_entities = m_entity_pc->TotalNumberOfParticles(LEV != 0);
324 m_ncells = CountCells(LEV);
325 if (Verbose() > 0) {
326 amrex::Print() << "Step " << m_istep << " Summary:" << std::endl;
327 amrex::Print() << " " << m_ntotal_messages << " total messages"
328 << std::endl;
329 for (int typ = 0; typ < particles::MessageTypes::NTYPES; typ++) {
330 amrex::Print() << " " << m_nmessages[typ] << " "
332 << "s (min: " << m_min_messages[typ]
333 << ", max: " << m_max_messages[typ] << ")"
334 << std::endl;
335 }
336 amrex::Print() << " " << m_nprocessed_messages
337 << " messages processed during this step" << std::endl;
338 amrex::Print() << " " << m_ntotal_entities << " total entities"
339 << std::endl;
340 amrex::Print() << " " << m_nentities << " entities" << std::endl;
341 amrex::Print() << " " << m_ncells << " logical processes" << std::endl;
342 for (int n = 0; n < m_nrollbacks.size(); n++) {
343 const auto nrlbk = m_nrollbacks[n];
344 if (nrlbk > 0) {
345 amrex::Print() << " " << nrlbk << " logical processes doing "
346 << n << " rollbacks" << std::endl;
347 }
348 }
349 }
350 AMREX_ALWAYS_ASSERT(
352}
amrex::Long m_ntotal_entities
Total entity count.
Definition SPADES.H:309
amrex::Long m_ncells
Cell count.
Definition SPADES.H:291
amrex::Long m_nentities
Entity count.
Definition SPADES.H:312
amrex::Long m_ntotal_messages
Total message count.
Definition SPADES.H:294
Here is the caller graph for this function:

◆ time_step()

void spades::SPADES::time_step ( const amrex::Real time)
private

Advance by the time step.

Parameters
time[in] time
266{
267 BL_PROFILE("spades::SPADES::time_step()");
268
269 amrex::Print() << "Performing step " << m_istep + 1;
270 amrex::Print() << " at time = " << m_t_new << ", dt = " << m_dt
271 << ", gvt = " << m_gvt << ", lbts = " << m_lbts << std::endl;
272
273 advance(time, m_dt);
274
275 ++m_istep;
276
277 summary();
278}
void advance(const amrex::Real time, const amrex::Real dt)
Advance for a single time step.
Definition SPADES.cpp:280
Here is the call graph for this function:
Here is the caller graph for this function:

◆ update_gvt()

void spades::SPADES::update_gvt ( )

Update the global virtual time.

822{
823 BL_PROFILE("spades::SPADES::update_gvt()");
824 const amrex::Real gvt = m_message_pc->compute_gvt();
825 AMREX_ALWAYS_ASSERT(gvt >= m_gvt);
826 AMREX_ALWAYS_ASSERT(gvt >= m_state.min(constants::LVT_IDX, 0));
827 m_gvt = gvt;
828}
Here is the caller graph for this function:

◆ update_lbts()

void spades::SPADES::update_lbts ( )

Update the Lower Bound on Incoming Time Stamp.

831{
832 BL_PROFILE("spades::SPADES::update_lbts()");
833
834 amrex::MultiFab lbts(boxArray(LEV), DistributionMap(LEV), 1, 0);
835 lbts.setVal(constants::LARGE_NUM);
836
837#ifdef AMREX_USE_OMP
838#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
839#endif
840 for (amrex::MFIter mfi = m_message_pc->MakeMFIter(LEV); mfi.isValid();
841 ++mfi) {
842 const amrex::Box& box = mfi.tilebox();
843 const auto& msg_cnt_arr = m_message_pc->counts().const_array(mfi);
844 const auto& msg_offsets_arr = m_message_pc->offsets().const_array(mfi);
845 const auto& lbts_arr = lbts.array(mfi);
846 const auto index = std::make_pair(mfi.index(), mfi.LocalTileIndex());
847 auto& msg_pti = m_message_pc->GetParticles(LEV)[index];
848 const auto msg_parrs = m_message_pc->particle_arrays(msg_pti);
849
850 amrex::ParallelFor(
851 box, [=] AMREX_GPU_DEVICE(
852 int i, int j, int AMREX_D_PICK(, , k)) noexcept {
853 const amrex::IntVect iv(AMREX_D_DECL(i, j, k));
854
855 if (msg_cnt_arr(iv, particles::MessageTypes::MESSAGE) > 0) {
856 const auto msg_getter = particles::Get(
857 iv, msg_cnt_arr, msg_offsets_arr, msg_parrs);
858 const auto prcv_soa =
860 lbts_arr(iv) =
861 msg_parrs.m_rdata[particles::CommonRealData::timestamp]
862 [prcv_soa];
863 }
864 });
865 }
866
867 m_lbts = lbts.min(0, 0);
868 AMREX_ALWAYS_ASSERT(m_lbts >= m_gvt);
869}
Here is the caller graph for this function:

◆ write_checkpoint_file()

void spades::SPADES::write_checkpoint_file ( ) const
private

Write checkpoint file to disk.

1108{
1109 BL_PROFILE("spades::SPADES::write_checkpoint_file()");
1110 const auto& varnames = m_state_varnames;
1111
1112 // chk00010 write a checkpoint file with this root directory
1113 // chk00010/Header this contains information you need to save (e.g.,
1114 // finest_level, t_new, etc.) and also
1115 // the BoxArrays at each level
1116 // chk00010/Level_0/
1117 // chk00010/Level_1/
1118 // etc. these subdirectories will hold the MultiFab data
1119 // at each level of refinement
1120
1121 const std::string& checkpointname = chk_file_name(m_istep);
1122
1123 amrex::Print() << "Writing checkpoint file " << checkpointname
1124 << " at time " << m_t_new << std::endl;
1125
1126 const int nlevels = finest_level + 1;
1127
1128 // ---- prebuild a hierarchy of directories
1129 // ---- dirName is built first. if dirName exists, it is renamed. then
1130 // build
1131 // ---- dirName/subDirPrefix_0 .. dirName/subDirPrefix_nlevels-1
1132 // ---- if callBarrier is true, call ParallelDescriptor::Barrier()
1133 // ---- after all directories are built
1134 // ---- ParallelDescriptor::IOProcessor() creates the directories
1135 amrex::PreBuildDirectorHierarchy(checkpointname, "Level_", nlevels, true);
1136
1137 // write Header file
1138 if (amrex::ParallelDescriptor::IOProcessor()) {
1139
1140 const std::string header_file_name(checkpointname + "/Header");
1141 amrex::VisMF::IO_Buffer io_buffer(amrex::VisMF::IO_Buffer_Size);
1142 std::ofstream header_file;
1143 header_file.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
1144 header_file.open(
1145 header_file_name.c_str(),
1146 std::ofstream::out | std::ofstream::trunc | std::ofstream::binary);
1147
1148 if (!header_file.good()) {
1149 amrex::FileOpenFailed(header_file_name);
1150 }
1151
1152 header_file.precision(constants::HEADER_FILE_PRECISION);
1153
1154 // write out title line
1155 header_file << "Checkpoint file for SPADES\n";
1156
1157 // write out finest_level
1158 header_file << finest_level << "\n";
1159
1160 // write out istep
1161 header_file << m_istep << " ";
1162 header_file << "\n";
1163
1164 // write out dt
1165 header_file << m_dt << " ";
1166 header_file << "\n";
1167
1168 // write out t_new
1169 header_file << m_t_new << " ";
1170 header_file << "\n";
1171
1172 // write the BoxArray at each level
1173
1174 boxArray(LEV).writeOn(header_file);
1175 header_file << '\n';
1176 header_file.close();
1177 }
1178
1179 // write the MultiFab data to, e.g., chk00010/Level_0/
1180
1181 amrex::VisMF::Write(
1182 m_state, amrex::MultiFabFileFullPrefix(
1183 LEV, checkpointname, "Level_", varnames[0]));
1184
1185 m_message_pc->Checkpoint(checkpointname, m_message_pc->identifier());
1186
1187 m_entity_pc->Checkpoint(checkpointname, m_entity_pc->identifier());
1188
1189 write_info_file(checkpointname);
1190
1191 write_rng_file(checkpointname);
1192}
void write_rng_file(const std::string &path) const
Write random number generator seed info.
Definition SPADES.cpp:1343
std::string chk_file_name(const int step) const
Get checkpoint file name.
Definition SPADES.cpp:1052
void write_info_file(const std::string &path) const
Write job info to disk.
Definition SPADES.cpp:1288
static constexpr int HEADER_FILE_PRECISION
Precision to use in header files.
Definition Constants.H:48
Here is the call graph for this function:
Here is the caller graph for this function:

◆ write_data_file()

void spades::SPADES::write_data_file ( const bool is_init) const
private

Write simulation information.

Parameters
is_init[in] boolean indicating if this is the initializing step
1382{
1383 BL_PROFILE("spades::SPADES::write_data_file()");
1384
1385 if (amrex::ParallelDescriptor::IOProcessor()) {
1386
1387 amrex::Print() << "Writing simulation data to " << m_data_fname
1388 << " at time " << m_t_new << std::endl;
1389
1390 if (is_init) {
1391 std::ofstream fh(m_data_fname.c_str(), std::ios::out);
1392 fh.precision(m_data_precision);
1393 fh << "step,gvt,lbts,lps,total_messages,";
1394 for (int typ = 0; typ < particles::MessageTypes::NTYPES; typ++) {
1395 fh << m_message_counts_varnames[typ] + "s,"
1396 << "min_" + m_message_counts_varnames[typ] + "s,"
1397 << "max_" + m_message_counts_varnames[typ] + "s,";
1398 }
1399 fh << "total_entities,"
1400 "entities,"
1401 "processed_messages_per_step,"
1402 "min_time,avg_time,max_time,min_rate,avg_"
1403 "rate,"
1404 "max_rate";
1405 for (int i = 0; i < m_nrollbacks.size(); i++) {
1406 fh << ",rollback_" << i;
1407 }
1408 fh << "\n";
1409 fh.close();
1410 }
1411
1412 std::ofstream fh(m_data_fname.c_str(), std::ios::app);
1413 fh.precision(m_data_precision);
1414 if (!fh.good()) {
1415 amrex::FileOpenFailed(m_data_fname);
1416 }
1417
1418 const auto n_processed_messages = m_nprocessed_messages;
1419 fh << m_t_new << "," << m_gvt << "," << m_lbts << "," << m_ncells << ","
1420 << m_ntotal_messages << ",";
1421 for (int typ = 0; typ < particles::MessageTypes::NTYPES; typ++) {
1422 fh << m_nmessages[typ] << "," << m_min_messages[typ] << ","
1423 << m_max_messages[typ] << ",";
1424 }
1425 fh << n_processed_messages << "," << m_ntotal_entities << ","
1426 << m_nentities << "," << m_min_timings[0] << "," << m_avg_timings[0]
1427 << "," << m_max_timings[0] << "," << m_min_timings[1] << ","
1428 << m_avg_timings[1] << "," << m_max_timings[1];
1429
1430 for (const auto& rlbk : m_nrollbacks) {
1431 fh << "," << rlbk;
1432 }
1433 fh << "\n";
1434 fh.close();
1435 }
1436}
const int m_data_precision
Data precision for data output.
Definition SPADES.H:387
Here is the caller graph for this function:

◆ write_info_file()

void spades::SPADES::write_info_file ( const std::string & path) const
private

Write job info to disk.

Parameters
path[in] path for file
1289{
1290 BL_PROFILE("spades::SPADES::write_info_file()");
1291 if (!amrex::ParallelDescriptor::IOProcessor()) {
1292 return;
1293 }
1294
1295 const std::string dash_line = "\n" + std::string(78, '-') + "\n";
1296 const std::string fname(path + "/spades_info");
1297 std::ofstream fh(fname.c_str(), std::ios::out);
1298 if (!fh.good()) {
1299 amrex::FileOpenFailed(fname);
1300 }
1301
1302 fh << dash_line << "Grid information: " << std::endl;
1303
1304 fh << " Level: " << LEV << "\n"
1305 << " num. boxes = " << boxArray().size() << "\n"
1306 << " maximum zones = ";
1307
1308 for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
1309 fh << Geom(LEV).Domain().length(dir) << " ";
1310 }
1311 fh << "\n";
1312
1313 fh << dash_line << "Input file parameters: " << std::endl;
1314 amrex::ParmParse::dumpTable(fh, true);
1315 fh.close();
1316}
Here is the caller graph for this function:

◆ write_plot_file()

void spades::SPADES::write_plot_file ( )
private

Write plotfile to disk.

1085{
1086 BL_PROFILE("spades::SPADES::write_plot_file()");
1087 const std::string& plotfilename = plot_file_name(m_istep);
1088 plot_file_mf();
1089 const auto& varnames = plot_file_var_names();
1090
1091 amrex::Print() << "Writing plot file " << plotfilename << " at time "
1092 << m_t_new << std::endl;
1093
1094 amrex::VisMF::SetNOutFiles(m_nfiles);
1095 amrex::WriteSingleLevelPlotfile(
1096 plotfilename, m_plt_mf, varnames, Geom(LEV), m_t_new, m_istep);
1097 if (m_write_messages) {
1098 m_message_pc->write_plot_file(plotfilename);
1099 }
1100 if (m_write_entities) {
1101 m_entity_pc->write_plot_file(plotfilename);
1102 }
1103
1104 write_info_file(plotfilename);
1105}
std::string plot_file_name(const int step) const
Get plotfile name.
Definition SPADES.cpp:1047
void plot_file_mf()
Put together the MultiFab for output.
Definition SPADES.cpp:1057
Here is the call graph for this function:
Here is the caller graph for this function:

◆ write_rng_file()

void spades::SPADES::write_rng_file ( const std::string & path) const
private

Write random number generator seed info.

Parameters
path[in] path for file
1344{
1345 BL_PROFILE("spades::SPADES::write_rng_file()");
1346
1347 const std::string base(path + "/rng");
1348 const std::string& fname = amrex::Concatenate(
1349 base, amrex::ParallelDescriptor::MyProc(), m_rng_file_name_digits);
1350 std::ofstream fh(fname.c_str(), std::ios::out);
1351 if (!fh.good()) {
1352 amrex::FileOpenFailed(fname);
1353 }
1354 amrex::SaveRandomState(fh);
1355 fh.close();
1356}
Here is the caller graph for this function:

Member Data Documentation

◆ LEV

int spades::SPADES::LEV {0}
staticconstexpr

Level index.

166{0};

◆ m_avg_timings

amrex::Vector<amrex::Real> spades::SPADES::m_avg_timings
private

Average timings for each step.

◆ m_chk_file

std::string spades::SPADES::m_chk_file {"chk"}
private

Checkpoint prefix (optional user input)

363{"chk"};

◆ m_chk_int

int spades::SPADES::m_chk_int = -1
private

Checkpoint frequency (optional user input)

◆ m_data_fname

std::string spades::SPADES::m_data_fname {"data.csv"}
private

Filename for simulation data.

384{"data.csv"};

◆ m_data_precision

const int spades::SPADES::m_data_precision {6}
private

Data precision for data output.

387{6};

◆ m_dt

amrex::Real spades::SPADES::m_dt {constants::LARGE_NUM}
private

Time step.

◆ m_entities_per_lp

int spades::SPADES::m_entities_per_lp {1}
private

Number of entities per logical process (optional user input)

399{1};

◆ m_entity_counts_varnames

amrex::Vector<std::string> spades::SPADES::m_entity_counts_varnames
private

Entity count names for output.

◆ m_entity_pc

std::unique_ptr<particles::EntityParticleContainer> spades::SPADES::m_entity_pc
private

Entity article container.

◆ m_file_name_digits

int spades::SPADES::m_file_name_digits {constants::FILE_NAME_DIGITS}
private

Digits used in the plotfile and checkpoint file names.

static constexpr int FILE_NAME_DIGITS
Number of digits to use in output file names.
Definition Constants.H:42

◆ m_gvt

amrex::Real spades::SPADES::m_gvt {constants::LOW_NUM}
private

Global virtual time.

◆ m_ic_op

std::unique_ptr<ic::InitializerBase> spades::SPADES::m_ic_op
private

Initial condition operator.

◆ m_ic_type

std::string spades::SPADES::m_ic_type
private

Initial condition type (optional user input)

◆ m_istep

int spades::SPADES::m_istep {0}
private

Current step.

273{0};

◆ m_lambda

amrex::Real spades::SPADES::m_lambda {1.0}
private

Width of exponential distribution (optional user input) smaller values of lambda = larger variance in random values larger values of lambda = smaller variance in random values

410{1.0};

◆ m_lbts

amrex::Real spades::SPADES::m_lbts {constants::LOW_NUM}
private

Lower bound on incoming time stamp.

◆ m_lookahead

amrex::Real spades::SPADES::m_lookahead {1.0}
private

Lookahead value (optional user input)

390{1.0};

◆ m_max_messages

amrex::Vector<amrex::Long> spades::SPADES::m_max_messages
private

Max of message counts.

◆ m_max_step

int spades::SPADES::m_max_step = std::numeric_limits<int>::max()
private

Maximum number of steps.

◆ m_max_timings

amrex::Vector<amrex::Real> spades::SPADES::m_max_timings
private

Max timings for each step.

◆ m_message_counts_varnames

amrex::Vector<std::string> spades::SPADES::m_message_counts_varnames
private

Message count names for output.

◆ m_message_pc

std::unique_ptr<particles::MessageParticleContainer> spades::SPADES::m_message_pc
private

Message particle container.

◆ m_messages_per_lp

int spades::SPADES::m_messages_per_lp {1}
private

Number of messages per logical process (optional user input)

402{1};

◆ m_messages_per_step

int spades::SPADES::m_messages_per_step {1}
private

Number of messages to process in each step (optional user input)

396{1};

◆ m_min_messages

amrex::Vector<amrex::Long> spades::SPADES::m_min_messages
private

Min of message counts.

◆ m_min_timings

amrex::Vector<amrex::Real> spades::SPADES::m_min_timings
private

Min timings for each step.

◆ m_ncells

amrex::Long spades::SPADES::m_ncells {0}
private

Cell count.

291{0};

◆ m_nentities

amrex::Long spades::SPADES::m_nentities {0}
private

Entity count.

312{0};

◆ m_nfiles

int spades::SPADES::m_nfiles {constants::TWO_TO_EIGHT}
private

Number of plot and checkpoint data files per write.

static constexpr int TWO_TO_EIGHT
Definition Constants.H:39

◆ m_nmessages

amrex::Vector<amrex::Long> spades::SPADES::m_nmessages
private

Message counts.

◆ m_nprocessed_messages

amrex::Long spades::SPADES::m_nprocessed_messages {0}
private

Count of processed messages.

306{0};

◆ m_nrollbacks

amrex::Vector<int> spades::SPADES::m_nrollbacks
private

Number of rollbacks.

◆ m_ntotal_entities

amrex::Long spades::SPADES::m_ntotal_entities {0}
private

Total entity count.

309{0};

◆ m_ntotal_messages

amrex::Long spades::SPADES::m_ntotal_messages {0}
private

Total message count.

294{0};

◆ m_plot_file

std::string spades::SPADES::m_plot_file {"plt"}
private

Plotfile prefix (optional user input)

357{"plt"};

◆ m_plot_int

int spades::SPADES::m_plot_int = -1
private

Plotfile frequency (optional user input)

◆ m_plt_mf

amrex::MultiFab spades::SPADES::m_plt_mf
private

Multifabs for output.

◆ m_restart_chkfile

std::string spades::SPADES::m_restart_chkfile
private

Restart file name, restart from this checkpoint if it is not empty.

◆ m_rng_file_name_digits

int spades::SPADES::m_rng_file_name_digits {constants::RNG_FILE_NAME_DIGITS}
private

Digits used in the rng seed file names.

static constexpr int RNG_FILE_NAME_DIGITS
Number of digits to use in random seed file names.
Definition Constants.H:45

◆ m_seed

int spades::SPADES::m_seed {0}
private

Random number generator seed (optional user input)

405{0};

◆ m_spades_varnames

amrex::Vector<std::string> spades::SPADES::m_spades_varnames
private

All variable names for output.

◆ m_state

amrex::MultiFab spades::SPADES::m_state
private

Multifabs to store the solution.

◆ m_state_ngrow

const int spades::SPADES::m_state_ngrow = 0
private

Number of ghost cells.

◆ m_state_varnames

amrex::Vector<std::string> spades::SPADES::m_state_varnames
private

State variable names for output.

◆ m_stop_time

amrex::Real spades::SPADES::m_stop_time = std::numeric_limits<amrex::Real>::max()
private

Stop time.

◆ m_t_new

amrex::Real spades::SPADES::m_t_new {0.0}
private

New time.

279{0.0};

◆ m_t_old

amrex::Real spades::SPADES::m_t_old {constants::LOW_NUM}
private

Old time.

◆ m_window_size

amrex::Real spades::SPADES::m_window_size {constants::LARGE_NUM}
private

Window size for processing messages (optional user input)

◆ m_write_entities

bool spades::SPADES::m_write_entities {false}
private

Boolean to output entities (optional user input)

381{false};

◆ m_write_messages

bool spades::SPADES::m_write_messages {false}
private

Boolean to output messages (optional user input)

378{false};

The documentation for this class was generated from the following files:
  • /home/runner/work/spades/spades/Source/SPADES.H
  • /home/runner/work/spades/spades/Source/SPADES.cpp