Flow123d  jenkins-Flow123d-windows-release-multijob-285
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 
45 #include "io/output_time.hh"
46 #include "system/system.hh"
47 #include "system/sys_profiler.hh"
48 #include "system/xio.h"
49 
50 #include "fem/dofhandler.hh"
51 #include "fields/field_fe.hh"
52 #include "fields/generic_field.hh"
53 
54 #include "mesh/partitioning.hh"
55 
57 
58 
59 namespace it = Input::Type;
60 
62  = DarcyFlowMH::EqData().make_output_field_selection("DarcyMHOutput_Selection", "Selection of fields available for output.")
63  .copy_values(OutputFields().make_output_field_selection("").close())
64  .close();
65 
67  = it::Record("DarcyMHOutput", "Parameters of MH output.")
69  "Parameters of output stream.")
70  .declare_key("output_fields", it::Array(OutputFields::output_selection),
71  it::Default::obligatory(), "List of fields to write to output file.")
72 // .declare_key("balance_output", it::FileName::output(), it::Default("water_balance.txt"),
73 // "Output file for water balance table.")
74  .declare_key("compute_errors", it::Bool(), it::Default("false"),
75  "SPECIAL PURPOSE. Computing errors pro non-compatible coupling.")
77  "Output file with raw data form MH module.");
78 
79 
81 {
82  *this += field_ele_pressure.name("pressure_p0").units(UnitSI().m());
83  *this += field_node_pressure.name("pressure_p1").units(UnitSI().m());
84  *this += field_ele_piezo_head.name("piezo_head_p0").units(UnitSI().m());
85  *this += field_ele_flux.name("velocity_p0").units(UnitSI().m().s(-1));
86  *this += subdomain.name("subdomain")
89  *this += region_id.name("region_id")
92 
93  fields_for_output += *this;
94 
95  *this += pressure_diff.name("pressure_diff").units(UnitSI().m());
96  *this += velocity_diff.name("velocity_diff").units(UnitSI().m().s(-1));
97  *this += div_diff.name("div_diff").units(UnitSI().s(-1));
98 
99 }
100 
101 
103 : darcy_flow(flow),
104  mesh_(&darcy_flow->mesh()),
105  in_rec_(in_rec),
106  balance_output_file(NULL),
107  raw_output_file(NULL)
108 {
109 
110 
111 
112  // we need to add data from the flow equation at this point, not in constructor of OutputFields
115 
116  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
118  auto ele_pressure_ptr=ele_pressure.create_field<3, FieldValue<3>::Scalar>(1);
120 
121  dh = new DOFHandlerMultiDim(*mesh_);
123  corner_pressure.resize(dh->n_global_dofs());
124  VecCreateSeqWithArray(PETSC_COMM_SELF, 1, dh->n_global_dofs(), &(corner_pressure[0]), &vec_corner_pressure);
125 
126  auto corner_ptr = make_shared< FieldFE<3, FieldValue<3>::Scalar> >();
127  corner_ptr->set_fe_data(dh, &map1, &map2, &map3, &vec_corner_pressure);
128 
131 
133  auto ele_piezo_head_ptr=ele_piezo_head.create_field<3, FieldValue<3>::Scalar>(1);
135 
137  auto ele_flux_ptr=ele_flux.create_field<3, FieldValue<3>::VectorFixed>(3);
139 
142 
144 
145 
147  output_stream->add_admissible_field_names(in_rec.val<Input::Array>("output_fields"));
149 
150  int rank;
152 
153  if (rank == 0) {
154  // temporary solution for balance output
155 // balance_output_file = xfopen( in_rec.val<FilePath>("balance_output"), "wt");
156 
157  // optionally open raw output file
158  FilePath raw_output_file_path;
159  if (in_rec.opt_val("raw_flow_output", raw_output_file_path)) {
160  xprintf(Msg, "Opening raw output: %s\n", string(raw_output_file_path).c_str() );
161  raw_output_file = xfopen(raw_output_file_path, "wt");
162  }
163 
164  }
165 
166  flow->balance_->units(
167  get_output_fields().field_ele_pressure.units()
168  *flow->data_.cross_section.units()*UnitSI().md(1)
169  *flow->data_.storativity.units());
170 }
171 
172 
173 
175 
178  VecDestroy(&vec_corner_pressure);
179  delete output_stream;
180 
181  delete dh;
182 };
183 
184 
185 
186 
187 
188 //=============================================================================
189 // CONVERT SOLUTION, CALCULATE BALANCES, ETC...
190 //=============================================================================
191 
192 
194 {
195  START_TIMER("Darcy output");
196 
197  if (darcy_flow->time().is_current( TimeGovernor::marks().type_output() )) {
198 
201 
203 
204 // water_balance();
205 
206  if (in_rec_.val<bool>("compute_errors")) compute_l2_difference();
207 
211 
213 
214  }
215 
216 }
217 
218 
219 //=============================================================================
220 // FILL TH "SCALAR" FIELD FOR ALL ELEMENTS IN THE MESH
221 //=============================================================================
222 
224  unsigned int sol_size;
225  double *sol;
226 
227  darcy_flow->get_solution_vector(sol, sol_size);
228  unsigned int soi = mesh_->n_sides();
229  unsigned int i = 0;
230  FOR_ELEMENTS(mesh_,ele) {
231  ele_pressure[i] = sol[ soi];
232  ele_piezo_head[i] = sol[soi ] + ele->centre()[Mesh::z_coord];
233  i++; soi++;
234  }
235 }
236 
237 
238 
239 /****
240  * compute Darcian velocity in centre of elements
241  *
242  */
245  MHFEValues fe_values;
246 
247  int i_side=0;
248  FOR_ELEMENTS(mesh_, ele) {
249  arma::vec3 flux_in_centre;
250  flux_in_centre.zeros();
251 
252  fe_values.update(ele,
256 
257  for (unsigned int li = 0; li < ele->n_sides(); li++) {
258  flux_in_centre += dh.side_flux( *(ele->side( li ) ) )
259  * fe_values.RT0_value( ele, ele->centre(), li )
260  / darcy_flow->data_.cross_section.value(ele->centre(), ele->element_accessor() );
261  }
262 
263  for(unsigned int j=0; j<3; j++)
264  ele_flux[3*i_side+j]=flux_in_centre[j];
265 
266  i_side++;
267  }
268 
269 }
270 
271 
273 {
274  unsigned int ndofs = max(dh->fe<1>()->n_dofs(), max(dh->fe<2>()->n_dofs(), dh->fe<3>()->n_dofs()));
275  unsigned int indices[ndofs];
276  unsigned int i_node;
277  FOR_ELEMENTS(mesh_, ele)
278  {
279  dh->get_dof_indices(ele, indices);
280  FOR_ELEMENT_NODES(ele, i_node)
281  {
282  corner_pressure[indices[i_node]] = node_scalar[mesh_->node_vector.index(ele->node[i_node])];
283  }
284  }
285 }
286 
287 
288 //=============================================================================
289 // pressure interpolation
290 //
291 // some sort of weighted average over elements and sides/edges of an node
292 //
293 // TODO:
294 // questions:
295 // - explain details of averaging and motivation for this type
296 // - edge scalars are stronger since thay are computed twice by each side
297 // - why division by (nod.aux-1)
298 // - ?implement without TED, TSD global arrays with single loop over elements, sides and one loop over nodes for normalization
299 //
300 //=============================================================================
301 
303 
304  vector<double> scalars(mesh_->n_nodes());
305 
306  double dist; //!< tmp variable for storing particular distance node --> element, node --> side*/
307 
308  /** Iterators */
309  const Node * node;
310 
311  int n_nodes = mesh_->node_vector.size(); //!< number of nodes in the mesh */
312  int node_index = 0; //!< index of each node */
313 
314  int* sum_elements = new int [n_nodes]; //!< sum elements joined to node */
315  int* sum_sides = new int [n_nodes]; //!< sum sides joined to node */
316  double* sum_ele_dist = new double [n_nodes]; //!< sum distances to all joined elements */
317  double* sum_side_dist = new double [n_nodes]; //!< Sum distances to all joined sides */
318 
319  /** tmp variables, will be replaced by ini keys
320  * TODO include them into ini file*/
321  bool count_elements = true; //!< scalar is counted via elements*/
322  bool count_sides = true; //!< scalar is counted via sides */
323 
324 
326 
327  /** init arrays */
328  for (int i = 0; i < n_nodes; i++){
329  sum_elements[i] = 0;
330  sum_sides[i] = 0;
331  sum_ele_dist[i] = 0.0;
332  sum_side_dist[i] = 0.0;
333  scalars[i] = 0.0;
334  };
335 
336  /**first pass - calculate sums (weights)*/
337  if (count_elements){
338  FOR_ELEMENTS(mesh_, ele)
339  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
340  node = ele->node[li]; //!< get Node pointer from element */
341  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
342 
343  dist = sqrt(
344  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
345  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
346  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
347  );
348  sum_ele_dist[node_index] += dist;
349  sum_elements[node_index]++;
350  }
351  }
352  if (count_sides){
353  FOR_SIDES(mesh_, side) {
354  for (unsigned int li = 0; li < side->n_nodes(); li++) {
355  node = side->node(li);//!< get Node pointer from element */
356  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
357  dist = sqrt(
358  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
359  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
360  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
361  );
362 
363  sum_side_dist[node_index] += dist;
364  sum_sides[node_index]++;
365  }
366  }
367  }
368 
369  /**second pass - calculate scalar */
370  if (count_elements){
371  FOR_ELEMENTS(mesh_, ele)
372  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
373  node = ele->node[li];//!< get Node pointer from element */
374  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
375 
376  /**TODO - calculate it again or store it in prior pass*/
377  dist = sqrt(
378  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
379  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
380  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
381  );
382  scalars[node_index] += ele_pressure[ele.index()] *
383  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
384  (sum_elements[node_index] + sum_sides[node_index] - 1);
385  }
386  }
387  if (count_sides) {
388  FOR_SIDES(mesh_, side) {
389  for (unsigned int li = 0; li < side->n_nodes(); li++) {
390  node = side->node(li);//!< get Node pointer from element */
391  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
392 
393  /**TODO - calculate it again or store it in prior pass*/
394  dist = sqrt(
395  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
396  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
397  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
398  );
399 
400 
401  scalars[node_index] += dh.side_scalar( *side ) *
402  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
403  (sum_sides[node_index] + sum_elements[node_index] - 1);
404  }
405  }
406  }
407 
408  /** free memory */
409  delete [] sum_elements;
410  delete [] sum_sides;
411  delete [] sum_ele_dist;
412  delete [] sum_side_dist;
413 
414  make_corner_scalar(scalars);
415 }
416 
417 
418 //=============================================================================
419 //
420 //=============================================================================
421 
424  if (balance_output_file == NULL) return;
425 
426  //BOUNDARY
427  //struct Boundary *bcd;
428  std::vector<double> bcd_balance( mesh_->region_db().boundary_size(), 0.0 );
429  std::vector<double> bcd_plus_balance( mesh_->region_db().boundary_size(), 0.0 );
430  std::vector<double> bcd_minus_balance( mesh_->region_db().boundary_size(), 0.0 );
431 
432  using namespace std;
433  //printing the head of water balance file
434  unsigned int c = 5; //column number without label
435  unsigned int w = 14; //column width
436  unsigned int wl = 2*(w-5)+7; //label column width
437  stringstream s; //helpful stringstream
438  string bc_head_format = "# %-*s%-*s%-*s%-*s%-*s%-*s\n",
439  bc_format = "%*s%-*d%-*s %-*g%-*g%-*g%-*g\n",
440  bc_total_format = "# %-*s%-*g%-*g%-*g\n\n\n";
441  s << setw((w*c+wl-15)/2) << setfill('-') << "-"; //drawing half line
442  fprintf(balance_output_file,"# %s WATER BALANCE %s\n",s.str().c_str(), s.str().c_str());
443  fprintf(balance_output_file,"# Time of computed water balance: %f\n\n\n",darcy_flow->time().t());
444 
445  fprintf(balance_output_file,"# Boundary water balance:\n");
446  fprintf(balance_output_file,bc_head_format.c_str(),w,"[boundary_id]",wl,"[label]",
447  w,"[total_balance]",w,"[total_outflow]",w,"[total_inflow]",w,"[time]");
448  s.clear();
449  s.str(std::string());
450  s << setw(w*c+wl) << setfill('-') << "-";
451  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //drawing long line
452 
453  //computing water balance over boundaries
454  FOR_BOUNDARIES(mesh_, bcd) {
455  // !! there can be more sides per one boundary
456  double flux = dh.side_flux( *(bcd->side()) );
457 
458  Region r = bcd->region();
459  if (! r.is_valid()) xprintf(Msg, "Invalid region, ele % d, edg: % d\n", bcd->bc_ele_idx_, bcd->edge_idx_);
460  unsigned int bc_region_idx = r.boundary_idx();
461  bcd_balance[bc_region_idx] += flux;
462 
463  if (flux > 0) bcd_plus_balance[bc_region_idx] += flux;
464  else bcd_minus_balance[bc_region_idx] += flux;
465  }
466  //printing water balance over boundaries
467  const RegionSet & b_set = mesh_->region_db().get_region_set("BOUNDARY");
468  double total_balance = 0, // for computing total balance on boundary
469  total_inflow = 0,
470  total_outflow = 0;
471  for( RegionSet::const_iterator reg = b_set.begin(); reg != b_set.end(); ++reg) {
472  total_balance += bcd_balance[reg->boundary_idx()];
473  total_outflow += bcd_plus_balance[reg->boundary_idx()];
474  total_inflow += bcd_minus_balance[reg->boundary_idx()];
475  fprintf(balance_output_file, bc_format.c_str(),2,"",w,reg->id(),wl,reg->label().c_str(),
476  w, bcd_balance[reg->boundary_idx()], w, bcd_plus_balance[reg->boundary_idx()],
477  w, bcd_minus_balance[reg->boundary_idx()], w, darcy_flow->time().t());
478  }
479  //total boundary balance
480  fprintf(balance_output_file,"# %s\n",s.str().c_str()); // drawing long line
481  fprintf(balance_output_file, bc_total_format.c_str(),w+wl+2,"total boundary balance",
482  w,total_balance, w, total_outflow, w, total_inflow);
483 
484  //SOURCES
485  string src_head_format = "# %-*s%-*s%-*s%-*s%-*s\n",
486  src_format = "%*s%-*d%-*s %-*g%-*s%-*g\n",
487  src_total_format = "# %-*s%-*g\n\n\n";
488  //computing water balance of sources
489  fprintf(balance_output_file,"# Source fluxes over material subdomains:\n"); //head
490  fprintf(balance_output_file,src_head_format.c_str(),w,"[region_id]",wl,"[label]",
491  w,"[total_balance]",2*w,"",w,"[time]");
492  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //long line
493  std::vector<double> src_balance( mesh_->region_db().bulk_size(), 0.0 ); // initialize by zero
494  FOR_ELEMENTS(mesh_, elm) {
495  src_balance[elm->element_accessor().region().bulk_idx()] += elm->measure() *
496  darcy_flow->data_.cross_section.value(elm->centre(), elm->element_accessor()) *
497  darcy_flow->data_.water_source_density.value(elm->centre(), elm->element_accessor());
498  }
499 
500  total_balance = 0;
501  //printing water balance of sources
502  const RegionSet & bulk_set = mesh_->region_db().get_region_set("BULK");
503  for( RegionSet::const_iterator reg = bulk_set.begin(); reg != bulk_set.end(); ++reg)
504  {
505  total_balance += src_balance[reg->bulk_idx()];
506  //"%*s%-*d%-*s %-*g%-*s%-*g\n";
507  fprintf(balance_output_file, src_format.c_str(), 2,"", w, reg->id(), wl,
508  reg->label().c_str(), w, src_balance[reg->bulk_idx()],2*w,"", w,darcy_flow->time().t());
509  }
510  //total sources balance
511  fprintf(balance_output_file,"# %s\n",s.str().c_str()); //drawing long line
512  fprintf(balance_output_file, src_total_format.c_str(),w+wl+2,"total sources balance",
513  w,total_balance);
514 }
515 
516 
517 /*
518  * Output of internal flow data.
519  */
521 {
523 
524  if (raw_output_file == NULL) return;
525 
526  char dbl_fmt[ 16 ]= "%.8g ";
527  // header
528  xfprintf( raw_output_file, "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n]\n");
529  xfprintf( raw_output_file, "$FlowField\nT=");
530  xfprintf( raw_output_file, dbl_fmt, darcy_flow->time().t());
531  xfprintf( raw_output_file, "\n%d\n", mesh_->n_elements() );
532 
533  unsigned int i;
534  int cit = 0;
535  FOR_ELEMENTS( mesh_, ele ) {
536  xfprintf( raw_output_file, "%d ", ele.id());
537  xfprintf( raw_output_file, dbl_fmt, ele_pressure[cit]);
538  for (i = 0; i < 3; i++)
539  xfprintf( raw_output_file, dbl_fmt, ele_flux[3*cit+i]);
540 
541  xfprintf( raw_output_file, " %d ", ele->n_sides());
542  for (i = 0; i < ele->n_sides(); i++)
543  xfprintf( raw_output_file, dbl_fmt, dh.side_scalar( *(ele->side(i) ) ) );
544  for (i = 0; i < ele->n_sides(); i++)
545  xfprintf( raw_output_file, dbl_fmt, dh.side_flux( *(ele->side(i) ) ) );
546 
547  xfprintf( raw_output_file, "\n" );
548  cit ++;
549  }
550  xfprintf( raw_output_file, "$EndFlowField\n\n" );
551 }
552 
553 
555 #include "fem/fe_p.hh"
556 #include "fem/fe_values.hh"
557 #include "fem/mapping_p1.hh"
558 #include "fields/field_python.hh"
559 #include "fields/field_values.hh"
560 
562 
563 /*
564 * Calculate approximation of L2 norm for:
565  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
566  * 2) difference between RT velocities and analytical solution
567  * 3) difference of divergence
568  * */
569 
570 struct DiffData {
571  double pressure_error[2], velocity_error[2], div_error[2];
576 
577 
578  double * solution;
579  const MH_DofHandler * dh;
581 
582  //std::vector< std::vector<double> > *ele_flux;
586 };
587 
588 template <int dim>
589 void l2_diff_local(ElementFullIter &ele, FEValues<dim,3> &fe_values, ExactSolution &anal_sol, DiffData &result) {
590 
591  fe_values.reinit(ele);
592  result.fe_values.update(ele,
593  result.data_->anisotropy,
594  result.data_->cross_section,
595  result.data_->conductivity);
596  double conductivity = result.data_->conductivity.value(ele->centre(), ele->element_accessor() );
597  double cross = result.data_->cross_section.value(ele->centre(), ele->element_accessor() );
598 
599 
600  // get coefficients on the current element
601  vector<double> fluxes(dim+1);
602  vector<double> pressure_traces(dim+1);
603 
604  for (unsigned int li = 0; li < ele->n_sides(); li++) {
605  fluxes[li] = result.dh->side_flux( *(ele->side( li ) ) );
606  pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
607  }
608  double pressure_mean = result.dh->element_scalar(ele);
609 
610  arma::vec analytical(5);
611  arma::vec3 flux_in_q_point;
612  arma::vec3 anal_flux;
613 
614  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
615 
616  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
617  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
618  double mean_x_squared=0;
619  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
620  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
621  {
622  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
623  * arma::dot( ele->node[i_node]->point(), ele->node[j_node]->point());
624  }
625 
626  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
627  arma::vec3 q_point = fe_values.point(i_point);
628 
629 
630  analytical = anal_sol.value(q_point, ele->element_accessor() );
631  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
632 
633  // compute postprocesed pressure
634  diff = 0;
635  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
636  unsigned int oposite_node = (i_shape + dim) % (dim + 1);
637 
638  diff += fluxes[ i_shape ] *
639  ( arma::dot( q_point, q_point )/ 2
640  - mean_x_squared / 2
641  - arma::dot( q_point, ele->node[oposite_node]->point() )
642  + arma::dot( ele->centre(), ele->node[oposite_node]->point() )
643  );
644  }
645 
646  diff = - (1.0 / conductivity) * diff / dim / ele->measure() / cross + pressure_mean ;
647  diff = ( diff - analytical[0]);
648  pressure_diff += diff * diff * fe_values.JxW(i_point);
649 
650 
651  // velocity difference
652  flux_in_q_point.zeros();
653  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
654  flux_in_q_point += fluxes[ i_shape ]
655  * result.fe_values.RT0_value( ele, q_point, i_shape )
656  / cross;
657  }
658 
659  flux_in_q_point -= anal_flux;
660  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
661 
662  // divergence diff
663  diff = 0;
664  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
665  diff = ( diff / ele->measure() / cross - analytical[4]);
666  divergence_diff += diff * diff * fe_values.JxW(i_point);
667 
668  }
669 
670 
671  result.velocity_diff[ele.index()] = velocity_diff;
672  result.velocity_error[dim-1] += velocity_diff;
673  if (dim == 2) {
674  result.mask_vel_error += (result.velocity_mask[ ele.index() ])? 0 : velocity_diff;
675  }
676 
677  result.pressure_diff[ele.index()] = pressure_diff;
678  result.pressure_error[dim-1] += pressure_diff;
679 
680  result.div_diff[ele.index()] = divergence_diff;
681  result.div_error[dim-1] += divergence_diff;
682 
683 }
684 
685 
686 
687 
688 
689 
691  DBGMSG("l2 norm output\n");
692  ofstream os( FilePath("solution_error", FilePath::output_file) );
693 
694  const unsigned int order = 4; // order of Gauss quadrature
695 
696  // we create trivial Dofhandler , for P0 elements, to get access to, FEValues on individual elements
697  // this we use to integrate our own functions - difference of postprocessed pressure and analytical solution
698  FE_P_disc<0,1,3> fe_1d;
699  FE_P_disc<0,2,3> fe_2d;
700 
701  QGauss<1> quad_1d( order );
702  QGauss<2> quad_2d( order );
703 
704  MappingP1<1,3> mapp_1d;
705  MappingP1<2,3> mapp_2d;
706 
707  FEValues<1,3> fe_values_1d(mapp_1d, quad_1d, fe_1d, update_JxW_values | update_quadrature_points);
708  FEValues<2,3> fe_values_2d(mapp_2d, quad_2d, fe_2d, update_JxW_values | update_quadrature_points);
709 
710  FilePath source_file( "analytical_module.py", FilePath::input_file);
711  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
712  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
713 
714  ExactSolution anal_sol_2d(5);
715  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
716 
717 
718  DiffData result;
719 
720  // mask 2d elements crossing 1d
721  result.velocity_mask.resize(mesh_->n_elements(),0);
722  for(Intersection & isec : mesh_->intersections) {
723  result.velocity_mask[ mesh_->element.index( isec.slave_iter() ) ]++;
724  }
725 
726  result.pressure_diff.resize( mesh_->n_elements() );
727  result.velocity_diff.resize( mesh_->n_elements() );
728  result.div_diff.resize( mesh_->n_elements() );
729 
730  result.pressure_error[0] = 0;
731  result.velocity_error[0] = 0;
732  result.div_error[0] = 0;
733  result.pressure_error[1] = 0;
734  result.velocity_error[1] = 0;
735  result.div_error[1] = 0;
736  result.mask_vel_error=0;
737 
738  //result.ele_flux = &( ele_flux );
739 
740  auto vel_diff_ptr = result.velocity_diff.create_field<3, FieldValue<3>::Scalar>(1);
742  auto pressure_diff_ptr = result.pressure_diff.create_field<3, FieldValue<3>::Scalar>(1);
743  output_fields.pressure_diff.set_field(mesh_->region_db().get_region_set("ALL"), pressure_diff_ptr, 0);
744  auto div_diff_ptr = result.div_diff.create_field<3, FieldValue<3>::Scalar>(1);
745  output_fields.div_diff.set_field(mesh_->region_db().get_region_set("ALL"), div_diff_ptr, 0);
746 
748 
749  unsigned int solution_size;
750  darcy_flow->get_solution_vector(result.solution, solution_size);
751 
752  result.dh = &( darcy_flow->get_mh_dofhandler());
753  result.darcy = darcy_flow;
754  result.data_ = &(darcy_flow->data_);
755 
756  FOR_ELEMENTS( mesh_, ele) {
757 
758  switch (ele->dim()) {
759  case 1:
760 
761  l2_diff_local<1>( ele, fe_values_1d, anal_sol_1d, result);
762  break;
763  case 2:
764  l2_diff_local<2>( ele, fe_values_2d, anal_sol_2d, result);
765  break;
766  }
767  }
768 
769  os << "l2 norm output\n\n"
770  << "pressure error 1d: " << sqrt(result.pressure_error[0]) << endl
771  << "pressure error 2d: " << sqrt(result.pressure_error[1]) << endl
772  << "velocity error 1d: " << sqrt(result.velocity_error[0]) << endl
773  << "velocity error 2d: " << sqrt(result.velocity_error[1]) << endl
774  << "masked velocity error 2d: " << sqrt(result.mask_vel_error) <<endl
775  << "div error 1d: " << sqrt(result.div_error[0]) << endl
776  << "div error 2d: " << sqrt(result.div_error[1]);
777 }
778 
FieldSet & data()
Definition: equation.hh:188
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
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
VectorSeqDouble ele_pressure
static auto subdomain(Mesh &mesh) -> IndexField
#define DBGMSG(...)
Definition: global_defs.h:196
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
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
Accessor to input data conforming to declared Array.
Definition: accessors.hh:558
Definition: system.hh:72
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:252
static Input::Type::Record input_type
void add_admissible_field_names(const Input::Array &in_array)
Registers names of output fields that can be written using this stream.
Definition: output_time.cc:165
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:364
Class for declaration of the input of type Bool.
Definition: type_base.hh:332
MappingP1< 1, 3 > map1
VectorSeqDouble velocity_diff
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:89
Selection & copy_values(const Selection &sel)
int xfprintf(FILE *out, const char *fmt,...)
FPRINTF WITH ERROR HANDLING.
Definition: xio.cc:395
void set_time(const TimeStep &time)
Definition: field_set.hh:186
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:361
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:257
Fields computed from the mesh data.
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
VectorSeqDouble pressure_diff
bool is_current(const TimeMark::Type &mask) const
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
const RegionDB & region_db() const
Definition: mesh.h:155
const MH_DofHandler & get_mh_dofhandler()
Field< 3, FieldValue< 3 >::Integer > subdomain
const TimeStep & step(int index=-1) const
Field< 3, FieldValue< 3 >::Scalar > storativity
double t() const
I/O functions with filename storing, able to track current line in opened file. All standard stdio fu...
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
unsigned int n_sides()
Definition: mesh.cc:172
static TimeMarks & marks()
Class for declaration of inputs sequences.
Definition: type_base.hh:239
double div_error[2]
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:75
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:401
static Default optional()
Definition: type_record.hh:102
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:141
FE_P_disc< 1, 2, 3 > fe2
static OutputTime * create_output_stream(const Input::Record &in_rec)
This method delete all object instances of class OutputTime stored in output_streams vector...
Definition: output_time.cc:141
boost::shared_ptr< Balance > balance_
object for calculation and writing the water balance to file.
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
VectorSeqDouble ele_piezo_head
MappingP1< 2, 3 > map2
Definitions of basic Lagrangean finite elements with polynomial shape functions.
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:42
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:327
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:100
static auto region_id(Mesh &mesh) -> IndexField
void resize(unsigned int size)
Create shared pointer and PETSC vector with given size.
const Element * slave_iter() const
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:383
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:312
FILE * balance_output_file
Temporary solution for writing balance into separate file.
Input::Record in_rec_
Accessor to the input record for the DarcyFlow output.
VectorSeqDouble ele_flux
void mark_output_times(const TimeGovernor &tg)
Definition: output_time.cc:182
VectorSeqDouble div_diff
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:32
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
void write_time_frame()
Definition: output_time.cc:215
Field< 3, FieldValue< 3 >::Scalar > conductivity
Field< 3, FieldValue< 3 >::Integer > region_id
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:75
Field< 3, FieldValue< 3 >::Scalar > pressure_diff
const OutputFields & get_output_fields()
void output(OutputTime *stream)
Definition: field_set.cc:143
const Selection & close() const
unsigned int boundary_size() const
Definition: region.cc:245
static FileName output()
Definition: type_base.hh:470
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:346
TimeGovernor const & time()
Definition: equation.hh:143
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:210
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:244
std::vector< int > velocity_mask
vector< Intersection > intersections
Definition: mesh.h:224
#define FOR_BOUNDARIES(_mesh_, i)
Definition: mesh.h:382
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:76
std::shared_ptr< FieldElementwise< spacedim, Value > > create_field(unsigned int n_comp)
Create and return shared pointer to FieldElementwise object.
Definitions of particular quadrature rules on simplices.
FieldCommon & name(const string &name)
Definition: field_common.hh:73
#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:137
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:169
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:67
FieldCommon & flags(FieldFlag::Flags::Mask mask)
double side_scalar(const Side &side) const
temporary replacement for DofHandler accessor, scalar (pressure) on edge of the side ...
Class for representation SI units of Fields.
Definition: unit_si.hh:31
UnitSI & md(int exp=-1)
Method sets value of exponent for m^{-d}, where d is dimension of region.
Definition: unit_si.cc:91
FE_P_disc< 1, 1, 3 > fe1
void l2_diff_local(ElementFullIter &ele, FEValues< dim, 3 > &fe_values, ExactSolution &anal_sol, DiffData &result)
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:44
#define FOR_SIDES(_mesh_, it)
Definition: mesh.h:408
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)
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:203
FILE * xfopen(const std::string &fname, const char *mode)
Definition: xio.cc:246
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:205
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:430