Skip to content

Commit 957262b

Browse files
Close BTD Snapshot when last slice is filled in. (#2516)
* Add logic to close snapshot when the last slice, with index 0, is filled * fix eol * call flush as soon as last valid Zslice is full * clean up unwanted print statements * add comments in code and doxygen comments * remove commented line * fix eol * add doxygen comments and describe the new parameter introduced in the function * fix typo from Axel Co-authored-by: Axel Huebl <[email protected]> * Update Source/Diagnostics/BTDiagnostics.cpp * Update Source/Diagnostics/BTDiagnostics.cpp * snake_style * fix eol * modified code so we dont use level info when setting snapshot full info Co-authored-by: Axel Huebl <[email protected]>
1 parent e459e3a commit 957262b

File tree

5 files changed

+80
-27
lines changed

5 files changed

+80
-27
lines changed

Source/Diagnostics/BTDiagnostics.H

Lines changed: 15 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -181,6 +181,15 @@ private:
181181
* will be equal to the predicted m_max_buffer_multifabs for each snapshot.
182182
*/
183183
amrex::Vector<int> m_max_buffer_multifabs;
184+
/** Vector of integers to indicate if the snapshot is full (0 not full, 1 for full).
185+
* If the snapshot is full, then the snapshot files are closed.
186+
*/
187+
amrex::Vector<int> m_snapshot_full;
188+
/** Vector of integers to store if the last valid slice in the lab-frame
189+
* is being filled. When the value is 1, then the buffer is flushed, and
190+
* m_snapshot_full is set to 1 for that snapshot index.
191+
*/
192+
amrex::Vector<int> m_lastValidZSlice;
184193
/** Vector of counters tracking number of times the buffer of multifab is
185194
* flushed out and emptied before being refilled again for each snapshot */
186195
amrex::Vector<int> m_buffer_flush_counter;
@@ -276,6 +285,12 @@ private:
276285
void IncrementBufferFlushCounter(int i_buffer) {
277286
m_buffer_flush_counter[i_buffer]++;
278287
}
288+
/** Set Snapshot full status to 1 if the last valid zslice, m_lastValidZSlice,
289+
* for the ith snapshot index, given by, i_buffer, is filled.
290+
*
291+
* \param[in] i_buffer snapshot index
292+
*/
293+
void SetSnapshotFullStatus (const int i_buffer);
279294
/** Vector of field-data stored in the cell-centered multifab, m_cell_centered_data.
280295
* All the fields are stored regardless of the specific fields to plot selected
281296
* by the user.

Source/Diagnostics/BTDiagnostics.cpp

Lines changed: 28 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -83,8 +83,13 @@ void BTDiagnostics::DerivedInitData ()
8383
m_buffer_flush_counter.resize(m_num_buffers);
8484
// allocate vector of geometry objects corresponding to each snapshot
8585
m_geom_snapshot.resize( m_num_buffers );
86+
m_snapshot_full.resize( m_num_buffers );
87+
m_lastValidZSlice.resize( m_num_buffers );
8688
for (int i = 0; i < m_num_buffers; ++i) {
8789
m_geom_snapshot[i].resize(nmax_lev);
90+
// initialize snapshot full boolean to false
91+
m_snapshot_full[i] = 0;
92+
m_lastValidZSlice[i] = 0;
8893
}
8994

9095
for (int lev = 0; lev < nmax_lev; ++lev) {
@@ -156,9 +161,14 @@ BTDiagnostics::DoDump (int step, int i_buffer, bool force_flush)
156161
// timestep < 0, i.e., at initialization time when step == -1
157162
if (step < 0 )
158163
return false;
159-
// buffer for this lab snapshot is full, time to dump it and continue
160-
// to collect more slices afterwards
161-
else if (buffer_full(i_buffer))
164+
// Do not call dump if the snapshot is already full and the files are closed.
165+
else if (m_snapshot_full[i_buffer] == 1)
166+
return false;
167+
// If buffer for this lab snapshot is full then dump it and continue to collect
168+
// slices afterwards; or
169+
// If last z-slice in the lab-frame snapshot is filled, call dump to
170+
// write the buffer and close the file.
171+
else if (buffer_full(i_buffer) || m_lastValidZSlice[i_buffer] == 1)
162172
return true;
163173
// forced: at the end of the simulation
164174
// empty: either lab snapshot was already fully written and buffer was reset
@@ -442,9 +452,12 @@ BTDiagnostics::PrepareFieldDataForOutput ()
442452
i_buffer, ZSliceInDomain,
443453
m_current_z_boost[i_buffer],
444454
m_buffer_box[i_buffer],
445-
k_index_zlab(i_buffer, lev), m_max_box_size );
455+
k_index_zlab(i_buffer, lev), m_max_box_size,
456+
m_snapshot_full[i_buffer] );
446457

447458
if (ZSliceInDomain) ++m_buffer_counter[i_buffer];
459+
// when the 0th z-index is filled, then set lastValidZSlice to 1
460+
if (k_index_zlab(i_buffer, lev) == 0) m_lastValidZSlice[i_buffer] = 1;
448461
}
449462
}
450463
}
@@ -465,15 +478,23 @@ BTDiagnostics::k_index_zlab (int i_buffer, int lev)
465478
amrex::Real prob_domain_zmin_lab = m_prob_domain_lab[i_buffer].lo( m_moving_window_dir );
466479
amrex::IntVect ref_ratio = amrex::IntVect(1);
467480
if (lev > 0 ) ref_ratio = WarpX::RefRatio(lev-1);
468-
int k_lab = static_cast<unsigned>( (
481+
int k_lab = static_cast<int>(floor (
469482
( m_current_z_lab[i_buffer]
470483
- (prob_domain_zmin_lab + 0.5*dz_lab(warpx.getdt(lev), ref_ratio[m_moving_window_dir]) ) )
471484
/ dz_lab( warpx.getdt(lev), ref_ratio[m_moving_window_dir] )
472485
) );
473486
return k_lab;
474487
}
475488

489+
void
490+
BTDiagnostics::SetSnapshotFullStatus (const int i_buffer)
491+
{
492+
if (m_snapshot_full[i_buffer] == 1) return;
493+
// if the last valid z-index of the snapshot, which is 0, is filled, then
494+
// set the snapshot full integer to 1
495+
if (m_lastValidZSlice[i_buffer] == 1) m_snapshot_full[i_buffer] = 1;
476496

497+
}
477498

478499
void
479500
BTDiagnostics::DefineFieldBufferMultiFab (const int i_buffer, const int lev)
@@ -610,8 +631,8 @@ BTDiagnostics::Flush (int i_buffer)
610631
file_name = amrex::Concatenate(m_file_prefix,i_buffer,5);
611632
file_name = file_name+"/buffer";
612633
}
613-
bool isLastBTDFlush = ( ( m_max_buffer_multifabs[i_buffer]
614-
- m_buffer_flush_counter[i_buffer]) == 1) ? true : false;
634+
SetSnapshotFullStatus(i_buffer);
635+
bool isLastBTDFlush = ( m_snapshot_full[i_buffer] == 1 ) ? true : false;
615636
bool const isBTD = true;
616637
double const labtime = m_t_lab[i_buffer];
617638
m_flush_format->WriteToFile(

Source/Diagnostics/ComputeDiagFunctors/BackTransformFunctor.H

Lines changed: 14 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -60,21 +60,22 @@ public:
6060

6161
/** \brief Prepare data required to back-transform fields for lab-frame snapshot, i_buffer
6262
*
63-
* \param[in] i_buffer index of the snapshot
64-
* \param[in] ZSliceInDomain if the z-slice at current_z_boost is within
65-
* the boosted-frame and lab-frame domain.
66-
* The fields are sliced and back-transformed only
67-
* if this value is true.
68-
* \param[in] current_z_boost current z-coordinate of the slice in boosted-frame
69-
* \param[in] buffer_box Box with index-space in lab-frame for the ith buffer
70-
* \param[in] k_index_zlab k-index in the lab-frame corresponding to the
71-
* current z co-ordinate in the lab-frame for the ith buffer.
72-
* \param[in] max_box_size maximum box size for the multifab to generate box arrays
63+
* \param[in] i_buffer index of the snapshot
64+
* \param[in] z_slice_in_domain if the z-slice at current_z_boost is within
65+
* the boosted-frame and lab-frame domain.
66+
* The fields are sliced and back-transformed only if this value is true.
67+
* \param[in] buffer_box Box with index-space in lab-frame for the ith buffer
68+
* \param[in] k_index_zlab k-index in the lab-frame corresponding to the
69+
* current z co-ordinate in the lab-frame for the ith buffer.
70+
* \param[in] max_box_size maximum box size for the multifab to generate box arrays
71+
* \param[in] snapshot_full if the current snapshot, with index, i_buffer, is full (1)
72+
or not (0). If it is full, then Lorentz-transform is not performed
73+
by setting m_perform_backtransform to 0.
7374
*/
74-
void PrepareFunctorData ( int i_buffer, bool ZSliceInDomain,
75+
void PrepareFunctorData ( int i_buffer, bool z_slice_in_domain,
7576
amrex::Real current_z_boost,
7677
amrex::Box buffer_box, const int k_index_zlab,
77-
const int max_box_size ) override;
78+
const int max_box_size, const int snapshot_full ) override;
7879
/** Allocate and initialize member variables and arrays required to back-transform
7980
* field-data from boosted-frame to lab-frame.
8081
*/
@@ -103,7 +104,7 @@ private:
103104
amrex::Vector<amrex::Real> m_current_z_boost;
104105
/** Vector of 0s and 1s stored to check if back-transformation is to be performed
105106
* for the ith buffer. The value is set 0 (false) or 1 (true) using the boolean
106-
* ZSliceInDomain in PrepareFunctorData().
107+
* z_slice_in_domain in PrepareFunctorData().
107108
*/
108109
amrex::Vector<int> m_perform_backtransform;
109110
/** Vector of k-index correspoding to the current lab-frame z co-ordinate for each buffer */

Source/Diagnostics/ComputeDiagFunctors/BackTransformFunctor.cpp

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -119,15 +119,15 @@ BackTransformFunctor::operator ()(amrex::MultiFab& mf_dst, int /*dcomp*/, const
119119

120120
void
121121
BackTransformFunctor::PrepareFunctorData (int i_buffer,
122-
bool ZSliceInDomain, amrex::Real current_z_boost,
122+
bool z_slice_in_domain, amrex::Real current_z_boost,
123123
amrex::Box buffer_box, const int k_index_zlab,
124-
const int max_box_size )
124+
const int max_box_size, const int snapshot_full)
125125
{
126126
m_buffer_box[i_buffer] = buffer_box;
127127
m_current_z_boost[i_buffer] = current_z_boost;
128128
m_k_index_zlab[i_buffer] = k_index_zlab;
129129
m_perform_backtransform[i_buffer] = 0;
130-
if (ZSliceInDomain) m_perform_backtransform[i_buffer] = 1;
130+
if (z_slice_in_domain == true and snapshot_full == 0) m_perform_backtransform[i_buffer] = 1;
131131
m_max_box_size = max_box_size;
132132
}
133133

Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H

Lines changed: 20 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -30,14 +30,30 @@ public:
3030
* multifab */
3131
int nComp () const { return m_ncomp; }
3232

33-
virtual void PrepareFunctorData ( int i_buffer, bool ZSliceInDomain,
33+
/** \brief Prepare data required to process fields in the operator()
34+
* Note that this function has parameters that are specific to
35+
* back-transformed diagnostics, that are unused for regular diagnostics.
36+
*
37+
* \param[in] i_buffer index of the back-transform snapshot
38+
* \param[in] z_slice_in_domain if the z-slice at current_z_boost is within
39+
* the boosted-frame and lab-frame domain.
40+
* The fields are sliced and back-transformed only if this value is true.
41+
* \param[in] buffer_box Box with index-space in lab-frame for the ith buffer
42+
* \param[in] k_index_zlab k-index in the lab-frame corresponding to the
43+
* current z co-ordinate in the lab-frame for the ith buffer.
44+
* \param[in] max_box_size maximum box size for the multifab to generate box arrays
45+
* \param[in] snapshot_full if the current snapshot, with index, i_buffer, is full (1)
46+
or not (0). If it is full, then Lorentz-transform
47+
is not performed by setting m_perform_backtransform to 0;
48+
*/
49+
virtual void PrepareFunctorData ( int i_buffer, bool z_slice_in_domain,
3450
amrex::Real current_z_boost,
3551
amrex::Box buffer_box, const int k_index_zlab,
36-
const int max_box_size) {
52+
const int max_box_size, const int snapshot_full) {
3753
amrex::ignore_unused(i_buffer,
38-
ZSliceInDomain,
54+
z_slice_in_domain,
3955
current_z_boost, buffer_box,
40-
k_index_zlab, max_box_size);
56+
k_index_zlab, max_box_size, snapshot_full);
4157
}
4258
virtual void InitData() {}
4359
private:

0 commit comments

Comments
 (0)