Flow123d  release_3.0.0-1191-g7740876
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"
40 #include "fem/fe_values_views.hh"
42 #include "fields/field_fe.hh"
44 #include "fields/generic_field.hh"
45 
46 #include "mesh/long_idx.hh"
47 #include "mesh/mesh.h"
48 #include "mesh/partitioning.hh"
49 #include "mesh/accessors.hh"
50 #include "mesh/node_accessor.hh"
51 #include "mesh/range_wrapper.hh"
52 
53 // #include "coupling/balance.hh"
56 
57 namespace it = Input::Type;
58 
59 
62  DarcyMH::EqData eq_data;
63  output_fields += eq_data;
64  return output_fields.make_output_type("Flow_Darcy_MH", "");
65 }
66 
67 
69 
70  static it::Record& rec = it::Record("Output_DarcyMHSpecific", "Specific Darcy flow MH output.")
72  .declare_key("compute_errors", it::Bool(), it::Default("false"),
73  "SPECIAL PURPOSE. Computes error norms of the solution, particulary suited for non-compatible coupling models.")
75  "Output file with raw data from MH module.")
76  .close();
77 
79  return output_fields.make_output_type_from_record(rec,
80  "Flow_Darcy_MH_specific",
81  "");
82 }
83 
84 
87 {
88 
89  *this += field_ele_pressure.name("pressure_p0").units(UnitSI().m())
91  .description("Pressure solution - P0 interpolation.");
92  *this += field_node_pressure.name("pressure_p1").units(UnitSI().m())
94  .description("Pressure solution - P1 interpolation.");
95  *this += field_ele_piezo_head.name("piezo_head_p0").units(UnitSI().m())
97  .description("Piezo head solution - P0 interpolation.");
98  *this += field_ele_flux.name("velocity_p0").units(UnitSI().m().s(-1))
100  .description("Velocity solution - P0 interpolation.");
101  *this += subdomain.name("subdomain")
104  .description("Subdomain ids of the domain decomposition.");
105  *this += region_id.name("region_id")
108  .description("Region ids.");
109 }
110 
111 
113 : EquationOutput()
114 {
115  *this += pressure_diff.name("pressure_diff").units(UnitSI().m())
117  .description("Error norm of the pressure solution. [Experimental]");
118  *this += velocity_diff.name("velocity_diff").units(UnitSI().m().s(-1))
120  .description("Error norm of the velocity solution. [Experimental]");
121  *this += div_diff.name("div_diff").units(UnitSI().s(-1))
123  .description("Error norm of the divergence of the velocity solution. [Experimental]");
124 }
125 
127 : darcy_flow(flow),
128  mesh_(&darcy_flow->mesh()),
129  compute_errors_(false),
130  fe0(1),
132 {
133  Input::Record in_rec_output = main_mh_in_rec.val<Input::Record>("output");
134 
136  main_mh_in_rec.val<Input::Record>("output_stream"),
138  prepare_output(in_rec_output);
139 
140  auto in_rec_specific = main_mh_in_rec.find<Input::Record>("output_specific");
141  if (in_rec_specific) {
142  in_rec_specific->opt_val("compute_errors", compute_errors_);
143 
144  // raw output
145  int rank;
147  if (rank == 0) {
148  // optionally open raw output file
149  FilePath raw_output_file_path;
150  if (in_rec_specific->opt_val("raw_flow_output", raw_output_file_path)) {
151  MessageOut() << "Opening raw flow output: " << raw_output_file_path << "\n";
152  try {
153  raw_output_file_path.open_stream(raw_output_file);
154  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, (*in_rec_specific))
155  }
156  }
157 
158  auto fields_array = in_rec_specific->val<Input::Array>("fields");
159  if(fields_array.size() > 0){
161  prepare_specific_output(*in_rec_specific);
162  }
163  }
164 }
165 
167 {
168  // we need to add data from the flow equation at this point, not in constructor of OutputFields
171 
172  all_element_idx_.resize(mesh_->n_elements());
173  for(unsigned int i=0; i<all_element_idx_.size(); i++) all_element_idx_[i] = i;
174 
175  // create shared pointer to a FieldFE and push this Field to output_field on all regions
177  ele_pressure_ptr=create_field<3, FieldValue<3>::Scalar>(ele_pressure, *mesh_, 1);
178  output_fields.field_ele_pressure.set_field(mesh_->region_db().get_region_set("ALL"), ele_pressure_ptr);
179 
180  ds = std::make_shared<EqualOrderDiscreteSpace>(mesh_, &fe0, &fe_data_1d.fe_p1, &fe_data_2d.fe_p1, &fe_data_3d.fe_p1);
181  DOFHandlerMultiDim dh_par(*mesh_);
182  dh_par.distribute_dofs(ds);
183  dh_ = dh_par.sequential();
184  corner_pressure.resize(dh_->n_global_dofs());
185 
186  auto corner_ptr = make_shared< FieldFE<3, FieldValue<3>::Scalar> >();
187  corner_ptr->set_fe_data(dh_, 0, corner_pressure);
188 
189  output_fields.field_node_pressure.set_field(mesh_->region_db().get_region_set("ALL"), corner_ptr);
191 
192  ele_piezo_head.resize(mesh_->n_elements());
193  ele_piezo_head_ptr=create_field<3, FieldValue<3>::Scalar>(ele_piezo_head, *mesh_, 1);
194  output_fields.field_ele_piezo_head.set_field(mesh_->region_db().get_region_set("ALL"), ele_piezo_head_ptr);
195 
196  ele_flux.resize(3*mesh_->n_elements());
197  ele_flux_ptr=create_field<3, FieldValue<3>::VectorFixed>(ele_flux, *mesh_, 3);
198  output_fields.field_ele_flux.set_field(mesh_->region_db().get_region_set("ALL"), ele_flux_ptr);
199 
202 
203  //output_stream->add_admissible_field_names(in_rec_output.val<Input::Array>("fields"));
204  //output_stream->mark_output_times(darcy_flow->time());
206 }
207 
209 {
211  diff_data.data_ = darcy_flow->data_.get();
212 
213  // mask 2d elements crossing 1d
217  diff_data.velocity_mask[ isec.bulk_ele_idx() ]++;
218  }
219  }
220 
224 
226 
227  diff_data.vel_diff_ptr = create_field<3, FieldValue<3>::Scalar>(diff_data.velocity_diff, *mesh_, 1);
229  diff_data.pressure_diff_ptr = create_field<3, FieldValue<3>::Scalar>(diff_data.pressure_diff, *mesh_, 1);
231  diff_data.div_diff_ptr = create_field<3, FieldValue<3>::Scalar>(diff_data.div_diff, *mesh_, 1);
233 
236 }
237 
239 {};
240 
241 
242 
243 
244 
245 //=============================================================================
246 // CONVERT SOLUTION, CALCULATE BALANCES, ETC...
247 //=============================================================================
248 
249 
251 {
252  START_TIMER("Darcy fields output");
253 
254  ElementSetRef observed_elements = output_stream->observe(mesh_)->observed_elements();
255  {
256  START_TIMER("post-process output fields");
257 
259 
263  else
264  make_element_scalar(observed_elements);
265 
268  else
269  make_element_vector(observed_elements);
270 
273  //else
274  // make_node_scalar_param(observed_elements);
275 
279 
280  // Internal output only if both ele_pressure and ele_flux are output.
283  {
285  }
286  }
287 
288  {
289  START_TIMER("evaluate output fields");
291  }
292 
293  if (compute_errors_)
294  {
295  START_TIMER("compute specific output fields");
297  }
298 
300  {
301  START_TIMER("evaluate output fields");
307  }
308 
309  {
310  START_TIMER("write time frame");
311  output_stream->write_time_frame();
312  }
313 
314 
315 }
316 
317 
318 //=============================================================================
319 // FILL TH "SCALAR" FIELD FOR ALL ELEMENTS IN THE MESH
320 //=============================================================================
321 
323 {
324  START_TIMER("DarcyFlowMHOutput::make_element_scalar");
325  unsigned int sol_size;
326  double *sol;
327 
328  darcy_flow->get_solution_vector(sol, sol_size);
329  unsigned int soi = mesh_->n_sides();
330  for(unsigned int i_ele : element_indices) {
332  ele_pressure[i_ele] = sol[ soi + i_ele];
333  ele_piezo_head[i_ele] = sol[soi + i_ele ]
334  - (darcy_flow->data_->gravity_[3] + arma::dot(darcy_flow->data_->gravity_vec_,ele.centre()));
335  }
336 }
337 
338 
339 
340 /****
341  * compute Darcian velocity in centre of elements
342  *
343  */
345  START_TIMER("DarcyFlowMHOutput::make_element_vector");
346  // need to call this to create mh solution vector
348 
349  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(darcy_flow->data_);
350  arma::vec3 flux_in_center;
351  for(unsigned int i_ele : element_indices) {
353 
354  flux_in_center = multidim_assembler[ele->dim() -1]->make_element_vector(ele);
355 
356  // place it in the sequential vector
357  for(unsigned int j=0; j<3; j++) ele_flux[3*i_ele + j]=flux_in_center[j];
358  }
359 }
360 
361 
363 {
364  START_TIMER("DarcyFlowMHOutput::make_corner_scalar");
365  unsigned int ndofs = dh_->max_elem_dofs();
366  std::vector<LongIdx> indices(ndofs);
367  unsigned int i_node;
368  /*for (DHCellAccessor cell : dh_->local_range()) {
369  cell.get_dof_indices(indices);
370  for (i_node=0; i_node<cell.elm()->n_nodes(); i_node++)
371  {
372  corner_pressure[indices[i_node]] = node_scalar[ cell.elm().node_accessor(i_node).idx() ];
373  }
374  }*/
375  for (auto ele : mesh_->elements_range()) {
376  dh_->cell_accessor_from_element(ele.idx()).get_dof_indices(indices);
377  for (i_node=0; i_node<ele->n_nodes(); i_node++)
378  {
379  corner_pressure[indices[i_node]] = node_scalar[ ele.node_accessor(i_node).idx() ];
380  }
381  }
382 }
383 
384 
385 //=============================================================================
386 // pressure interpolation
387 //
388 // some sort of weighted average over elements and sides/edges of an node
389 //
390 // TODO:
391 // questions:
392 // - explain details of averaging and motivation for this type
393 // - edge scalars are stronger since thay are computed twice by each side
394 // - why division by (nod.aux-1)
395 // - ?implement without TED, TSD global arrays with single loop over elements, sides and one loop over nodes for normalization
396 //
397 //=============================================================================
398 
400  START_TIMER("DarcyFlowMHOutput::make_node_scalar_param");
401 
402  vector<double> scalars(mesh_->n_nodes());
403 
404  double dist; //!< tmp variable for storing particular distance node --> element, node --> side*/
405 
406  /** Accessors */
407  NodeAccessor<3> node;
408 
409  int n_nodes = mesh_->n_nodes(); //!< number of nodes in the mesh */
410  int node_index = 0; //!< index of each node */
411 
412  int* sum_elements = new int [n_nodes]; //!< sum elements joined to node */
413  int* sum_sides = new int [n_nodes]; //!< sum sides joined to node */
414  double* sum_ele_dist = new double [n_nodes]; //!< sum distances to all joined elements */
415  double* sum_side_dist = new double [n_nodes]; //!< Sum distances to all joined sides */
416 
417  /** tmp variables, will be replaced by ini keys
418  * TODO include them into ini file*/
419  bool count_elements = true; //!< scalar is counted via elements*/
420  bool count_sides = true; //!< scalar is counted via sides */
421 
422 
424 
425  /** init arrays */
426  for (int i = 0; i < n_nodes; i++){
427  sum_elements[i] = 0;
428  sum_sides[i] = 0;
429  sum_ele_dist[i] = 0.0;
430  sum_side_dist[i] = 0.0;
431  scalars[i] = 0.0;
432  };
433 
434  /**first pass - calculate sums (weights)*/
435  if (count_elements){
436  for (auto ele : mesh_->elements_range())
437  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
438  node = ele.node_accessor(li); //!< get NodeAccessor from element */
439  node_index = node.idx(); //!< get nod index from mesh */
440 
441  dist = sqrt(
442  ((node->getX() - ele.centre()[ 0 ])*(node->getX() - ele.centre()[ 0 ])) +
443  ((node->getY() - ele.centre()[ 1 ])*(node->getY() - ele.centre()[ 1 ])) +
444  ((node->getZ() - ele.centre()[ 2 ])*(node->getZ() - ele.centre()[ 2 ]))
445  );
446  sum_ele_dist[node_index] += dist;
447  sum_elements[node_index]++;
448  }
449  }
450  if (count_sides){
451  for (auto ele : mesh_->elements_range())
452  for(SideIter side = ele.side(0); side->side_idx() < ele->n_sides(); ++side) {
453  for (unsigned int li = 0; li < side->n_nodes(); li++) {
454  node = side->node(li);//!< get NodeAccessor from element */
455  node_index = node.idx(); //!< get nod index from mesh */
456  dist = sqrt(
457  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
458  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
459  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
460  );
461 
462  sum_side_dist[node_index] += dist;
463  sum_sides[node_index]++;
464  }
465  }
466  }
467 
468  /**second pass - calculate scalar */
469  if (count_elements){
470  for (auto ele : mesh_->elements_range())
471  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
472  node = ele.node_accessor(li);//!< get NodeAccessor from element */
473  node_index = ele.node_accessor(li).idx(); //!< get nod index from mesh */
474 
475  /**TODO - calculate it again or store it in prior pass*/
476  dist = sqrt(
477  ((node->getX() - ele.centre()[ 0 ])*(node->getX() - ele.centre()[ 0 ])) +
478  ((node->getY() - ele.centre()[ 1 ])*(node->getY() - ele.centre()[ 1 ])) +
479  ((node->getZ() - ele.centre()[ 2 ])*(node->getZ() - ele.centre()[ 2 ]))
480  );
481  scalars[node_index] += ele_pressure[ ele.idx() ] *
482  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
483  (sum_elements[node_index] + sum_sides[node_index] - 1);
484  }
485  }
486  if (count_sides) {
487  for (auto ele : mesh_->elements_range())
488  for(SideIter side = ele.side(0); side->side_idx() < ele->n_sides(); ++side) {
489  for (unsigned int li = 0; li < side->n_nodes(); li++) {
490  node = side->node(li);//!< get NodeAccessor from element */
491  node_index = node.idx(); //!< get nod index from mesh */
492 
493  /**TODO - calculate it again or store it in prior pass*/
494  dist = sqrt(
495  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
496  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
497  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
498  );
499 
500 
501  scalars[node_index] += dh.side_scalar( *side ) *
502  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
503  (sum_sides[node_index] + sum_elements[node_index] - 1);
504  }
505  }
506  }
507 
508  /** free memory */
509  delete [] sum_elements;
510  delete [] sum_sides;
511  delete [] sum_ele_dist;
512  delete [] sum_side_dist;
513 
514  make_corner_scalar(scalars);
515 }
516 
517 
518 /*
519  * Output of internal flow data.
520  */
522 {
523  START_TIMER("DarcyFlowMHOutput::output_internal_flow_data");
525 
526  if (! raw_output_file.is_open()) return;
527 
528  //char dbl_fmt[ 16 ]= "%.8g ";
529  // header
530  raw_output_file << "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n]\n";
531  raw_output_file << fmt::format("$FlowField\nT={}\n", darcy_flow->time().t());
532  raw_output_file << fmt::format("{}\n" , mesh_->n_elements() );
533 
534  int cit = 0;
535  for (auto ele : mesh_->elements_range()) {
536  raw_output_file << fmt::format("{} {} ", ele.index(), ele_pressure[cit]);
537  for (unsigned int i = 0; i < 3; i++)
538  raw_output_file << ele_flux[3*cit+i] << " ";
539 
540  raw_output_file << ele->n_sides() << " ";
541 
542  for (unsigned int i = 0; i < ele->n_sides(); i++) {
543  raw_output_file << dh.side_scalar( *(ele.side(i) ) ) << " ";
544  }
545  for (unsigned int i = 0; i < ele->n_sides(); i++) {
546  raw_output_file << dh.side_flux( *(ele.side(i) ) ) << " ";
547  }
548 
549  raw_output_file << endl;
550  cit ++;
551  }
552  raw_output_file << "$EndFlowField\n" << endl;
553 }
554 
555 
557 #include "fem/fe_p.hh"
558 #include "fem/fe_values.hh"
559 #include "fem/mapping_p1.hh"
560 #include "fields/field_python.hh"
561 #include "fields/field_values.hh"
562 
564 
565 /*
566 * Calculate approximation of L2 norm for:
567  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
568  * 2) difference between RT velocities and analytical solution
569  * 3) difference of divergence
570  * */
571 
572 template <int dim>
574  FEValues<dim,3> &fe_values, FEValues<dim,3> &fv_rt,
575  ExactSolution &anal_sol, DarcyFlowMHOutput::DiffData &result) {
576 
577  fv_rt.reinit(ele);
578  fe_values.reinit(ele);
579 
580  double conductivity = result.data_->conductivity.value(ele.centre(), ele );
581  double cross = result.data_->cross_section.value(ele.centre(), ele );
582 
583 
584  // get coefficients on the current element
585  vector<double> fluxes(dim+1);
586 // vector<double> pressure_traces(dim+1);
587 
588  for (unsigned int li = 0; li < ele->n_sides(); li++) {
589  fluxes[li] = result.dh->side_flux( *(ele.side( li ) ) );
590 // pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
591  }
592  double pressure_mean = result.dh->element_scalar(ele);
593 
594  arma::vec analytical(5);
595  arma::vec3 flux_in_q_point;
596  arma::vec3 anal_flux;
597 
598  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
599 
600  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
601  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
602  double mean_x_squared=0;
603  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
604  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
605  {
606  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
607  * arma::dot( ele.node(i_node)->point(), ele.node(j_node)->point());
608  }
609 
610  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
611  arma::vec3 q_point = fe_values.point(i_point);
612 
613 
614  analytical = anal_sol.value(q_point, ele );
615  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
616 
617  // compute postprocesed pressure
618  diff = 0;
619  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
620  unsigned int oposite_node = RefElement<dim>::oposite_node(i_shape);
621 
622  diff += fluxes[ i_shape ] *
623  ( arma::dot( q_point, q_point )/ 2
624  - mean_x_squared / 2
625  - arma::dot( q_point, ele.node(oposite_node)->point() )
626  + arma::dot( ele.centre(), ele.node(oposite_node)->point() )
627  );
628  }
629 
630  diff = - (1.0 / conductivity) * diff / dim / ele.measure() / cross + pressure_mean ;
631  diff = ( diff - analytical[0]);
632  pressure_diff += diff * diff * fe_values.JxW(i_point);
633 
634 
635  // velocity difference
636  flux_in_q_point.zeros();
637  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
638  flux_in_q_point += fluxes[ i_shape ]
639  * fv_rt.vector_view(0).value(i_shape, i_point)
640  / cross;
641  }
642 
643  flux_in_q_point -= anal_flux;
644  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
645 
646  // divergence diff
647  diff = 0;
648  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
649  diff = ( diff / ele.measure() / cross - analytical[4]);
650  divergence_diff += diff * diff * fe_values.JxW(i_point);
651 
652  }
653 
654 
655  result.velocity_diff[ ele.idx() ] = sqrt(velocity_diff);
656  result.velocity_error[dim-1] += velocity_diff;
657  if (dim == 2 && result.velocity_mask.size() != 0 ) {
658  result.mask_vel_error += (result.velocity_mask[ ele.idx() ])? 0 : velocity_diff;
659  }
660 
661  result.pressure_diff[ ele.idx() ] = sqrt(pressure_diff);
662  result.pressure_error[dim-1] += pressure_diff;
663 
664  result.div_diff[ ele.idx() ] = sqrt(divergence_diff);
665  result.div_error[dim-1] += divergence_diff;
666 
667 }
668 
670 : fe_p0(0), fe_p1(1), order(4), quad(dim, order),
671  fe_values(mapp,quad,fe_p0,update_JxW_values | update_quadrature_points),
672  fv_rt(mapp,quad,fe_rt,update_values | update_quadrature_points)
673 {}
674 
675 
677  DebugOut() << "l2 norm output\n";
678  ofstream os( FilePath("solution_error", FilePath::output_file) );
679 
680  FilePath source_file( "analytical_module.py", FilePath::input_file);
681  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
682  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
683 
684  ExactSolution anal_sol_2d(5);
685  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
686 
687  ExactSolution anal_sol_3d(5);
688  anal_sol_3d.set_python_field_from_file( source_file, "all_values_3d");
689 
692  for(unsigned int j=0; j<3; j++){
693  diff_data.pressure_error[j] = 0;
694  diff_data.velocity_error[j] = 0;
695  diff_data.div_error[j] = 0;
696  }
697 
698  //diff_data.ele_flux = &( ele_flux );
699 
700  unsigned int solution_size;
702 
703 
704  for (auto ele : mesh_->elements_range()) {
705 
706  switch (ele->dim()) {
707  case 1:
708 
709  l2_diff_local<1>( ele, fe_data_1d.fe_values, fe_data_1d.fv_rt, anal_sol_1d, diff_data);
710  break;
711  case 2:
712  l2_diff_local<2>( ele, fe_data_2d.fe_values, fe_data_2d.fv_rt, anal_sol_2d, diff_data);
713  break;
714  case 3:
715  l2_diff_local<3>( ele, fe_data_3d.fe_values, fe_data_3d.fv_rt, anal_sol_3d, diff_data);
716  break;
717  }
718  }
719 
720  // square root for L2 norm
721  for(unsigned int j=0; j<3; j++){
724  diff_data.div_error[j] = sqrt(diff_data.div_error[j]);
725  }
727 
728  os << "l2 norm output\n\n"
729  << "pressure error 1d: " << diff_data.pressure_error[0] << endl
730  << "pressure error 2d: " << diff_data.pressure_error[1] << endl
731  << "pressure error 3d: " << diff_data.pressure_error[2] << endl
732  << "velocity error 1d: " << diff_data.velocity_error[0] << endl
733  << "velocity error 2d: " << diff_data.velocity_error[1] << endl
734  << "velocity error 3d: " << diff_data.velocity_error[2] << endl
735  << "masked velocity error 2d: " << diff_data.mask_vel_error <<endl
736  << "div error 1d: " << diff_data.div_error[0] << endl
737  << "div error 2d: " << diff_data.div_error[1] << endl
738  << "div error 3d: " << diff_data.div_error[2];
739 }
740 
741 
TimeGovernor & time()
Definition: equation.hh:148
FieldSet & data()
Definition: equation.hh:193
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...
Shape function values.
Definition: update_flags.hh:87
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > vel_diff_ptr
static auto subdomain(Mesh &mesh) -> IndexField
struct DarcyFlowMHOutput::DiffData diff_data
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
Output class for darcy_flow_mh model.
unsigned int n_nodes() const
Definition: elements.h:129
virtual void prepare_specific_output(Input::Record in_rec)
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:329
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
Classes with algorithms for computation of intersections of meshes.
double getY() const
Definition: nodes.hh:58
virtual void prepare_output(Input::Record in_rec)
unsigned int side_idx() const
Definition: side_impl.hh:82
OutputSpecificFields output_specific_fields
Specific quantities for output in DarcyFlowMH - error estimates etc.
Field< 3, FieldValue< 3 >::Scalar > field_node_pressure
MixedMeshIntersections & mixed_intersections()
Definition: mesh.cc:710
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Class for declaration of the input of type Bool.
Definition: type_base.hh:459
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:243
void output(TimeStep step)
const Input::Type::Instance & make_output_type_from_record(Input::Type::Record &in_rec, const string &equation_name, const string &aditional_description="")
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:64
static std::shared_ptr< OutputTime > create_output_stream(const std::string &equation_name, const Input::Record &in_rec, std::string unit_str)
This method delete all object instances of class OutputTime stored in output_streams vector...
Definition: output_time.cc:185
bool is_output_specific_fields
Output specific field stuff.
void l2_diff_local(ElementAccessor< 3 > &ele, FEValues< dim, 3 > &fe_values, FEValues< dim, 3 > &fv_rt, FieldPython< 3, FieldValue< 3 >::Vector > &anal_sol, DiffData &result)
Computes L2 error on an element.
Standard quantities for output in DarcyFlowMH.
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:89
Class FEValues calculates finite element data on the actual cells such as shape function values...
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
Mat< double, N, 1 > vec
Definition: armor.hh:211
SideIter side(const unsigned int loc_index)
Definition: accessors.hh:137
void fill_output_data(VectorMPI &vec_seq, std::shared_ptr< FieldFE< spacedim, Value > > field_ptr)
Definition: field_fe.hh:309
const RegionDB & region_db() const
Definition: mesh.h:143
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > ele_pressure_ptr
Field of pressure head in barycenters of elements.
virtual unsigned int n_nodes() const
Definition: mesh.h:129
const TimeStep & step(int index=-1) const
std::shared_ptr< DOFHandlerMultiDim > dh_
double t() const
Field< 3, FieldValue< 3 >::Scalar > subdomain
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
std::vector< unsigned int > all_element_idx_
std::shared_ptr< OutputTime > output_stream
std::shared_ptr< EqData > data_
Field< 3, FieldValue< 3 >::Scalar > pressure_diff
unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.hh:407
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
std::shared_ptr< DiscreteSpace > ds
virtual ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
Definition: mesh.cc:727
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:318
double getZ() const
Definition: nodes.hh:60
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Definition: accessors.hh:285
unsigned int dim() const
Definition: elements.h:124
Field< 3, FieldValue< 3 >::Scalar > velocity_diff
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
double element_scalar(ElementAccessor< 3 > &ele) const
temporary replacement for DofHandler accessor, scalar (pressure) on element
void open_stream(Stream &stream) const
Definition: file_path.cc:211
bool opt_val(const string &key, Ret &value) const
Field< 3, FieldValue< 3 >::Scalar > field_ele_pressure
Transformed quadrature points.
Global macros to enhance readability and debugging, general constants.
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
static const Input::Type::Instance & get_input_type_specific()
double JxW(const unsigned int point_no) override
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.hh:336
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:292
const Ret val(const string &key) const
unsigned int n_sides() const
Definition: elements.h:135
static auto region_id(Mesh &mesh) -> IndexField
void resize(unsigned int local_size)
Definition: vector_mpi.hh:76
Type dot(const Mat< Type, nRows, nCols > &a, const Mat< Type, nRows, nCols > &b)
Definition: armor.hh:176
double getX() const
Definition: nodes.hh:56
#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:152
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:389
static const Input::Type::Instance & get_input_type()
Field< 3, FieldValue< 3 >::Scalar > div_diff
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:503
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > pressure_diff_ptr
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1041
ofstream raw_output_file
Raw data output file.
double measure() const
Computes the measure of the element.
Definition: accessors.hh:254
unsigned int n_sides() const
Definition: mesh.cc:227
const FEValuesViews::Vector< dim, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:388
#define MPI_Comm_rank
Definition: mpi.h:236
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
double side_flux(const Side &side) const
temporary replacement for DofHandler accessor, flux through given side
FieldCommon & description(const string &description)
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > ele_piezo_head_ptr
Field of piezo-metric head in barycenters of elements.
std::string get_unit_string() const
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:346
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
Field< 3, FieldValue< 3 >::Scalar > region_id
bool compute_errors_
Specific experimental error computing.
static Input::Type::Record & get_input_type()
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:198
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
Definitions of particular quadrature rules on simplices.
FieldCommon & name(const string &name)
#define MPI_COMM_WORLD
Definition: mpi.h:123
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > div_diff_ptr
mixed-hybrid model of linear Darcy flow, possibly unsteady.
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:183
Record type proxy class.
Definition: type_record.hh:182
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > ele_flux_ptr
Field of flux in barycenter of every element.
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)
arma::vec3 & point()
Definition: nodes.hh:67
static FileName output()
The factory function for declaring type FileName for input files.
Definition: type_base.cc:533
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:252
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:111
void make_element_scalar(ElementSetRef element_indices)
std::vector< IntersectionLocal< 1, 2 > > intersection_storage12_
Stores 1D-2D intersections.
Implementation of range helper class.
unsigned int idx() const
Definitions of Raviart-Thomas finite elements.
MortarMethod mortar_method_
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
FieldCommon & output_type(OutputTime::DiscreteSpace rt)
void output()
Calculate values for output.
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:541
Transformed quadrature weights.
const MH_DofHandler & get_mh_dofhandler() override
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
const Node * node(unsigned int ni) const
Definition: accessors.hh:145
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:347