Flow123d  JS_before_hm-915-gc632be9
balance.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file balance.cc
15  * @ingroup transport
16  * @brief Mass balance
17  */
18 
19 #include <iostream>
20 #include <iomanip>
21 #include <unordered_map>
22 
23 #include "system/system.hh"
24 #include "system/sys_profiler.hh"
25 #include "system/index_types.hh"
26 
27 #include <petscmat.h>
28 #include "mesh/mesh.h"
29 #include "mesh/accessors.hh"
30 #include "io/output_time_set.hh"
31 #include "coupling/balance.hh"
32 #include "tools/unit_si.hh"
33 #include "tools/time_governor.hh"
34 #include "la/distribution.hh"
35 #include "fem/dofhandler.hh"
36 #include "fem/dh_cell_accessor.hh"
37 
38 using namespace Input::Type;
39 
40 bool Balance::do_yaml_output_ = true;
41 
43  return Selection("Balance_output_format", "Format of output file for balance.")
44  .add_value(Balance::legacy, "legacy", "Legacy format used by previous program versions.")
45  .add_value(Balance::txt, "txt", "Excel format with tab delimiter.")
46  .add_value(Balance::gnuplot, "gnuplot", "Format compatible with GnuPlot datafile with fixed column width.")
47  .close();
48 }
49 
51  return Record("Balance", "Balance of a conservative quantity, boundary fluxes and sources.")
52  .declare_key("times", OutputTimeSet::get_input_type(), Default("[]"), "" )
53  .declare_key("add_output_times", Bool(), Default("true"), "Add all output times of the balanced equation to the balance output times set. "
54  "Note that this is not the time set of the output stream.")
55  //.declare_key("balance_on", Bool(), Default("true"), "Balance is computed if the value is true.")
56  .declare_key("format", Balance::get_format_selection_input_type(), Default("\"txt\""), "Format of output file.")
57  .declare_key("cumulative", Bool(), Default("false"), "Compute cumulative balance over time. "
58  "If true, then balance is calculated at each computational time step, which can slow down the program.")
59  .declare_key("file", FileName::output(), Default::read_time("File name generated from the balanced quantity: <quantity_name>_balance.*"), "File name for output of balance.")
60  .close();
61 }
62 
65 }
66 
67 /*
68 std::shared_ptr<Balance> Balance::make_balance(
69  const std::string &file_prefix,
70  const Mesh *mesh,
71  const Input::Record &in_rec,
72  TimeGovernor &tg)
73 {
74  auto ptr = make_shared<Balance>(file_prefix, mesh, in_rec, tg);
75  auto &marks = TimeGovernor::marks();
76  auto balance_output_type = tg.equation_mark_type() | TimeGovernor::marks().type_balance_output();
77  if (marks.begin(balance_output_type) == marks.end(balance_output_type) ) {
78  // no balance output time => force balance off
79  ptr.reset();
80  }
81  return ptr;
82 }*/
83 
84 
85 
86 Balance::Balance(const std::string &file_prefix, const Mesh *mesh)
87  : file_prefix_(file_prefix),
88  mesh_(mesh),
89  last_time_(),
90  initial_(true),
91  allocation_done_(false),
92  balance_on_(true),
93  output_line_counter_(0),
94  output_yaml_header_(false)
95 
96 {
97  MPI_Comm_rank(PETSC_COMM_WORLD, &rank_);
99 }
100 
101 
102 
104 {
105  if (rank_ == 0) {
106  output_.close();
107  if (do_yaml_output_) output_yaml_.close();
108  }
109  if (! allocation_done_) return;
110 
111  for (unsigned int c=0; c<quantities_.size(); ++c)
112  {
113  chkerr(MatDestroy(&(region_mass_matrix_[c])));
114  chkerr(MatDestroy(&(be_flux_matrix_[c])));
115  chkerr(MatDestroy(&(region_source_matrix_[c])));
116  chkerr(MatDestroy(&(region_source_rhs_[c])));
117  chkerr(VecDestroy(&(be_flux_vec_[c])));
118  chkerr(VecDestroy(&(region_mass_vec_[c])));
119  }
120  delete[] region_mass_matrix_;
121  delete[] be_flux_matrix_;
122  delete[] be_flux_vec_;
123  delete[] region_source_matrix_;
124  delete[] region_source_rhs_;
125  delete[] region_mass_vec_;
126 }
127 
128 
130  const Input::Record &in_rec,
131  TimeGovernor &tg)
132 {
134 
135  time_ = &tg;
136 
137  auto &marks = TimeGovernor::marks();
138 
139  output_mark_type_ = tg.equation_mark_type() | marks.type_output(),
140  balance_output_type_ = tg.equation_fixed_mark_type() | marks.type_balance_output();
141 
142  cumulative_ = in_rec.val<bool>("cumulative");
143  output_format_ = in_rec.val<OutputFormat>("format");
144 
145  OutputTimeSet time_set;
146  time_set.read_from_input( in_rec.val<Input::Array>("times"), tg, balance_output_type_);
147  add_output_times_ = in_rec.val<bool>("add_output_times");
148 
149  input_record_ = in_rec;
150 
151 }
152 
153 
154 void Balance::units(const UnitSI &unit)
155 {
157 
158  units_ = unit;
159  units_.undef(false);
160 }
161 
162 unsigned int Balance::add_quantity(const string &name)
163 {
165 
166  Quantity q(quantities_.size(), name);
167  quantities_.push_back(q);
168 
169  return q.index_;
170 }
171 
172 
174 {
175  vector<unsigned int> indices;
176  for (auto name : names)
177  indices.push_back(add_quantity(name));
178  return indices;
179 }
180 
181 
182 void Balance::allocate(unsigned int n_loc_dofs,
183  unsigned int max_dofs_per_boundary)
184 {
186  n_loc_dofs_seq_ = n_loc_dofs;
187  n_loc_dofs_par_ = n_loc_dofs;
188  max_dofs_per_boundary_ = max_dofs_per_boundary;
189 }
190 
191 void Balance::allocate(const std::shared_ptr<DOFHandlerMultiDim>& dh,
192  unsigned int max_dofs_per_boundary)
193 {
195  // for sequential matrices, we need to include ghost values
196  n_loc_dofs_seq_ = dh->get_local_to_global_map().size();
197  // for parallel matrices, we use the local size from dof distribution
198  n_loc_dofs_par_ = dh->distr()->lsize();
199  max_dofs_per_boundary_ = max_dofs_per_boundary;
200 }
201 
203 {
204  if (allocation_done_) return;
205 
206  auto &marks = TimeGovernor::marks();
207  if (add_output_times_)
208  marks.add_to_type_all(output_mark_type_, balance_output_type_);
209  // if there are no balance marks turn balance off
210  if (marks.begin(balance_output_type_) == marks.end(balance_output_type_) )
211  {
212  balance_on_ = false;
213  cumulative_ = false;
214  return;
215  }
216 
217  // Max. number of regions to which a single dof can contribute.
218  // TODO: estimate or compute this number directly (from mesh or dof handler).
219  const int n_bulk_regs_per_dof = min(10, (int)mesh_->region_db().bulk_size());
220  const unsigned int n_quant = quantities_.size();
221  const unsigned int n_bdr_reg = mesh_->region_db().boundary_size();
222  const unsigned int n_blk_reg = mesh_->region_db().bulk_size();
223 
224 
225  // enumerate boundary edges by unique id and
226  // create map: (local bulk ele idx, side idx) -> boundary edge id
227  // create vector that maps: boundary edge id -> region boundary index
228  unsigned int be_id = 0;
229  for (unsigned int loc_el = 0; loc_el < mesh_->get_el_ds()->lsize(); loc_el++)
230  {
232  if (elm->boundary_idx_ != nullptr)
233  {
234  for(unsigned int si=0; si<elm->n_sides(); si++)
235  {
236  if (elm.side(si)->is_boundary()){
237  Boundary bcd = elm.side(si)->cond();
238  LongIdx ele_side_uid = get_boundary_edge_uid(elm.side(si));
239  be_id_map_[ele_side_uid] = be_id;
240  be_regions_.push_back(bcd.region().boundary_idx());
241  be_id++;
242  }
243  }
244  }
245  }
246 
247 
248  fluxes_ .resize(n_quant, vector<double>(n_bdr_reg, 0));
249  fluxes_in_ .resize(n_quant, vector<double>(n_bdr_reg, 0));
250  fluxes_out_.resize(n_quant, vector<double>(n_bdr_reg, 0));
251 
252  masses_ .resize(n_quant, vector<double>(n_blk_reg, 0));
253  sources_in_ .resize(n_quant, vector<double>(n_blk_reg, 0));
254  sources_out_.resize(n_quant, vector<double>(n_blk_reg, 0));
255 
256  sum_fluxes_ .resize(n_quant, 0);
257  sum_fluxes_in_ .resize(n_quant, 0);
258  sum_fluxes_out_ .resize(n_quant, 0);
259  sum_masses_ .resize(n_quant, 0);
260  sum_sources_ .resize(n_quant, 0);
261  sum_sources_in_ .resize(n_quant, 0);
262  sum_sources_out_.resize(n_quant, 0);
263 
264  if (cumulative_)
265  {
266  initial_mass_ .resize(n_quant, 0);
267  integrated_fluxes_ .resize(n_quant, 0);
268  integrated_sources_.resize(n_quant, 0);
269  increment_sources_.resize(n_quant, 0);
270  increment_fluxes_.resize(n_quant, 0);
271  }
272 
273 
274 
275  region_mass_matrix_ = new Mat[n_quant];
276  be_flux_matrix_ = new Mat[n_quant];
277  region_source_matrix_ = new Mat[n_quant];
278  region_source_rhs_ = new Mat[n_quant];
279  region_mass_vec_ = new Vec[n_quant];
280  be_flux_vec_ = new Vec[n_quant];
281 
282 
283  for (unsigned int c=0; c<n_quant; ++c)
284  {
285  chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
287  (rank_==0)?mesh_->region_db().bulk_size():0,
288  PETSC_DECIDE,
289  PETSC_DECIDE,
290  (rank_==0)?n_bulk_regs_per_dof:0,
291  0,
292  (rank_==0)?0:n_bulk_regs_per_dof,
293  0,
294  &(region_mass_matrix_[c])));
295 
296  chkerr(MatCreateSeqAIJ(PETSC_COMM_SELF,
299  n_bulk_regs_per_dof,
300  NULL,
301  &(region_source_matrix_[c])));
302 
303  chkerr(MatCreateSeqAIJ(PETSC_COMM_SELF,
306  n_bulk_regs_per_dof,
307  NULL,
308  &(region_source_rhs_[c])));
309 
310  chkerr(MatCreateAIJ(PETSC_COMM_WORLD,
311  be_regions_.size(), // n local rows, number of local boundary edges
312  n_loc_dofs_par_, // n local cols (local rows of multiplying vector)
313  PETSC_DECIDE, // n global rows
314  PETSC_DECIDE, // n global cols
315  max_dofs_per_boundary_, // allocation, local poriton
316  0,
317  0,
318  0,
319  &(be_flux_matrix_[c])));
320 
321  chkerr(VecCreateMPI(PETSC_COMM_WORLD,
322  (rank_==0)?mesh_->region_db().bulk_size():0,
323  PETSC_DECIDE,
324  &(region_mass_vec_[c])));
325 
326  chkerr(VecCreateMPI(PETSC_COMM_WORLD,
327  be_regions_.size(),
328  PETSC_DECIDE,
329  &(be_flux_vec_[c])));
330  }
331 
332  // set be_offset_, used in add_flux_matrix_values()
333  chkerr(VecGetOwnershipRange(be_flux_vec_[0], &be_offset_, NULL));
334 
335  if (rank_ == 0) {
336  // set default value by output_format_
337  std::string default_file_name;
338  switch (output_format_)
339  {
340  case txt:
341  default_file_name = file_prefix_ + "_balance.txt";
342  break;
343  case gnuplot:
344  default_file_name = file_prefix_ + "_balance.dat";
345  break;
346  case legacy:
347  default_file_name = file_prefix_ + "_balance.txt";
348  break;
349  }
350 
352  try {
354  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, input_record_)
355 
356 
357  // set file name of YAML output
358  if (do_yaml_output_) {
359  string yaml_file_name = file_prefix_ + "_balance.yaml";
361  }
362  }
363 
364  allocation_done_ = true;
365 }
366 
367 
368 
370 {
371  if (! balance_on_) return false;
372 
373  auto &marks = TimeGovernor::marks();
374  bool res = marks.current(time_->step(), balance_output_type_) != marks.end(balance_output_type_);
375 
376  //cout << "flag: " << res << " time: " << step.end() << " type:" << hex << balance_output_type_ << endl;
377  return res;
378 }
379 
380 void Balance::start_mass_assembly(unsigned int quantity_idx)
381 {
382  lazy_initialize();
383  if (! balance_on_) return;
384  chkerr(MatZeroEntries(region_mass_matrix_[quantity_idx]));
385  chkerr(VecZeroEntries(region_mass_vec_[quantity_idx]));
386 }
387 
388 
389 void Balance::start_flux_assembly(unsigned int quantity_idx)
390 {
391  lazy_initialize();
392  if (! balance_on_) return;
393  chkerr(MatZeroEntries(be_flux_matrix_[quantity_idx]));
394  chkerr(VecZeroEntries(be_flux_vec_[quantity_idx]));
395 }
396 
397 
398 void Balance::start_source_assembly(unsigned int quantity_idx)
399 {
400  lazy_initialize();
401  if (! balance_on_) return;
402  chkerr(MatZeroEntries(region_source_matrix_[quantity_idx]));
403  chkerr(MatZeroEntries(region_source_rhs_[quantity_idx]));
404 }
405 
406 
407 void Balance::finish_mass_assembly(unsigned int quantity_idx)
408 {
410  if (! balance_on_) return;
411 
412  chkerr(MatAssemblyBegin(region_mass_matrix_[quantity_idx], MAT_FINAL_ASSEMBLY));
413  chkerr(MatAssemblyEnd(region_mass_matrix_[quantity_idx], MAT_FINAL_ASSEMBLY));
414  chkerr(VecAssemblyBegin(region_mass_vec_[quantity_idx]));
415  chkerr(VecAssemblyEnd(region_mass_vec_[quantity_idx]));
416 }
417 
418 void Balance::finish_flux_assembly(unsigned int quantity_idx)
419 {
421  if (! balance_on_) return;
422 
423  chkerr(MatAssemblyBegin(be_flux_matrix_[quantity_idx], MAT_FINAL_ASSEMBLY));
424  chkerr(MatAssemblyEnd(be_flux_matrix_[quantity_idx], MAT_FINAL_ASSEMBLY));
425  chkerr(VecAssemblyBegin(be_flux_vec_[quantity_idx]));
426  chkerr(VecAssemblyEnd(be_flux_vec_[quantity_idx]));
427 }
428 
429 void Balance::finish_source_assembly(unsigned int quantity_idx)
430 {
432  if (! balance_on_) return;
433 
434  chkerr(MatAssemblyBegin(region_source_matrix_[quantity_idx], MAT_FINAL_ASSEMBLY));
435  chkerr(MatAssemblyEnd(region_source_matrix_[quantity_idx], MAT_FINAL_ASSEMBLY));
436  chkerr(MatAssemblyBegin(region_source_rhs_[quantity_idx], MAT_FINAL_ASSEMBLY));
437  chkerr(MatAssemblyEnd(region_source_rhs_[quantity_idx], MAT_FINAL_ASSEMBLY));
438 }
439 
440 
441 void Balance::add_mass_values(unsigned int quantity_idx,
442  const DHCellAccessor &dh_cell,
443  const LocDofVec &loc_dof_indices,
444  const std::vector<double> &mat_values,
445  double vec_value)
446 {
448  if (! balance_on_) return;
449 
450  // map local dof indices to global
451  uint m = mat_values.size();
452  int row_dofs[m];
453  for (uint i=0; i<m; i++)
454  row_dofs[i]= dh_cell.dh()->get_local_to_global_map()[loc_dof_indices[i]];
455 
456  PetscInt reg_array[1] = { (int)dh_cell.elm().region_idx().bulk_idx() };
457 
458  chkerr_assert(MatSetValues(region_mass_matrix_[quantity_idx],
459  m,
460  row_dofs,
461  1,
462  reg_array,
463  &(mat_values[0]),
464  ADD_VALUES));
465 
466  chkerr_assert(VecSetValue(region_mass_vec_[quantity_idx],
467  dh_cell.elm().region_idx().bulk_idx(),
468  vec_value,
469  ADD_VALUES));
470 }
471 
472 void Balance::add_flux_values(unsigned int quantity_idx,
473  const DHCellSide &side,
474  const LocDofVec &loc_dof_indices,
475  const std::vector<double> &mat_values,
476  double vec_value)
477 {
479  if (! balance_on_) return;
480 
481  // filling row elements corresponding to a boundary edge
482 
483  // map local dof indices to global
484  uint m = mat_values.size();
485  int col_dofs[m];
486  for (uint i=0; i<m; i++)
487  col_dofs[i]= side.cell().dh()->get_local_to_global_map()[loc_dof_indices[i]];
488 
489  SideIter s = SideIter(side.side());
490  PetscInt glob_be_idx[1] = { int(be_offset_ + be_id_map_[get_boundary_edge_uid(s)]) };
491 
492  chkerr_assert(MatSetValues(be_flux_matrix_[quantity_idx],
493  1,
494  glob_be_idx,
495  m,
496  col_dofs,
497  &(mat_values[0]),
498  ADD_VALUES));
499 
500  chkerr_assert(VecSetValue(be_flux_vec_[quantity_idx],
501  glob_be_idx[0],
502  vec_value,
503  ADD_VALUES));
504 }
505 
506 void Balance::add_source_values(unsigned int quantity_idx,
507  unsigned int region_idx,
508  const LocDofVec &loc_dof_indices,
509  const vector<double> &mult_mat_values,
510  const vector<double> &add_mat_values)
511 {
513  if (! balance_on_) return;
514 
515  PetscInt reg_array[1] = { (int)region_idx };
516 
517  chkerr_assert(MatSetValues(region_source_matrix_[quantity_idx],
518  loc_dof_indices.size(),
519  loc_dof_indices.memptr(),
520  1,
521  reg_array,
522  &(mult_mat_values[0]),
523  ADD_VALUES));
524 
525  chkerr_assert(MatSetValues(region_source_rhs_[quantity_idx],
526  loc_dof_indices.size(),
527  loc_dof_indices.memptr(),
528  1,
529  reg_array,
530  &(add_mat_values[0]),
531  ADD_VALUES));
532 }
533 
534 
535 void Balance::add_cumulative_source(unsigned int quantity_idx, double source)
536 {
538  if (!cumulative_) return;
539 
540  if (rank_ == 0)
541  increment_sources_[quantity_idx] += source;
542 }
543 
544 
545 void Balance::calculate_cumulative(unsigned int quantity_idx,
546  const Vec &solution)
547 {
549  if (!cumulative_) return;
550  if (time_->tlevel() <= 0) return;
551 
552  // sources
553  double temp_source = 0;
554  int lsize, n_cols_mat, n_cols_rhs;
555  //const int *cols;
556  const double *vals_mat, *vals_rhs, *sol_array;
557  // chkerr(VecGetLocalSize(solution, &lsize));
558  // chkerr(ISLocalToGlobalMappingGetSize(solution->mapping, &lsize); // cannot do for const
559  lsize = n_loc_dofs_seq_;
560  chkerr(VecGetArrayRead(solution, &sol_array));
561 
562  // computes transpose multiplication and sums region_source_rhs_ over dofs
563  // resulting in a vector of sources for each region
564  // transpose(region_source_matrix_) * solution + region_source_rhs_*ones(n_blk_reg)
565  // the region vector is then summed up to temp_source
566  for (int i=0; i<lsize; ++i){
567  chkerr(MatGetRow(region_source_matrix_[quantity_idx], i, &n_cols_mat, NULL, &vals_mat));
568  chkerr(MatGetRow(region_source_rhs_[quantity_idx], i, &n_cols_rhs, NULL, &vals_rhs));
569 
570  ASSERT_DBG(n_cols_mat == n_cols_rhs);
571 
572  for (int j=0; j<n_cols_mat; ++j)
573  temp_source += vals_mat[j]*sol_array[i] + vals_rhs[j];
574 
575  chkerr(MatRestoreRow(region_source_matrix_[quantity_idx], i, &n_cols_mat, NULL, &vals_mat));
576  chkerr(MatRestoreRow(region_source_rhs_[quantity_idx], i, &n_cols_rhs, NULL, &vals_rhs));
577  }
578  chkerr(VecRestoreArrayRead(solution, &sol_array));
579 
580  increment_sources_[quantity_idx] += temp_source*time_->dt();
581 
582 
583  // fluxes
584  Vec temp;
585  chkerr(VecCreateMPI(PETSC_COMM_WORLD,
586  be_regions_.size(),
587  PETSC_DECIDE,
588  &temp));
589 
590  chkerr(MatMultAdd(be_flux_matrix_[quantity_idx], solution, be_flux_vec_[quantity_idx], temp));
591 
592  double sum_fluxes;
593  chkerr(VecSum(temp, &sum_fluxes));
594  chkerr(VecDestroy(&temp));
595 
596  if (rank_ == 0)
597  // sum fluxes in one step
598  // Since internally we keep outgoing fluxes, we change sign
599  // to write to output _incoming_ fluxes.
600  increment_fluxes_[quantity_idx] += -1.0 * sum_fluxes*time_->dt();
601 }
602 
603 
604 void Balance::calculate_mass(unsigned int quantity_idx,
605  const Vec &solution,
606  vector<double> &output_array)
607 {
609  if (! balance_on_) return;
610 
611  Vec bulk_vec;
612 
613  chkerr(VecCreateMPIWithArray(PETSC_COMM_WORLD,
614  1,
615  (rank_==0)?mesh_->region_db().bulk_size():0,
616  PETSC_DECIDE,
617  &(output_array[0]),
618  &bulk_vec));
619 
620  // compute mass on regions: M'.u
621  chkerr(VecZeroEntries(bulk_vec));
622  chkerr(MatMultTransposeAdd(region_mass_matrix_[quantity_idx],
623  solution,
624  region_mass_vec_[quantity_idx],
625  bulk_vec));
626  chkerr(VecDestroy(&bulk_vec));
627 }
628 
629 void Balance::calculate_instant(unsigned int quantity_idx, const Vec& solution)
630 {
631  if ( !is_current() ) return;
632 
633  calculate_mass(quantity_idx, solution, masses_[quantity_idx]);
634 
635  // compute positive/negative sources
636  for (unsigned int r=0; r<mesh_->region_db().bulk_size(); ++r)
637  {
638  sources_in_[quantity_idx][r] = 0;
639  sources_out_[quantity_idx][r] = 0;
640  }
641 
642  int lsize, n_cols_mat, n_cols_rhs;
643  const int *cols; // the columns must be same - matrices created and filled in the same way
644  const double *vals_mat, *vals_rhs, *sol_array;
645  // chkerr(VecGetLocalSize(solution, &lsize));
646  // chkerr(ISLocalToGlobalMappingGetSize(solution->mapping, &lsize); // cannot do for const
647  lsize = n_loc_dofs_seq_;
648  chkerr(VecGetArrayRead(solution, &sol_array));
649 
650  // computes transpose multiplication and sums region_source_rhs_ over dofs
651  // resulting in a vector of sources for each region, one positive, one negative
652  // transpose(region_source_matrix_) * solution + region_source_rhs_*ones(n_blk_reg)
653  for (int i=0; i<lsize; ++i)
654  {
655  int row = i;
656  chkerr(MatGetRow(region_source_matrix_[quantity_idx], row, &n_cols_mat, &cols, &vals_mat));
657  chkerr(MatGetRow(region_source_rhs_[quantity_idx], row, &n_cols_rhs, NULL, &vals_rhs));
658 
659  ASSERT_DBG(n_cols_mat == n_cols_rhs);
660 
661  for (int j=0; j<n_cols_mat; ++j)
662  {
663  int col = cols[j];
664 
665  double f = vals_mat[j]*sol_array[i] + vals_rhs[j];
666  if (f > 0) sources_in_[quantity_idx][col] += f;
667  else sources_out_[quantity_idx][col] += f;
668  }
669  }
670  chkerr(VecRestoreArrayRead(solution, &sol_array));
671 
672  // calculate flux
673  Vec temp;
674  chkerr(VecCreateMPI(PETSC_COMM_WORLD,
675  be_regions_.size(),
676  PETSC_DECIDE,
677  &temp));
678 
679  chkerr(MatMultAdd(be_flux_matrix_[quantity_idx], solution, be_flux_vec_[quantity_idx], temp));
680  // Since internally we keep outgoing fluxes, we change sign
681  // to write to output _incoming_ fluxes.
682  chkerr(VecScale(temp, -1));
683 
684  // compute positive/negative fluxes
685  fluxes_in_[quantity_idx].assign(mesh_->region_db().boundary_size(), 0);
686  fluxes_out_[quantity_idx].assign(mesh_->region_db().boundary_size(), 0);
687  const double *flux_array;
688 // int lsize;
689  chkerr(VecGetArrayRead(temp, &flux_array));
690  chkerr(VecGetLocalSize(temp, &lsize));
691  for (int e=0; e<lsize; ++e)
692  {
693  if (flux_array[e] < 0)
694  fluxes_out_[quantity_idx][be_regions_[e]] += flux_array[e];
695  else
696  fluxes_in_[quantity_idx][be_regions_[e]] += flux_array[e];
697  }
698  chkerr(VecRestoreArrayRead(temp, &flux_array));
699  chkerr(VecDestroy(&temp));
700 }
701 
702 
703 
704 
705 
707 {
709  if (! balance_on_) return;
710  if (! is_current() ) return;
711 
712  // gather results from processes and sum them up
713  const unsigned int n_quant = quantities_.size();
714  const unsigned int n_blk_reg = mesh_->region_db().bulk_size();
715  const unsigned int n_bdr_reg = mesh_->region_db().boundary_size();
716  const int buf_size = n_quant*2*n_blk_reg + n_quant*2*n_bdr_reg + n_quant;
717  double sendbuffer[buf_size], recvbuffer[buf_size];
718  for (unsigned int qi=0; qi<n_quant; qi++)
719  {
720  for (unsigned int ri=0; ri<n_blk_reg; ri++)
721  {
722  sendbuffer[qi*2*n_blk_reg + + ri] = sources_in_[qi][ri];
723  sendbuffer[qi*2*n_blk_reg + n_blk_reg + ri] = sources_out_[qi][ri];
724  }
725  for (unsigned int ri=0; ri<n_bdr_reg; ri++)
726  {
727  sendbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + + ri] = fluxes_in_[qi][ri];
728  sendbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + n_bdr_reg + ri] = fluxes_out_[qi][ri];
729  }
730  if (cumulative_)
731  {
732  sendbuffer[n_quant*2*n_blk_reg + n_quant*2*n_bdr_reg + qi] = increment_sources_[qi];
733  }
734  }
735 
736  MPI_Reduce(&sendbuffer,recvbuffer,buf_size,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);
737  // for other than 0th process update last_time and finish,
738  // on process #0 sum balances over all regions and calculate
739  // cumulative balance over time.
740  if (rank_ == 0)
741  {
742  // update balance vectors
743  for (unsigned int qi=0; qi<n_quant; qi++)
744  {
745  for (unsigned int ri=0; ri<n_blk_reg; ri++)
746  {
747  sources_in_[qi][ri] = recvbuffer[qi*2*n_blk_reg + + ri];
748  sources_out_[qi][ri] = recvbuffer[qi*2*n_blk_reg + n_blk_reg + ri];
749  }
750  for (unsigned int ri=0; ri<n_bdr_reg; ri++)
751  {
752  fluxes_in_[qi][ri] = recvbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + + ri];
753  fluxes_out_[qi][ri] = recvbuffer[n_quant*2*n_blk_reg + qi*2*n_bdr_reg + n_bdr_reg + ri];
754  }
755  if (cumulative_)
756  {
757  increment_sources_[qi] = recvbuffer[n_quant*2*n_blk_reg + n_quant*2*n_bdr_reg + qi];
758  }
759  }
760  }
761 
762 
763  // The convention for input/output of fluxes is that positive means inward.
764  // Therefore in the following code we switch sign of fluxes.
765  if (rank_ == 0)
766  {
767  sum_fluxes_.assign(n_quant, 0);
768  sum_fluxes_in_.assign(n_quant, 0);
769  sum_fluxes_out_.assign(n_quant, 0);
770  sum_masses_.assign(n_quant, 0);
771  sum_sources_.assign(n_quant, 0);
772  sum_sources_in_.assign(n_quant, 0);
773  sum_sources_out_.assign(n_quant, 0);
774 
775  // sum all boundary fluxes
776  const RegionSet & b_set = mesh_->region_db().get_region_set(".BOUNDARY");
777  for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
778  {
779  for (unsigned int qi=0; qi<n_quant; qi++)
780  {
781  sum_fluxes_[qi] += fluxes_in_ [qi][reg->boundary_idx()] + fluxes_out_[qi][reg->boundary_idx()];
782  sum_fluxes_in_[qi] += fluxes_in_ [qi][reg->boundary_idx()];
783  sum_fluxes_out_[qi] += fluxes_out_[qi][reg->boundary_idx()];
784  }
785  }
786 
787  // sum all volume sources
788  const RegionSet & bulk_set = mesh_->region_db().get_region_set("BULK");
789  for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
790  {
791  for (unsigned int qi=0; qi<n_quant; qi++)
792  {
793  sum_masses_[qi] += masses_[qi][reg->bulk_idx()];
794  sum_sources_[qi] += sources_in_[qi][reg->bulk_idx()] + sources_out_[qi][reg->bulk_idx()];
795  sum_sources_in_[qi] += sources_in_[qi][reg->bulk_idx()];
796  sum_sources_out_[qi] += sources_out_[qi][reg->bulk_idx()];
797  }
798  }
799 
800  // cumulative balance over time
801  if (cumulative_)
802  {
803  // save initial time and mass
804  if (initial_)
805  {
807  for (unsigned int qi=0; qi<n_quant; qi++)
808  initial_mass_[qi] = sum_masses_[qi];
809  initial_ = false;
810  }
811 
812  for (unsigned int qi=0; qi<n_quant; qi++)
813  {
816  }
817  }
818  }
819 
820  last_time_ = time_->t();
821 
822 
823  // perform actual output
824  switch (output_format_)
825  {
826  case txt:
827  output_csv(time_->t(), '\t', "");
828  break;
829  case gnuplot:
830  output_csv(time_->t(), ' ', "#", 30);
831  break;
832  case legacy:
833  output_legacy(time_->t());
834  break;
835  }
836  // output in YAML format
837  output_yaml(time_->t());
838 
839  if (rank_ == 0)
840  {
841  sum_fluxes_.assign(n_quant, 0);
842  sum_sources_.assign(n_quant, 0);
843  increment_fluxes_.assign(n_quant, 0);
844  }
845  increment_sources_.assign(n_quant, 0);
846 }
847 
848 
849 void Balance::output_legacy(double time)
850 {
851  // write output only on process #0
852  if (rank_ != 0) return;
853 
854  const unsigned int n_quant = quantities_.size();
855 
856  // print the head of mass balance file
857  unsigned int c = 6; //column number without label
858  unsigned int w = 14; //column width
859  unsigned int wl = 2*(w-5)+7; //label column width
860  string bc_head_format = "# %-*s%-*s%-*s%-*s%-*s%-*s\n",
861  bc_format = "%*s%-*d%-*s%-*s%-*g%-*g%-*g\n",
862  bc_total_format = "# %-*s%-*s%-*g%-*g%-*g\n";
863 
864  output_ << "# " << setw((w*c+wl-14)/2) << setfill('-') << "--"
865  << " MASS BALANCE "
866  << setw((w*c+wl-14)/2) << setfill('-') << "" << endl
867  << "# Time: " << (time / time_->get_coef()) << "[" << time_->get_unit_string() << "]\n\n\n";
868 
869  // header for table of boundary fluxes
870  output_ << "# Mass flux through boundary [M/T]:\n# "
871  << setiosflags(ios::left) << setfill(' ')
872  << setw(w) << "[boundary_id]"
873  << setw(wl) << "[label]"
874  << setw(w) << "[substance]"
875  << setw(w) << "[total flux]"
876  << setw(w) << "[outward flux]"
877  << setw(w) << "[inward flux]"
878  << endl;
879 
880  // draw long line
881  output_ << "# " << setw(w*c+wl) << setfill('-') << "" << setfill(' ') << endl;
882 
883  // print mass fluxes over boundaries
884  const RegionSet & b_set = mesh_->region_db().get_region_set(".BOUNDARY");
885  for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg) {
886  for (unsigned int qi=0; qi<n_quant; qi++) {
887  output_ << setw(2) << ""
888  << setw(w) << (int)reg->id()
889  << setw(wl) << reg->label().c_str()
890  << setw(w) << quantities_[qi].name_.c_str()
891  << setw(w) << fluxes_in_[qi][reg->boundary_idx()] + fluxes_out_[qi][reg->boundary_idx()]
892  << setw(w) << fluxes_out_[qi][reg->boundary_idx()]
893  << setw(w) << fluxes_in_[qi][reg->boundary_idx()]
894  << endl;
895  }
896  }
897 
898  // draw long line
899  output_ << "# " << setw(w*c+wl) << setfill('-') << "" << setfill(' ') << endl;
900 
901  // total boundary balance
902  for (unsigned int qi=0; qi<n_quant; qi++)
903  output_ << "# " << setiosflags(ios::left)
904  << setw(w+wl) << "Total mass flux of substance [M/T]"
905  << setw(w) << quantities_[qi].name_.c_str()
906  << setw(w) << sum_fluxes_[qi]
907  << setw(w) << sum_fluxes_out_[qi]
908  << setw(w) << sum_fluxes_in_[qi]
909  << endl;
910  output_ << "\n\n";
911 
912 
913  // header for table of volume sources and masses
914  string src_head_format = "# %-*s%-*s%-*s%-*s%-*s\n",
915  src_format = "%*s%-*d%-*s%-*s%-*g%-*g\n",
916  src_total_format = "# %-*s%-*s%-*g%-*g\n";
917  output_ << "# Mass [M] and sources [M/T] on regions:\n"
918  << "# " << setiosflags(ios::left)
919  << setw(w) << "[region_id]"
920  << setw(wl) << "[label]"
921  << setw(w) << "[substance]"
922  << setw(w) << "[total_mass]"
923  << setw(w) << "[total_source]"
924  << endl;
925 
926  // draw long line
927  output_ << "# " << setw(w*c+wl) << setfill('-') << "" << setfill(' ') << endl;
928 
929  // print balance of volume sources and masses
930  const RegionSet & bulk_set = mesh_->region_db().get_region_set("BULK");
931  for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
932  {
933  for (unsigned int qi=0; qi<n_quant; qi++)
934  {
935  output_ << setw(2) << ""
936  << setw(w) << (int)reg->id()
937  << setw(wl) << reg->label().c_str()
938  << setw(w) << quantities_[qi].name_.c_str()
939  << setw(w) << masses_[qi][reg->bulk_idx()]
940  << setw(w) << sources_in_[qi][reg->bulk_idx()] + sources_out_[qi][reg->bulk_idx()]
941  << endl;
942  }
943  }
944 
945  // draw long line
946  output_ << "# " << setw(w*c+wl) << setfill('-') << "" << setfill(' ') << endl;
947 
948  // total sources balance
949  for (unsigned int qi=0; qi<n_quant; qi++)
950  output_ << "# " << setiosflags(ios::left) << setw(w+wl) << "Total mass [M] and sources [M/T]"
951  << setw(w) << quantities_[qi].name_.c_str()
952  << setw(w) << sum_masses_[qi]
953  << setw(w) << sum_sources_[qi]
954  << endl;
955 
956  if (cumulative_)
957  {
958  // Print cumulative sources
959  output_ << "\n\n# Cumulative mass balance on time interval ["
960  << setiosflags(ios::left) << time_->init_time() << ","
961  << setiosflags(ios::left) << time << "]\n"
962  << "# Initial mass [M] + sources integrated over time [M] - flux integrated over time [M] = current mass [M]\n"
963  << "# " << setiosflags(ios::left)
964  << setw(w) << "[substance]"
965  << setw(w) << "[A=init. mass]"
966  << setw(w) << "[B=source]"
967  << setw(w) << "[C=flux]"
968  << setw(w) << "[A+B-C]"
969  << setw(w) << "[D=curr. mass]"
970  << setw(w) << "[A+B-C-D=err.]"
971  << setw(w) << "[rel. error]"
972  << endl;
973 
974  for (unsigned int qi=0; qi<n_quant; qi++)
975  {
976  double denominator = max(fabs(initial_mass_[qi]+integrated_sources_[qi]-integrated_fluxes_[qi]),fabs(sum_masses_[qi]));
977  output_ << " " << setiosflags(ios::left)
978  << setw(w) << quantities_[qi].name_.c_str()
979  << setw(w) << initial_mass_[qi]
980  << setw(w) << integrated_sources_[qi]
981  << setw(w) << integrated_fluxes_[qi]
982  << setw(w) << initial_mass_[qi]+integrated_sources_[qi]-integrated_fluxes_[qi]
983  << setw(w) << sum_masses_[qi]
985  << setw(w) << fabs(initial_mass_[qi]+integrated_sources_[qi]+integrated_fluxes_[qi]-sum_masses_[qi])/(denominator==0?1:denominator)
986  << endl;
987  }
988  }
989 
990  output_ << endl << endl;
991 }
992 
993 
994 std::string Balance::csv_zero_vals(unsigned int cnt, char delimiter)
995 {
996  std::stringstream ss;
997  for (unsigned int i=0; i<cnt; i++) ss << format_csv_val(0, delimiter);
998  return ss.str();
999 }
1000 
1001 
1002 void Balance::output_csv(double time, char delimiter, const std::string& comment_string, unsigned int repeat)
1003 {
1004  // write output only on process #0
1005  if (rank_ != 0) return;
1006 
1007  const unsigned int n_quant = quantities_.size();
1008 
1009  // print data header only on first line
1010  if (repeat==0 && output_line_counter_==0) format_csv_output_header(delimiter, comment_string);
1011 
1012  // print sources and masses over bulk regions
1013  const RegionSet & bulk_set = mesh_->region_db().get_region_set("BULK");
1014  for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
1015  {
1016  for (unsigned int qi=0; qi<n_quant; qi++)
1017  {
1018  // print data header (repeat header after every "repeat" lines)
1019  if (repeat && (output_line_counter_%repeat == 0)) format_csv_output_header(delimiter, comment_string);
1020 
1021  output_ << format_csv_val(time / time_->get_coef(), delimiter, true)
1022  << format_csv_val(reg->label(), delimiter)
1023  << format_csv_val(quantities_[qi].name_, delimiter)
1024  << csv_zero_vals(3, delimiter)
1025  << format_csv_val(masses_[qi][reg->bulk_idx()], delimiter)
1026  << format_csv_val(sources_in_[qi][reg->bulk_idx()] + sources_out_[qi][reg->bulk_idx()], delimiter)
1027  << format_csv_val(sources_in_[qi][reg->bulk_idx()], delimiter)
1028  << format_csv_val(sources_out_[qi][reg->bulk_idx()], delimiter)
1029  << csv_zero_vals(5, delimiter) << endl;
1031  }
1032  }
1033 
1034  // print mass fluxes over boundaries
1035  const RegionSet & b_set = mesh_->region_db().get_region_set(".BOUNDARY");
1036  for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
1037  {
1038  for (unsigned int qi=0; qi<n_quant; qi++) {
1039  // print data header (repeat header after every "repeat" lines)
1040  if (repeat && (output_line_counter_%repeat == 0)) format_csv_output_header(delimiter, comment_string);
1041 
1042  output_ << format_csv_val(time / time_->get_coef(), delimiter, true)
1043  << format_csv_val(reg->label(), delimiter)
1044  << format_csv_val(quantities_[qi].name_, delimiter)
1045  << format_csv_val(fluxes_in_[qi][reg->boundary_idx()] + fluxes_out_[qi][reg->boundary_idx()], delimiter)
1046  << format_csv_val(fluxes_in_[qi][reg->boundary_idx()], delimiter)
1047  << format_csv_val(fluxes_out_[qi][reg->boundary_idx()], delimiter)
1048  << csv_zero_vals(9, delimiter) << endl;
1050  }
1051  }
1052 
1053  if (cumulative_)
1054  {
1055  for (unsigned int qi=0; qi<n_quant; qi++)
1056  {
1057  // print data header (repeat header after every "repeat" lines)
1058  if (repeat && (output_line_counter_%repeat == 0)) format_csv_output_header(delimiter, comment_string);
1059 
1060  double error = sum_masses_[qi] - (initial_mass_[qi] + integrated_sources_[qi] + integrated_fluxes_[qi]);
1061  output_ << format_csv_val(time / time_->get_coef(), delimiter, true)
1062  << format_csv_val("ALL", delimiter)
1063  << format_csv_val(quantities_[qi].name_, delimiter)
1064  << format_csv_val(sum_fluxes_[qi], delimiter)
1065  << format_csv_val(sum_fluxes_in_[qi], delimiter)
1066  << format_csv_val(sum_fluxes_out_[qi], delimiter)
1067  << format_csv_val(sum_masses_[qi], delimiter)
1068  << format_csv_val(sum_sources_[qi], delimiter)
1069  << format_csv_val(sum_sources_in_[qi], delimiter)
1070  << format_csv_val(sum_sources_out_[qi], delimiter)
1071  << format_csv_val(increment_fluxes_[qi], delimiter)
1072  << format_csv_val(increment_sources_[qi], delimiter)
1073  << format_csv_val(integrated_fluxes_[qi], delimiter)
1074  << format_csv_val(integrated_sources_[qi], delimiter)
1075  << format_csv_val(error, delimiter) << endl;
1077  }
1078  }
1079 
1080 }
1081 
1082 
1083 void Balance::format_csv_output_header(char delimiter, const std::string& comment_string)
1084 {
1085  std::stringstream ss;
1086  if (delimiter == ' ') {
1087  ss << setw(output_column_width-comment_string.size()) << "\"time [" << time_->get_unit_string() << "]\"";
1088  } else {
1089  ss << "\"time [" << time_->get_unit_string() << "]\"";
1090  }
1091 
1092  output_ << comment_string << ss.str()
1093  << format_csv_val("region", delimiter)
1094  << format_csv_val("quantity [" + units_.format_text() + "]", delimiter)
1095  << format_csv_val("flux", delimiter)
1096  << format_csv_val("flux_in", delimiter)
1097  << format_csv_val("flux_out", delimiter)
1098  << format_csv_val("mass", delimiter)
1099  << format_csv_val("source", delimiter)
1100  << format_csv_val("source_in", delimiter)
1101  << format_csv_val("source_out", delimiter)
1102  << format_csv_val("flux_increment", delimiter)
1103  << format_csv_val("source_increment", delimiter)
1104  << format_csv_val("flux_cumulative", delimiter)
1105  << format_csv_val("source_cumulative", delimiter)
1106  << format_csv_val("error", delimiter)
1107  << endl;
1108 }
1109 
1110 std::string Balance::format_csv_val(std::string val, char delimiter, bool initial)
1111 {
1112  std::stringstream ss;
1113  std::replace( val.begin(), val.end(), '\"', '\'');
1114 
1115  if (!initial) ss << delimiter;
1116  if (delimiter == ' ') {
1117  std::stringstream sval;
1118  sval << "\"" << val << "\"";
1119  ss << " " << setw(output_column_width-1) << sval.str();
1120  } else {
1121  ss << "\"" << val << "\"";
1122  }
1123 
1124  return ss.str();
1125 }
1126 
1127 std::string Balance::format_csv_val(double val, char delimiter, bool initial)
1128 {
1129  std::stringstream ss;
1130 
1131  if (!initial) ss << delimiter;
1132  if (delimiter == ' ') {
1133  ss << " " << setw(output_column_width-1) << val;
1134  } else {
1135  ss << val;
1136  }
1137  return ss.str();
1138 }
1139 
1140 void Balance::output_yaml(double time)
1141 {
1142 
1143  // write output only on process #0
1144  if (!do_yaml_output_ || rank_ != 0) return;
1145 
1146  const unsigned int n_quant = quantities_.size();
1147 
1148  // print data header only once
1149  if (!output_yaml_header_) {
1150  output_yaml_ << "column_names: [ flux, flux_in, flux_out, mass, source, source_in, source_out, flux_increment, "
1151  << "source_increment, flux_cumulative, source_cumulative, error ]" << endl;
1152  output_yaml_ << "data:" << endl;
1153  output_yaml_header_ = true;
1154  }
1155 
1156  output_yaml_ << setfill(' ');
1157 
1158  // print sources and masses over bulk regions
1159  const RegionSet & bulk_set = mesh_->region_db().get_region_set("BULK");
1160  for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
1161  {
1162  for (unsigned int qi=0; qi<n_quant; qi++)
1163  {
1164  output_yaml_ << " - time: " << (time / time_->get_coef()) << endl;
1165  output_yaml_ << setw(4) << "" << "region: " << reg->label() << endl;
1166  output_yaml_ << setw(4) << "" << "quantity: " << quantities_[qi].name_ << endl;
1167  output_yaml_ << setw(4) << "" << "data: " << "[ 0, 0, 0, " << masses_[qi][reg->bulk_idx()] << ", "
1168  << sources_in_[qi][reg->bulk_idx()] + sources_out_[qi][reg->bulk_idx()] << ", "
1169  << sources_in_[qi][reg->bulk_idx()] << ", " << sources_out_[qi][reg->bulk_idx()]
1170  << ", 0, 0, 0, 0, 0 ]" << endl;
1171  }
1172  }
1173 
1174  // print mass fluxes over boundaries
1175  const RegionSet & b_set = mesh_->region_db().get_region_set(".BOUNDARY");
1176  for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg)
1177  {
1178  for (unsigned int qi=0; qi<n_quant; qi++) {
1179  output_yaml_ << " - time: " << (time / time_->get_coef()) << endl;
1180  output_yaml_ << setw(4) << "" << "region: " << reg->label() << endl;
1181  output_yaml_ << setw(4) << "" << "quantity: " << quantities_[qi].name_ << endl;
1182  output_yaml_ << setw(4) << "" << "data: " << "[ "
1183  << fluxes_in_[qi][reg->boundary_idx()] + fluxes_out_[qi][reg->boundary_idx()] << ", "
1184  << fluxes_in_[qi][reg->boundary_idx()] << ", " << fluxes_out_[qi][reg->boundary_idx()]
1185  << ", 0, 0, 0, 0, 0, 0, 0, 0, 0 ]" << endl;
1186  }
1187  }
1188 
1189  if (cumulative_)
1190  {
1191  for (unsigned int qi=0; qi<n_quant; qi++)
1192  {
1193  double error = sum_masses_[qi] - (initial_mass_[qi] + integrated_sources_[qi] + integrated_fluxes_[qi]);
1194  output_yaml_ << " - time: " << (time / time_->get_coef()) << endl;
1195  output_yaml_ << setw(4) << "" << "region: ALL" << endl;
1196  output_yaml_ << setw(4) << "" << "quantity: " << quantities_[qi].name_ << endl;
1197  output_yaml_ << setw(4) << "" << "data: " << "[ " << sum_fluxes_[qi] << ", "
1198  << sum_fluxes_in_[qi] << ", " << sum_fluxes_out_[qi] << ", "
1199  << sum_masses_[qi] << ", " << sum_sources_[qi] << ", "
1200  << sum_sources_in_[qi] << ", " << sum_sources_out_[qi] << ", "
1201  << increment_fluxes_[qi] << ", " << increment_sources_[qi] << ", "
1202  << integrated_fluxes_[qi] << ", " << integrated_sources_[qi] << ", "
1203  << error << " ]" << endl;
1204  }
1205  }
1206 }
1207 
1208 
1209 
1210 
1211 
1212 
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)
Definition: balance.cc:506
UnitSI units_
Units of conserved quantities.
Definition: balance.hh:460
void lazy_initialize()
Definition: balance.cc:202
static const Input::Type::Record & get_input_type()
Main balance input record type.
Definition: balance.cc:50
unsigned int n_loc_dofs_par_
Allocation parameters. Set by the allocate method used in the lazy_initialize.
Definition: balance.hh:434
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
std::vector< double > integrated_sources_
Definition: balance.hh:515
unsigned int add_quantity(const string &name)
Definition: balance.cc:162
std::vector< double > sum_fluxes_out_
Definition: balance.hh:507
void calculate_mass(unsigned int quantity_idx, const Vec &solution, vector< double > &output_array)
Definition: balance.cc:604
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:329
LongIdx get_boundary_edge_uid(SideIter side)
Definition: balance.hh:426
TimeMark::Type output_mark_type_
TimeMark type for output of particular equation.
Definition: balance.hh:527
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
TimeMark::Type balance_output_type_
TimeMark type for balance output of particular equation.
Definition: balance.hh:524
unsigned int * boundary_idx_
Definition: elements.h:80
unsigned int uint
std::string format_text() const
Definition: unit_si.cc:127
int tlevel() const
void finish_source_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix and vectors for computing source.
Definition: balance.cc:429
unsigned int bulk_size() const
Definition: region.cc:269
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Class for declaration of the input of type Bool.
Definition: type_base.hh:459
gnuplot
Definition: balance.hh:150
bool initial_
true before calculating the mass at initial time, otherwise false
Definition: balance.hh:530
void allocate(const std::shared_ptr< DOFHandlerMultiDim > &dh, unsigned int max_dofs_per_boundary)
Definition: balance.cc:191
std::vector< std::vector< double > > sources_in_
Definition: balance.hh:501
void add_cumulative_source(unsigned int quantity_idx, double source)
Definition: balance.cc:535
std::vector< std::vector< double > > fluxes_in_
Definition: balance.hh:498
std::vector< double > sum_fluxes_in_
Definition: balance.hh:506
std::vector< std::vector< double > > fluxes_out_
Definition: balance.hh:499
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:64
const TimeGovernor * time_
Definition: balance.hh:557
#define MPI_SUM
Definition: mpi.h:196
Mat * region_source_rhs_
Matrices for calculation of signed source (n_dofs x n_bulk_regions).
Definition: balance.hh:473
Definition: mesh.h:78
Input::Record input_record_
Record for current balance.
Definition: balance.hh:555
Cell accessor allow iterate over DOF handler cells.
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.
Definition: system.hh:148
void output_legacy(double time)
Perform output in old format (for compatibility)
Definition: balance.cc:849
SideIter side(const unsigned int loc_index)
const RegionDB & region_db() const
Definition: mesh.h:143
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
Region region()
const TimeStep & step(int index=-1) const
bool allocation_done_
true before allocating necessary internal structures (Petsc matrices etc.)
Definition: balance.hh:536
double t() const
void calculate_cumulative(unsigned int quantity_idx, const Vec &solution)
Definition: balance.cc:545
Basic time management functionality for unsteady (and steady) solvers (class Equation).
Side side() const
Return Side of given cell and side_idx.
static TimeMarks & marks()
#define MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm)
Definition: mpi.h:608
void calculate_instant(unsigned int quantity_idx, const Vec &solution)
Definition: balance.cc:629
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Basic time management class.
std::vector< double > sum_fluxes_
Definition: balance.hh:505
virtual ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
Definition: mesh.cc:839
bool balance_on_
If the balance is on. Balance is off in the case of no balance output time marks. ...
Definition: balance.hh:539
std::vector< double > increment_sources_
Definition: balance.hh:518
Vec * be_flux_vec_
Vectors for calculation of flux (n_boundary_edges).
Definition: balance.hh:476
static void set_yaml_output()
Set global variable to output balance files into YAML format (in addition to the table format)...
Definition: balance.cc:63
static const Input::Type::Selection & get_format_selection_input_type()
Input selection for file format.
Definition: balance.cc:42
std::vector< double > integrated_fluxes_
Definition: balance.hh:516
unsigned int boundary_idx() const
Returns index of the region in the boundary set.
Definition: region.hh:86
std::vector< Quantity > quantities_
Names of conserved quantities.
Definition: balance.hh:455
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...
Definition: balance.cc:994
std::vector< std::vector< double > > sources_out_
Definition: balance.hh:502
Mat * be_flux_matrix_
Matrices for calculation of flux (n_boundary_edges x n_dofs).
Definition: balance.hh:467
double last_time_
time of last calculated balance
Definition: balance.hh:521
void open_stream(Stream &stream) const
Definition: file_path.cc:211
double init_time() const
std::vector< std::vector< double > > fluxes_
Definition: balance.hh:497
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
bool is_boundary() const
Returns true for side on the boundary.
std::ofstream output_
Handle for file for output in given OutputFormat of balance and total fluxes over individual regions ...
Definition: balance.hh:446
void finish_mass_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix for computing mass.
Definition: balance.cc:407
TimeMark::Type equation_fixed_mark_type() const
void format_csv_output_header(char delimiter, const std::string &comment_string)
Print output header.
Definition: balance.cc:1083
int rank_
MPI rank.
Definition: balance.hh:546
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
TimeMark::Type equation_mark_type() const
const Ret val(const string &key) const
unsigned int n_sides() const
Definition: elements.h:132
void init_from_input(const Input::Record &in_rec, TimeGovernor &tg)
Definition: balance.cc:129
Mat * region_mass_matrix_
Matrices for calculation of mass (n_dofs x n_bulk_regions).
Definition: balance.hh:464
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
std::vector< double > sum_sources_out_
Definition: balance.hh:511
unsigned int max_dofs_per_boundary_
Definition: balance.hh:436
bool cumulative_
if true then cumulative balance is computed
Definition: balance.hh:533
unsigned int n_loc_dofs_seq_
Definition: balance.hh:435
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:503
Mat * region_source_matrix_
Matrices for calculation of source (n_dofs x n_bulk_regions).
Definition: balance.hh:470
std::vector< double > increment_fluxes_
Definition: balance.hh:517
void chkerr_assert(unsigned int ierr)
Definition: system.hh:159
bool is_current()
Returns true if the current time step is marked for the balance output.
Definition: balance.cc:369
Balance(const std::string &file_prefix, const Mesh *mesh)
Definition: balance.cc:86
void output_csv(double time, char delimiter, const std::string &comment_string, unsigned int repeat=0)
Perform output in csv format.
Definition: balance.cc:1002
const DOFHandlerMultiDim * dh() const
Return DOF handler.
Distribution * get_el_ds() const
Definition: mesh.h:154
~Balance()
Definition: balance.cc:103
static bool do_yaml_output_
Definition: balance.hh:431
#define MPI_DOUBLE
Definition: mpi.h:156
#define MPI_Comm_rank
Definition: mpi.h:236
std::vector< unsigned int > add_quantities(const std::vector< string > &names)
Definition: balance.cc:173
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
const Mesh * mesh_
Definition: balance.hh:457
static Default read_time(const std::string &description)
The factory function to make an default value that will be specified at the time when a key will be r...
Definition: type_record.hh:97
Support classes for parallel programing.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
void output_yaml(double time)
Perform output in yaml format.
Definition: balance.cc:1140
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
void output()
Perform output to file for given time instant.
Definition: balance.cc:706
std::vector< double > sum_sources_
Definition: balance.hh:509
std::string get_unit_string() const
void start_flux_assembly(unsigned int quantity_idx)
Definition: balance.cc:389
std::vector< double > initial_mass_
Definition: balance.hh:512
static const unsigned int output_column_width
Size of column in output (used if delimiter is space)
Definition: balance.hh:391
std::ofstream output_yaml_
Definition: balance.hh:449
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
const Selection & close() const
Close the Selection, no more values can be added.
unsigned int boundary_size() const
Definition: region.cc:262
double dt() const
void finish_flux_assembly(unsigned int quantity_idx)
This method must be called after assembling the matrix and vector for computing flux.
Definition: balance.cc:418
std::vector< unsigned int > be_regions_
Maps local boundary edge to its region boundary index.
Definition: balance.hh:489
RegionIdx region_idx() const
Definition: accessors.hh:168
OutputFormat output_format_
Format of output file.
Definition: balance.hh:452
#define ASSERT_DBG(expr)
std::unordered_map< LongIdx, unsigned int > be_id_map_
Definition: balance.hh:486
static const Input::Type::Array get_input_type()
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:91
void units(const UnitSI &unit)
Setter for units of conserved quantities.
Definition: balance.cc:154
int be_offset_
Offset for local part of vector of boundary edges.
Definition: balance.hh:492
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)
Definition: balance.cc:472
std::vector< std::vector< double > > masses_
Definition: balance.hh:500
const std::vector< LongIdx > & get_local_to_global_map() const
Get the map between local dof indices and the global ones.
Definition: dofhandler.hh:258
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)
Definition: balance.cc:441
void start_mass_assembly(unsigned int quantity_idx)
Definition: balance.cc:380
Record type proxy class.
Definition: type_record.hh:182
std::vector< double > sum_masses_
Definition: balance.hh:508
Vec * region_mass_vec_
Vectors for calculation of mass (n_bulk_regions).
Definition: balance.hh:479
unsigned int output_line_counter_
hold count of line printed into output_
Definition: balance.hh:549
bool add_output_times_
Add output time marks to balance output time marks.
Definition: balance.hh:542
Class for representation SI units of Fields.
Definition: unit_si.hh:40
Boundary cond() const
double get_coef() const
static FileName output()
The factory function for declaring type FileName for input files.
Definition: type_base.cc:533
OutputFormat
Definition: balance.hh:146
LongIdx * get_el_4_loc() const
Definition: mesh.h:160
void undef(bool val=true)
Definition: unit_si.cc:197
void start_source_assembly(unsigned int quantity_idx)
Definition: balance.cc:398
Template for classes storing finite set of named values.
Side accessor allows to iterate over sides of DOF handler cell.
FilePath balance_output_file_
File path for output_ stream.
Definition: balance.hh:443
std::string file_prefix_
Save prefix passed in in constructor.
Definition: balance.hh:440
bool output_yaml_header_
marks whether YAML output has printed header
Definition: balance.hh:552
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.
Definition: balance.cc:1110
std::vector< double > sum_sources_in_
Definition: balance.hh:510
unsigned int lsize(int proc) const
get local size