Flow123d  release_2.1.0-84-g6a13a75
darcy_flow_mh_output.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file darcy_flow_mh_output.cc
15  * @ingroup flow
16  * @brief Output class for darcy_flow_mh model.
17  * @author Jan Brezina
18  */
19 
20 #include <vector>
21 #include <iostream>
22 #include <sstream>
23 #include <string>
24 
25 #include <system/global_defs.h>
26 
27 #include "flow/darcy_flow_mh.hh"
30 
31 #include "io/output_time.hh"
32 #include "io/observe.hh"
33 #include "system/system.hh"
34 #include "system/sys_profiler.hh"
35 
36 #include "fields/field_set.hh"
37 #include "fem/dofhandler.hh"
38 #include "fem/fe_values.hh"
39 #include "fem/fe_rt.hh"
41 #include "fields/field_fe.hh"
42 #include "fields/generic_field.hh"
43 
44 #include "mesh/partitioning.hh"
45 
46 #include "coupling/balance.hh"
47 
48 
49 namespace it = Input::Type;
50 
51 
54  DarcyMH::EqData eq_data;
55  output_fields += eq_data;
56  output_fields += output_fields.error_fields_for_output;
57  return output_fields.make_output_type("Flow_Darcy_MH", "");
58 }
59 
61  return it::Record("Output_DarcyMHSpecific", "Specific Darcy flow MH output.")
62  .declare_key("compute_errors", it::Bool(), it::Default("false"),
63  "SPECIAL PURPOSE. Computing errors pro non-compatible coupling.")
65  "Output file with raw data form MH module.")
66  .close();
67 }
68 
69 
72 {
73 
74  *this += field_ele_pressure.name("pressure_p0").units(UnitSI().m());
75  *this += field_node_pressure.name("pressure_p1").units(UnitSI().m());
76  *this += field_ele_piezo_head.name("piezo_head_p0").units(UnitSI().m());
77  *this += field_ele_flux.name("velocity_p0").units(UnitSI().m().s(-1));
78  *this += subdomain.name("subdomain")
81  *this += region_id.name("region_id")
84 
85  error_fields_for_output += pressure_diff.name("pressure_diff").units(UnitSI().m());
86  error_fields_for_output += velocity_diff.name("velocity_diff").units(UnitSI().m().s(-1));
87  error_fields_for_output += div_diff.name("div_diff").units(UnitSI().s(-1));
88 
89 }
90 
91 
93 : darcy_flow(flow),
94  mesh_(&darcy_flow->mesh()),
95  compute_errors_(false)
96 {
97  Input::Record in_rec_output = main_mh_in_rec.val<Input::Record>("output");
98 
99 
100  // we need to add data from the flow equation at this point, not in constructor of OutputFields
102  output_fields.set_mesh(*mesh_);
103 
104  all_element_idx_.resize(mesh_->n_elements());
105  for(unsigned int i=0; i<all_element_idx_.size(); i++) all_element_idx_[i] = i;
106 
107  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
109  auto ele_pressure_ptr=ele_pressure.create_field<3, FieldValue<3>::Scalar>(1);
110  output_fields.field_ele_pressure.set_field(mesh_->region_db().get_region_set("ALL"), ele_pressure_ptr);
111 
112  dh = new DOFHandlerMultiDim(*mesh_);
114  corner_pressure.resize(dh->n_global_dofs());
115  VecCreateSeqWithArray(PETSC_COMM_SELF, 1, dh->n_global_dofs(), &(corner_pressure[0]), &vec_corner_pressure);
116 
117  auto corner_ptr = make_shared< FieldFE<3, FieldValue<3>::Scalar> >();
118  corner_ptr->set_fe_data(dh, &map1, &map2, &map3, &vec_corner_pressure);
119 
120  output_fields.field_node_pressure.set_field(mesh_->region_db().get_region_set("ALL"), corner_ptr);
121  output_fields.field_node_pressure.output_type(OutputTime::NODE_DATA);
122 
124  auto ele_piezo_head_ptr=ele_piezo_head.create_field<3, FieldValue<3>::Scalar>(1);
125  output_fields.field_ele_piezo_head.set_field(mesh_->region_db().get_region_set("ALL"), ele_piezo_head_ptr);
126 
128  auto ele_flux_ptr=ele_flux.create_field<3, FieldValue<3>::VectorFixed>(3);
129  output_fields.field_ele_flux.set_field(mesh_->region_db().get_region_set("ALL"), ele_flux_ptr);
130 
131  output_fields.subdomain = GenericField<3>::subdomain(*mesh_);
132  output_fields.region_id = GenericField<3>::region_id(*mesh_);
133 
134  output_stream = OutputTime::create_output_stream("flow", *mesh_, main_mh_in_rec.val<Input::Record>("output_stream"));
135  //output_stream->add_admissible_field_names(in_rec_output.val<Input::Array>("fields"));
136  //output_stream->mark_output_times(darcy_flow->time());
137  output_fields.initialize(output_stream, in_rec_output, darcy_flow->time() );
138 
139  int rank;
141 
142  auto in_rec_specific = main_mh_in_rec.find<Input::Record>("output_specific");
143  if (in_rec_specific) {
144  in_rec_specific->opt_val("compute_errors", compute_errors_);
145  if (rank == 0) {
146  // optionally open raw output file
147  FilePath raw_output_file_path;
148  if (in_rec_specific->opt_val("raw_flow_output", raw_output_file_path)) {
149  MessageOut() << "Opening raw output: " << raw_output_file_path << "\n";
150  try {
151  raw_output_file_path.open_stream(raw_output_file);
152  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, (*in_rec_specific))
153  }
154 
155  }
156  }
157 }
158 
159 
160 
162 {
163  if (dh) {
164  chkerr(VecDestroy(&vec_corner_pressure));
165  delete dh;
166  }
167 };
168 
169 
170 
171 
172 
173 //=============================================================================
174 // CONVERT SOLUTION, CALCULATE BALANCES, ETC...
175 //=============================================================================
176 
177 
179 {
180  START_TIMER("Darcy fields output");
181 
182  ElementSetRef observed_elements = output_stream->observe()->observed_elements();
183  {
184  START_TIMER("post-process output fields");
185 
187 
191  else
192  make_element_scalar(observed_elements);
193 
196  else
197  make_element_vector(observed_elements);
198 
201  //else
202  // make_node_scalar_param(observed_elements);
203 
204  // Internal output only if both ele_pressure and ele_flux are output.
208 
210  }
211 
212  {
213  START_TIMER("evaluate output fields");
215  }
216 
217  {
218  START_TIMER("write time frame");
219  output_stream->write_time_frame();
220  }
221 
222 
223 }
224 
225 
226 //=============================================================================
227 // FILL TH "SCALAR" FIELD FOR ALL ELEMENTS IN THE MESH
228 //=============================================================================
229 
231 {
232  START_TIMER("DarcyFlowMHOutput::make_element_scalar");
233  unsigned int sol_size;
234  double *sol;
235 
236  darcy_flow->get_solution_vector(sol, sol_size);
237  unsigned int soi = mesh_->n_sides();
238  for(unsigned int i_ele : element_indices) {
239  ElementFullIter ele = mesh_->element(i_ele);
240  ele_pressure[i_ele] = sol[ soi + i_ele];
241  ele_piezo_head[i_ele] = sol[soi + i_ele ] + ele->centre()[Mesh::z_coord];
242  }
243 }
244 
245 
246 
247 /****
248  * compute Darcian velocity in centre of elements
249  *
250  */
252  START_TIMER("DarcyFlowMHOutput::make_element_vector");
253  // need to call this to create mh solution vector
255 
256  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(darcy_flow->data_);
257  arma::vec3 flux_in_center;
258  for(unsigned int i_ele : element_indices) {
259  ElementFullIter ele = mesh_->element(i_ele);
260 
261  flux_in_center = multidim_assembler[ele->dim() -1]->make_element_vector(ele);
262 
263  // place it in the sequential vector
264  for(unsigned int j=0; j<3; j++) ele_flux[3*i_ele + j]=flux_in_center[j];
265  }
266 }
267 
268 
270 {
271  START_TIMER("DarcyFlowMHOutput::make_corner_scalar");
272  unsigned int ndofs = max(dh->fe<1>()->n_dofs(), max(dh->fe<2>()->n_dofs(), dh->fe<3>()->n_dofs()));
273  unsigned int indices[ndofs];
274  unsigned int i_node;
275  FOR_ELEMENTS(mesh_, ele)
276  {
277  dh->get_dof_indices(ele, indices);
278  FOR_ELEMENT_NODES(ele, i_node)
279  {
280  corner_pressure[indices[i_node]] = node_scalar[mesh_->node_vector.index(ele->node[i_node])];
281  }
282  }
283 }
284 
285 
286 //=============================================================================
287 // pressure interpolation
288 //
289 // some sort of weighted average over elements and sides/edges of an node
290 //
291 // TODO:
292 // questions:
293 // - explain details of averaging and motivation for this type
294 // - edge scalars are stronger since thay are computed twice by each side
295 // - why division by (nod.aux-1)
296 // - ?implement without TED, TSD global arrays with single loop over elements, sides and one loop over nodes for normalization
297 //
298 //=============================================================================
299 
301  START_TIMER("DarcyFlowMHOutput::make_node_scalar_param");
302 
303  vector<double> scalars(mesh_->n_nodes());
304 
305  double dist; //!< tmp variable for storing particular distance node --> element, node --> side*/
306 
307  /** Iterators */
308  const Node * node;
309 
310  int n_nodes = mesh_->node_vector.size(); //!< number of nodes in the mesh */
311  int node_index = 0; //!< index of each node */
312 
313  int* sum_elements = new int [n_nodes]; //!< sum elements joined to node */
314  int* sum_sides = new int [n_nodes]; //!< sum sides joined to node */
315  double* sum_ele_dist = new double [n_nodes]; //!< sum distances to all joined elements */
316  double* sum_side_dist = new double [n_nodes]; //!< Sum distances to all joined sides */
317 
318  /** tmp variables, will be replaced by ini keys
319  * TODO include them into ini file*/
320  bool count_elements = true; //!< scalar is counted via elements*/
321  bool count_sides = true; //!< scalar is counted via sides */
322 
323 
325 
326  /** init arrays */
327  for (int i = 0; i < n_nodes; i++){
328  sum_elements[i] = 0;
329  sum_sides[i] = 0;
330  sum_ele_dist[i] = 0.0;
331  sum_side_dist[i] = 0.0;
332  scalars[i] = 0.0;
333  };
334 
335  /**first pass - calculate sums (weights)*/
336  if (count_elements){
337  FOR_ELEMENTS(mesh_, ele)
338  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
339  node = ele->node[li]; //!< get Node pointer from element */
340  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
341 
342  dist = sqrt(
343  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
344  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
345  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
346  );
347  sum_ele_dist[node_index] += dist;
348  sum_elements[node_index]++;
349  }
350  }
351  if (count_sides){
352  FOR_SIDES(mesh_, side) {
353  for (unsigned int li = 0; li < side->n_nodes(); li++) {
354  node = side->node(li);//!< get Node pointer from element */
355  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
356  dist = sqrt(
357  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
358  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
359  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
360  );
361 
362  sum_side_dist[node_index] += dist;
363  sum_sides[node_index]++;
364  }
365  }
366  }
367 
368  /**second pass - calculate scalar */
369  if (count_elements){
370  FOR_ELEMENTS(mesh_, ele)
371  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
372  node = ele->node[li];//!< get Node pointer from element */
373  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
374 
375  /**TODO - calculate it again or store it in prior pass*/
376  dist = sqrt(
377  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
378  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
379  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
380  );
381  scalars[node_index] += ele_pressure[ele.index()] *
382  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
383  (sum_elements[node_index] + sum_sides[node_index] - 1);
384  }
385  }
386  if (count_sides) {
387  FOR_SIDES(mesh_, side) {
388  for (unsigned int li = 0; li < side->n_nodes(); li++) {
389  node = side->node(li);//!< get Node pointer from element */
390  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
391 
392  /**TODO - calculate it again or store it in prior pass*/
393  dist = sqrt(
394  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
395  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
396  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
397  );
398 
399 
400  scalars[node_index] += dh.side_scalar( *side ) *
401  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
402  (sum_sides[node_index] + sum_elements[node_index] - 1);
403  }
404  }
405  }
406 
407  /** free memory */
408  delete [] sum_elements;
409  delete [] sum_sides;
410  delete [] sum_ele_dist;
411  delete [] sum_side_dist;
412 
413  make_corner_scalar(scalars);
414 }
415 
416 
417 /*
418  * Output of internal flow data.
419  */
421 {
422  START_TIMER("DarcyFlowMHOutput::output_internal_flow_data");
424 
425  if (! raw_output_file.is_open()) return;
426 
427  //char dbl_fmt[ 16 ]= "%.8g ";
428  // header
429  raw_output_file << "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n]\n";
430  raw_output_file << fmt::format("$FlowField\nT={}\n", darcy_flow->time().t());
431  raw_output_file << fmt::format("{}\n" , mesh_->n_elements() );
432 
433  ;
434  int cit = 0;
435  FOR_ELEMENTS( mesh_, ele ) {
436  raw_output_file << fmt::format("{} {} ", ele.id(), ele_pressure[cit]);
437  for (unsigned int i = 0; i < 3; i++)
438  raw_output_file << ele_flux[3*cit+i] << " ";
439 
440  raw_output_file << ele->n_sides() << " ";
441  std::vector< std::vector<unsigned int > > old_to_new_side =
442  { {0, 1},
443  {0, 1, 2},
444  {0, 1, 2, 3}
445  };
446  for (unsigned int i = 0; i < ele->n_sides(); i++) {
447  unsigned int i_new_side = old_to_new_side[ele->dim()-1][i];
448  raw_output_file << dh.side_scalar( *(ele->side(i_new_side) ) ) << " ";
449  }
450  for (unsigned int i = 0; i < ele->n_sides(); i++) {
451  unsigned int i_new_side = old_to_new_side[ele->dim()-1][i];
452  raw_output_file << dh.side_flux( *(ele->side(i_new_side) ) ) << " ";
453  }
454 
455  raw_output_file << endl;
456  cit ++;
457  }
458  raw_output_file << "$EndFlowField\n" << endl;
459 }
460 
461 
463 #include "fem/fe_p.hh"
464 #include "fem/fe_values.hh"
465 #include "fem/mapping_p1.hh"
466 #include "fields/field_python.hh"
467 #include "fields/field_values.hh"
468 
470 
471 /*
472 * Calculate approximation of L2 norm for:
473  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
474  * 2) difference between RT velocities and analytical solution
475  * 3) difference of divergence
476  * */
477 
478 struct DiffData {
479  double pressure_error[2], velocity_error[2], div_error[2];
484 
485 
486  double * solution;
487  const MH_DofHandler * dh;
488 
489  //std::vector< std::vector<double> > *ele_flux;
493 };
494 
495 template <int dim>
497  FEValues<dim,3> &fe_values, FEValues<dim,3> &fv_rt,
498  ExactSolution &anal_sol, DiffData &result) {
499 
500  fv_rt.reinit(ele);
501  fe_values.reinit(ele);
502 
503  double conductivity = result.data_->conductivity.value(ele->centre(), ele->element_accessor() );
504  double cross = result.data_->cross_section.value(ele->centre(), ele->element_accessor() );
505 
506 
507  // get coefficients on the current element
508  vector<double> fluxes(dim+1);
509 // vector<double> pressure_traces(dim+1);
510 
511  for (unsigned int li = 0; li < ele->n_sides(); li++) {
512  fluxes[li] = result.dh->side_flux( *(ele->side( li ) ) );
513 // pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
514  }
515  double pressure_mean = result.dh->element_scalar(ele);
516 
517  arma::vec analytical(5);
518  arma::vec3 flux_in_q_point;
519  arma::vec3 anal_flux;
520 
521  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
522 
523  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
524  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
525  double mean_x_squared=0;
526  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
527  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
528  {
529  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
530  * arma::dot( ele->node[i_node]->point(), ele->node[j_node]->point());
531  }
532 
533  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
534  arma::vec3 q_point = fe_values.point(i_point);
535 
536 
537  analytical = anal_sol.value(q_point, ele->element_accessor() );
538  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
539 
540  // compute postprocesed pressure
541  diff = 0;
542  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
543  unsigned int oposite_node = RefElement<dim>::oposite_node(i_shape);
544 
545  diff += fluxes[ i_shape ] *
546  ( arma::dot( q_point, q_point )/ 2
547  - mean_x_squared / 2
548  - arma::dot( q_point, ele->node[oposite_node]->point() )
549  + arma::dot( ele->centre(), ele->node[oposite_node]->point() )
550  );
551  }
552 
553  diff = - (1.0 / conductivity) * diff / dim / ele->measure() / cross + pressure_mean ;
554  diff = ( diff - analytical[0]);
555  pressure_diff += diff * diff * fe_values.JxW(i_point);
556 
557 
558  // velocity difference
559  flux_in_q_point.zeros();
560  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
561  flux_in_q_point += fluxes[ i_shape ]
562  * fv_rt.shape_vector(i_shape, i_point)
563  / cross;
564  }
565 
566  flux_in_q_point -= anal_flux;
567  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
568 
569  // divergence diff
570  diff = 0;
571  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
572  diff = ( diff / ele->measure() / cross - analytical[4]);
573  divergence_diff += diff * diff * fe_values.JxW(i_point);
574 
575  }
576 
577 
578  result.velocity_diff[ele.index()] = velocity_diff;
579  result.velocity_error[dim-1] += velocity_diff;
580  if (dim == 2) {
581  result.mask_vel_error += (result.velocity_mask[ ele.index() ])? 0 : velocity_diff;
582  }
583 
584  result.pressure_diff[ele.index()] = pressure_diff;
585  result.pressure_error[dim-1] += pressure_diff;
586 
587  result.div_diff[ele.index()] = divergence_diff;
588  result.div_error[dim-1] += divergence_diff;
589 
590 }
591 
592 
593 
594 
595 
596 
598  DebugOut() << "l2 norm output\n";
599  ofstream os( FilePath("solution_error", FilePath::output_file) );
600 
601  const unsigned int order = 4; // order of Gauss quadrature
602 
603  // we create trivial Dofhandler , for P0 elements, to get access to, FEValues on individual elements
604  // this we use to integrate our own functions - difference of postprocessed pressure and analytical solution
605  FE_P_disc<0,1,3> fe_1d;
606  FE_P_disc<0,2,3> fe_2d;
607 
608  QGauss<1> quad_1d( order );
609  QGauss<2> quad_2d( order );
610 
611  MappingP1<1,3> mapp_1d;
612  MappingP1<2,3> mapp_2d;
613 
614  FEValues<1,3> fe_values_1d(mapp_1d, quad_1d, fe_1d, update_JxW_values | update_quadrature_points);
615  FEValues<2,3> fe_values_2d(mapp_2d, quad_2d, fe_2d, update_JxW_values | update_quadrature_points);
616 
617  // FEValues for velocity.
618  FE_RT0<1,3> fe_rt1d;
619  FE_RT0<2,3> fe_rt2d;
620  FEValues<1,3> fv_rt1d(mapp_1d,quad_1d, fe_rt1d, update_values | update_quadrature_points);
621  FEValues<2,3> fv_rt2d(mapp_2d,quad_2d, fe_rt2d, update_values | update_quadrature_points);
622 
623  FilePath source_file( "analytical_module.py", FilePath::input_file);
624  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
625  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
626 
627  ExactSolution anal_sol_2d(5);
628  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
629 
630 
631  DiffData result;
632 
633  // mask 2d elements crossing 1d
634  result.velocity_mask.resize(mesh_->n_elements(),0);
635  for(Intersection & isec : mesh_->intersections) {
636  result.velocity_mask[ mesh_->element.index( isec.slave_iter() ) ]++;
637  }
638 
639  result.pressure_diff.resize( mesh_->n_elements() );
640  result.velocity_diff.resize( mesh_->n_elements() );
641  result.div_diff.resize( mesh_->n_elements() );
642 
643  result.pressure_error[0] = 0;
644  result.velocity_error[0] = 0;
645  result.div_error[0] = 0;
646  result.pressure_error[1] = 0;
647  result.velocity_error[1] = 0;
648  result.div_error[1] = 0;
649  result.mask_vel_error=0;
650 
651  //result.ele_flux = &( ele_flux );
652 
654 
655  auto vel_diff_ptr = result.velocity_diff.create_field<3, FieldValue<3>::Scalar>(1);
657  auto pressure_diff_ptr = result.pressure_diff.create_field<3, FieldValue<3>::Scalar>(1);
658  output_fields.pressure_diff.set_field(mesh_->region_db().get_region_set("ALL"), pressure_diff_ptr, 0);
659  auto div_diff_ptr = result.div_diff.create_field<3, FieldValue<3>::Scalar>(1);
660  output_fields.div_diff.set_field(mesh_->region_db().get_region_set("ALL"), div_diff_ptr, 0);
661 
664 
665 
666  unsigned int solution_size;
667  darcy_flow->get_solution_vector(result.solution, solution_size);
668 
669  result.dh = &( darcy_flow->get_mh_dofhandler());
670  result.darcy = darcy_flow;
671  result.data_ = darcy_flow->data_.get();
672 
673  FOR_ELEMENTS( mesh_, ele) {
674 
675  switch (ele->dim()) {
676  case 1:
677 
678  l2_diff_local<1>( ele, fe_values_1d, fv_rt1d, anal_sol_1d, result);
679  break;
680  case 2:
681  l2_diff_local<2>( ele, fe_values_2d, fv_rt2d, anal_sol_2d, result);
682  break;
683  }
684  }
685 
686  os << "l2 norm output\n\n"
687  << "pressure error 1d: " << sqrt(result.pressure_error[0]) << endl
688  << "pressure error 2d: " << sqrt(result.pressure_error[1]) << endl
689  << "velocity error 1d: " << sqrt(result.velocity_error[0]) << endl
690  << "velocity error 2d: " << sqrt(result.velocity_error[1]) << endl
691  << "masked velocity error 2d: " << sqrt(result.mask_vel_error) <<endl
692  << "div error 1d: " << sqrt(result.div_error[0]) << endl
693  << "div error 2d: " << sqrt(result.div_error[1]);
694 }
695 
TimeGovernor & time()
Definition: equation.hh:148
FieldSet & data()
Definition: equation.hh:193
const MH_DofHandler * dh
void set_python_field_from_file(const FilePath &file_name, const string &func_name)
void get_solution_vector(double *&vec, unsigned int &vec_size) override
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
Shape function values.
Definition: update_flags.hh:87
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.hh:292
VectorSeqDouble ele_pressure
static auto subdomain(Mesh &mesh) -> IndexField
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
Output class for darcy_flow_mh model.
#define FOR_ELEMENT_NODES(i, j)
Definition: elements.h:188
double getY() const
Definition: nodes.hh:59
Definition: nodes.hh:32
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:56
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:405
Class for declaration of the input of type Bool.
Definition: type_base.hh:452
MappingP1< 1, 3 > map1
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:231
VectorSeqDouble velocity_diff
void output(TimeStep step)
void make_element_vector(ElementSetRef element_indices)
std::string format(CStringRef format_str, ArgList args)
Definition: format.h:3141
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:57
DOFHandlerMultiDim * dh
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
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:245
Fields computed from the mesh data.
Iterator< Ret > find(const string &key) const
void make_corner_scalar(vector< double > &node_scalar)
Helper class that stores data of generic types.
Definition: type_generic.hh:88
Class FEValues calculates finite element data on the actual cells such as shape function values...
int index() const
Definition: sys_vector.hh:78
VectorSeqDouble pressure_diff
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:142
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
const RegionDB & region_db() const
Definition: mesh.h:147
Field< 3, FieldValue< 3 >::Integer > subdomain
const TimeStep & step(int index=-1) const
double t() const
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
std::vector< unsigned int > all_element_idx_
unsigned int n_sides()
Definition: mesh.cc:170
std::shared_ptr< OutputTime > output_stream
std::shared_ptr< EqData > data_
unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.hh:330
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:301
double div_error[2]
FiniteElement< dim, 3 > * fe() const
Returns finite element object for given space dimension.
void make_node_scalar_param(ElementSetRef element_indices)
Calculate nodes scalar, store it in double* node_scalars instead of node->scalar. ...
static unsigned int oposite_node(unsigned int sid)
Definition: ref_element.hh:251
Symmetric Gauss-Legendre quadrature formulae on simplices.
double getZ() const
Definition: nodes.hh:61
unsigned int size() const
Returns size of the container. This is independent of the allocated space.
Definition: sys_vector.hh:391
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:119
void open_stream(Stream &stream) const
Definition: file_path.cc:211
bool opt_val(const string &key, Ret &value) const
Discontinuous Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:202
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:133
FE_P_disc< 1, 2, 3 > fe2
double element_scalar(ElementFullIter &ele) const
temporary replacement for DofHandler accessor, scalar (pressure) on element
double pressure_error[2]
VectorSeqDouble ele_piezo_head
MappingP1< 2, 3 > map2
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:33
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:58
Accessor to the data with type Type::Record.
Definition: accessors.hh:286
const Ret val(const string &key) const
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:57
#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:248
double velocity_error[2]
unsigned int index(const T *pointer) const
Definition: sys_vector.hh:373
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:350
static const Input::Type::Instance & get_input_type()
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:488
VectorSeqDouble ele_flux
ofstream raw_output_file
Raw data output file.
static const Input::Type::Record & get_input_type_specific()
arma::vec::fixed< spacedim > shape_vector(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.hh:254
VectorSeqDouble div_diff
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:48
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
static std::shared_ptr< OutputTime > create_output_stream(const std::string &equation_name, Mesh &mesh, const Input::Record &in_rec)
This method delete all object instances of class OutputTime stored in output_streams vector...
Definition: output_time.cc:191
Field< 3, FieldValue< 3 >::Integer > region_id
double side_flux(const Side &side) const
temporary replacement for DofHandler accessor, flux through given side
Field< 3, FieldValue< 3 >::Scalar > pressure_diff
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:334
bool compute_errors_
Specific experimental error computing.
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:196
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:157
std::vector< int > velocity_mask
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:149
vector< Intersection > intersections
Definition: mesh.h:235
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:86
#define MPI_COMM_WORLD
Definition: mpi.h:123
Field< 3, FieldValue< 3 >::Scalar > div_diff
FE_P_disc< 1, 3, 3 > fe3
unsigned int n_nodes() const
Definition: mesh.h:129
vector< double > corner_pressure
mixed-hybrid model of linear Darcy flow, possibly unsteady.
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:173
Record type proxy class.
Definition: type_record.hh:177
DarcyMH::EqData * data_
bool is_field_output_time(const FieldCommon &field, TimeStep step) const
const Input::Type::Instance & make_output_type(const string &equation_name, const string &aditional_description="")
Field< 3, FieldValue< 3 >::Scalar > cross_section
Field< 3, FieldValue< 3 >::Scalar > conductivity
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:40
DarcyFlowMHOutput(DarcyMH *flow, Input::Record in_rec)
FE_P_disc< 1, 1, 3 > fe1
static FileName output()
The factory function for declaring type FileName for input files.
Definition: type_base.cc:508
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
#define FOR_SIDES(_mesh_, it)
Definition: mesh.h:449
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:240
void make_element_scalar(ElementSetRef element_indices)
void l2_diff_local(ElementFullIter &ele, FEValues< dim, 3 > &fe_values, FEValues< dim, 3 > &fv_rt, ExactSolution &anal_sol, DiffData &result)
Definitions of Raviart-Thomas finite elements.
const unsigned int n_global_dofs() const
Getter for the number of all mesh dofs required by the given finite element.
Definition: dofhandler.hh:61
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:214
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:216
void output()
Calculate values for output.
Transformed quadrature weights.
const MH_DofHandler & get_mh_dofhandler() override
Mixed-hybrid of steady Darcy flow with sources and variable density.
arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
Definition: fe_values.hh:302