Flow123d  jenkins-Flow123d-windows32-release-multijob-28
darcy_flow_mh_output.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  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id: darcy_flow_mh.hh 877 2011-02-04 13:13:25Z jakub.sistek $
21  * $Revision: 877 $
22  * $LastChangedBy: jakub.sistek $
23  * $LastChangedDate: 2011-02-04 14:13:25 +0100 (Fri, 04 Feb 2011) $
24  *
25  * @file
26  * @brief Output class for darcy_flow_mh model.
27  * @ingroup flow
28  *
29  * @author Jan Brezina
30  *
31  */
32 
33 
34 #include <vector>
35 #include <iostream>
36 #include <sstream>
37 #include <string>
38 
39 #include <system/global_defs.h>
40 
41 #include "flow/mh_fe_values.hh"
42 #include "flow/darcy_flow_mh.hh"
44 #include "system/system.hh"
45 #include "system/sys_profiler.hh"
46 #include "system/xio.h"
47 
48 #include "fem/dofhandler.hh"
49 #include "fields/field_fe.hh"
50 
51 #include "io/output.h"
52 #include "mesh/partitioning.hh"
53 
54 
55 
56 namespace it = Input::Type;
57 
59  = DarcyFlowMH::EqData().make_output_field_selection("DarcyMHOutput_Selection", "Selection of fields available for output.")
60  .copy_values(OutputFields().make_output_field_selection("").close())
61  .close();
62 
64  = it::Record("DarcyMHOutput", "Parameters of MH output.")
66  "Parameters of output stream.")
67  .declare_key("output_fields", it::Array(OutputFields::output_selection),
68  it::Default::obligatory(), "List of fields to write to output file.")
69  .declare_key("balance_output", it::FileName::output(), it::Default("water_balance.txt"),
70  "Output file for water balance table.")
71  .declare_key("compute_errors", it::Bool(), it::Default("false"),
72  "SPECIAL PURPOSE. Computing errors pro non-compatible coupling.")
74  "Output file with raw data form MH module.");
75 
76 
78 {
79 // *this += darcy_flow->get_data();
80  *this += field_ele_pressure.name("pressure_p0").units("L");
81  *this += field_node_pressure.name("pressure_p1").units("L");
82  *this += field_ele_piezo_head.name("piezo_head_p0").units("L");
83  *this += field_ele_flux.name("velocity_p0").units("L/T");
84  *this += subdomain.name("subdomain").units("0");
85 
86  fields_for_output += *this;
87 
88  *this += pressure_diff.name("pressure_diff").units("0");
89  *this += velocity_diff.name("velocity_diff").units("0");
90  *this += div_diff.name("div_diff").units("0");
91 
92 }
93 
94 
96 : darcy_flow(flow),
97  mesh_(&darcy_flow->mesh()),
98  in_rec_(in_rec),
99  balance_output_file(NULL),
100  raw_output_file(NULL)
101 {
102 
103 
104 
105  // we need to add data from the flow equation at this point, not in constructor of OutputFields
108 
109  //VecCreateSeqWithArray(PETSC_COMM_SELF, 1, dh->n_global_dofs(), corner_pressure, &vec_corner_pressure);
110 
111  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
112  ele_pressure.resize(mesh_->n_elements());
113  auto ele_pressure_ptr=make_shared< FieldElementwise<3, FieldValue<3>::Scalar> >(ele_pressure, 1);
115 
116  dh = new DOFHandlerMultiDim(*mesh_);
118  corner_pressure.resize(dh->n_global_dofs());
119  VecCreateSeqWithArray(PETSC_COMM_SELF, 1, dh->n_global_dofs(), &(corner_pressure[0]), &vec_corner_pressure);
120 
121  auto corner_ptr = make_shared< FieldFE<3, FieldValue<3>::Scalar> >();
122  corner_ptr->set_fe_data(dh, &map1, &map2, &map3, &vec_corner_pressure);
123 
126 
127  ele_piezo_head.resize(mesh_->n_elements());
129  make_shared< FieldElementwise<3, FieldValue<3>::Scalar> >(ele_piezo_head, 1));
130 
131  ele_flux.resize(3*mesh_->n_elements());
133  make_shared< FieldElementwise<3, FieldValue<3>::VectorFixed> >(ele_flux, 3));
134 
135  auto &vec_int_sub = mesh_->get_part()->seq_output_partition();
136  subdomains.resize(vec_int_sub.size());
137  for(unsigned int i=0; i<subdomains.size();i++)
138  subdomains[i]=vec_int_sub[i];
140  make_shared< FieldElementwise<3, FieldValue<3>::Integer> >(subdomains, 1));
141 
143 
144 
148  //DBGMSG("output stream: %p\n", output_stream);
149 
150  int rank;
151  // int ierr
153 
154  if (rank == 0) {
155  // temporary solution for balance output
156  balance_output_file = xfopen( in_rec.val<FilePath>("balance_output"), "wt");
157 
158  // optionally open raw output file
159  FilePath raw_output_file_path;
160  if (in_rec.opt_val("raw_flow_output", raw_output_file_path)) {
161  xprintf(Msg, "Opening raw output: %s\n", string(raw_output_file_path).c_str() );
162  raw_output_file = xfopen(raw_output_file_path, "wt");
163  }
164 
165  }
166 }
167 
168 
169 
171  //if (output_writer != NULL) delete output_writer;
172 
175  //delete [] ele_pressure;
176  //delete [] node_pressure;
177  //delete [] ele_flux;
178  VecDestroy(&vec_corner_pressure);
179  //delete [] corner_pressure;
180  //delete [] ele_piezo_head;
181  delete output_stream;
182 
183  delete dh;
184 };
185 
186 
187 
188 
189 
190 //=============================================================================
191 // CONVERT SOLUTION, CALCULATE BALANCES, ETC...
192 //=============================================================================
193 
194 
196 {
197  START_TIMER("Darcy output");
198 
199 // std::string eleVectorName = "velocity_elements";
200 // std::string eleVectorUnit = "L/T";
201 
202  //cout << "DMHO_output: rank: " << rank << "\t output_writer: " << output_writer << endl;
203 
204  // skip initial output for steady solver
205  //darcy_flow->time().view("Darcy output output");
206 
207  //if (darcy_flow->time().is_steady() && darcy_flow->time().tlevel() ==0) return;
208 
209  //DBGVAR( darcy_flow->time().is_current( TimeGovernor::marks().type_output() ) );
210  if (darcy_flow->time().is_current( TimeGovernor::marks().type_output() )) {
211 
214 
216 
217  water_balance();
218 
219  if (in_rec_.val<bool>("compute_errors")) compute_l2_difference();
220 
224 
226 
227  }
228 
229 }
230 
231 
232 //=============================================================================
233 // FILL TH "SCALAR" FIELD FOR ALL ELEMENTS IN THE MESH
234 //=============================================================================
235 
237  unsigned int sol_size;
238  double *sol;
239 
240  darcy_flow->get_solution_vector(sol, sol_size);
241  unsigned int soi = mesh_->n_sides();
242  unsigned int i = 0;
243  FOR_ELEMENTS(mesh_,ele) {
244  ele_pressure[i] = sol[ soi];
245  ele_piezo_head[i] = sol[soi ] + ele->centre()[Mesh::z_coord];
246  i++; soi++;
247  }
248 }
249 
250 
251 
252 /****
253  * compute Darcian velocity in centre of elements
254  *
255  */
258  MHFEValues fe_values;
259 
260  int i_side=0;
261  FOR_ELEMENTS(mesh_, ele) {
262  arma::vec3 flux_in_centre;
263  flux_in_centre.zeros();
264 
265  fe_values.update(ele,
269 
270  for (unsigned int li = 0; li < ele->n_sides(); li++) {
271  flux_in_centre += dh.side_flux( *(ele->side( li ) ) )
272  * fe_values.RT0_value( ele, ele->centre(), li )
273  / darcy_flow->data_.cross_section.value(ele->centre(), ele->element_accessor() );
274  }
275 
276  for(unsigned int j=0; j<3; j++)
277  ele_flux[3*i_side+j]=flux_in_centre[j];
278 
279  i_side++;
280  }
281 
282 }
283 
284 
286 {
287  unsigned int ndofs = max(dh->fe<1>()->n_dofs(), max(dh->fe<2>()->n_dofs(), dh->fe<3>()->n_dofs()));
288  unsigned int indices[ndofs];
289  unsigned int i_node;
290  FOR_ELEMENTS(mesh_, ele)
291  {
292  dh->get_dof_indices(ele, indices);
293  FOR_ELEMENT_NODES(ele, i_node)
294  {
295  corner_pressure[indices[i_node]] = node_scalar[mesh_->node_vector.index(ele->node[i_node])];
296  }
297  }
298 }
299 
300 
301 //=============================================================================
302 // pressure interpolation
303 //
304 // some sort of weighted average over elements and sides/edges of an node
305 //
306 // TODO:
307 // questions:
308 // - explain details of averaging and motivation for this type
309 // - edge scalars are stronger since thay are computed twice by each side
310 // - why division by (nod.aux-1)
311 // - ?implement without TED, TSD global arrays with single loop over elements, sides and one loop over nodes for normalization
312 //
313 //=============================================================================
314 
316 
317  vector<double> scalars(mesh_->n_nodes());
318 
319  double dist; //!< tmp variable for storing particular distance node --> element, node --> side*/
320 
321  /** Iterators */
322  const Node * node;
323  //ElementIter ele;
324  //struct Side* side;
325 
326  int n_nodes = mesh_->node_vector.size(); //!< number of nodes in the mesh */
327  int node_index = 0; //!< index of each node */
328 
329  int* sum_elements = new int [n_nodes]; //!< sum elements joined to node */
330  int* sum_sides = new int [n_nodes]; //!< sum sides joined to node */
331  double* sum_ele_dist = new double [n_nodes]; //!< sum distances to all joined elements */
332  double* sum_side_dist = new double [n_nodes]; //!< Sum distances to all joined sides */
333 
334  /** tmp variables, will be replaced by ini keys
335  * TODO include them into ini file*/
336  bool count_elements = true; //!< scalar is counted via elements*/
337  bool count_sides = true; //!< scalar is counted via sides */
338 
339 
341 
342  /** init arrays */
343  for (int i = 0; i < n_nodes; i++){
344  sum_elements[i] = 0;
345  sum_sides[i] = 0;
346  sum_ele_dist[i] = 0.0;
347  sum_side_dist[i] = 0.0;
348  scalars[i] = 0.0;
349  };
350 
351  /**first pass - calculate sums (weights)*/
352  if (count_elements){
353  FOR_ELEMENTS(mesh_, ele)
354  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
355  node = ele->node[li]; //!< get Node pointer from element */
356  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
357 
358  dist = sqrt(
359  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
360  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
361  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
362  );
363  sum_ele_dist[node_index] += dist;
364  sum_elements[node_index]++;
365  }
366  }
367  if (count_sides){
368  FOR_SIDES(mesh_, side) {
369  for (unsigned int li = 0; li < side->n_nodes(); li++) {
370  node = side->node(li);//!< get Node pointer from element */
371  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
372  dist = sqrt(
373  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
374  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
375  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
376  );
377 
378  sum_side_dist[node_index] += dist;
379  sum_sides[node_index]++;
380  }
381  }
382  }
383 
384  /**second pass - calculate scalar */
385  if (count_elements){
386  FOR_ELEMENTS(mesh_, ele)
387  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
388  node = ele->node[li];//!< get Node pointer from element */
389  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
390 
391  /**TODO - calculate it again or store it in prior pass*/
392  dist = sqrt(
393  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
394  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
395  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
396  );
397  scalars[node_index] += ele_pressure[ele.index()] *
398  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
399  (sum_elements[node_index] + sum_sides[node_index] - 1);
400  }
401  }
402  if (count_sides) {
403  FOR_SIDES(mesh_, side) {
404  for (unsigned int li = 0; li < side->n_nodes(); li++) {
405  node = side->node(li);//!< get Node pointer from element */
406  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
407 
408  /**TODO - calculate it again or store it in prior pass*/
409  dist = sqrt(
410  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
411  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
412  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
413  );
414 
415 
416  scalars[node_index] += dh.side_scalar( *side ) *
417  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
418  (sum_sides[node_index] + sum_elements[node_index] - 1);
419  }
420  }
421  }
422 
423 // xprintf(Msg, "**********************************************************************************************\n");
424 // for (int i =0; i<n_nodes; i++){
425 // xprintf(Msg, "make_node_scalar_param id: %i, %f\n", i, scalars[i]);
426 // };
427 // xprintf(Msg, "**********************************************************************************************\n");
428 
429  /** free memory */
430  delete [] sum_elements;
431  delete [] sum_sides;
432  delete [] sum_ele_dist;
433  delete [] sum_side_dist;
434 
435  make_corner_scalar(scalars);
436 }
437 
438 
439 /*
440 void DarcyFlowMHOutput::make_node_scalar() {
441  int li;
442  NodeIter nod;
443  struct Side *sde;
444  double dist;
445  int max_side_id = 0;
446 
447  double **TED;
448  double **TSD;
449 
450 
451  TED = (double **) xmalloc((mesh_->element.size() + 1) * sizeof (double *));
452 
453  FOR_SIDES(mesh_, sde)
454  if (max_side_id <= sde->id)
455  max_side_id = sde->id;
456 
457  TSD = (double **) xmalloc((max_side_id + 1) * sizeof (double *));
458 
459  FOR_ELEMENTS(mesh_, ele)
460  TED[ele.index()] = (double*) xmalloc(ele->n_nodes * sizeof (double));
461  FOR_SIDES(mesh_, sde)
462  TSD[sde->id] = (double*) xmalloc(sde->n_nodes * sizeof (double));
463 
464  FOR_NODES(mesh_, nod ) {
465  nod->scalar = 0.0;
466  nod->faux = 0.0;
467  nod->aux = 0;
468  }
469  FOR_ELEMENTS(mesh_, ele)
470  for (li = 0; li < ele->n_nodes; li++) {
471  nod = ele->node[li];
472 
473  dist = sqrt(
474  ((nod->getX() - ele->centre()[ 0 ])*(nod->getX() - ele->centre()[ 0 ])) +
475  ((nod->getY() - ele->centre()[ 1 ])*(nod->getY() - ele->centre()[ 1 ])) +
476  ((nod->getZ() - ele->centre()[ 2 ])*(nod->getZ() - ele->centre()[ 2 ]))
477  );
478 
479  TED[ele.index()][li] = dist;
480  nod->faux += dist; // nod->faux += 1 / dist;
481  nod->aux++;
482  }
483 
484  FOR_SIDES(mesh_, sde) {
485  for (li = 0; li < sde->n_nodes; li++) {
486  nod = sde->node[li];
487 
488  dist = sqrt(
489  ((nod->getX() - sde->centre()[ 0 ])*(nod->getX() - sde->centre()[ 0 ])) +
490  ((nod->getY() - sde->centre()[ 1 ])*(nod->getY() - sde->centre()[ 1 ])) +
491  ((nod->getZ() - sde->centre()[ 2 ])*(nod->getZ() - sde->centre()[ 2 ]))
492  );
493 
494  TSD[sde->id][li] = dist;
495  nod->faux += dist; // nod->faux += 1 / dist;
496  nod->aux++;
497  }
498  }
499  FOR_ELEMENTS(mesh_, ele)
500  for (li = 0; li < ele->n_nodes; li++) {
501  nod = ele->node[li];
502  nod->scalar += ele->scalar * (1 - TED[ele.index()][li] / nod->faux)
503  / (nod->aux - 1); // 1 / (dist * nod->faux);
504  }
505  FOR_SIDES(mesh_, sde)
506  for (li = 0; li < sde->n_nodes; li++) {
507  nod = sde->node[li];
508  nod->scalar += sde->scalar * (1 - TSD[sde->id][li] / nod->faux)
509  / (nod->aux - 1); // 1 / (dist * nod->faux);
510  }
511  xfree(TED);
512  xfree(TSD);
513 }
514 */
515 /*
516 //=============================================================================
517 //
518 //=============================================================================
519 void make_node_vector(Mesh* mesh)
520 {
521  int ni;
522  int nei;
523  TNode* nod;
524  ElementIter ele;
525  double cnt;
526 
527  for( ni = 0; ni < mesh->n_nodes; ni++ ) {
528  nod = mesh->node + ni;
529  nod->vector[ 0 ] = 0.0;
530  nod->vector[ 1 ] = 0.0;
531  cnt = 0.0;
532  for( nei = 0; nei < nod->n_elements(); nei++ ) {
533  ele = nod->element[ nei ];
534  nod->vector[ 0 ] += ele->vector[ 0 ];
535  nod->vector[ 1 ] += ele->vector[ 1 ];
536  cnt += 1.0;
537  }
538  nod->vector[ 0 ] /= cnt;
539  nod->vector[ 1 ] /= cnt;
540 
541  }
542 }
543  */
544 //=============================================================================
545 // FILL TH "FLUX" FIELD FOR ALL VV NEIGHBOURS IN THE MESH
546 //=============================================================================
547 /*
548 void DarcyFlowMHOutput::make_neighbour_flux() {
549  struct Neighbour *ngh;
550 
551  FOR_NEIGHBOURS(mesh_, ngh) {
552  if (ngh->type != VV_2E)
553  continue;
554 
555  ngh->flux = ngh->sigma * ngh->geom_factor * (ngh->element[1]->scalar - ngh->element[0]->scalar);
556  continue;
557  }
558 }*/
559 
560 //=============================================================================
561 //
562 //=============================================================================
563 
566  if (balance_output_file == NULL) return;
567 
568  //BOUNDARY
569  //struct Boundary *bcd;
570  std::vector<double> bcd_balance( mesh_->region_db().boundary_size(), 0.0 );
571  std::vector<double> bcd_plus_balance( mesh_->region_db().boundary_size(), 0.0 );
572  std::vector<double> bcd_minus_balance( mesh_->region_db().boundary_size(), 0.0 );
573 
574  using namespace std;
575  //printing the head of water balance file
576  unsigned int c = 5; //column number without label
577  unsigned int w = 14; //column width
578  unsigned int wl = 2*(w-5)+7; //label column width
579  stringstream s; //helpful stringstream
580  string bc_head_format = "# %-*s%-*s%-*s%-*s%-*s%-*s\n",
581  bc_format = "%*s%-*d%-*s %-*g%-*g%-*g%-*g\n",
582  bc_total_format = "# %-*s%-*g%-*g%-*g\n\n\n";
583  s << setw((w*c+wl-15)/2) << setfill('-') << "-"; //drawing half line
584  fprintf(balance_output_file,"# %s WATER BALANCE %s\n",s.str().c_str(), s.str().c_str());
585  fprintf(balance_output_file,"# Time of computed water balance: %f\n\n\n",darcy_flow->time().t());
586 
587  fprintf(balance_output_file,"# Boundary water balance:\n");
588  fprintf(balance_output_file,bc_head_format.c_str(),w,"[boundary_id]",wl,"[label]",
589  w,"[total_balance]",w,"[total_outflow]",w,"[total_inflow]",w,"[time]");
590  s.clear();
591  s.str(std::string());
592  s << setw(w*c+wl) << setfill('-') << "-";
593  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //drawing long line
594 
595  //computing water balance over boundaries
596  FOR_BOUNDARIES(mesh_, bcd) {
597  // !! there can be more sides per one boundary
598  double flux = dh.side_flux( *(bcd->side()) );
599 
600  Region r = bcd->region();
601  //DBGMSG("flux: %f side: %d %d reg: %s\n", flux, bcd->side()->element()->index(), bcd->side()->el_idx(), r.label().c_str() );
602  if (! r.is_valid()) xprintf(Msg, "Invalid region, ele % d, edg: % d\n", bcd->bc_ele_idx_, bcd->edge_idx_);
603  unsigned int bc_region_idx = r.boundary_idx();
604  bcd_balance[bc_region_idx] += flux;
605 
606  if (flux > 0) bcd_plus_balance[bc_region_idx] += flux;
607  else bcd_minus_balance[bc_region_idx] += flux;
608  }
609  //printing water balance over boundaries
610  //DBGMSG("DB[boundary] size: %u\n", mesh_->region_db().boundary_size());
611  const RegionSet & b_set = mesh_->region_db().get_region_set("BOUNDARY");
612  double total_balance = 0, // for computing total balance on boundary
613  total_inflow = 0,
614  total_outflow = 0;
615  for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg) {
616  //DBGMSG("writing reg->idx() and id() and boundary_idx(): %d\t%d\t%d\n", reg->idx(), reg->id(), reg->boundary_idx());
617  total_balance += bcd_balance[reg->boundary_idx()];
618  total_outflow += bcd_plus_balance[reg->boundary_idx()];
619  total_inflow += bcd_minus_balance[reg->boundary_idx()];
620  fprintf(balance_output_file, bc_format.c_str(),2,"",w,reg->id(),wl,reg->label().c_str(),
621  w, bcd_balance[reg->boundary_idx()], w, bcd_plus_balance[reg->boundary_idx()],
622  w, bcd_minus_balance[reg->boundary_idx()], w, darcy_flow->time().t());
623  }
624  //total boundary balance
625  fprintf(balance_output_file,"# %s\n",s.str().c_str()); // drawing long line
626  fprintf(balance_output_file, bc_total_format.c_str(),w+wl+2,"total boundary balance",
627  w,total_balance, w, total_outflow, w, total_inflow);
628 
629  //SOURCES
630  string src_head_format = "# %-*s%-*s%-*s%-*s%-*s\n",
631  src_format = "%*s%-*d%-*s %-*g%-*s%-*g\n",
632  src_total_format = "# %-*s%-*g\n\n\n";
633  //computing water balance of sources
634  fprintf(balance_output_file,"# Source fluxes over material subdomains:\n"); //head
635  fprintf(balance_output_file,src_head_format.c_str(),w,"[region_id]",wl,"[label]",
636  w,"[total_balance]",2*w,"",w,"[time]");
637  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //long line
638  std::vector<double> src_balance( mesh_->region_db().bulk_size(), 0.0 ); // initialize by zero
639  FOR_ELEMENTS(mesh_, elm) {
640  //DBGMSG("writing reg->idx() and id() and bulk_idx(): %d\t%d\t%d\n",
641  // elm->element_accessor().region().idx(),
642  // elm->element_accessor().region().id(),
643  // elm->element_accessor().region().bulk_idx());
644  src_balance[elm->element_accessor().region().bulk_idx()] += elm->measure() *
645  darcy_flow->data_.cross_section.value(elm->centre(), elm->element_accessor()) *
646  darcy_flow->data_.water_source_density.value(elm->centre(), elm->element_accessor());
647  }
648 
649  total_balance = 0;
650  //printing water balance of sources
651  //DBGMSG("DB[bulk] size: %u\n", mesh_->region_db().bulk_size());
652  const RegionSet & bulk_set = mesh_->region_db().get_region_set("BULK");
653  for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
654  {
655  total_balance += src_balance[reg->bulk_idx()];
656  //"%*s%-*d%-*s %-*g%-*s%-*g\n";
657  fprintf(balance_output_file, src_format.c_str(), 2,"", w, reg->id(), wl,
658  reg->label().c_str(), w, src_balance[reg->bulk_idx()],2*w,"", w,darcy_flow->time().t());
659  }
660  //total sources balance
661  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //drawing long line
662  fprintf(balance_output_file, src_total_format.c_str(),w+wl+2,"total sources balance",
663  w,total_balance);
664 }
665 
666 //=============================================================================
667 // UNUSED FUNCTION
668 //=============================================================================
669 /*
670 double calc_water_balance(Mesh* mesh, int c_water) {
671  double rc;
672  struct Boundary *bcd;
673 
674  rc = 0.0;
675 
676  return rc;
677 }
678 */
679 
680 /*
681  * Output of internal flow data.
682  */
683 
685 {
687 
688  if (raw_output_file == NULL) return;
689 
690  char dbl_fmt[ 16 ]= "%.8g ";
691  // header
692  xfprintf( raw_output_file, "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n]\n");
693  xfprintf( raw_output_file, "$FlowField\nT=");
694  xfprintf( raw_output_file, dbl_fmt, darcy_flow->time().t());
695  xfprintf( raw_output_file, "\n%d\n", mesh_->n_elements() );
696 
697  unsigned int i;
698  int cit = 0;
699  FOR_ELEMENTS( mesh_, ele ) {
700  //xfprintf( raw_output_file, "%d ", cit);
701  xfprintf( raw_output_file, "%d ", ele.id());
702  xfprintf( raw_output_file, dbl_fmt, ele_pressure[cit]);
703  for (i = 0; i < 3; i++)
704  xfprintf( raw_output_file, dbl_fmt, ele_flux[3*cit+i]);
705 
706  xfprintf( raw_output_file, " %d ", ele->n_sides());
707  for (i = 0; i < ele->n_sides(); i++)
708  xfprintf( raw_output_file, dbl_fmt, dh.side_scalar( *(ele->side(i) ) ) );
709  for (i = 0; i < ele->n_sides(); i++)
710  xfprintf( raw_output_file, dbl_fmt, dh.side_flux( *(ele->side(i) ) ) );
711 
712  //xfprintf( raw_output_file, "%d ", ele->n_neighs_vv);
713  //for (i = 0; i < ele->n_neighs_vv; i++)
714  // xfprintf( raw_output_file, "%d ", ele->neigh_vv[i]->id);
715 
716  xfprintf( raw_output_file, "\n" );
717  cit ++;
718  }
719  xfprintf( raw_output_file, "$EndFlowField\n\n" );
720 }
721 
722 
724 #include "fem/fe_p.hh"
725 #include "fem/fe_values.hh"
726 #include "fem/mapping_p1.hh"
727 //#include "functions/function_python.hh"
728 #include "fields/field_python.hh"
729 #include "fields/field_values.hh"
730 
732 
733 /*
734 * Calculate approximation of L2 norm for:
735  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
736  * 2) difference between RT velocities and analytical solution
737  * 3) difference of divergence
738  * */
739 
740 struct DiffData {
741  double pressure_error[2], velocity_error[2], div_error[2];
745 
746  double * solution;
747  const MH_DofHandler * dh;
749 
750  //double *ele_flux;
751 
754 };
755 
756 template <int dim>
757 void l2_diff_local(ElementFullIter &ele, FEValues<dim,3> &fe_values, ExactSolution &anal_sol, DiffData &result) {
758 
759  fe_values.reinit(ele);
760  result.fe_values.update(ele,
761  result.data_->anisotropy,
762  result.data_->cross_section,
763  result.data_->conductivity);
764  double conductivity = result.data_->conductivity.value(ele->centre(), ele->element_accessor() );
765  double cross = result.data_->cross_section.value(ele->centre(), ele->element_accessor() );
766 
767 
768  // get coefficients on the current element
769  vector<double> fluxes(dim+1);
770  vector<double> pressure_traces(dim+1);
771 
772  for (unsigned int li = 0; li < ele->n_sides(); li++) {
773  fluxes[li] = result.dh->side_flux( *(ele->side( li ) ) );
774  pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
775  }
776  double pressure_mean = result.dh->element_scalar(ele);
777 
778  arma::vec analytical(5);
779  arma::vec3 flux_in_q_point;
780  arma::vec3 anal_flux;
781 
782  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
783 
784  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
785  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
786  double mean_x_squared=0;
787  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
788  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
789  {
790  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
791  * arma::dot( ele->node[i_node]->point(), ele->node[j_node]->point());
792  }
793 
794  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
795  arma::vec3 q_point = fe_values.point(i_point);
796 
797 
798  analytical = anal_sol.value(q_point, ele->element_accessor() );
799  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
800 
801  // compute postprocesed pressure
802  diff = 0;
803  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
804  unsigned int oposite_node = (i_shape + dim) % (dim + 1);
805 
806  diff += fluxes[ i_shape ] *
807  ( arma::dot( q_point, q_point )/ 2
808  - mean_x_squared / 2
809  - arma::dot( q_point, ele->node[oposite_node]->point() )
810  + arma::dot( ele->centre(), ele->node[oposite_node]->point() )
811  );
812  }
813 
814  diff = - (1.0 / conductivity) * diff / dim / ele->measure() / cross + pressure_mean ;
815  diff = ( diff - analytical[0]);
816  pressure_diff += diff * diff * fe_values.JxW(i_point);
817 
818 
819  // velocity difference
820  flux_in_q_point.zeros();
821  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
822  flux_in_q_point += fluxes[ i_shape ]
823  * result.fe_values.RT0_value( ele, q_point, i_shape )
824  / cross;
825  }
826 
827  flux_in_q_point -= anal_flux;
828  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
829 
830  // divergence diff
831  diff = 0;
832  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
833  diff = ( diff / ele->measure() / cross - analytical[4]);
834  divergence_diff += diff * diff * fe_values.JxW(i_point);
835 
836  }
837 
838 
839  result.velocity_diff[ele.index()] = velocity_diff;
840  result.velocity_error[dim-1] += velocity_diff;
841 
842  result.pressure_diff[ele.index()] = pressure_diff;
843  result.pressure_error[dim-1] += pressure_diff;
844 
845  result.div_diff[ele.index()] = divergence_diff;
846  result.div_error[dim-1] += divergence_diff;
847 
848 }
849 
850 
851 
852 
853 
854 
856  DBGMSG("l2 norm output\n");
857  ofstream os( FilePath("solution_error", FilePath::output_file) );
858 
859  const unsigned int order = 4; // order of Gauss quadrature
860 
861  // we create trivial Dofhandler , for P0 elements, to get access to, FEValues on individual elements
862  // this we use to integrate our own functions - difference of postprocessed pressure and analytical solution
863  FE_P_disc<0,1,3> fe_1d;
864  FE_P_disc<0,2,3> fe_2d;
865 
866  QGauss<1> quad_1d( order );
867  QGauss<2> quad_2d( order );
868 
869  MappingP1<1,3> mapp_1d;
870  MappingP1<2,3> mapp_2d;
871 
872  FEValues<1,3> fe_values_1d(mapp_1d, quad_1d, fe_1d, update_JxW_values | update_quadrature_points);
873  FEValues<2,3> fe_values_2d(mapp_2d, quad_2d, fe_2d, update_JxW_values | update_quadrature_points);
874 
875  FilePath source_file( "analytical_module.py", FilePath::input_file);
876  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
877  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
878 
879  ExactSolution anal_sol_2d(5);
880  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
881 
882 
883  static DiffData result;
884 
885  result.pressure_diff.resize( mesh_->n_elements() );
886  result.velocity_diff.resize( mesh_->n_elements() );
887  result.div_diff.resize( mesh_->n_elements() );
888 
889  result.pressure_error[0] = 0;
890  result.velocity_error[0] = 0;
891  result.div_error[0] = 0;
892  result.pressure_error[1] = 0;
893  result.velocity_error[1] = 0;
894  result.div_error[1] = 0;
895 
896  //result.ele_flux = &(ele_flux[0]);
897 
898  //output_writer->register_elem_data(mesh_, "pressure_diff", "0", in_rec_.val<Input::Record>("output_stream") ,result.pressure_diff);
899  //output_writer->register_elem_data(mesh_, "velocity_diff", "0", in_rec_.val<Input::Record>("output_stream"),result.pressure_diff);
900  //output_writer->register_elem_data(mesh_, "div_diff", "0", in_rec_.val<Input::Record>("output_stream"),result.pressure_diff);
901 
902 
903 
904  auto vel_diff_ptr = std::make_shared< FieldElementwise<3, FieldValue<3>::Scalar> >(&(result.velocity_diff[0]), 1, mesh_->n_elements());
906  auto pressure_diff_ptr = std::make_shared< FieldElementwise<3, FieldValue<3>::Scalar> >(&(result.pressure_diff[0]), 1, mesh_->n_elements());
907  output_fields.pressure_diff.set_field(mesh_->region_db().get_region_set("ALL"), pressure_diff_ptr, 0);
908  auto div_diff_ptr = std::make_shared< FieldElementwise<3, FieldValue<3>::Scalar> >(&(result.div_diff[0]), 1, mesh_->n_elements());
909  output_fields.div_diff.set_field(mesh_->region_db().get_region_set("ALL"), div_diff_ptr, 0);
910 
912 
913  unsigned int solution_size;
914  darcy_flow->get_solution_vector(result.solution, solution_size);
915 
916  result.dh = &( darcy_flow->get_mh_dofhandler());
917  result.darcy = darcy_flow;
918  result.data_ = &(darcy_flow->data_);
919 
920  FOR_ELEMENTS( mesh_, ele) {
921 
922  switch (ele->dim()) {
923  case 1:
924 
925  l2_diff_local<1>( ele, fe_values_1d, anal_sol_1d, result);
926  break;
927  case 2:
928  l2_diff_local<2>( ele, fe_values_2d, anal_sol_2d, result);
929  break;
930  }
931  }
932 
933  os << "l2 norm output\n\n"
934  << "pressure error 1d: " << sqrt(result.pressure_error[0]) << endl
935  << "pressure error 2d: " << sqrt(result.pressure_error[1]) << endl
936  << "velocity error 1d: " << sqrt(result.velocity_error[0]) << endl
937  << "velocity error 2d: " << sqrt(result.velocity_error[1]) << endl
938  << "div error 1d: " << sqrt(result.div_error[0]) << endl
939  << "div error 2d: " << sqrt(result.div_error[1]);
940 }
941 
FieldSet & data()
Definition: equation.hh:197
const MH_DofHandler * dh
void set_python_field_from_file(const FilePath &file_name, const string &func_name)
Field< 3, FieldValue< 3 >::VectorFixed > field_ele_flux
vector< double > ele_piezo_head
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Field< 3, FieldValue< 3 >::Scalar > velocity_diff
void set_limit_side(LimitSide side)
Definition: field_set.hh:171
#define DBGMSG(...)
Definition: global_defs.h:195
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
vector< double > div_diff
Output class for darcy_flow_mh model.
DarcyFlowMH_Steady * darcy_flow
#define FOR_ELEMENT_NODES(i, j)
Definition: elements.h:150
virtual void get_solution_vector(double *&vec, unsigned int &vec_size)
const arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
Definition: fe_values.cc:184
Header: The functions for all outputs.
Accessor to input data conforming to declared Array.
Definition: accessors.hh:521
Definition: system.hh:75
double getY() const
Definition: nodes.hh:71
static Input::Type::Record input_type
The specification of output stream.
Definition: output_time.hh:57
Definition: nodes.hh:44
unsigned int bulk_size() const
Definition: region.cc:256
static Input::Type::Record input_type
Field< 3, FieldValue< 3 >::Scalar > field_node_pressure
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:39
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:359
Class for declaration of the input of type Bool.
Definition: type_base.hh:321
MappingP1< 1, 3 > map1
int xfclose(FILE *stream)
FCLOSE WITH ERROR HANDLING.
Definition: xio.cc:309
DarcyFlowMHOutput(DarcyFlowMH_Steady *flow, Input::Record in_rec)
const double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.cc:178
DOFHandlerMultiDim * dh
static Default obligatory()
Definition: type_record.hh:87
Selection & copy_values(const Selection &sel)
int xfprintf(FILE *out, const char *fmt,...)
FPRINTF WITH ERROR HANDLING.
Definition: xio.cc:395
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:365
void distribute_dofs(FiniteElement< 1, 3 > &fe1d, FiniteElement< 2, 3 > &fe2d, FiniteElement< 3, 3 > &fe3d, const unsigned int offset=0)
Distributes degrees of freedom on the mesh needed for the given finite elements.
Definition: dofhandler.cc:260
void make_corner_scalar(vector< double > &node_scalar)
Class FEValues calculates finite element data on the actual cells such as shape function values...
int index() const
Definition: sys_vector.hh:88
bool is_current(const TimeMark::Type &mask) const
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
const RegionDB & region_db() const
Definition: mesh.h:148
const MH_DofHandler & get_mh_dofhandler()
Field< 3, FieldValue< 3 >::Integer > subdomain
double t() const
I/O functions with filename storing, able to track current line in opened file. All standard stdio fu...
Partitioning * get_part()
Definition: mesh.cc:163
unsigned int n_sides()
Definition: mesh.cc:152
static TimeMarks & marks()
Class for declaration of inputs sequences.
Definition: type_base.hh:230
double div_error[2]
FieldCommon & units(const string &units)
Set basic units of the field.
FiniteElement< dim, 3 > * fe() const
Symmetric Gauss-Legendre quadrature formulae on simplices.
const unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.cc:202
DarcyFlowMH_Steady * darcy
unsigned int boundary_idx() const
Returns index of the region in the boundary set.
Definition: region.hh:78
double getZ() const
Definition: nodes.hh:73
unsigned int size() const
Returns size of the container. This is independent of the allocated space.
Definition: sys_vector.hh:412
static Default optional()
Definition: type_record.hh:100
Data for Darcy flow equation.
bool opt_val(const string &key, Ret &value) const
Discontinuous Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:213
Field< 3, FieldValue< 3 >::Scalar > field_ele_pressure
Transformed quadrature points.
Global macros to enhance readability and debugging, general constants.
unsigned int n_elements() const
Definition: mesh.h:137
FE_P_disc< 1, 2, 3 > fe2
static OutputTime * create_output_stream(const Input::Record &in_rec)
This method write all registered data to output streams.
Definition: output.cc:171
double element_scalar(ElementFullIter &ele) const
temporary replacement for DofHandler accessor, scalar (pressure) on element
double pressure_error[2]
void update(ElementFullIter ele, FieldType &cond_anisothropy, FieldType_Scalar &cross_section, FieldType_Scalar &conductivity)
Definition: mh_fe_values.cc:34
MappingP1< 2, 3 > map2
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Field< 3, FieldValue< 3 >::TensorFixed > anisotropy
FILE * raw_output_file
Raw data output file.
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:104
double getX() const
Definition: nodes.hh:69
#define START_TIMER(tag)
Starts a timer with specified tag.
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:251
double velocity_error[2]
MHFEValues fe_values
unsigned int index(const T *pointer) const
Definition: sys_vector.hh:394
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:399
FILE * balance_output_file
Temporary solution for writing balance into separate file.
vector< double > ele_pressure
vector< double > velocity_diff
Input::Record in_rec_
Accessor to the input record for the DarcyFlow output.
vector< double > subdomains
void mark_output_times(const TimeGovernor &tg)
Definition: output.cc:331
void set_time(const TimeGovernor &time)
Definition: field_set.hh:186
static Input::Type::Selection output_selection
MappingP1< 3, 3 > map3
#define MPI_Comm_rank
Definition: mpi.h:236
Dedicated class for storing path to input and output files.
Definition: file_path.hh:34
vector< int > & seq_output_partition()
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
void write_time_frame()
Definition: output.cc:362
Field< 3, FieldValue< 3 >::Scalar > conductivity
double side_flux(const Side &side) const
temporary replacement for DofHandler accessor, flux through given side
arma::vec3 RT0_value(ElementFullIter ele, arma::vec3 point, unsigned int face)
Definition: mh_fe_values.cc:76
Field< 3, FieldValue< 3 >::Scalar > pressure_diff
void output(OutputTime *stream)
Definition: field_set.cc:137
const Selection & close() const
unsigned int boundary_size() const
Definition: region.cc:249
static FileName output()
Definition: type_base.hh:449
void get_dof_indices(const CellIterator &cell, unsigned int indices[]) const override
Returns the global indices of dofs associated to the cell.
Definition: dofhandler.cc:349
void add_admissible_field_names(const Input::Array &in_array, const Input::Type::Selection &in_sel)
Registers names of output fields that can be written using this stream.
Definition: output.cc:231
TimeGovernor const & time()
Definition: equation.hh:152
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:205
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:244
#define FOR_BOUNDARIES(_mesh_, i)
Definition: mesh.h:377
Mixed-hybrid of steady Darcy flow with sources and variable density.
Input::Type::Selection make_output_field_selection(const string &name, const string &desc="")
Definition: field_set.cc:70
Definitions of particular quadrature rules on simplices.
FieldCommon & name(const string &name)
Definition: field_common.hh:71
#define MPI_COMM_WORLD
Definition: mpi.h:123
Calculates finite element data on the actual cell.
Definition: fe_values.hh:370
Field< 3, FieldValue< 3 >::Scalar > div_diff
FE_P_disc< 1, 3, 3 > fe3
unsigned int n_nodes() const
Definition: mesh.h:133
vector< double > corner_pressure
mixed-hybrid model of linear Darcy flow, possibly unsteady.
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:157
Record type proxy class.
Definition: type_record.hh:161
Field< 3, FieldValue< 3 >::Scalar > cross_section
DarcyFlowMH_Steady::EqData * data_
bool is_valid() const
Returns false if the region has undefined/invalid value.
Definition: region.hh:70
double side_scalar(const Side &side) const
temporary replacement for DofHandler accessor, scalar (pressure) on edge of the side ...
FE_P_disc< 1, 1, 3 > fe1
void l2_diff_local(ElementFullIter &ele, FEValues< dim, 3 > &fe_values, ExactSolution &anal_sol, DiffData &result)
vector< double > ele_flux
#define FOR_SIDES(_mesh_, it)
Definition: mesh.h:403
void make_node_scalar_param()
Calculate nodes scalar, store it in double* node_scalars instead of node->scalar. ...
Template for classes storing finite set of named values.
const unsigned int n_global_dofs() const
Getter for the number of all mesh dofs required by the given finite element.
Definition: dofhandler.hh:72
Field< 3, FieldValue< 3 >::Scalar > water_source_density
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
FieldCommon & output_type(OutputTime::DiscreteSpace rt)
vector< double > pressure_diff
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:196
FILE * xfopen(const std::string &fname, const char *mode)
Definition: xio.cc:246
void output()
Calculate values for output.
Transformed quadrature weights.
Record & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:390