WarpX
Public Member Functions | Public Attributes | Private Attributes | List of all members
BTDSpeciesHeaderImpl Class Reference

Class to read, modify, and write species header when back-transformed diag format is selected as plotfile. This class enables multiple BTD particle buffers to be interweaved and stitched into a single plotfile with a single Header. More...

#include <BTD_Plotfile_Header_Impl.H>

Public Member Functions

 BTDSpeciesHeaderImpl (std::string const &Headerfile_path, std::string const &species_name)
 
void ReadHeader ()
 
void WriteHeader ()
 
void set_DataIndex (int lev, int box_id, int data_index)
 
void AddTotalParticles (const int new_particles)
 
void IncrementParticleBoxArraySize ()
 
void AppendParticleInfoForNewBox (int data_index, int particles_per_box, int offset)
 

Public Attributes

amrex::Vector< amrex::Vector< int > > m_which_data
 
amrex::Vector< amrex::Vector< int > > m_particles_per_box
 
amrex::Vector< amrex::Vector< int > > m_offset_per_box
 
amrex::Vector< intm_particleBoxArray_size
 
int m_total_particles
 
int m_nextid
 
int m_finestLevel
 

Private Attributes

std::string m_Header_path
 
std::string m_species_name
 
std::string m_file_version
 
int m_spacedim
 
int m_num_output_real
 
int m_num_output_int
 
amrex::Vector< std::string > m_real_comp_names
 
amrex::Vector< std::string > m_int_comp_names
 
bool m_is_checkpoint
 

Detailed Description

Class to read, modify, and write species header when back-transformed diag format is selected as plotfile. This class enables multiple BTD particle buffers to be interweaved and stitched into a single plotfile with a single Header.

Constructor & Destructor Documentation

◆ BTDSpeciesHeaderImpl()

BTDSpeciesHeaderImpl::BTDSpeciesHeaderImpl ( std::string const &  Headerfile_path,
std::string const &  species_name 
)

Constructor.

Parameters
[in]Headerfile_pathstring containing path of Headerfile
[in]species_namestring containing species name

Member Function Documentation

◆ AddTotalParticles()

void BTDSpeciesHeaderImpl::AddTotalParticles ( const int  new_particles)
inline

Add new_particles to existing to obtain the total number of particles of the species.

Parameters
[in]new_particlestotal particles in the new buffer

◆ AppendParticleInfoForNewBox()

void BTDSpeciesHeaderImpl::AppendParticleInfoForNewBox ( int  data_index,
int  particles_per_box,
int  offset 
)

Append particle info of the newly added box, namely, m_which_data, m_particles_per_box, m_offset_per_box.

◆ IncrementParticleBoxArraySize()

void BTDSpeciesHeaderImpl::IncrementParticleBoxArraySize ( )
inline

Increment number of boxes in a box array by 1, with every flush.

◆ ReadHeader()

void BTDSpeciesHeaderImpl::ReadHeader ( )

Reads the Header file for BTD species

◆ set_DataIndex()

void BTDSpeciesHeaderImpl::set_DataIndex ( int  lev,
int  box_id,
int  data_index 
)

Set data Index of the data-file, DATAXXXXX, that the particles belong to

◆ WriteHeader()

void BTDSpeciesHeaderImpl::WriteHeader ( )

Writes the meta-data of species Header file, with path, m_Header_path

Member Data Documentation

◆ m_file_version

std::string BTDSpeciesHeaderImpl::m_file_version
private

String containing file version of the species output

◆ m_finestLevel

int BTDSpeciesHeaderImpl::m_finestLevel

finest level of the particle output

◆ m_Header_path

std::string BTDSpeciesHeaderImpl::m_Header_path
private

Path of the headerfile

◆ m_int_comp_names

amrex::Vector<std::string> BTDSpeciesHeaderImpl::m_int_comp_names
private

The real integer names

◆ m_is_checkpoint

bool BTDSpeciesHeaderImpl::m_is_checkpoint
private

if the output is a checkpoint

◆ m_nextid

int BTDSpeciesHeaderImpl::m_nextid

Id of the next particle

◆ m_num_output_int

int BTDSpeciesHeaderImpl::m_num_output_int
private

Number of integer attributes

◆ m_num_output_real

int BTDSpeciesHeaderImpl::m_num_output_real
private

Number of real attributes

◆ m_offset_per_box

amrex::Vector< amrex::Vector<int> > BTDSpeciesHeaderImpl::m_offset_per_box

Vector of particle offset per box of each particle data file for each level

◆ m_particleBoxArray_size

amrex::Vector<int> BTDSpeciesHeaderImpl::m_particleBoxArray_size

Vector of box array size per for each level

◆ m_particles_per_box

amrex::Vector< amrex::Vector<int> > BTDSpeciesHeaderImpl::m_particles_per_box

Vector of particles per box for each level

◆ m_real_comp_names

amrex::Vector<std::string> BTDSpeciesHeaderImpl::m_real_comp_names
private

The real attribute names

◆ m_spacedim

int BTDSpeciesHeaderImpl::m_spacedim
private

Number of dimensions stored in the plotfile, should be same as AMREX_SPACEDIM

◆ m_species_name

std::string BTDSpeciesHeaderImpl::m_species_name
private

Species name of the particles flushed out

◆ m_total_particles

int BTDSpeciesHeaderImpl::m_total_particles

Total number of particles of species, m_species_name, in the particle output

◆ m_which_data

amrex::Vector< amrex::Vector<int> > BTDSpeciesHeaderImpl::m_which_data

Vector of data indices of all particle data for each level


The documentation for this class was generated from the following files: