Flow123d
release_3.0.0-968-gc87a28e79
|
Go to the documentation of this file.
21 #include <unordered_map>
42 return Selection(
"Balance_output_format",
"Format of output file for balance.")
50 return Record(
"Balance",
"Balance of a conservative quantity, boundary fluxes and sources.")
52 .
declare_key(
"add_output_times",
Bool(),
Default(
"true"),
"Add all output times of the balanced equation to the balance output times set. "
53 "Note that this is not the time set of the output stream.")
57 "If true, then balance is calculated at each computational time step, which can slow down the program.")
86 : file_prefix_(file_prefix),
90 allocation_done_(false),
92 output_line_counter_(0),
93 output_yaml_header_(false)
181 for (
auto name : names)
188 unsigned int max_dofs_per_boundary)
221 unsigned int be_id = 0;
222 for (
unsigned int loc_el = 0; loc_el <
mesh_->
get_el_ds()->lsize(); loc_el++)
227 for(
unsigned int si=0; si<elm->
n_sides(); si++)
277 for (
unsigned int c=0; c<n_quant; ++c)
279 chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
284 (
rank_==0)?n_bulk_regs_per_dof:0,
286 (
rank_==0)?0:n_bulk_regs_per_dof,
290 chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
295 (
rank_==0)?n_bulk_regs_per_dof:0,
297 (
rank_==0)?0:n_bulk_regs_per_dof,
301 chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
306 (
rank_==0)?n_bulk_regs_per_dof:0,
308 (
rank_==0)?0:n_bulk_regs_per_dof,
312 chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
323 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
328 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
333 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
339 chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
351 for (
unsigned int loc_el=0; loc_el<
be_regions_.size(); ++loc_el)
363 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
371 chkerr(VecCreateMPI(PETSC_COMM_WORLD,
382 std::string default_file_name;
491 unsigned int region_idx,
498 PetscInt reg_array[1] = { (int)region_idx };
530 unsigned int region_idx,
537 PetscInt reg_array[1] = { (int)region_idx };
550 unsigned int region_idx,
575 unsigned int region_idx,
582 PetscInt reg_array[1] = { (int)region_idx };
614 chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD,
622 chkerr(VecZeroEntries(bulk_vec));
627 chkerr(VecSum(bulk_vec, &sum_sources));
628 chkerr(VecDestroy(&bulk_vec));
638 chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD,
646 chkerr(VecZeroEntries(boundary_vec));
652 chkerr(VecScale(temp, -1));
654 chkerr(VecDestroy(&temp));
657 chkerr(VecSum(boundary_vec, &sum_fluxes));
658 chkerr(VecDestroy(&boundary_vec));
675 chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD,
683 chkerr(VecZeroEntries(bulk_vec));
688 chkerr(VecDestroy(&bulk_vec));
700 chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD,
708 chkerr(VecZeroEntries(bulk_vec));
717 const double *sol_array, *mat_array, *rhs_array;
718 chkerr(VecGetLocalSize(solution, &lsize));
719 chkerr(VecDuplicate(solution, &mat_r));
720 chkerr(VecDuplicate(solution, &rhs_r));
721 chkerr(VecGetArrayRead(solution, &sol_array));
727 VecGetArrayRead(mat_r, &mat_array);
728 VecGetArrayRead(rhs_r, &rhs_array);
732 for (
int i=0; i<lsize; ++i)
734 double f = mat_array[i]*sol_array[i] + rhs_array[i];
739 VecRestoreArrayRead(mat_r, &mat_array);
740 VecRestoreArrayRead(rhs_r, &rhs_array);
742 chkerr(VecRestoreArrayRead(solution, &sol_array));
743 chkerr(VecDestroy(&rhs_r));
744 chkerr(VecDestroy(&mat_r));
745 chkerr(VecDestroy(&bulk_vec));
751 chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1,
753 PETSC_DECIDE, &(
fluxes_[quantity_idx][0]), &boundary_vec));
756 chkerr(VecZeroEntries(boundary_vec));
762 chkerr(VecScale(temp, -1));
768 const double *flux_array;
770 chkerr(VecGetArrayRead(temp, &flux_array));
771 chkerr(VecGetLocalSize(temp, &lsize));
772 for (
int e=0; e<lsize; ++e)
774 if (flux_array[e] < 0)
779 chkerr(VecRestoreArrayRead(temp, &flux_array));
780 chkerr(VecDestroy(&temp));
781 chkerr(VecDestroy(&boundary_vec));
798 const int buf_size = n_quant*2*n_blk_reg + n_quant*2*n_bdr_reg;
799 double sendbuffer[buf_size], recvbuffer[buf_size];
800 for (
unsigned int qi=0; qi<n_quant; qi++)
802 for (
unsigned int ri=0; ri<n_blk_reg; ri++)
804 sendbuffer[qi*2*n_blk_reg + + ri] =
sources_in_[qi][ri];
805 sendbuffer[qi*2*n_blk_reg + n_blk_reg + ri] =
sources_out_[qi][ri];
807 for (
unsigned int ri=0; ri<n_bdr_reg; ri++)
809 sendbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + + ri] =
fluxes_in_[qi][ri];
810 sendbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + n_bdr_reg + ri] =
fluxes_out_[qi][ri];
820 for (
unsigned int qi=0; qi<n_quant; qi++)
822 for (
unsigned int ri=0; ri<n_blk_reg; ri++)
824 sources_in_[qi][ri] = recvbuffer[qi*2*n_blk_reg + + ri];
825 sources_out_[qi][ri] = recvbuffer[qi*2*n_blk_reg + n_blk_reg + ri];
827 for (
unsigned int ri=0; ri<n_bdr_reg; ri++)
829 fluxes_in_[qi][ri] = recvbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + + ri];
830 fluxes_out_[qi][ri] = recvbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + n_bdr_reg + ri];
850 for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
852 for (
unsigned int qi=0; qi<n_quant; qi++)
862 for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
864 for (
unsigned int qi=0; qi<n_quant; qi++)
880 for (
unsigned int qi=0; qi<n_quant; qi++)
885 for (
unsigned int qi=0; qi<n_quant; qi++)
925 if (
rank_ != 0)
return;
932 unsigned int wl = 2*(w-5)+7;
933 string bc_head_format =
"# %-*s%-*s%-*s%-*s%-*s%-*s\n",
934 bc_format =
"%*s%-*d%-*s%-*s%-*g%-*g%-*g\n",
935 bc_total_format =
"# %-*s%-*s%-*g%-*g%-*g\n";
937 output_ <<
"# " << setw((w*c+wl-14)/2) << setfill(
'-') <<
"--"
939 << setw((w*c+wl-14)/2) << setfill(
'-') <<
"" << endl
943 output_ <<
"# Mass flux through boundary [M/T]:\n# "
944 << setiosflags(ios::left) << setfill(
' ')
945 << setw(w) <<
"[boundary_id]"
946 << setw(wl) <<
"[label]"
947 << setw(w) <<
"[substance]"
948 << setw(w) <<
"[total flux]"
949 << setw(w) <<
"[outward flux]"
950 << setw(w) <<
"[inward flux]"
954 output_ <<
"# " << setw(w*c+wl) << setfill(
'-') <<
"" << setfill(
' ') << endl;
958 for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg) {
959 for (
unsigned int qi=0; qi<n_quant; qi++) {
961 << setw(w) << (int)reg->id()
962 << setw(wl) << reg->label().c_str()
964 << setw(w) <<
fluxes_[qi][reg->boundary_idx()]
966 << setw(w) <<
fluxes_in_[qi][reg->boundary_idx()]
972 output_ <<
"# " << setw(w*c+wl) << setfill(
'-') <<
"" << setfill(
' ') << endl;
975 for (
unsigned int qi=0; qi<n_quant; qi++)
976 output_ <<
"# " << setiosflags(ios::left)
977 << setw(w+wl) <<
"Total mass flux of substance [M/T]"
987 string src_head_format =
"# %-*s%-*s%-*s%-*s%-*s\n",
988 src_format =
"%*s%-*d%-*s%-*s%-*g%-*g\n",
989 src_total_format =
"# %-*s%-*s%-*g%-*g\n";
990 output_ <<
"# Mass [M] and sources [M/T] on regions:\n"
991 <<
"# " << setiosflags(ios::left)
992 << setw(w) <<
"[region_id]"
993 << setw(wl) <<
"[label]"
994 << setw(w) <<
"[substance]"
995 << setw(w) <<
"[total_mass]"
996 << setw(w) <<
"[total_source]"
1000 output_ <<
"# " << setw(w*c+wl) << setfill(
'-') <<
"" << setfill(
' ') << endl;
1004 for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
1006 for (
unsigned int qi=0; qi<n_quant; qi++)
1009 << setw(w) << (int)reg->id()
1010 << setw(wl) << reg->label().c_str()
1012 << setw(w) <<
masses_[qi][reg->bulk_idx()]
1013 << setw(w) <<
sources_[qi][reg->bulk_idx()]
1019 output_ <<
"# " << setw(w*c+wl) << setfill(
'-') <<
"" << setfill(
' ') << endl;
1022 for (
unsigned int qi=0; qi<n_quant; qi++)
1023 output_ <<
"# " << setiosflags(ios::left) << setw(w+wl) <<
"Total mass [M] and sources [M/T]"
1032 output_ <<
"\n\n# Cumulative mass balance on time interval ["
1034 << setiosflags(ios::left) << time <<
"]\n"
1035 <<
"# Initial mass [M] + sources integrated over time [M] - flux integrated over time [M] = current mass [M]\n"
1036 <<
"# " << setiosflags(ios::left)
1037 << setw(w) <<
"[substance]"
1038 << setw(w) <<
"[A=init. mass]"
1039 << setw(w) <<
"[B=source]"
1040 << setw(w) <<
"[C=flux]"
1041 << setw(w) <<
"[A+B-C]"
1042 << setw(w) <<
"[D=curr. mass]"
1043 << setw(w) <<
"[A+B-C-D=err.]"
1044 << setw(w) <<
"[rel. error]"
1047 for (
unsigned int qi=0; qi<n_quant; qi++)
1050 output_ <<
" " << setiosflags(ios::left)
1069 std::stringstream ss;
1070 for (
unsigned int i=0; i<cnt; i++) ss <<
format_csv_val(0, delimiter);
1078 if (
rank_ != 0)
return;
1087 for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
1089 for (
unsigned int qi=0; qi<n_quant; qi++)
1109 for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
1111 for (
unsigned int qi=0; qi<n_quant; qi++) {
1128 for (
unsigned int qi=0; qi<n_quant; qi++)
1158 std::stringstream ss;
1159 if (delimiter ==
' ') {
1165 output_ << comment_string << ss.str()
1185 std::stringstream ss;
1186 std::replace( val.begin(), val.end(),
'\"',
'\'');
1188 if (!initial) ss << delimiter;
1189 if (delimiter ==
' ') {
1190 std::stringstream sval;
1191 sval <<
"\"" << val <<
"\"";
1194 ss <<
"\"" << val <<
"\"";
1202 std::stringstream ss;
1204 if (!initial) ss << delimiter;
1205 if (delimiter ==
' ') {
1223 output_yaml_ <<
"column_names: [ flux, flux_in, flux_out, mass, source, source_in, source_out, flux_increment, "
1224 <<
"source_increment, flux_cumulative, source_cumulative, error ]" << endl;
1233 for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
1235 for (
unsigned int qi=0; qi<n_quant; qi++)
1238 output_yaml_ << setw(4) <<
"" <<
"region: " << reg->label() << endl;
1240 output_yaml_ << setw(4) <<
"" <<
"data: " <<
"[ 0, 0, 0, " <<
masses_[qi][reg->bulk_idx()] <<
", "
1242 <<
sources_out_[qi][reg->bulk_idx()] <<
", 0, 0, 0, 0, 0 ]" << endl;
1248 for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
1250 for (
unsigned int qi=0; qi<n_quant; qi++) {
1252 output_yaml_ << setw(4) <<
"" <<
"region: " << reg->label() << endl;
1254 output_yaml_ << setw(4) <<
"" <<
"data: " <<
"[ " <<
fluxes_[qi][reg->boundary_idx()] <<
", "
1256 <<
", 0, 0, 0, 0, 0, 0, 0, 0, 0 ]" << endl;
1262 for (
unsigned int qi=0; qi<n_quant; qi++)
1266 output_yaml_ << setw(4) <<
"" <<
"region: ALL" << endl;
1274 << error <<
" ]" << endl;
Distribution * get_el_ds() const
void calculate_instant(unsigned int quantity_idx, const Vec &solution)
void add_mass_vec_value(unsigned int quantity_idx, unsigned int region_idx, double value)
static const unsigned int output_column_width
Size of column in output (used if delimiter is space)
std::vector< unsigned int > be_regions_
Maps local boundary edge to its region boundary index.
std::vector< double > sum_fluxes_out_
Balance(const std::string &file_prefix, const Mesh *mesh)
Basic time management class.
UnitSI units_
Units of conserved quantities.
std::string format_text() const
const TimeGovernor * time_
ofstream output_
Handle for file for output in given OutputFormat of balance and total fluxes over individual regions ...
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
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.
const unsigned int index_
Internal index within list of quantities.
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
unsigned int add_quantity(const string &name)
TimeMark::Type balance_output_type_
TimeMark type for balance output of particular equation.
TimeMark::Type equation_fixed_mark_type() const
Vec * region_source_vec_
Vectors for calculation of source (n_bulk_regions).
Support classes for parallel programing.
std::vector< double > integrated_sources_
void format_csv_output_header(char delimiter, const std::string &comment_string)
Print output header.
std::vector< double > sum_sources_
void finish_source_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix and vectors for computing source.
static constexpr bool value
void add_flux_matrix_values(unsigned int quantity_idx, SideIter side, const std::vector< LongIdx > &dof_indices, const std::vector< double > &values)
Dedicated class for storing path to input and output files.
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
bool output_yaml_header_
marks whether YAML output has printed header
void start_flux_assembly(unsigned int quantity_idx)
std::vector< double > increment_fluxes_
void init_from_input(const Input::Record &in_rec, TimeGovernor &tg)
static bool do_yaml_output_
std::vector< Quantity > quantities_
Names of conserved quantities.
void open_stream(Stream &stream) const
void units(const UnitSI &unit)
Setter for units of conserved quantities.
void add_mass_matrix_values(unsigned int quantity_idx, unsigned int region_idx, const std::vector< LongIdx > &dof_indices, const std::vector< double > &values)
unsigned int boundary_size() const
void start_mass_assembly(unsigned int quantity_idx)
Vec * region_mass_vec_
Vectors for calculation of mass (n_bulk_regions).
Mat * be_flux_matrix_
Matrices for calculation of flux (n_boundary_edges x n_dofs).
void finish_mass_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix for computing mass.
void add_flux_vec_value(unsigned int quantity_idx, SideIter side, double value)
void output_yaml(double time)
Perform output in yaml format.
Mat * region_mass_matrix_
Matrices for calculation of mass (n_dofs x n_bulk_regions).
std::string get_unit_string() const
bool balance_on_
If the balance is on. Balance is off in the case of no balance output time marks.
unsigned int bulk_size() const
bool cumulative_
if true then cumulative balance is computed
std::vector< std::vector< double > > masses_
Basic time management functionality for unsteady (and steady) solvers (class Equation).
const TimeStep & step(int index=-1) const
Vec * be_flux_vec_
Vectors for calculation of flux (n_boundary_edges).
void chkerr_assert(unsigned int ierr)
std::vector< double > sum_fluxes_in_
Class for representation SI units of Fields.
const RegionDB & region_db() const
std::vector< double > sum_masses_
std::vector< std::vector< double > > fluxes_
RegionSet get_region_set(const std::string &set_name) const
#define MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm)
unsigned int output_line_counter_
hold count of line printed into output_
void output()
Perform output to file for given time instant.
std::vector< double > initial_mass_
bool add_output_times_
Add output time marks to balance output time marks.
void allocate(unsigned int n_loc_dofs, unsigned int max_dofs_per_boundary)
void calculate_mass(unsigned int quantity_idx, const Vec &solution, vector< double > &output_array)
void read_from_input(Input::Array in_array, const TimeGovernor &tg)
Mat * region_source_matrix_
Matrices for calculation of source (n_dofs x n_bulk_regions).
bool is_current()
Returns true if the current time step is marked for the balance output.
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,...
void calculate_cumulative(unsigned int quantity_idx, const Vec &solution)
OutputFormat output_format_
Format of output file.
std::vector< std::vector< double > > fluxes_in_
int be_offset_
Offset for local part of vector of boundary edges.
std::vector< std::vector< double > > sources_in_
void undef(bool val=true)
std::vector< double > sum_fluxes_
std::vector< double > increment_sources_
std::string file_prefix_
Save prefix passed in in constructor.
unsigned int n_sides() const
double last_time_
time of last calculated balance
static const Input::Type::Array get_input_type()
static TimeMarks & marks()
static void set_yaml_output()
Set global variable to output balance files into YAML format (in addition to the table format).
void output_csv(double time, char delimiter, const std::string &comment_string, unsigned int repeat=0)
Perform output in csv format.
std::vector< std::vector< double > > fluxes_out_
static const Input::Type::Selection & get_format_selection_input_type()
Input selection for file format.
void add_source_vec_values(unsigned int quantity_idx, unsigned int region_idx, const std::vector< LongIdx > &dof_values, const std::vector< double > &values)
LongIdx get_boundary_edge_uid(SideIter side)
Vec ones_
auxiliary vectors for summation of matrix columns
void add_source_matrix_values(unsigned int quantity_idx, unsigned int region_idx, const std::vector< LongIdx > &dof_indices, const std::vector< double > &values)
unsigned int boundary_idx() const
Returns index of the region in the boundary set.
Input::Record input_record_
Record for current balance.
TimeMark::Type output_mark_type_
TimeMark type for output of particular equation.
TimeMark::Type equation_mark_type() const
void start_source_assembly(unsigned int quantity_idx)
unsigned int n_loc_dofs_
Allocation parameters. Set by the allocate method used in the lazy_initialize.
std::vector< double > sum_sources_in_
std::vector< std::vector< double > > sources_
bool allocation_done_
true before allocating necessary internal structures (Petsc matrices etc.)
std::vector< double > sum_sources_out_
FilePath balance_output_file_
File path for output_ stream.
bool initial_
true before calculating the mass at initial time, otherwise false
unsigned int max_dofs_per_boundary_
void add_cumulative_source(unsigned int quantity_idx, double source)
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
SideIter side(const unsigned int loc_index)
void finish_flux_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix and vector for computing flux.
virtual ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
std::vector< double > integrated_fluxes_
Mat * region_source_rhs_
Matrices for calculation of signed source (n_dofs x n_bulk_regions).
std::vector< std::vector< double > > sources_out_
unsigned int * boundary_idx_
std::unordered_map< LongIdx, unsigned int > be_id_map_
static const Input::Type::Record & get_input_type()
Main balance input record type.
void output_legacy(double time)
Perform output in old format (for compatibility)
std::vector< unsigned int > add_quantities(const std::vector< string > &names)
LongIdx * get_el_4_loc() const