Flow123d  intersections_paper-476-gbe68821
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"
49 
50 namespace it = Input::Type;
51 
52 
55  DarcyMH::EqData eq_data;
56  output_fields += eq_data;
57  output_fields += output_fields.error_fields_for_output;
58  return output_fields.make_output_type("Flow_Darcy_MH", "");
59 }
60 
62  return it::Record("Output_DarcyMHSpecific", "Specific Darcy flow MH output.")
63  .declare_key("compute_errors", it::Bool(), it::Default("false"),
64  "SPECIAL PURPOSE. Computing errors pro non-compatible coupling.")
66  "Output file with raw data form MH module.")
67  .close();
68 }
69 
70 
73 {
74 
75  *this += field_ele_pressure.name("pressure_p0").units(UnitSI().m());
76  *this += field_node_pressure.name("pressure_p1").units(UnitSI().m());
77  *this += field_ele_piezo_head.name("piezo_head_p0").units(UnitSI().m());
78  *this += field_ele_flux.name("velocity_p0").units(UnitSI().m().s(-1));
79  *this += subdomain.name("subdomain")
82  *this += region_id.name("region_id")
85 
86  error_fields_for_output += pressure_diff.name("pressure_diff").units(UnitSI().m());
87  error_fields_for_output += velocity_diff.name("velocity_diff").units(UnitSI().m().s(-1));
88  error_fields_for_output += div_diff.name("div_diff").units(UnitSI().s(-1));
89 
90 }
91 
92 
94 : darcy_flow(flow),
95  mesh_(&darcy_flow->mesh()),
96  compute_errors_(false)
97 {
98  Input::Record in_rec_output = main_mh_in_rec.val<Input::Record>("output");
99 
100 
101  // we need to add data from the flow equation at this point, not in constructor of OutputFields
103  output_fields.set_mesh(*mesh_);
104 
105  all_element_idx_.resize(mesh_->n_elements());
106  for(unsigned int i=0; i<all_element_idx_.size(); i++) all_element_idx_[i] = i;
107 
108  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
110  auto ele_pressure_ptr=ele_pressure.create_field<3, FieldValue<3>::Scalar>(1);
111  output_fields.field_ele_pressure.set_field(mesh_->region_db().get_region_set("ALL"), ele_pressure_ptr);
112 
113  dh_ = make_shared<DOFHandlerMultiDim>(*mesh_);
114  dh_->distribute_dofs(fe1, fe2, fe3);
115  corner_pressure.resize(dh_->n_global_dofs());
116 
117  auto corner_ptr = make_shared< FieldFE<3, FieldValue<3>::Scalar> >();
118  corner_ptr->set_fe_data(dh_, &map1, &map2, &map3, &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 
123  ele_piezo_head.resize(mesh_->n_elements());
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 
127  ele_flux.resize(3*mesh_->n_elements());
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 
164 
165 
166 
167 
168 //=============================================================================
169 // CONVERT SOLUTION, CALCULATE BALANCES, ETC...
170 //=============================================================================
171 
172 
174 {
175  START_TIMER("Darcy fields output");
176 
177  ElementSetRef observed_elements = output_stream->observe()->observed_elements();
178  {
179  START_TIMER("post-process output fields");
180 
182 
186  else
187  make_element_scalar(observed_elements);
188 
191  else
192  make_element_vector(observed_elements);
193 
196  //else
197  // make_node_scalar_param(observed_elements);
198 
199  // Internal output only if both ele_pressure and ele_flux are output.
203 
205  }
206 
207  {
208  START_TIMER("evaluate output fields");
210  }
211 
212  {
213  START_TIMER("write time frame");
214  output_stream->write_time_frame();
215  }
216 
217 
218 }
219 
220 
221 //=============================================================================
222 // FILL TH "SCALAR" FIELD FOR ALL ELEMENTS IN THE MESH
223 //=============================================================================
224 
226 {
227  START_TIMER("DarcyFlowMHOutput::make_element_scalar");
228  unsigned int sol_size;
229  double *sol;
230 
231  darcy_flow->get_solution_vector(sol, sol_size);
232  unsigned int soi = mesh_->n_sides();
233  for(unsigned int i_ele : element_indices) {
234  ElementFullIter ele = mesh_->element(i_ele);
235  ele_pressure[i_ele] = sol[ soi + i_ele];
236  ele_piezo_head[i_ele] = sol[soi + i_ele ]
237  - (darcy_flow->data_->gravity_[3] + arma::dot(darcy_flow->data_->gravity_vec_,ele->centre()));
238  }
239 }
240 
241 
242 
243 /****
244  * compute Darcian velocity in centre of elements
245  *
246  */
248  START_TIMER("DarcyFlowMHOutput::make_element_vector");
249  // need to call this to create mh solution vector
251 
252  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(darcy_flow->data_);
253  arma::vec3 flux_in_center;
254  for(unsigned int i_ele : element_indices) {
255  ElementFullIter ele = mesh_->element(i_ele);
256 
257  flux_in_center = multidim_assembler[ele->dim() -1]->make_element_vector(ele);
258 
259  // place it in the sequential vector
260  for(unsigned int j=0; j<3; j++) ele_flux[3*i_ele + j]=flux_in_center[j];
261  }
262 }
263 
264 
266 {
267  START_TIMER("DarcyFlowMHOutput::make_corner_scalar");
268  unsigned int ndofs = max(dh_->fe<1>()->n_dofs(), max(dh_->fe<2>()->n_dofs(), dh_->fe<3>()->n_dofs()));
269  unsigned int indices[ndofs];
270  unsigned int i_node;
271  FOR_ELEMENTS(mesh_, ele)
272  {
273  dh_->get_dof_indices(ele, indices);
274  FOR_ELEMENT_NODES(ele, i_node)
275  {
276  corner_pressure[indices[i_node]] = node_scalar[mesh_->node_vector.index(ele->node[i_node])];
277  }
278  }
279 }
280 
281 
282 //=============================================================================
283 // pressure interpolation
284 //
285 // some sort of weighted average over elements and sides/edges of an node
286 //
287 // TODO:
288 // questions:
289 // - explain details of averaging and motivation for this type
290 // - edge scalars are stronger since thay are computed twice by each side
291 // - why division by (nod.aux-1)
292 // - ?implement without TED, TSD global arrays with single loop over elements, sides and one loop over nodes for normalization
293 //
294 //=============================================================================
295 
297  START_TIMER("DarcyFlowMHOutput::make_node_scalar_param");
298 
299  vector<double> scalars(mesh_->n_nodes());
300 
301  double dist; //!< tmp variable for storing particular distance node --> element, node --> side*/
302 
303  /** Iterators */
304  const Node * node;
305 
306  int n_nodes = mesh_->node_vector.size(); //!< number of nodes in the mesh */
307  int node_index = 0; //!< index of each node */
308 
309  int* sum_elements = new int [n_nodes]; //!< sum elements joined to node */
310  int* sum_sides = new int [n_nodes]; //!< sum sides joined to node */
311  double* sum_ele_dist = new double [n_nodes]; //!< sum distances to all joined elements */
312  double* sum_side_dist = new double [n_nodes]; //!< Sum distances to all joined sides */
313 
314  /** tmp variables, will be replaced by ini keys
315  * TODO include them into ini file*/
316  bool count_elements = true; //!< scalar is counted via elements*/
317  bool count_sides = true; //!< scalar is counted via sides */
318 
319 
321 
322  /** init arrays */
323  for (int i = 0; i < n_nodes; i++){
324  sum_elements[i] = 0;
325  sum_sides[i] = 0;
326  sum_ele_dist[i] = 0.0;
327  sum_side_dist[i] = 0.0;
328  scalars[i] = 0.0;
329  };
330 
331  /**first pass - calculate sums (weights)*/
332  if (count_elements){
333  FOR_ELEMENTS(mesh_, ele)
334  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
335  node = ele->node[li]; //!< get Node pointer from element */
336  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
337 
338  dist = sqrt(
339  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
340  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
341  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
342  );
343  sum_ele_dist[node_index] += dist;
344  sum_elements[node_index]++;
345  }
346  }
347  if (count_sides){
348  FOR_SIDES(mesh_, side) {
349  for (unsigned int li = 0; li < side->n_nodes(); li++) {
350  node = side->node(li);//!< get Node pointer from element */
351  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
352  dist = sqrt(
353  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
354  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
355  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
356  );
357 
358  sum_side_dist[node_index] += dist;
359  sum_sides[node_index]++;
360  }
361  }
362  }
363 
364  /**second pass - calculate scalar */
365  if (count_elements){
366  FOR_ELEMENTS(mesh_, ele)
367  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
368  node = ele->node[li];//!< get Node pointer from element */
369  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
370 
371  /**TODO - calculate it again or store it in prior pass*/
372  dist = sqrt(
373  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
374  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
375  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
376  );
377  scalars[node_index] += ele_pressure[ele.index()] *
378  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
379  (sum_elements[node_index] + sum_sides[node_index] - 1);
380  }
381  }
382  if (count_sides) {
383  FOR_SIDES(mesh_, side) {
384  for (unsigned int li = 0; li < side->n_nodes(); li++) {
385  node = side->node(li);//!< get Node pointer from element */
386  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
387 
388  /**TODO - calculate it again or store it in prior pass*/
389  dist = sqrt(
390  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
391  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
392  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
393  );
394 
395 
396  scalars[node_index] += dh.side_scalar( *side ) *
397  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
398  (sum_sides[node_index] + sum_elements[node_index] - 1);
399  }
400  }
401  }
402 
403  /** free memory */
404  delete [] sum_elements;
405  delete [] sum_sides;
406  delete [] sum_ele_dist;
407  delete [] sum_side_dist;
408 
409  make_corner_scalar(scalars);
410 }
411 
412 
413 /*
414  * Output of internal flow data.
415  */
417 {
418  START_TIMER("DarcyFlowMHOutput::output_internal_flow_data");
420 
421  if (! raw_output_file.is_open()) return;
422 
423  //char dbl_fmt[ 16 ]= "%.8g ";
424  // header
425  raw_output_file << "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n]\n";
426  raw_output_file << fmt::format("$FlowField\nT={}\n", darcy_flow->time().t());
427  raw_output_file << fmt::format("{}\n" , mesh_->n_elements() );
428 
429  ;
430  int cit = 0;
431  FOR_ELEMENTS( mesh_, ele ) {
432  raw_output_file << fmt::format("{} {} ", ele.id(), ele_pressure[cit]);
433  for (unsigned int i = 0; i < 3; i++)
434  raw_output_file << ele_flux[3*cit+i] << " ";
435 
436  raw_output_file << ele->n_sides() << " ";
437  std::vector< std::vector<unsigned int > > old_to_new_side =
438  { {0, 1},
439  {0, 1, 2},
440  {0, 1, 2, 3}
441  };
442  for (unsigned int i = 0; i < ele->n_sides(); i++) {
443  unsigned int i_new_side = old_to_new_side[ele->dim()-1][i];
444  raw_output_file << dh.side_scalar( *(ele->side(i_new_side) ) ) << " ";
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_flux( *(ele->side(i_new_side) ) ) << " ";
449  }
450 
451  raw_output_file << endl;
452  cit ++;
453  }
454  raw_output_file << "$EndFlowField\n" << endl;
455 }
456 
457 
459 #include "fem/fe_p.hh"
460 #include "fem/fe_values.hh"
461 #include "fem/mapping_p1.hh"
462 #include "fields/field_python.hh"
463 #include "fields/field_values.hh"
464 
466 
467 /*
468 * Calculate approximation of L2 norm for:
469  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
470  * 2) difference between RT velocities and analytical solution
471  * 3) difference of divergence
472  * */
473 
474 struct DiffData {
475  double pressure_error[2], velocity_error[2], div_error[2];
480 
481 
482  double * solution;
483  const MH_DofHandler * dh;
484 
485  //std::vector< std::vector<double> > *ele_flux;
489 };
490 
491 template <int dim>
493  FEValues<dim,3> &fe_values, FEValues<dim,3> &fv_rt,
494  ExactSolution &anal_sol, DiffData &result) {
495 
496  fv_rt.reinit(ele);
497  fe_values.reinit(ele);
498 
499  double conductivity = result.data_->conductivity.value(ele->centre(), ele->element_accessor() );
500  double cross = result.data_->cross_section.value(ele->centre(), ele->element_accessor() );
501 
502 
503  // get coefficients on the current element
504  vector<double> fluxes(dim+1);
505 // vector<double> pressure_traces(dim+1);
506 
507  for (unsigned int li = 0; li < ele->n_sides(); li++) {
508  fluxes[li] = result.dh->side_flux( *(ele->side( li ) ) );
509 // pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
510  }
511  double pressure_mean = result.dh->element_scalar(ele);
512 
513  arma::vec analytical(5);
514  arma::vec3 flux_in_q_point;
515  arma::vec3 anal_flux;
516 
517  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
518 
519  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
520  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
521  double mean_x_squared=0;
522  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
523  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
524  {
525  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
526  * arma::dot( ele->node[i_node]->point(), ele->node[j_node]->point());
527  }
528 
529  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
530  arma::vec3 q_point = fe_values.point(i_point);
531 
532 
533  analytical = anal_sol.value(q_point, ele->element_accessor() );
534  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
535 
536  // compute postprocesed pressure
537  diff = 0;
538  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
539  unsigned int oposite_node = RefElement<dim>::oposite_node(i_shape);
540 
541  diff += fluxes[ i_shape ] *
542  ( arma::dot( q_point, q_point )/ 2
543  - mean_x_squared / 2
544  - arma::dot( q_point, ele->node[oposite_node]->point() )
545  + arma::dot( ele->centre(), ele->node[oposite_node]->point() )
546  );
547  }
548 
549  diff = - (1.0 / conductivity) * diff / dim / ele->measure() / cross + pressure_mean ;
550  diff = ( diff - analytical[0]);
551  pressure_diff += diff * diff * fe_values.JxW(i_point);
552 
553 
554  // velocity difference
555  flux_in_q_point.zeros();
556  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
557  flux_in_q_point += fluxes[ i_shape ]
558  * fv_rt.shape_vector(i_shape, i_point)
559  / cross;
560  }
561 
562  flux_in_q_point -= anal_flux;
563  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
564 
565  // divergence diff
566  diff = 0;
567  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
568  diff = ( diff / ele->measure() / cross - analytical[4]);
569  divergence_diff += diff * diff * fe_values.JxW(i_point);
570 
571  }
572 
573 
574  result.velocity_diff[ele.index()] = velocity_diff;
575  result.velocity_error[dim-1] += velocity_diff;
576  if (dim == 2 && result.velocity_mask.size() != 0 ) {
577  result.mask_vel_error += (result.velocity_mask[ ele.index() ])? 0 : velocity_diff;
578  }
579 
580  result.pressure_diff[ele.index()] = pressure_diff;
581  result.pressure_error[dim-1] += pressure_diff;
582 
583  result.div_diff[ele.index()] = divergence_diff;
584  result.div_error[dim-1] += divergence_diff;
585 
586 }
587 
588 
589 
590 
591 
592 
594  DebugOut() << "l2 norm output\n";
595  ofstream os( FilePath("solution_error", FilePath::output_file) );
596 
597  const unsigned int order = 4; // order of Gauss quadrature
598 
599  // we create trivial Dofhandler , for P0 elements, to get access to, FEValues on individual elements
600  // this we use to integrate our own functions - difference of postprocessed pressure and analytical solution
601  FE_P_disc<0,1,3> fe_1d;
602  FE_P_disc<0,2,3> fe_2d;
603 
604  QGauss<1> quad_1d( order );
605  QGauss<2> quad_2d( order );
606 
607  MappingP1<1,3> mapp_1d;
608  MappingP1<2,3> mapp_2d;
609 
610  FEValues<1,3> fe_values_1d(mapp_1d, quad_1d, fe_1d, update_JxW_values | update_quadrature_points);
611  FEValues<2,3> fe_values_2d(mapp_2d, quad_2d, fe_2d, update_JxW_values | update_quadrature_points);
612 
613  // FEValues for velocity.
614  FE_RT0<1,3> fe_rt1d;
615  FE_RT0<2,3> fe_rt2d;
616  FEValues<1,3> fv_rt1d(mapp_1d,quad_1d, fe_rt1d, update_values | update_quadrature_points);
617  FEValues<2,3> fv_rt2d(mapp_2d,quad_2d, fe_rt2d, update_values | update_quadrature_points);
618 
619  FilePath source_file( "analytical_module.py", FilePath::input_file);
620  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
621  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
622 
623  ExactSolution anal_sol_2d(5);
624  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
625 
626 
627  DiffData result;
628  result.dh = &( darcy_flow->get_mh_dofhandler());
629  result.darcy = darcy_flow;
630  result.data_ = darcy_flow->data_.get();
631 
632  // mask 2d elements crossing 1d
633  if (result.data_->mortar_method_ != DarcyMH::NoMortar) {
634  result.velocity_mask.resize(mesh_->n_elements(),0);
636  result.velocity_mask[ isec.bulk_ele_idx() ]++;
637  }
638  }
639 
640  result.pressure_diff.resize( mesh_->n_elements() );
641  result.velocity_diff.resize( mesh_->n_elements() );
642  result.div_diff.resize( mesh_->n_elements() );
643 
644  result.pressure_error[0] = 0;
645  result.velocity_error[0] = 0;
646  result.div_error[0] = 0;
647  result.pressure_error[1] = 0;
648  result.velocity_error[1] = 0;
649  result.div_error[1] = 0;
650  result.mask_vel_error=0;
651 
652  //result.ele_flux = &( ele_flux );
653 
655 
656  auto vel_diff_ptr = result.velocity_diff.create_field<3, FieldValue<3>::Scalar>(1);
658  auto pressure_diff_ptr = result.pressure_diff.create_field<3, FieldValue<3>::Scalar>(1);
659  output_fields.pressure_diff.set_field(mesh_->region_db().get_region_set("ALL"), pressure_diff_ptr, 0);
660  auto div_diff_ptr = result.div_diff.create_field<3, FieldValue<3>::Scalar>(1);
661  output_fields.div_diff.set_field(mesh_->region_db().get_region_set("ALL"), div_diff_ptr, 0);
662 
665 
666 
667  unsigned int solution_size;
668  darcy_flow->get_solution_vector(result.solution, solution_size);
669 
670 
671  FOR_ELEMENTS( mesh_, ele) {
672 
673  switch (ele->dim()) {
674  case 1:
675 
676  l2_diff_local<1>( ele, fe_values_1d, fv_rt1d, anal_sol_1d, result);
677  break;
678  case 2:
679  l2_diff_local<2>( ele, fe_values_2d, fv_rt2d, anal_sol_2d, result);
680  break;
681  }
682  }
683 
684  os << "l2 norm output\n\n"
685  << "pressure error 1d: " << sqrt(result.pressure_error[0]) << endl
686  << "pressure error 2d: " << sqrt(result.pressure_error[1]) << endl
687  << "velocity error 1d: " << sqrt(result.velocity_error[0]) << endl
688  << "velocity error 2d: " << sqrt(result.velocity_error[1]) << endl
689  << "masked velocity error 2d: " << sqrt(result.mask_vel_error) <<endl
690  << "div error 1d: " << sqrt(result.div_error[0]) << endl
691  << "div error 2d: " << sqrt(result.div_error[1]);
692 }
693 
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:155
Classes with algorithms for computation of intersections of meshes.
double getY() const
Definition: nodes.hh:59
Definition: nodes.hh:32
Field< 3, FieldValue< 3 >::Scalar > field_node_pressure
MixedMeshIntersections & mixed_intersections()
Definition: mesh.cc:689
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:444
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:241
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
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
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
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
const RegionDB & region_db() const
Definition: mesh.h:160
Field< 3, FieldValue< 3 >::Integer > subdomain
const TimeStep & step(int index=-1) const
std::shared_ptr< DOFHandlerMultiDim > dh_
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:198
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]
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.cc:286
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:146
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.
double getX() const
Definition: nodes.hh:57
#define START_TIMER(tag)
Starts a timer with specified tag.
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
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
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:97
#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:142
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:488
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:250
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)
std::vector< IntersectionLocal< 1, 2 > > intersection_storage12_
Stores 1D-2D intersections.
VectorSeqDouble corner_pressure
Definitions of Raviart-Thomas finite elements.
MortarMethod mortar_method_
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:233
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:235
void output()
Calculate values for output.
Transformed quadrature weights.
const MH_DofHandler & get_mh_dofhandler() override
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
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