Flow123d
mass_balance.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  *
17  * You should have received a copy of the GNU General Public License along with this program; if not,
18  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
19  *
20  *
21  * $Id$
22  * $Revision$
23  * $LastChangedBy$
24  * $LastChangedDate$
25  *
26  * @file
27  * @ingroup transport
28  * @brief Mass balance
29  *
30  *
31  */
32 
33 #include <iostream>
34 #include <iomanip>
35 
36 #include "system/system.hh"
37 #include "system/sys_profiler.hh"
38 #include "system/xio.h"
39 
40 #include <petscmat.h>
41 #include "mesh/mesh.h"
42 
44 
45 using namespace Input::Type;
46 
48  = Record("MassBalance", "Balance of mass, boundary fluxes and sources for transport of substances.")
49  .declare_key("cumulative", Bool(), Default("false"), "Compute cumulative balance over time. "
50  "If true, then balance is calculated at each computational time step, which can slow down the program.")
51  .declare_key("file", FileName::output(), Default("mass_balance.txt"), "File name for output of mass balance.")
52 ;
53 
54 
55 
57  : equation_(eq), initial_time(0), last_time(-1), initial(true)
58 {
59  balance_output_file = xfopen( in_rec.val<FilePath>("file"), "wt");
60 
61  cumulative = in_rec.val<bool>("cumulative");
62 
63  const int n_bcd_reg_ = equation_->region_db()->boundary_size();
64  const int n_blk_reg_ = equation_->region_db()->bulk_size();
65 
66  bcd_balance.resize(eq->n_substances(), vector<double>(n_bcd_reg_, 0));
67  bcd_plus_balance.resize(eq->n_substances(), vector<double>(n_bcd_reg_, 0));
68  bcd_minus_balance.resize(eq->n_substances(), vector<double>(n_bcd_reg_, 0));
69  mass.resize(eq->n_substances(), vector<double>(n_blk_reg_, 0));
70  src_balance.resize(eq->n_substances(), vector<double>(n_blk_reg_, 0));
71 
72  bcd_total_balance.resize(eq->n_substances(), 0.);
73  bcd_total_inflow.resize(eq->n_substances(), 0.);
74  bcd_total_outflow.resize(eq->n_substances(), 0.);
75  mass_total.resize(eq->n_substances(), 0.);
76  src_total_balance.resize(eq->n_substances(), 0.);
77 
78  initial_mass.resize(eq->n_substances(), 0.);
79  integrated_sources.resize(eq->n_substances(), 0.);
80  integrated_fluxes.resize(eq->n_substances(), 0.);
81 }
82 
84 {
86 }
87 
88 
89 void MassBalance::calculate(double time) {
90  F_ENTRY;
91 
92  // return if we already calculated at the given time
93  if (last_time == time) return;
94 
95  const int n_bcd_reg_ = equation_->region_db()->boundary_size();
96  const int n_blk_reg_ = equation_->region_db()->bulk_size();
97  const int n_subst = equation_->n_substances();
98 
99  for (int i=0; i<n_subst; i++)
100  {
101  for (int j=0; j<n_bcd_reg_; j++)
102  {
103  bcd_balance[i][j] = 0;
104  bcd_plus_balance[i][j] = 0;
105  bcd_minus_balance[i][j] = 0;
106  }
107  for (int j=0; j<n_blk_reg_; j++)
108  {
109  mass[i][j] = 0;
110  src_balance[i][j] = 0;
111  }
112  bcd_total_balance[i] = 0;
113  bcd_total_outflow[i] = 0;
114  bcd_total_inflow[i] = 0;
115  mass_total[i] = 0;
116  src_total_balance[i] = 0;
117  }
118 
119  // compute fluxes, mass and volume sources
122 
123 
124  // gather results from processes and sum them up
125  int buf_size = n_subst*(3*n_bcd_reg_ + 2*n_blk_reg_);
126  double sendbuffer[buf_size], recvbuffer[buf_size];
127  for (int i=0; i<n_subst; i++)
128  {
129  for (int j=0; j<n_bcd_reg_; j++)
130  {
131  sendbuffer[i*3*n_bcd_reg_+ j] = bcd_balance[i][j];
132  sendbuffer[i*3*n_bcd_reg_+ n_bcd_reg_+j] = bcd_plus_balance[i][j];
133  sendbuffer[i*3*n_bcd_reg_+2*n_bcd_reg_+j] = bcd_minus_balance[i][j];
134  }
135  for (int j=0; j<n_blk_reg_; j++)
136  {
137  sendbuffer[n_subst*3*n_bcd_reg_+i*2*n_blk_reg_ +j] = mass[i][j];
138  sendbuffer[n_subst*3*n_bcd_reg_+i*2*n_blk_reg_+n_blk_reg_+j] = src_balance[i][j];
139  }
140  }
141  MPI_Reduce(&sendbuffer,recvbuffer,buf_size,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);
142 
143  int rank;
144  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
145 
146  // for other than 0th process update last_time and finish,
147  // on process #0 sum balances over all regions and calculate
148  // cumulative balance over time.
149  if (rank != 0)
150  {
151  last_time = time;
152  return;
153  }
154 
155 
156  // update balance vectors
157  for (int i=0; i<n_subst; i++)
158  {
159  for (int j=0; j<n_bcd_reg_; j++)
160  {
161  bcd_balance[i][j] = recvbuffer[i*3*n_bcd_reg_+ j];
162  bcd_plus_balance[i][j] = recvbuffer[i*3*n_bcd_reg_+ n_bcd_reg_+j];
163  bcd_minus_balance[i][j] = recvbuffer[i*3*n_bcd_reg_+2*n_bcd_reg_+j];
164  }
165  for (int j=0; j<n_blk_reg_; j++)
166  {
167  mass[i][j] = recvbuffer[n_subst*3*n_bcd_reg_+i*2*n_blk_reg_ +j];
168  src_balance[i][j] = recvbuffer[n_subst*3*n_bcd_reg_+i*2*n_blk_reg_+n_blk_reg_+j];
169  }
170  }
171 
172 
173  // sum all boundary fluxes
174  const RegionSet & b_set = equation_->region_db()->get_region_set("BOUNDARY");
175  for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg) {
176  //DBGMSG("writing reg->idx() and id() and boundary_idx(): %d\t%d\t%d\n", reg->idx(), reg->id(), reg->boundary_idx());
177  for (int sbi=0; sbi<n_subst; sbi++) {
178  bcd_total_balance[sbi] += bcd_balance[sbi][reg->boundary_idx()];
179  bcd_total_outflow[sbi] += bcd_plus_balance[sbi][reg->boundary_idx()];
180  bcd_total_inflow[sbi] += bcd_minus_balance[sbi][reg->boundary_idx()];
181  }
182  }
183 
184  // sum all volume sources
185  const RegionSet & bulk_set = equation_->region_db()->get_region_set("BULK");
186  for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
187  {
188  for (int sbi=0; sbi<n_subst; sbi++)
189  {
190  mass_total[sbi] += mass[sbi][reg->bulk_idx()];
191  src_total_balance[sbi] += src_balance[sbi][reg->bulk_idx()];
192  }
193  }
194 
195  if (!cumulative) return;
196 
197  // cumulative balance over time
198 
199  // quantities need for the balance over time interval
200  static vector<double> last_sources(n_subst, 0.),
201  last_fluxes(n_subst, 0.);
202 
203  // save initial time and mass
204  if (initial)
205  {
206  initial_time = time;
208  for (int i=0; i<n_subst; i++)
209  initial_mass[i] = mass_total[i];
210  initial = false;
211  }
212 
213 
214  // sum sources and fluxes according to the time integration scheme
215  // TODO: Think if we really need TimeIntegrationScheme or the cummulative
216  // quantities can be calculated in the same way for both explicit and implicit
217  // methods.
218  switch (equation_->time_scheme())
219  {
221 // for (int i=0; i<n_subst; i++)
222 // {
223 // integrated_sources[i] += last_sources[i]*(time-last_time);
224 // integrated_fluxes[i] += last_fluxes[i]*(time-last_time);
225 // }
226 // break;
228  for (int i=0; i<n_subst; i++)
229  {
232  }
233  break;
235  for (int i=0; i<n_subst; i++)
236  {
237  integrated_sources[i] += (last_sources[i]+src_total_balance[i])*0.5*(time-last_time);
238  integrated_fluxes[i] += (last_fluxes[i]+bcd_total_balance[i])*0.5*(time-last_time);
239  }
240  break;
241  default:
242  break;
243  }
244 
245  last_time = time;
246  for (int i=0; i<n_subst; i++)
247  {
248  last_sources[i] = src_total_balance[i];
249  last_fluxes[i] = bcd_total_balance[i];
250  }
251 }
252 
253 
254 void MassBalance::output(double time)
255 {
256  // calculate balances for the given time
257  if (last_time != time) calculate(time);
258 
259  int rank;
260  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
261 
262  // write output only on process #0
263  if (rank != 0) return;
264 
265  const int n_subst = equation_->n_substances();
266 
267  // print the head of mass balance file
268  unsigned int c = 6; //column number without label
269  unsigned int w = 14; //column width
270  unsigned int wl = 2*(w-5)+7; //label column width
271  stringstream s; //helpful stringstream
272  string bc_head_format = "# %-*s%-*s%-*s%-*s%-*s%-*s\n",
273  bc_format = "%*s%-*d%-*s%-*s%-*g%-*g%-*g\n",
274  bc_total_format = "# %-*s%-*s%-*g%-*g%-*g\n";
275  s << setw((w*c+wl-14)/2) << setfill('-') << "--"; //drawing half line
276  fprintf(balance_output_file,"# %s MASS BALANCE %s\n",s.str().c_str(), s.str().c_str());
277  fprintf(balance_output_file,"# Time: %f\n\n\n",time);
278 
279  // header for table of boundary fluxes
280  fprintf(balance_output_file,"# Mass flux through boundary [M/T]:\n");
281  fprintf(balance_output_file,bc_head_format.c_str(),w,"[boundary_id]",wl,"[label]",
282  w,"[substance]",w,"[total flux]",w,"[outward flux]",w,"[inward flux]");
283  s.clear();
284  s.str(std::string());
285  s << setw(w*c+wl) << setfill('-') << "-";
286  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //drawing long line
287 
288  // print mass fluxes over boundaries
289  const RegionSet & b_set = equation_->region_db()->get_region_set("BOUNDARY");
290  for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg) {
291  for (int sbi=0; sbi<n_subst; sbi++) {
292  fprintf(balance_output_file, bc_format.c_str(),2,"",w,reg->id(),wl,reg->label().c_str(),
293  w, equation_->substance_names()[sbi].c_str(),w, bcd_balance[sbi][reg->boundary_idx()],
294  w, bcd_plus_balance[sbi][reg->boundary_idx()],
295  w, bcd_minus_balance[sbi][reg->boundary_idx()]);
296  }
297  }
298  // total boundary balance
299  fprintf(balance_output_file,"# %s\n",s.str().c_str()); // drawing long line
300  for (int sbi=0; sbi<n_subst; sbi++)
301  fprintf(balance_output_file, bc_total_format.c_str(),w+wl,"Total mass flux of substance [M/T]",
302  w,equation_->substance_names()[sbi].c_str(),w,bcd_total_balance[sbi], w, bcd_total_outflow[sbi], w, bcd_total_inflow[sbi]);
303  fprintf(balance_output_file, "\n\n");
304 
305 
306  // header for table of volume sources and masses
307  string src_head_format = "# %-*s%-*s%-*s%-*s%-*s\n",
308  src_format = "%*s%-*d%-*s%-*s%-*g%-*g\n",
309  src_total_format = "# %-*s%-*s%-*g%-*g\n";
310  fprintf(balance_output_file,"# Mass [M] and sources [M/T] on regions:\n"); //head
311  fprintf(balance_output_file,src_head_format.c_str(),w,"[region_id]",wl,"[label]",
312  w,"[substance]",w,"[total_mass]",w,"[total_source]");
313  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //long line
314 
315  // print balance of volume sources and masses
316  const RegionSet & bulk_set = equation_->region_db()->get_region_set("BULK");
317  for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
318  {
319  for (int sbi=0; sbi<n_subst; sbi++)
320  {
321  fprintf(balance_output_file, src_format.c_str(), 2,"", w, reg->id(), wl,
322  reg->label().c_str(), w, equation_->substance_names()[sbi].c_str(),
323  w,mass[sbi][reg->bulk_idx()],
324  w,src_balance[sbi][reg->bulk_idx()]);
325  }
326  }
327  // total sources balance
328  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //drawing long line
329  for (int sbi=0; sbi<n_subst; sbi++)
330  fprintf(balance_output_file, src_total_format.c_str(),w+wl,"Total mass [M] and sources [M/T]",
331  w,equation_->substance_names()[sbi].c_str(),
332  w,mass_total[sbi],
333  w,src_total_balance[sbi]);
334 
335  if (cumulative)
336  {
337  // Print cumulative sources
338  fprintf(balance_output_file, "\n\n# Cumulative mass balance on time interval [%-g,%-g]\n"
339  "# Initial mass [M] + sources integrated over time [M] - flux integrated over time [M] = current mass [M]\n"
340  "# %-*s%-*s%-*s%-*s%-*s%-*s%-*s%-*s\n",
341  initial_time, time,
342  w,"[substance]",
343  w,"[A=init. mass]",
344  w,"[B=source]",
345  w,"[C=flux]",
346  w,"[A+B-C]",
347  w,"[D=curr. mass]",
348  w,"[A+B-C-D=err.]",
349  w,"[rel. error]");
350 
351  for (int i=0; i<n_subst; i++)
352  {
353  double denominator = max(fabs(initial_mass[i]+integrated_sources[i]-integrated_fluxes[i]),fabs(mass_total[i]));
354  fprintf(balance_output_file, " %-*s%-*g%-*g%-*g%-*g%-*g%-*g%-*g\n",
355  w,equation_->substance_names()[i].c_str(),
356  w,initial_mass[i],
357  w,integrated_sources[i],
358  w,integrated_fluxes[i],
359  w,initial_mass[i]+integrated_sources[i]-integrated_fluxes[i],
360  w,mass_total[i],
361  w,initial_mass[i]+integrated_sources[i]-integrated_fluxes[i]-mass_total[i],
362  w,fabs(initial_mass[i]+integrated_sources[i]-integrated_fluxes[i]-mass_total[i])/(denominator==0?1:denominator));
363  }
364  }
365 
366  fprintf(balance_output_file, "\n\n");
367 }
368 
369