WarpX
Loading...
Searching...
No Matches
BTDiagnostics.H
Go to the documentation of this file.
1/* Copyright 2021 Revathi Jambunathan
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef WARPX_BTDIAGNOSTICS_H_
8#define WARPX_BTDIAGNOSTICS_H_
9
10#include "Diagnostics.H"
13#include "Utils/WarpXConst.H"
15
16#include <AMReX_Box.H>
17#include <AMReX_Geometry.H>
18#include <AMReX_IntVect.H>
19#include <AMReX_REAL.H>
20#include <AMReX_RealBox.H>
21#include <AMReX_Vector.H>
22
23#include <limits>
24#include <memory>
25#include <string>
26
27class BTDiagnostics final : public Diagnostics
28{
29public:
30
31 BTDiagnostics (int i, const std::string& name, DiagTypes diag_type);
32
33private:
35 bool m_plot_raw_fields = false;
39 void ReadParameters ();
49 void Flush (int i_buffer, bool force_flush) override;
58 bool DoDump (int step, int i_buffer, bool force_flush=false) override;
66 bool DoComputeAndPack (int step, bool force_flush=false) override;
68 void DerivedInitData () override;
75 void InitializeFieldFunctors (int lev) override;
82 void InitializeFieldFunctorsRZopenPMD (int lev) override;
94 void AddRZModesToOutputNames (const std::string& field, int ncomp, bool cellcenter_data);
95
99 void InitializeParticleBuffer (const MultiParticleContainer& mpc) override;
100
104 void PrepareBufferData () override;
106 void UpdateBufferData () override;
114 void PrepareFieldDataForOutput () override;
116 void PrepareParticleDataForOutput() override;
117
126 bool GetZSliceInDomainFlag (int i_buffer, int lev);
127
134 bool GetKIndexInSnapshotBoxFlag (int i_buffer, int lev);
135
142 void InitializeBufferData ( int i_buffer , int lev, bool restart=false) override;
151
159 amrex::Real m_gamma_boost;
160 amrex::Real m_beta_boost;
167
169 int m_num_snapshots_lab = std::numeric_limits<int>::lowest();
171 amrex::Real m_dt_snapshots_lab = std::numeric_limits<amrex::Real>::lowest();
175 amrex::Real m_dz_snapshots_lab = 0.0;
176
178 int m_buffer_size = 256;
179
257 void DefineCellCenteredMultiFab(int lev);
264 void DefineFieldBufferMultiFab (int i_buffer, int lev);
265
272 void DefineSnapshotGeometry (int i_buffer, int lev);
273
278 [[nodiscard]] amrex::Real UpdateCurrentZBoostCoordinate(amrex::Real t_lab, amrex::Real t_boost) const
279 {
280 const amrex::Real current_z_boost = (t_lab / m_gamma_boost - t_boost) * PhysConst::c / m_beta_boost;
281 return current_z_boost;
282 }
283
287 [[nodiscard]] amrex::Real UpdateCurrentZLabCoordinate(amrex::Real t_lab, amrex::Real t_boost) const
288 {
289 const amrex::Real current_z_lab = (t_lab - t_boost / m_gamma_boost ) * PhysConst::c / m_beta_boost;
290 return current_z_lab;
291 }
292
297 [[nodiscard]] amrex::Real dz_lab (amrex::Real dt, amrex::Real ref_ratio) const;
303 [[nodiscard]] int k_index_zlab (int i_buffer, int lev) const;
309 bool buffer_full (int i_buffer) {
310 return (k_index_zlab(i_buffer,0) == m_buffer_box[i_buffer].smallEnd(m_moving_window_dir));
311 }
312
317 bool buffer_empty (int i_buffer) {
318 return ( m_buffer_counter[i_buffer] == 0) ;
319 }
320
323
327 void NullifyFirstFlush(int i_buffer) {m_first_flush_after_restart[i_buffer] = 0;}
328
333
337 void ResetBufferCounter(int i_buffer) {
338 m_buffer_counter[i_buffer] = 0;
339 }
340
343 void IncrementBufferFlushCounter(int i_buffer) {
344 m_buffer_flush_counter[i_buffer]++;
345 }
346
351 void SetSnapshotFullStatus (int i_buffer);
356#ifdef WARPX_DIM_RZ
358 "Br", "Bt", "Bz",
359 "jr", "jt", "jz", "rho"};
361 "Br", "Bt", "Bz",
362 "jr", "jt", "jz",
363 "rho"};
364#else
366 "Bx", "By", "Bz",
367 "jx", "jy", "jz", "rho"};
368#endif
369
370
374 void MergeBuffersForPlotfile (int i_snapshot);
378 void InterleaveBufferAndSnapshotHeader (const std::string& buffer_Header,
379 const std::string& snapshot_Header);
383 void InterleaveFabArrayHeader (const std::string& Buffer_FabHeader_path,
384 const std::string& snapshot_FabHeader_path,
385 const std::string& newsnapshot_FabFilename);
389 void InterleaveSpeciesHeader(const std::string& buffer_species_Header_path,
390 const std::string& snapshot_species_Header_path,
391 const std::string& species_name, int new_data_index);
392
396 void InterleaveParticleDataHeader(const std::string& buffer_ParticleHdrFilename,
397 const std::string& snapshot_ParticleHdrFilename);
400 void InitializeParticleFunctors () override;
401
403 void ResetTotalParticlesInBuffer(int i_buffer);
405 void ClearParticleBuffer(int i_buffer);
407 void RedistributeParticleBuffer (int i_buffer);
409
410 amrex::Real gettlab (int i_buffer) override {return m_t_lab[i_buffer];}
411 void settlab (int i_buffer, amrex::Real tlab) override {m_t_lab[i_buffer] = tlab; }
412 int get_buffer_k_index_hi (int i_buffer) override {return m_buffer_k_index_hi[i_buffer]; }
413 void set_buffer_k_index_hi (int i_buffer, int kindex) override {m_buffer_k_index_hi[i_buffer] = kindex;}
414 amrex::Real get_snapshot_domain_lo (int i_buffer, int idim) override {return m_snapshot_domain_lab[i_buffer].lo(idim); }
415 amrex::Real get_snapshot_domain_hi (int i_buffer, int idim) override {return m_snapshot_domain_lab[i_buffer].hi(idim); }
416 int get_flush_counter (int i_buffer) override {return m_buffer_flush_counter[i_buffer]; }
417 void set_flush_counter ([[maybe_unused]] int i_buffer, int flush_counter) override { m_buffer_flush_counter[i_buffer] = flush_counter; }
418 int get_last_valid_Zslice ( int i_buffer) override { return m_lastValidZSlice[i_buffer]; }
419 void set_last_valid_Zslice ( int i_buffer, int last_valid_Zslice) override { m_lastValidZSlice[i_buffer] = last_valid_Zslice; }
420 int get_snapshot_full_flag ( int i_buffer) override { return m_snapshot_full[i_buffer]; }
421 void set_snapshot_full ( int i_buffer, int snapshot_full) override { m_snapshot_full[i_buffer] = snapshot_full; }
422};
423#endif // WARPX_BTDIAGNOSTICS_H_
DiagTypes
Definition Diagnostics.H:25
void set_flush_counter(int i_buffer, int flush_counter) override
Definition BTDiagnostics.H:417
amrex::Real m_beta_boost
Definition BTDiagnostics.H:160
int m_moving_window_dir
Definition BTDiagnostics.H:165
amrex::Real m_moving_window_beta
Definition BTDiagnostics.H:166
void DerivedInitData() override
Definition BTDiagnostics.cpp:70
amrex::Vector< std::string > m_cellcenter_varnames
Definition BTDiagnostics.H:357
void Flush(int i_buffer, bool force_flush) override
Flush m_mf_output and particles to file. Currently, a temporary customized output format for the buff...
Definition BTDiagnostics.cpp:1038
void DefineCellCenteredMultiFab(int lev)
Definition BTDiagnostics.cpp:513
void InitializeFieldFunctorsRZopenPMD(int lev) override
Definition BTDiagnostics.cpp:673
amrex::Real dz_lab(amrex::Real dt, amrex::Real ref_ratio) const
Definition BTDiagnostics.cpp:891
void IncrementBufferFlushCounter(int i_buffer)
Definition BTDiagnostics.H:343
void InterleaveFabArrayHeader(const std::string &Buffer_FabHeader_path, const std::string &snapshot_FabHeader_path, const std::string &newsnapshot_FabFilename)
Definition BTDiagnostics.cpp:1372
amrex::Real m_dz_snapshots_lab
Definition BTDiagnostics.H:175
bool DoDump(int step, int i_buffer, bool force_flush=false) override
Definition BTDiagnostics.cpp:300
bool GetKIndexInSnapshotBoxFlag(int i_buffer, int lev)
Definition BTDiagnostics.cpp:1031
amrex::Vector< std::string > m_cellcenter_varnames_fields
Definition BTDiagnostics.H:360
void AddRZModesToOutputNames(const std::string &field, int ncomp, bool cellcenter_data)
Definition BTDiagnostics.cpp:733
amrex::Real m_dt_snapshots_lab
Definition BTDiagnostics.H:171
amrex::Vector< int > m_lastValidZSlice
Definition BTDiagnostics.H:233
void set_last_valid_Zslice(int i_buffer, int last_valid_Zslice) override
Definition BTDiagnostics.H:419
BTDiagnostics(int i, const std::string &name, DiagTypes diag_type)
Definition BTDiagnostics.cpp:58
void PrepareFieldDataForOutput() override
Definition BTDiagnostics.cpp:806
bool m_do_back_transformed_particles
Definition BTDiagnostics.H:150
void DefineFieldBufferMultiFab(int i_buffer, int lev)
Definition BTDiagnostics.cpp:922
amrex::Vector< amrex::RealBox > m_buffer_domain_lab
Definition BTDiagnostics.H:184
void ReadParameters()
Definition BTDiagnostics.cpp:210
amrex::Vector< int > m_buffer_counter
Definition BTDiagnostics.H:219
amrex::Vector< amrex::Real > m_old_z_boost
Definition BTDiagnostics.H:201
void set_snapshot_full(int i_buffer, int snapshot_full) override
Definition BTDiagnostics.H:421
amrex::Vector< int > m_snapshot_full
Definition BTDiagnostics.H:228
int get_snapshot_full_flag(int i_buffer) override
Definition BTDiagnostics.H:420
void NullifyFirstFlush(int i_buffer)
Definition BTDiagnostics.H:327
amrex::Vector< amrex::Vector< std::unique_ptr< ComputeDiagFunctor const > > > m_cell_center_functors
Definition BTDiagnostics.H:251
amrex::Real gettlab(int i_buffer) override
Definition BTDiagnostics.H:410
utils::parser::BTDIntervalsParser m_intervals
Definition BTDiagnostics.H:41
void ClearParticleBuffer(int i_buffer)
Definition BTDiagnostics.cpp:1541
void MergeBuffersForPlotfile(int i_snapshot)
Definition BTDiagnostics.cpp:1156
void settlab(int i_buffer, amrex::Real tlab) override
Definition BTDiagnostics.H:411
void InitializeParticleBuffer(const MultiParticleContainer &mpc) override
Definition BTDiagnostics.cpp:1464
int m_buffer_size
Definition BTDiagnostics.H:178
amrex::Vector< amrex::Real > m_current_z_lab
Definition BTDiagnostics.H:195
amrex::Real get_snapshot_domain_hi(int i_buffer, int idim) override
Definition BTDiagnostics.H:415
bool buffer_full(int i_buffer)
Definition BTDiagnostics.H:309
amrex::Vector< int > m_field_buffer_multifab_defined
Definition BTDiagnostics.H:332
void InterleaveBufferAndSnapshotHeader(const std::string &buffer_Header, const std::string &snapshot_Header)
Definition BTDiagnostics.cpp:1327
amrex::Vector< int > m_max_buffer_multifabs
Definition BTDiagnostics.H:224
bool buffer_empty(int i_buffer)
Definition BTDiagnostics.H:317
void set_buffer_k_index_hi(int i_buffer, int kindex) override
Definition BTDiagnostics.H:413
void InterleaveSpeciesHeader(const std::string &buffer_species_Header_path, const std::string &snapshot_species_Header_path, const std::string &species_name, int new_data_index)
Definition BTDiagnostics.cpp:1402
void RedistributeParticleBuffer(int i_buffer)
Definition BTDiagnostics.cpp:1149
void SetSnapshotFullStatus(int i_buffer)
Definition BTDiagnostics.cpp:913
amrex::Vector< int > m_buffer_flush_counter
Definition BTDiagnostics.H:236
bool DoComputeAndPack(int step, bool force_flush=false) override
Definition BTDiagnostics.cpp:328
amrex::Vector< amrex::Box > m_buffer_box
Definition BTDiagnostics.H:192
void InterleaveParticleDataHeader(const std::string &buffer_ParticleHdrFilename, const std::string &snapshot_ParticleHdrFilename)
Definition BTDiagnostics.cpp:1426
bool m_plot_raw_fields_guards
Definition BTDiagnostics.H:37
amrex::Real UpdateCurrentZBoostCoordinate(amrex::Real t_lab, amrex::Real t_boost) const
Definition BTDiagnostics.H:278
void UpdateVarnamesForRZopenPMD()
Definition BTDiagnostics.cpp:612
void PrepareParticleDataForOutput() override
Definition BTDiagnostics.cpp:1494
amrex::Vector< amrex::Real > m_current_z_boost
Definition BTDiagnostics.H:198
void UpdateBufferData() override
Definition BTDiagnostics.cpp:784
void InitializeParticleFunctors() override
Definition BTDiagnostics.cpp:1447
void InitializeBufferData(int i_buffer, int lev, bool restart=false) override
Definition BTDiagnostics.cpp:339
void ResetTotalParticlesInBuffer(int i_buffer)
Definition BTDiagnostics.cpp:1532
std::string m_cell_centered_data_name
Definition BTDiagnostics.H:247
void DefineSnapshotGeometry(int i_buffer, int lev)
Definition BTDiagnostics.cpp:987
bool m_plot_raw_fields
Definition BTDiagnostics.H:35
amrex::Vector< int > m_buffer_k_index_hi
Definition BTDiagnostics.H:240
bool m_do_back_transformed_fields
Definition BTDiagnostics.H:146
amrex::Real m_gamma_boost
Definition BTDiagnostics.H:159
amrex::Vector< amrex::IntVect > m_snapshot_ncells_lab
Definition BTDiagnostics.H:186
amrex::Real UpdateCurrentZLabCoordinate(amrex::Real t_lab, amrex::Real t_boost) const
Definition BTDiagnostics.H:287
int k_index_zlab(int i_buffer, int lev) const
Definition BTDiagnostics.cpp:898
amrex::Vector< amrex::Vector< amrex::Geometry > > m_geom_snapshot
Definition BTDiagnostics.H:212
int get_last_valid_Zslice(int i_buffer) override
Definition BTDiagnostics.H:418
amrex::Real get_snapshot_domain_lo(int i_buffer, int idim) override
Definition BTDiagnostics.H:414
int get_flush_counter(int i_buffer) override
Definition BTDiagnostics.H:416
void ResetBufferCounter(int i_buffer)
Definition BTDiagnostics.H:337
amrex::Vector< int > m_first_flush_after_restart
Definition BTDiagnostics.H:322
amrex::Vector< amrex::Box > m_snapshot_box
Definition BTDiagnostics.H:189
void PrepareBufferData() override
Definition BTDiagnostics.cpp:761
void InitializeFieldFunctors(int lev) override
Definition BTDiagnostics.cpp:537
int m_num_snapshots_lab
Definition BTDiagnostics.H:169
amrex::Vector< amrex::Real > m_t_lab
Definition BTDiagnostics.H:181
int get_buffer_k_index_hi(int i_buffer) override
Definition BTDiagnostics.H:412
amrex::Vector< int > m_snapshot_geometry_defined
Definition BTDiagnostics.H:330
bool GetZSliceInDomainFlag(int i_buffer, int lev)
Definition BTDiagnostics.cpp:1010
amrex::Vector< amrex::RealBox > m_snapshot_domain_lab
Definition Diagnostics.H:341
Diagnostics(int i, std::string name, DiagTypes diag_type)
Definition Diagnostics.cpp:41
Definition MultiParticleContainer.H:68
This class is a parser for multiple slices of the form x,y,z,... where x, y and z are slices of the f...
Definition IntervalsParser.H:177
static constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:44