WarpX
BTDiagnostics.H
Go to the documentation of this file.
1 #ifndef WARPX_BTDIAGNOSTICS_H_
2 #define WARPX_BTDIAGNOSTICS_H_
3 
4 #include "Diagnostics.H"
6 #include "Utils/WarpXConst.H"
7 
8 #include <AMReX_Box.H>
9 #include <AMReX_Geometry.H>
10 #include <AMReX_IntVect.H>
11 #include <AMReX_MultiFab.H>
12 #include <AMReX_REAL.H>
13 #include <AMReX_RealBox.H>
14 #include <AMReX_Vector.H>
15 
16 #include <limits>
17 #include <memory>
18 #include <string>
19 
20 class
21 BTDiagnostics final : public Diagnostics
22 {
23 public:
24 
25  BTDiagnostics (int i, std::string name);
26 
27 private:
29  bool m_plot_raw_fields = false;
31  bool m_plot_raw_fields_guards = false;
33  bool m_plot_raw_rho = false;
35  bool m_plot_raw_F = false;
37  void ReadParameters ();
45  void Flush (int i_buffer) override;
53  bool DoDump (int step, int i_buffer, bool force_flush=false) override;
61  bool DoComputeAndPack (int step, bool force_flush=false) override;
63  void DerivedInitData () override;
70  void InitializeFieldFunctors (int lev) override;
75  void InitializeParticleBuffer () override {}
83  void PrepareFieldDataForOutput () override;
84 
92  bool GetZSliceInDomainFlag (const int i_buffer, const int lev);
93 
99  void InitializeFieldBufferData ( int i_buffer , int lev) override;
103  bool m_do_back_transformed_fields = true;
107  bool m_do_back_transformed_particles = true;
108 
116  amrex::Real m_gamma_boost;
117  amrex::Real m_beta_boost;
123 
125  int m_num_snapshots_lab = std::numeric_limits<int>::lowest();
127  amrex::Real m_dt_snapshots_lab = std::numeric_limits<amrex::Real>::lowest();
131  amrex::Real m_dz_snapshots_lab = 0.0;
132 
134  int m_buffer_size = 256;
136  int m_max_box_size = 256;
137 
139  amrex::Vector<amrex::Real> m_t_lab;
142  amrex::Vector<amrex::RealBox> m_prob_domain_lab;
145  amrex::Vector<amrex::RealBox> m_snapshot_domain_lab;
148  amrex::Vector<amrex::RealBox> m_buffer_domain_lab;
150  amrex::Vector<amrex::IntVect> m_snapshot_ncells_lab;
153  amrex::Vector<amrex::Box> m_snapshot_box;
156  amrex::Vector<amrex::Box> m_buffer_box;
159  amrex::Vector<amrex::Real> m_current_z_lab;
162  amrex::Vector<amrex::Real> m_current_z_boost;
173  amrex::Vector< amrex::Vector <amrex::Geometry> > m_geom_snapshot;
180  amrex::Vector<int> m_buffer_counter;
185  amrex::Vector<int> m_max_buffer_multifabs;
188  amrex::Vector<int> m_buffer_flush_counter;
195  amrex::Vector<std::unique_ptr<amrex::MultiFab> > m_cell_centered_data;
199  amrex::Vector< amrex::Vector <std::unique_ptr<ComputeDiagFunctor const> > > m_cell_center_functors;
205  void DefineCellCenteredMultiFab(int lev);
212  void DefineFieldBufferMultiFab (const int i_buffer, const int lev);
213 
220  void DefineSnapshotGeometry (const int i_buffer, const int lev);
221 
226  amrex::Real UpdateCurrentZBoostCoordinate(amrex::Real t_lab, amrex::Real t_boost)
227  {
228  amrex::Real current_z_boost = (t_lab / m_gamma_boost - t_boost) * PhysConst::c / m_beta_boost;
229  return current_z_boost;
230  }
235  amrex::Real UpdateCurrentZLabCoordinate(amrex::Real t_lab, amrex::Real t_boost)
236  {
237  amrex::Real current_z_lab = (t_lab - t_boost / m_gamma_boost ) * PhysConst::c / m_beta_boost;
238  return current_z_lab;
239  }
244  amrex::Real dz_lab (amrex::Real dt, amrex::Real ref_ratio);
250  int k_index_zlab (int i_buffer, int lev);
256  bool buffer_full (int i_buffer) {
257  return ( m_buffer_counter[i_buffer] == m_buffer_size );
258  }
259 
264  bool buffer_empty (int i_buffer) {
265  return ( m_buffer_counter[i_buffer] == 0) ;
266  }
267 
271  void ResetBufferCounter(int i_buffer) {
272  m_buffer_counter[i_buffer] = 0;
273  }
277  void IncrementBufferFlushCounter(int i_buffer) {
278  m_buffer_flush_counter[i_buffer]++;
279  }
284  amrex::Vector< std::string > m_cellcenter_varnames = {"Ex", "Ey", "Ez",
285  "Bx", "By", "Bz",
286  "jx", "jy", "jz", "rho"};
287 
289  void TMP_ClearSpeciesDataForBTD() override;
290 
294  void MergeBuffersForPlotfile (int i_snapshot);
298  void InterleaveBufferAndSnapshotHeader ( std::string buffer_Header,
299  std::string snapshot_Header);
303  void InterleaveFabArrayHeader( std::string Buffer_FabHeaderFilename,
304  std::string snapshot_FabHeaderFilename,
305  std::string newsnapshot_FabFilename);
306 };
307 
308 #endif // WARPX_BTDIAGNOSTICS_H_
void ResetBufferCounter(int i_buffer)
Definition: BTDiagnostics.H:271
amrex::Vector< amrex::RealBox > m_prob_domain_lab
Definition: BTDiagnostics.H:142
void InitializeParticleBuffer() override
Definition: BTDiagnostics.H:75
amrex::Vector< amrex::Vector< amrex::Geometry > > m_geom_snapshot
Definition: BTDiagnostics.H:173
amrex::Vector< amrex::RealBox > m_snapshot_domain_lab
Definition: BTDiagnostics.H:145
void IncrementBufferFlushCounter(int i_buffer)
Definition: BTDiagnostics.H:277
amrex::Vector< amrex::Box > m_snapshot_box
Definition: BTDiagnostics.H:153
amrex::Vector< int > m_buffer_counter
Definition: BTDiagnostics.H:180
int m_moving_window_dir
Definition: BTDiagnostics.H:122
amrex::Vector< amrex::Real > m_current_z_boost
Definition: BTDiagnostics.H:162
amrex::Vector< amrex::Real > m_t_lab
Definition: BTDiagnostics.H:139
amrex::Vector< amrex::RealBox > m_buffer_domain_lab
Definition: BTDiagnostics.H:148
amrex::Vector< amrex::IntVect > m_snapshot_ncells_lab
Definition: BTDiagnostics.H:150
amrex::Real UpdateCurrentZBoostCoordinate(amrex::Real t_lab, amrex::Real t_boost)
Definition: BTDiagnostics.H:226
amrex::Vector< amrex::Real > m_current_z_lab
Definition: BTDiagnostics.H:159
amrex::Vector< int > m_buffer_flush_counter
Definition: BTDiagnostics.H:188
i
Definition: check_interp_points_and_weights.py:171
amrex::Vector< amrex::Vector< std::unique_ptr< ComputeDiagFunctor const > > > m_cell_center_functors
Definition: BTDiagnostics.H:199
name
Definition: run_automated.py:204
amrex::Real m_gamma_boost
Definition: BTDiagnostics.H:116
amrex::Vector< int > m_max_buffer_multifabs
Definition: BTDiagnostics.H:185
amrex::Vector< amrex::Box > m_buffer_box
Definition: BTDiagnostics.H:156
bool buffer_full(int i_buffer)
Definition: BTDiagnostics.H:256
amrex::Real UpdateCurrentZLabCoordinate(amrex::Real t_lab, amrex::Real t_boost)
Definition: BTDiagnostics.H:235
bool buffer_empty(int i_buffer)
Definition: BTDiagnostics.H:264
amrex::Vector< std::unique_ptr< amrex::MultiFab > > m_cell_centered_data
Definition: BTDiagnostics.H:195
Definition: BTDiagnostics.H:20
amrex::Real m_beta_boost
Definition: BTDiagnostics.H:117
base class for diagnostics. Contains main routines to filter, compute and flush diagnostics.
Definition: Diagnostics.H:25