Flow123d  master-f44eb46
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
Balance Class Reference

#include <balance.hh>

Collaboration diagram for Balance:
Collaboration graph
[legend]

Classes

class  Quantity
 

Public Types

enum  OutputFormat { legacy , txt , gnuplot }
 

Public Member Functions

 Balance (const std::string &file_prefix, const Mesh *mesh)
 
 ~Balance ()
 
void init_from_input (const Input::Record &in_rec, TimeGovernor &tg)
 
void units (const UnitSI &unit)
 Setter for units of conserved quantities. More...
 
bool cumulative () const
 Getter for cumulative_. More...
 
unsigned int add_quantity (const string &name)
 
std::vector< unsigned int > add_quantities (const std::vector< string > &names)
 
void allocate (const std::shared_ptr< DOFHandlerMultiDim > &dh, unsigned int max_dofs_per_boundary)
 
void allocate (unsigned int n_loc_dofs, unsigned int max_dofs_per_boundary)
 
bool is_current ()
 Returns true if the current time step is marked for the balance output. More...
 
void start_mass_assembly (unsigned int quantity_idx)
 
void start_mass_assembly (std::vector< unsigned int > q_idx_vec)
 Variant of the start_mass_assembly() method for a set of quantities. More...
 
void start_flux_assembly (unsigned int quantity_idx)
 
void start_flux_assembly (std::vector< unsigned int > q_idx_vec)
 Variant of the start_flux_assembly() method for a set of quantities. More...
 
void start_source_assembly (unsigned int quantity_idx)
 
void start_source_assembly (std::vector< unsigned int > q_idx_vec)
 Variant of the start_source_assembly() method for a set of quantities. More...
 
void add_mass_values (unsigned int quantity_idx, const DHCellAccessor &dh_cell, const LocDofVec &loc_dof_indices, const std::vector< double > &mat_values, double vec_value)
 
void add_flux_values (unsigned int quantity_idx, const DHCellSide &side, const LocDofVec &loc_dof_indices, const std::vector< double > &mat_values, double vec_value)
 
void add_source_values (unsigned int quantity_idx, unsigned int region_idx, const LocDofVec &loc_dof_indices, const std::vector< double > &mult_mat_values, const std::vector< double > &add_mat_values)
 
void finish_mass_assembly (unsigned int quantity_idx)
 This method must be called after assembling the matrix for computing mass. More...
 
void finish_mass_assembly (std::vector< unsigned int > q_idx_vec)
 Variant of the finish_mass_assembly() method for a set of quantities. More...
 
void finish_flux_assembly (unsigned int quantity_idx)
 This method must be called after assembling the matrix and vector for computing flux. More...
 
void finish_flux_assembly (std::vector< unsigned int > q_idx_vec)
 Variant of the finish_flux_assembly() method for a set of quantities. More...
 
void finish_source_assembly (unsigned int quantity_idx)
 This method must be called after assembling the matrix and vectors for computing source. More...
 
void finish_source_assembly (std::vector< unsigned int > q_idx_vec)
 Variant of the finish_source_assembly() method for a set of quantities. More...
 
void calculate_cumulative (unsigned int quantity_idx, const Vec &solution)
 
void calculate_mass (unsigned int quantity_idx, const Vec &solution, vector< double > &output_array)
 
void calculate_instant (unsigned int quantity_idx, const Vec &solution)
 
void add_cumulative_source (unsigned int quantity_idx, double source)
 
void output ()
 Perform output to file for given time instant. More...
 

Static Public Member Functions

static const Input::Type::Recordget_input_type ()
 Main balance input record type. More...
 
static const Input::Type::Selectionget_format_selection_input_type ()
 Input selection for file format. More...
 
static void set_yaml_output ()
 Set global variable to output balance files into YAML format (in addition to the table format). More...
 

Private Member Functions

void lazy_initialize ()
 
void output_legacy (double time)
 Perform output in old format (for compatibility) More...
 
void output_csv (double time, char delimiter, const std::string &comment_string, unsigned int repeat=0)
 Perform output in csv format. More...
 
void output_yaml (double time)
 Perform output in yaml format. More...
 
std::string csv_zero_vals (unsigned int cnt, char delimiter)
 Return part of output represented by zero values. Count of zero values is given by cnt parameter. More...
 
void format_csv_output_header (char delimiter, const std::string &comment_string)
 Print output header. More...
 
std::string format_csv_val (std::string val, char delimiter, bool initial=false)
 Format string value of csv output. Wrap string into quotes and if delimiter is space, align text to column. More...
 
std::string format_csv_val (double val, char delimiter, bool initial=false)
 Format double value of csv output. If delimiter is space, align text to column. More...
 
LongIdx get_boundary_edge_uid (SideIter side)
 

Private Attributes

unsigned int n_loc_dofs_par_
 Allocation parameters. Set by the allocate method used in the lazy_initialize. More...
 
unsigned int n_loc_dofs_seq_
 
unsigned int max_dofs_per_boundary_
 
std::string file_prefix_
 Save prefix passed in in constructor. More...
 
FilePath balance_output_file_
 File path for output_ stream. More...
 
std::ofstream output_
 Handle for file for output in given OutputFormat of balance and total fluxes over individual regions and region sets. More...
 
std::ofstream output_yaml_
 
OutputFormat output_format_
 Format of output file. More...
 
std::vector< Quantityquantities_
 Names of conserved quantities. More...
 
const Meshmesh_
 
UnitSI units_
 Units of conserved quantities. More...
 
Mat * region_mass_matrix_
 Matrices for calculation of mass (n_dofs x n_bulk_regions). More...
 
Mat * be_flux_matrix_
 Matrices for calculation of flux (n_boundary_edges x n_dofs). More...
 
Mat * region_source_matrix_
 Matrices for calculation of source (n_dofs x n_bulk_regions). More...
 
Mat * region_source_rhs_
 Matrices for calculation of signed source (n_dofs x n_bulk_regions). More...
 
Vec * be_flux_vec_
 Vectors for calculation of flux (n_boundary_edges). More...
 
Vec * region_mass_vec_
 Vectors for calculation of mass (n_bulk_regions). More...
 
std::unordered_map< LongIdx, unsigned int > be_id_map_
 
std::vector< unsigned int > be_regions_
 Maps local boundary edge to its region boundary index. More...
 
int be_offset_
 Offset for local part of vector of boundary edges. More...
 
std::vector< std::vector< double > > fluxes_
 
std::vector< std::vector< double > > fluxes_in_
 
std::vector< std::vector< double > > fluxes_out_
 
std::vector< std::vector< double > > masses_
 
std::vector< std::vector< double > > sources_in_
 
std::vector< std::vector< double > > sources_out_
 
std::vector< double > sum_fluxes_
 
std::vector< double > sum_fluxes_in_
 
std::vector< double > sum_fluxes_out_
 
std::vector< double > sum_masses_
 
std::vector< double > sum_sources_
 
std::vector< double > sum_sources_in_
 
std::vector< double > sum_sources_out_
 
std::vector< double > initial_mass_
 
std::vector< double > integrated_sources_
 
std::vector< double > integrated_fluxes_
 
std::vector< double > increment_fluxes_
 
std::vector< double > increment_sources_
 
double last_time_
 time of last calculated balance More...
 
TimeMark::Type balance_output_type_
 TimeMark type for balance output of particular equation. More...
 
TimeMark::Type output_mark_type_
 TimeMark type for output of particular equation. More...
 
bool initial_
 true before calculating the mass at initial time, otherwise false More...
 
bool cumulative_
 if true then cumulative balance is computed More...
 
bool allocation_done_
 true before allocating necessary internal structures (Petsc matrices etc.) More...
 
bool balance_on_
 If the balance is on. Balance is off in the case of no balance output time marks. More...
 
bool add_output_times_
 Add output time marks to balance output time marks. More...
 
int rank_
 MPI rank. More...
 
unsigned int output_line_counter_
 hold count of line printed into output_ More...
 
bool output_yaml_header_
 marks whether YAML output has printed header More...
 
Input::Record input_record_
 Record for current balance. More...
 
const TimeGovernortime_
 

Static Private Attributes

static const unsigned int output_column_width = 20
 Size of column in output (used if delimiter is space) More...
 
static bool do_yaml_output_ = true
 

Detailed Description

Design of balance class - serves as storage and writer. Equations themselves call methods of Balance that add/modify mass, source and flux of various quantities and generate output.

One instance of Balance can handle several conservative quantities of the same type (e.g. mass of several substances or their phases).

The mass, flux and source are calculated as follows:

m(q,r) =  ( M'(q) * solution + mv(q) )[r]
f(q,r) = -( R' * ( F(q) * solution + fv(q) ) )[r]
s(q,r) =  ( S'(q) * solution + sv(q) )[r]

where M' stands for matrix transpose,

m(q,r)...mass of q-th substance in region r
f(q,r)...incoming flux of q-th substance in region r
s(q,r)...source of q-th substance in region r

and

M(q)...region_mass_matrix_      n_dofs x n_bulk_regions
F(q)...be_flux_matrix_          n_boundary_edges x n_dofs
S(q)...region_source_matrix_    n_dofs x n_bulk_regions
SV(q)..region_source_rhs_       n_dofs x n_bulk_regions

mv(q)..region_mass_vec_ n_bulk_regions fv(q)..be_flux_vec_ n_boundary_edges sv(q)..region_source_vec_ n_bulk_regions R......region_be_matrix_ n_boundary_edges x n_boundary_regions

Remark: Matrix F and the vector fv are such that F*solution+fv produces outcoming fluxes per boundary edge. However we write to output incoming flux due to users' convention and consistently with input interface.

Note that it holds:

sv(q) = column sum of SV(q)

Except for that, we also provide information on positive/negative flux and source:

fp(q,r) = ( R' * EFP(q) )[r],   EFP(q)[e] = max{ 0, ( F(q) * solution + fv(q) )[e] }
fn(q,r) = ( R' * EFN(q) )[r],   EFN(q)[e] = min{ 0, ( F(q) * solution + fv(q) )[e] }
sp(q,r) = sum_{i in DOFS } max{ 0, ( S(q)[i,r] * solution[i] + SV(q)[i,r] ) }
sn(q,r) = sum_{i in DOFS } min{ 0, ( S(q)[i,r] * solution[i] + SV(q)[i,r] ) }

where

fp(q,r)...positive (inward) flux of q-th quantity in region r
fn(q,r)...negative (outward) flux of q-th quantity in region r
sp(q,r)...positive source (spring) of q-th quantity in region r
sn(q,r)...negative source (sink) of q-th quantity in region r

Remark: The matrix R is needed only for calculation of signed fluxes. The reason is that to determine sign, we decompose flux to sum of local contributions per each boundary element and check its sign. It is not possible to decompose flux using shape functions, since their normal derivatives may have any sign.

Output values (if not relevant, zero is supplied):

#time bulk_region quantity 0 0 0 mass source source_in source_out 0 0 0 #time boundary_region quantity flux flux_in flux_out 0 0 0 0 0 0 0 #time ALL quantity flux flux_in flux_out mass source source_in source_out integrated_flux integrated_source error

error = current_mass - (initial_mass + integrated_source - integrated_flux)

Definition at line 119 of file balance.hh.

Member Enumeration Documentation

◆ OutputFormat

Possible formats of output file.

Enumerator
legacy 

legacy

txt 

csv

gnuplot 

gnuplot

Definition at line 146 of file balance.hh.

Constructor & Destructor Documentation

◆ Balance()

Balance::Balance ( const std::string &  file_prefix,
const Mesh mesh 
)

Constructor.

Parameters
file_prefixPrefix of output file name.
meshMesh.

Definition at line 89 of file balance.cc.

◆ ~Balance()

Balance::~Balance ( )

Destructor.

Definition at line 106 of file balance.cc.

Member Function Documentation

◆ add_cumulative_source()

void Balance::add_cumulative_source ( unsigned int  quantity_idx,
double  source 
)

Adds provided values to the cumulative sources.

Parameters
quantity_idxIndex of quantity.
sourcesSources per region.
dtActual time step.

Definition at line 542 of file balance.cc.

◆ add_flux_values()

void Balance::add_flux_values ( unsigned int  quantity_idx,
const DHCellSide side,
const LocDofVec loc_dof_indices,
const std::vector< double > &  mat_values,
double  vec_value 
)

Adds elements into matrix for computing (outgoing) flux. The flux matrix F is in format [n_boundary_edges x n_dofs] for each quantity. The flux vector fv is in format [n_boundary_edges] for each quantity. See class Balance description above for details.

Parameters
quantity_idxIndex of quantity.
sideDHCellSide iterator.
loc_dof_indicesLocal dof indices (to the solution vector) to be added.
mat_valuesValues to be added into matrix F.
vec_valueValue to be added into vector fv.

Definition at line 479 of file balance.cc.

◆ add_mass_values()

void Balance::add_mass_values ( unsigned int  quantity_idx,
const DHCellAccessor dh_cell,
const LocDofVec loc_dof_indices,
const std::vector< double > &  mat_values,
double  vec_value 
)

Adds elements into matrix for computing mass. The mass matrix M is in format [n_dofs x n_bulk_regions] for each quantity. The mass vector mv is in format [n_bulk_regions] for each quantity. The region mass is later computed using transpose multiplication (M'*u + mv)[r]. See class Balance description above for details.

Parameters
quantity_idxIndex of quantity.
dh_cellDofhandler cell accessor.
loc_dof_indicesLocal dof indices (to the solution vector) to be added.
mat_valuesValues to be added into matrix M.
vec_valueValue to be added into vector mv.

Definition at line 448 of file balance.cc.

◆ add_quantities()

std::vector< unsigned int > Balance::add_quantities ( const std::vector< string > &  names)

Define a set of conservative quantities.

Parameters
namesList of quantities' names.
Returns
List of quantities' indices.

Definition at line 176 of file balance.cc.

◆ add_quantity()

unsigned int Balance::add_quantity ( const string &  name)

Define a single conservative quantity.

Parameters
nameName of the quantity.
Returns
Quantity's internal index.

Definition at line 165 of file balance.cc.

Here is the caller graph for this function:

◆ add_source_values()

void Balance::add_source_values ( unsigned int  quantity_idx,
unsigned int  region_idx,
const LocDofVec loc_dof_indices,
const std::vector< double > &  mult_mat_values,
const std::vector< double > &  add_mat_values 
)

Adds elements into matrix and vector for computing source. The source region matrix S is in format [n_dofs x n_bulk_regions] for each quantity (multiplies solution). The source region rhs matrix SV is in format [n_dofs x n_bulk_regions] for each quantity (addition). The region source is then computed as (S'(q) * solution + sv(q))[r], while sv(q) being column sum of SV(q). The source term for a specific region and quantity is then

Parameters
quantity_idxIndex of quantity.
region_idxIndex of bulk region.
loc_dof_indicesLocal dof indices (to the solution vector) to be added.
mult_mat_valuesValues to be added into matrix S.
add_mat_valuesValues to be added into matrix SV.

Definition at line 513 of file balance.cc.

◆ allocate() [1/2]

void Balance::allocate ( const std::shared_ptr< DOFHandlerMultiDim > &  dh,
unsigned int  max_dofs_per_boundary 
)

Allocates matrices and vectors for balance based on DofHandler.

Parameters
dhDofHandler of the corresponding linear system (and solution).
max_dofs_per_boundaryNumber of dofs contributing to one boundary edge.

Definition at line 194 of file balance.cc.

◆ allocate() [2/2]

void Balance::allocate ( unsigned int  n_loc_dofs,
unsigned int  max_dofs_per_boundary 
)

Allocates matrices and vectors for balance.

Parameters
n_loc_dofsNumber of solution dofs on the local process.
max_dofs_per_boundaryNumber of dofs contributing to one boundary edge.

Definition at line 185 of file balance.cc.

◆ calculate_cumulative()

void Balance::calculate_cumulative ( unsigned int  quantity_idx,
const Vec &  solution 
)

Updates cumulative quantities for balance. This method can be called in substeps even if no output is generated. It calculates the sum of source and sum of (incoming) flux over time interval.

Parameters
quantity_idxIndex of quantity.
solutionSolution vector.

Definition at line 552 of file balance.cc.

◆ calculate_instant()

void Balance::calculate_instant ( unsigned int  quantity_idx,
const Vec &  solution 
)

Calculates actual mass, incoming flux and source.

Parameters
quantity_idxIndex of quantity.
solutionSolution vector.

Definition at line 636 of file balance.cc.

◆ calculate_mass()

void Balance::calculate_mass ( unsigned int  quantity_idx,
const Vec &  solution,
vector< double > &  output_array 
)

Calculates actual mass and save it to given vector.

Parameters
quantity_idxIndex of quantity.
solutionSolution vector.
output_arrayVector of output masses per region.

Definition at line 611 of file balance.cc.

Here is the caller graph for this function:

◆ csv_zero_vals()

std::string Balance::csv_zero_vals ( unsigned int  cnt,
char  delimiter 
)
private

Return part of output represented by zero values. Count of zero values is given by cnt parameter.

Definition at line 1002 of file balance.cc.

Here is the caller graph for this function:

◆ cumulative()

bool Balance::cumulative ( ) const
inline

Getter for cumulative_.

Definition at line 190 of file balance.hh.

◆ finish_flux_assembly() [1/2]

void Balance::finish_flux_assembly ( std::vector< unsigned int >  q_idx_vec)
inline

Variant of the finish_flux_assembly() method for a set of quantities.

Definition at line 334 of file balance.hh.

◆ finish_flux_assembly() [2/2]

void Balance::finish_flux_assembly ( unsigned int  quantity_idx)

This method must be called after assembling the matrix and vector for computing flux.

Definition at line 425 of file balance.cc.

Here is the caller graph for this function:

◆ finish_mass_assembly() [1/2]

void Balance::finish_mass_assembly ( std::vector< unsigned int >  q_idx_vec)
inline

Variant of the finish_mass_assembly() method for a set of quantities.

Definition at line 324 of file balance.hh.

◆ finish_mass_assembly() [2/2]

void Balance::finish_mass_assembly ( unsigned int  quantity_idx)

This method must be called after assembling the matrix for computing mass.

Definition at line 414 of file balance.cc.

Here is the caller graph for this function:

◆ finish_source_assembly() [1/2]

void Balance::finish_source_assembly ( std::vector< unsigned int >  q_idx_vec)
inline

Variant of the finish_source_assembly() method for a set of quantities.

Definition at line 344 of file balance.hh.

◆ finish_source_assembly() [2/2]

void Balance::finish_source_assembly ( unsigned int  quantity_idx)

This method must be called after assembling the matrix and vectors for computing source.

Definition at line 436 of file balance.cc.

Here is the caller graph for this function:

◆ format_csv_output_header()

void Balance::format_csv_output_header ( char  delimiter,
const std::string &  comment_string 
)
private

Print output header.

Definition at line 1096 of file balance.cc.

Here is the caller graph for this function:

◆ format_csv_val() [1/2]

std::string Balance::format_csv_val ( double  val,
char  delimiter,
bool  initial = false 
)
private

Format double value of csv output. If delimiter is space, align text to column.

Definition at line 1141 of file balance.cc.

◆ format_csv_val() [2/2]

std::string Balance::format_csv_val ( std::string  val,
char  delimiter,
bool  initial = false 
)
private

Format string value of csv output. Wrap string into quotes and if delimiter is space, align text to column.

Definition at line 1124 of file balance.cc.

Here is the caller graph for this function:

◆ get_boundary_edge_uid()

LongIdx Balance::get_boundary_edge_uid ( SideIter  side)
inlineprivate

Computes unique id of local boundary edge from local element index and element side index

Parameters
sideis a side of locally owned element

Definition at line 426 of file balance.hh.

Here is the caller graph for this function:

◆ get_format_selection_input_type()

const Selection & Balance::get_format_selection_input_type ( )
static

Input selection for file format.

Definition at line 42 of file balance.cc.

Here is the caller graph for this function:

◆ get_input_type()

const Record & Balance::get_input_type ( )
static

Main balance input record type.

Definition at line 53 of file balance.cc.

Here is the caller graph for this function:

◆ init_from_input()

void Balance::init_from_input ( const Input::Record in_rec,
TimeGovernor tg 
)

Initialize the balance object according to the input. The balance output time marks are set according to the already existing output time marks of the same equation. So, this method must be called after Output::

Parameters
in_recInput record of balance.
tgTimeGovernor of the equation. We need just equation mark type.

Definition at line 132 of file balance.cc.

◆ is_current()

bool Balance::is_current ( )

Returns true if the current time step is marked for the balance output.

Definition at line 376 of file balance.cc.

Here is the caller graph for this function:

◆ lazy_initialize()

void Balance::lazy_initialize ( )
private

Postponed allocation and initialization to allow calling setters in arbitrary order. In particular we need to perform adding of output times after the output time marks are set. On the other hand we need to read the input before we make the allocation.

The initialization is done during the first call of any start_*_assembly method.

Definition at line 205 of file balance.cc.

Here is the caller graph for this function:

◆ output()

void Balance::output ( )

Perform output to file for given time instant.

Definition at line 713 of file balance.cc.

◆ output_csv()

void Balance::output_csv ( double  time,
char  delimiter,
const std::string &  comment_string,
unsigned int  repeat = 0 
)
private

Perform output in csv format.

Definition at line 1010 of file balance.cc.

Here is the caller graph for this function:

◆ output_legacy()

void Balance::output_legacy ( double  time)
private

Perform output in old format (for compatibility)

Definition at line 856 of file balance.cc.

Here is the caller graph for this function:

◆ output_yaml()

void Balance::output_yaml ( double  time)
private

Perform output in yaml format.

Definition at line 1154 of file balance.cc.

Here is the caller graph for this function:

◆ set_yaml_output()

void Balance::set_yaml_output ( )
static

Set global variable to output balance files into YAML format (in addition to the table format).

Definition at line 66 of file balance.cc.

Here is the caller graph for this function:

◆ start_flux_assembly() [1/2]

void Balance::start_flux_assembly ( std::vector< unsigned int >  q_idx_vec)
inline

Variant of the start_flux_assembly() method for a set of quantities.

Definition at line 246 of file balance.hh.

◆ start_flux_assembly() [2/2]

void Balance::start_flux_assembly ( unsigned int  quantity_idx)

This method must be called before assembling the matrix and vector for fluxes. It actually erases the matrix and vector.

Definition at line 396 of file balance.cc.

Here is the caller graph for this function:

◆ start_mass_assembly() [1/2]

void Balance::start_mass_assembly ( std::vector< unsigned int >  q_idx_vec)
inline

Variant of the start_mass_assembly() method for a set of quantities.

Definition at line 233 of file balance.hh.

◆ start_mass_assembly() [2/2]

void Balance::start_mass_assembly ( unsigned int  quantity_idx)

This method must be called before assembling the matrix for computing mass. It actually erases the matrix.

Definition at line 387 of file balance.cc.

Here is the caller graph for this function:

◆ start_source_assembly() [1/2]

void Balance::start_source_assembly ( std::vector< unsigned int >  q_idx_vec)
inline

Variant of the start_source_assembly() method for a set of quantities.

Definition at line 259 of file balance.hh.

◆ start_source_assembly() [2/2]

void Balance::start_source_assembly ( unsigned int  quantity_idx)

This method must be called before assembling the matrix and vectors for sources. It actually erases the matrix and vectors.

Definition at line 405 of file balance.cc.

Here is the caller graph for this function:

◆ units()

void Balance::units ( const UnitSI unit)

Setter for units of conserved quantities.

Definition at line 157 of file balance.cc.

Member Data Documentation

◆ add_output_times_

bool Balance::add_output_times_
private

Add output time marks to balance output time marks.

Definition at line 542 of file balance.hh.

◆ allocation_done_

bool Balance::allocation_done_
private

true before allocating necessary internal structures (Petsc matrices etc.)

Definition at line 536 of file balance.hh.

◆ balance_on_

bool Balance::balance_on_
private

If the balance is on. Balance is off in the case of no balance output time marks.

Definition at line 539 of file balance.hh.

◆ balance_output_file_

FilePath Balance::balance_output_file_
private

File path for output_ stream.

Definition at line 443 of file balance.hh.

◆ balance_output_type_

TimeMark::Type Balance::balance_output_type_
private

TimeMark type for balance output of particular equation.

Definition at line 524 of file balance.hh.

◆ be_flux_matrix_

Mat* Balance::be_flux_matrix_
private

Matrices for calculation of flux (n_boundary_edges x n_dofs).

Definition at line 467 of file balance.hh.

◆ be_flux_vec_

Vec* Balance::be_flux_vec_
private

Vectors for calculation of flux (n_boundary_edges).

Definition at line 476 of file balance.hh.

◆ be_id_map_

std::unordered_map<LongIdx, unsigned int> Balance::be_id_map_
private

Maps unique identifier of (local bulk element idx, side idx) returned by get_boundary_edge_uid(side) to local boundary edge. Example usage: be_id = be_id_map_[get_boundary_edge_uid(side)]

Definition at line 486 of file balance.hh.

◆ be_offset_

int Balance::be_offset_
private

Offset for local part of vector of boundary edges.

Definition at line 492 of file balance.hh.

◆ be_regions_

std::vector<unsigned int> Balance::be_regions_
private

Maps local boundary edge to its region boundary index.

Definition at line 489 of file balance.hh.

◆ cumulative_

bool Balance::cumulative_
private

if true then cumulative balance is computed

Definition at line 533 of file balance.hh.

◆ do_yaml_output_

bool Balance::do_yaml_output_ = true
staticprivate

Definition at line 431 of file balance.hh.

◆ file_prefix_

std::string Balance::file_prefix_
private

Save prefix passed in in constructor.

Definition at line 440 of file balance.hh.

◆ fluxes_

std::vector<std::vector<double> > Balance::fluxes_
private

Definition at line 497 of file balance.hh.

◆ fluxes_in_

std::vector<std::vector<double> > Balance::fluxes_in_
private

Definition at line 498 of file balance.hh.

◆ fluxes_out_

std::vector<std::vector<double> > Balance::fluxes_out_
private

Definition at line 499 of file balance.hh.

◆ increment_fluxes_

std::vector<double> Balance::increment_fluxes_
private

Definition at line 517 of file balance.hh.

◆ increment_sources_

std::vector<double> Balance::increment_sources_
private

Definition at line 518 of file balance.hh.

◆ initial_

bool Balance::initial_
private

true before calculating the mass at initial time, otherwise false

Definition at line 530 of file balance.hh.

◆ initial_mass_

std::vector<double> Balance::initial_mass_
private

Definition at line 512 of file balance.hh.

◆ input_record_

Input::Record Balance::input_record_
private

Record for current balance.

Definition at line 555 of file balance.hh.

◆ integrated_fluxes_

std::vector<double> Balance::integrated_fluxes_
private

Definition at line 516 of file balance.hh.

◆ integrated_sources_

std::vector<double> Balance::integrated_sources_
private

Definition at line 515 of file balance.hh.

◆ last_time_

double Balance::last_time_
private

time of last calculated balance

Definition at line 521 of file balance.hh.

◆ masses_

std::vector<std::vector<double> > Balance::masses_
private

Definition at line 500 of file balance.hh.

◆ max_dofs_per_boundary_

unsigned int Balance::max_dofs_per_boundary_
private

Definition at line 436 of file balance.hh.

◆ mesh_

const Mesh* Balance::mesh_
private

Definition at line 457 of file balance.hh.

◆ n_loc_dofs_par_

unsigned int Balance::n_loc_dofs_par_
private

Allocation parameters. Set by the allocate method used in the lazy_initialize.

Definition at line 434 of file balance.hh.

◆ n_loc_dofs_seq_

unsigned int Balance::n_loc_dofs_seq_
private

Definition at line 435 of file balance.hh.

◆ output_

std::ofstream Balance::output_
private

Handle for file for output in given OutputFormat of balance and total fluxes over individual regions and region sets.

Definition at line 446 of file balance.hh.

◆ output_column_width

const unsigned int Balance::output_column_width = 20
staticprivate

Size of column in output (used if delimiter is space)

Definition at line 391 of file balance.hh.

◆ output_format_

OutputFormat Balance::output_format_
private

Format of output file.

Definition at line 452 of file balance.hh.

◆ output_line_counter_

unsigned int Balance::output_line_counter_
private

hold count of line printed into output_

Definition at line 549 of file balance.hh.

◆ output_mark_type_

TimeMark::Type Balance::output_mark_type_
private

TimeMark type for output of particular equation.

Definition at line 527 of file balance.hh.

◆ output_yaml_

std::ofstream Balance::output_yaml_
private

Definition at line 449 of file balance.hh.

◆ output_yaml_header_

bool Balance::output_yaml_header_
private

marks whether YAML output has printed header

Definition at line 552 of file balance.hh.

◆ quantities_

std::vector<Quantity> Balance::quantities_
private

Names of conserved quantities.

Definition at line 455 of file balance.hh.

◆ rank_

int Balance::rank_
private

MPI rank.

Definition at line 546 of file balance.hh.

◆ region_mass_matrix_

Mat* Balance::region_mass_matrix_
private

Matrices for calculation of mass (n_dofs x n_bulk_regions).

Definition at line 464 of file balance.hh.

◆ region_mass_vec_

Vec* Balance::region_mass_vec_
private

Vectors for calculation of mass (n_bulk_regions).

Definition at line 479 of file balance.hh.

◆ region_source_matrix_

Mat* Balance::region_source_matrix_
private

Matrices for calculation of source (n_dofs x n_bulk_regions).

Definition at line 470 of file balance.hh.

◆ region_source_rhs_

Mat* Balance::region_source_rhs_
private

Matrices for calculation of signed source (n_dofs x n_bulk_regions).

Definition at line 473 of file balance.hh.

◆ sources_in_

std::vector<std::vector<double> > Balance::sources_in_
private

Definition at line 501 of file balance.hh.

◆ sources_out_

std::vector<std::vector<double> > Balance::sources_out_
private

Definition at line 502 of file balance.hh.

◆ sum_fluxes_

std::vector<double> Balance::sum_fluxes_
private

Definition at line 505 of file balance.hh.

◆ sum_fluxes_in_

std::vector<double> Balance::sum_fluxes_in_
private

Definition at line 506 of file balance.hh.

◆ sum_fluxes_out_

std::vector<double> Balance::sum_fluxes_out_
private

Definition at line 507 of file balance.hh.

◆ sum_masses_

std::vector<double> Balance::sum_masses_
private

Definition at line 508 of file balance.hh.

◆ sum_sources_

std::vector<double> Balance::sum_sources_
private

Definition at line 509 of file balance.hh.

◆ sum_sources_in_

std::vector<double> Balance::sum_sources_in_
private

Definition at line 510 of file balance.hh.

◆ sum_sources_out_

std::vector<double> Balance::sum_sources_out_
private

Definition at line 511 of file balance.hh.

◆ time_

const TimeGovernor* Balance::time_
private

Definition at line 557 of file balance.hh.

◆ units_

UnitSI Balance::units_
private

Units of conserved quantities.

Definition at line 460 of file balance.hh.


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