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