21 #include <unordered_map> 41 return Selection(
"Balance_output_format",
"Format of output file for balance.")
49 return Record(
"Balance",
"Balance of a conservative quantity, boundary fluxes and sources.")
51 .
declare_key(
"add_output_times",
Bool(),
Default(
"true"),
"Add all output times of the balanced equation to the balance output times set. " 52 "Note that this is not the time set of the output stream.")
56 "If true, then balance is calculated at each computational time step, which can slow down the program.")
85 : file_prefix_(file_prefix),
89 allocation_done_(false),
91 output_line_counter_(0),
92 output_yaml_header_(false)
174 for (
auto name : names)
181 unsigned int max_dofs_per_boundary)
214 unsigned int be_id = 0;
220 for(
unsigned int si=0; si<elm->
n_sides(); si++)
268 for (
unsigned int c=0; c<n_quant; ++c)
270 chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
275 (
rank_==0)?n_bulk_regs_per_dof:0,
277 (rank_==0)?0:n_bulk_regs_per_dof,
281 chkerr(MatCreateSeqAIJ(PETSC_COMM_SELF,
288 chkerr(MatCreateSeqAIJ(PETSC_COMM_SELF,
295 chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
306 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
311 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
322 std::string default_file_name;
429 unsigned int region_idx,
436 PetscInt reg_array[1] = { (int)region_idx };
467 unsigned int region_idx,
475 PetscInt reg_array[1] = { (int)region_idx };
478 loc_dof_indices.size(),
479 &(loc_dof_indices[0]),
486 loc_dof_indices.size(),
487 &(loc_dof_indices[0]),
495 unsigned int region_idx,
537 double temp_source = 0;
538 int lsize, n_cols_mat, n_cols_rhs;
540 const double *vals_mat, *vals_rhs, *sol_array;
541 chkerr(VecGetLocalSize(solution, &lsize));
542 chkerr(VecGetArrayRead(solution, &sol_array));
548 for (
int i=0; i<lsize; ++i){
554 for (
int j=0; j<n_cols_mat; ++j)
555 temp_source += vals_mat[j]*sol_array[i] + vals_rhs[j];
560 chkerr(VecRestoreArrayRead(solution, &sol_array));
567 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
575 chkerr(VecSum(temp, &sum_fluxes));
576 chkerr(VecDestroy(&temp));
595 chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD,
603 chkerr(VecZeroEntries(bulk_vec));
608 chkerr(VecDestroy(&bulk_vec));
624 int lsize, n_cols_mat, n_cols_rhs;
626 const double *vals_mat, *vals_rhs, *sol_array;
627 chkerr(VecGetLocalSize(solution, &lsize));
628 chkerr(VecGetArrayRead(solution, &sol_array));
633 for (
int i=0; i<lsize; ++i)
641 for (
int j=0; j<n_cols_mat; ++j)
645 double f = vals_mat[j]*sol_array[i] + vals_rhs[j];
650 chkerr(VecRestoreArrayRead(solution, &sol_array));
654 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
662 chkerr(VecScale(temp, -1));
667 const double *flux_array;
669 chkerr(VecGetArrayRead(temp, &flux_array));
670 chkerr(VecGetLocalSize(temp, &lsize));
671 for (
int e=0; e<lsize; ++e)
673 if (flux_array[e] < 0)
676 fluxes_in_[quantity_idx][be_regions_[e]] += flux_array[e];
678 chkerr(VecRestoreArrayRead(temp, &flux_array));
679 chkerr(VecDestroy(&temp));
696 const int buf_size = n_quant*2*n_blk_reg + n_quant*2*n_bdr_reg + n_quant;
697 double sendbuffer[buf_size], recvbuffer[buf_size];
698 for (
unsigned int qi=0; qi<n_quant; qi++)
700 for (
unsigned int ri=0; ri<n_blk_reg; ri++)
702 sendbuffer[qi*2*n_blk_reg + + ri] =
sources_in_[qi][ri];
703 sendbuffer[qi*2*n_blk_reg + n_blk_reg + ri] =
sources_out_[qi][ri];
705 for (
unsigned int ri=0; ri<n_bdr_reg; ri++)
707 sendbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + + ri] =
fluxes_in_[qi][ri];
708 sendbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + n_bdr_reg + ri] =
fluxes_out_[qi][ri];
723 for (
unsigned int qi=0; qi<n_quant; qi++)
725 for (
unsigned int ri=0; ri<n_blk_reg; ri++)
727 sources_in_[qi][ri] = recvbuffer[qi*2*n_blk_reg + + ri];
728 sources_out_[qi][ri] = recvbuffer[qi*2*n_blk_reg + n_blk_reg + ri];
730 for (
unsigned int ri=0; ri<n_bdr_reg; ri++)
732 fluxes_in_[qi][ri] = recvbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + + ri];
733 fluxes_out_[qi][ri] = recvbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + n_bdr_reg + ri];
757 for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
759 for (
unsigned int qi=0; qi<n_quant; qi++)
769 for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
771 for (
unsigned int qi=0; qi<n_quant; qi++)
787 for (
unsigned int qi=0; qi<n_quant; qi++)
792 for (
unsigned int qi=0; qi<n_quant; qi++)
832 if (
rank_ != 0)
return;
839 unsigned int wl = 2*(w-5)+7;
840 string bc_head_format =
"# %-*s%-*s%-*s%-*s%-*s%-*s\n",
841 bc_format =
"%*s%-*d%-*s%-*s%-*g%-*g%-*g\n",
842 bc_total_format =
"# %-*s%-*s%-*g%-*g%-*g\n";
844 output_ <<
"# " << setw((w*c+wl-14)/2) << setfill(
'-') <<
"--" 846 << setw((w*c+wl-14)/2) << setfill(
'-') <<
"" << endl
850 output_ <<
"# Mass flux through boundary [M/T]:\n# " 851 << setiosflags(ios::left) << setfill(
' ')
852 << setw(w) <<
"[boundary_id]" 853 << setw(wl) <<
"[label]" 854 << setw(w) <<
"[substance]" 855 << setw(w) <<
"[total flux]" 856 << setw(w) <<
"[outward flux]" 857 << setw(w) <<
"[inward flux]" 861 output_ <<
"# " << setw(w*c+wl) << setfill(
'-') <<
"" << setfill(
' ') << endl;
865 for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg) {
866 for (
unsigned int qi=0; qi<n_quant; qi++) {
868 << setw(w) << (int)reg->id()
869 << setw(wl) << reg->label().c_str()
873 << setw(w) <<
fluxes_in_[qi][reg->boundary_idx()]
879 output_ <<
"# " << setw(w*c+wl) << setfill(
'-') <<
"" << setfill(
' ') << endl;
882 for (
unsigned int qi=0; qi<n_quant; qi++)
883 output_ <<
"# " << setiosflags(ios::left)
884 << setw(w+wl) <<
"Total mass flux of substance [M/T]" 894 string src_head_format =
"# %-*s%-*s%-*s%-*s%-*s\n",
895 src_format =
"%*s%-*d%-*s%-*s%-*g%-*g\n",
896 src_total_format =
"# %-*s%-*s%-*g%-*g\n";
897 output_ <<
"# Mass [M] and sources [M/T] on regions:\n" 898 <<
"# " << setiosflags(ios::left)
899 << setw(w) <<
"[region_id]" 900 << setw(wl) <<
"[label]" 901 << setw(w) <<
"[substance]" 902 << setw(w) <<
"[total_mass]" 903 << setw(w) <<
"[total_source]" 907 output_ <<
"# " << setw(w*c+wl) << setfill(
'-') <<
"" << setfill(
' ') << endl;
911 for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
913 for (
unsigned int qi=0; qi<n_quant; qi++)
916 << setw(w) << (int)reg->id()
917 << setw(wl) << reg->label().c_str()
919 << setw(w) <<
masses_[qi][reg->bulk_idx()]
926 output_ <<
"# " << setw(w*c+wl) << setfill(
'-') <<
"" << setfill(
' ') << endl;
929 for (
unsigned int qi=0; qi<n_quant; qi++)
930 output_ <<
"# " << setiosflags(ios::left) << setw(w+wl) <<
"Total mass [M] and sources [M/T]" 939 output_ <<
"\n\n# Cumulative mass balance on time interval [" 941 << setiosflags(ios::left) << time <<
"]\n" 942 <<
"# Initial mass [M] + sources integrated over time [M] - flux integrated over time [M] = current mass [M]\n" 943 <<
"# " << setiosflags(ios::left)
944 << setw(w) <<
"[substance]" 945 << setw(w) <<
"[A=init. mass]" 946 << setw(w) <<
"[B=source]" 947 << setw(w) <<
"[C=flux]" 948 << setw(w) <<
"[A+B-C]" 949 << setw(w) <<
"[D=curr. mass]" 950 << setw(w) <<
"[A+B-C-D=err.]" 951 << setw(w) <<
"[rel. error]" 954 for (
unsigned int qi=0; qi<n_quant; qi++)
957 output_ <<
" " << setiosflags(ios::left)
976 std::stringstream ss;
977 for (
unsigned int i=0; i<cnt; i++) ss <<
format_csv_val(0, delimiter);
982 void Balance::output_csv(
double time,
char delimiter,
const std::string& comment_string,
unsigned int repeat)
985 if (
rank_ != 0)
return;
994 for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
996 for (
unsigned int qi=0; qi<n_quant; qi++)
1016 for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
1018 for (
unsigned int qi=0; qi<n_quant; qi++) {
1035 for (
unsigned int qi=0; qi<n_quant; qi++)
1065 std::stringstream ss;
1066 if (delimiter ==
' ') {
1072 output_ << comment_string << ss.str()
1092 std::stringstream ss;
1093 std::replace( val.begin(), val.end(),
'\"',
'\'');
1095 if (!initial) ss << delimiter;
1096 if (delimiter ==
' ') {
1097 std::stringstream sval;
1098 sval <<
"\"" << val <<
"\"";
1101 ss <<
"\"" << val <<
"\"";
1109 std::stringstream ss;
1111 if (!initial) ss << delimiter;
1112 if (delimiter ==
' ') {
1130 output_yaml_ <<
"column_names: [ flux, flux_in, flux_out, mass, source, source_in, source_out, flux_increment, " 1131 <<
"source_increment, flux_cumulative, source_cumulative, error ]" << endl;
1140 for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
1142 for (
unsigned int qi=0; qi<n_quant; qi++)
1145 output_yaml_ << setw(4) <<
"" <<
"region: " << reg->label() << endl;
1147 output_yaml_ << setw(4) <<
"" <<
"data: " <<
"[ 0, 0, 0, " <<
masses_[qi][reg->bulk_idx()] <<
", " 1150 <<
", 0, 0, 0, 0, 0 ]" << endl;
1156 for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
1158 for (
unsigned int qi=0; qi<n_quant; qi++) {
1160 output_yaml_ << setw(4) <<
"" <<
"region: " << reg->label() << endl;
1165 <<
", 0, 0, 0, 0, 0, 0, 0, 0, 0 ]" << endl;
1171 for (
unsigned int qi=0; qi<n_quant; qi++)
1175 output_yaml_ << setw(4) <<
"" <<
"region: ALL" << endl;
1183 << error <<
" ]" << endl;
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
UnitSI units_
Units of conserved quantities.
static const Input::Type::Record & get_input_type()
Main balance input record type.
std::vector< double > integrated_sources_
unsigned int add_quantity(const string &name)
std::vector< double > sum_fluxes_out_
void calculate_mass(unsigned int quantity_idx, const Vec &solution, vector< double > &output_array)
RegionSet get_region_set(const std::string &set_name) const
LongIdx get_boundary_edge_uid(SideIter side)
TimeMark::Type output_mark_type_
TimeMark type for output of particular equation.
void allocate(unsigned int n_loc_dofs, unsigned int max_dofs_per_boundary)
TimeMark::Type balance_output_type_
TimeMark type for balance output of particular equation.
unsigned int * boundary_idx_
std::string format_text() const
void finish_source_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix and vectors for computing source.
unsigned int bulk_size() const
void add_mass_vec_value(unsigned int quantity_idx, unsigned int region_idx, double value)
bool initial_
true before calculating the mass at initial time, otherwise false
std::vector< std::vector< double > > sources_in_
void add_cumulative_source(unsigned int quantity_idx, double source)
std::vector< std::vector< double > > fluxes_in_
std::vector< double > sum_fluxes_in_
std::vector< std::vector< double > > fluxes_out_
const TimeGovernor * time_
Mat * region_source_rhs_
Matrices for calculation of signed source (n_dofs x n_bulk_regions).
Input::Record input_record_
Record for current balance.
ofstream output_
Handle for file for output in given OutputFormat of balance and total fluxes over individual regions ...
void read_from_input(Input::Array in_array, const TimeGovernor &tg)
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
void output_legacy(double time)
Perform output in old format (for compatibility)
SideIter side(const unsigned int loc_index)
void add_mass_matrix_values(unsigned int quantity_idx, unsigned int region_idx, const std::vector< LongIdx > &dof_indices, const std::vector< double > &values)
const RegionDB & region_db() const
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
const TimeStep & step(int index=-1) const
bool allocation_done_
true before allocating necessary internal structures (Petsc matrices etc.)
void calculate_cumulative(unsigned int quantity_idx, const Vec &solution)
Basic time management functionality for unsteady (and steady) solvers (class Equation).
static TimeMarks & marks()
#define MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm)
void calculate_instant(unsigned int quantity_idx, const Vec &solution)
Basic time management class.
std::vector< double > sum_fluxes_
virtual ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
bool balance_on_
If the balance is on. Balance is off in the case of no balance output time marks. ...
static constexpr bool value
std::vector< double > increment_sources_
Vec * be_flux_vec_
Vectors for calculation of flux (n_boundary_edges).
static void set_yaml_output()
Set global variable to output balance files into YAML format (in addition to the table format)...
static const Input::Type::Selection & get_format_selection_input_type()
Input selection for file format.
std::vector< double > integrated_fluxes_
unsigned int boundary_idx() const
Returns index of the region in the boundary set.
std::vector< Quantity > quantities_
Names of conserved quantities.
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...
std::vector< std::vector< double > > sources_out_
Mat * be_flux_matrix_
Matrices for calculation of flux (n_boundary_edges x n_dofs).
double last_time_
time of last calculated balance
void open_stream(Stream &stream) const
std::vector< std::vector< double > > fluxes_
bool is_boundary() const
Returns true for side on the boundary.
void finish_mass_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix for computing mass.
TimeMark::Type equation_fixed_mark_type() const
void format_csv_output_header(char delimiter, const std::string &comment_string)
Print output header.
TimeMark::Type equation_mark_type() const
unsigned int n_sides() const
void init_from_input(const Input::Record &in_rec, TimeGovernor &tg)
Mat * region_mass_matrix_
Matrices for calculation of mass (n_dofs x n_bulk_regions).
std::vector< double > sum_sources_out_
unsigned int max_dofs_per_boundary_
bool cumulative_
if true then cumulative balance is computed
Mat * region_source_matrix_
Matrices for calculation of source (n_dofs x n_bulk_regions).
std::vector< double > increment_fluxes_
void chkerr_assert(unsigned int ierr)
bool is_current()
Returns true if the current time step is marked for the balance output.
Balance(const std::string &file_prefix, const Mesh *mesh)
void output_csv(double time, char delimiter, const std::string &comment_string, unsigned int repeat=0)
Perform output in csv format.
Distribution * get_el_ds() const
static bool do_yaml_output_
std::vector< unsigned int > add_quantities(const std::vector< string > &names)
Dedicated class for storing path to input and output files.
Support classes for parallel programing.
void output_yaml(double time)
Perform output in yaml format.
void output()
Perform output to file for given time instant.
std::vector< double > sum_sources_
std::string get_unit_string() const
void start_flux_assembly(unsigned int quantity_idx)
std::vector< double > initial_mass_
static const unsigned int output_column_width
Size of column in output (used if delimiter is space)
void add_source_values(unsigned int quantity_idx, unsigned int region_idx, const std::vector< LongIdx > &loc_dof_indices, const std::vector< double > &mat_values, const std::vector< double > &vec_values)
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
unsigned int boundary_size() const
void finish_flux_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix and vector for computing flux.
std::vector< unsigned int > be_regions_
Maps local boundary edge to its region boundary index.
OutputFormat output_format_
Format of output file.
std::unordered_map< LongIdx, unsigned int > be_id_map_
static const Input::Type::Array get_input_type()
void units(const UnitSI &unit)
Setter for units of conserved quantities.
int be_offset_
Offset for local part of vector of boundary edges.
std::vector< std::vector< double > > masses_
void start_mass_assembly(unsigned int quantity_idx)
std::vector< double > sum_masses_
Vec * region_mass_vec_
Vectors for calculation of mass (n_bulk_regions).
unsigned int output_line_counter_
hold count of line printed into output_
bool add_output_times_
Add output time marks to balance output time marks.
Class for representation SI units of Fields.
unsigned int n_loc_dofs_
Allocation parameters. Set by the allocate method used in the lazy_initialize.
LongIdx * get_el_4_loc() const
void add_flux_matrix_values(unsigned int quantity_idx, SideIter side, const std::vector< LongIdx > &dof_indices, const std::vector< double > &values)
void undef(bool val=true)
void start_source_assembly(unsigned int quantity_idx)
void add_flux_vec_value(unsigned int quantity_idx, SideIter side, double value)
FilePath balance_output_file_
File path for output_ stream.
std::string file_prefix_
Save prefix passed in in constructor.
bool output_yaml_header_
marks whether YAML output has printed header
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.
std::vector< double > sum_sources_in_
unsigned int lsize(int proc) const
get local size