Flow123d  release_2.2.0-33-g759111d
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 ]
242  - (darcy_flow->data_->gravity_[3] + arma::dot(darcy_flow->data_->gravity_vec_,ele->centre()));
243  }
244 }
245 
246 
247 
248 /****
249  * compute Darcian velocity in centre of elements
250  *
251  */
253  START_TIMER("DarcyFlowMHOutput::make_element_vector");
254  // need to call this to create mh solution vector
256 
257  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(darcy_flow->data_);
258  arma::vec3 flux_in_center;
259  for(unsigned int i_ele : element_indices) {
260  ElementFullIter ele = mesh_->element(i_ele);
261 
262  flux_in_center = multidim_assembler[ele->dim() -1]->make_element_vector(ele);
263 
264  // place it in the sequential vector
265  for(unsigned int j=0; j<3; j++) ele_flux[3*i_ele + j]=flux_in_center[j];
266  }
267 }
268 
269 
271 {
272  START_TIMER("DarcyFlowMHOutput::make_corner_scalar");
273  unsigned int ndofs = max(dh->fe<1>()->n_dofs(), max(dh->fe<2>()->n_dofs(), dh->fe<3>()->n_dofs()));
274  unsigned int indices[ndofs];
275  unsigned int i_node;
276  FOR_ELEMENTS(mesh_, ele)
277  {
278  dh->get_dof_indices(ele, indices);
279  FOR_ELEMENT_NODES(ele, i_node)
280  {
281  corner_pressure[indices[i_node]] = node_scalar[mesh_->node_vector.index(ele->node[i_node])];
282  }
283  }
284 }
285 
286 
287 //=============================================================================
288 // pressure interpolation
289 //
290 // some sort of weighted average over elements and sides/edges of an node
291 //
292 // TODO:
293 // questions:
294 // - explain details of averaging and motivation for this type
295 // - edge scalars are stronger since thay are computed twice by each side
296 // - why division by (nod.aux-1)
297 // - ?implement without TED, TSD global arrays with single loop over elements, sides and one loop over nodes for normalization
298 //
299 //=============================================================================
300 
302  START_TIMER("DarcyFlowMHOutput::make_node_scalar_param");
303 
304  vector<double> scalars(mesh_->n_nodes());
305 
306  double dist; //!< tmp variable for storing particular distance node --> element, node --> side*/
307 
308  /** Iterators */
309  const Node * node;
310 
311  int n_nodes = mesh_->node_vector.size(); //!< number of nodes in the mesh */
312  int node_index = 0; //!< index of each node */
313 
314  int* sum_elements = new int [n_nodes]; //!< sum elements joined to node */
315  int* sum_sides = new int [n_nodes]; //!< sum sides joined to node */
316  double* sum_ele_dist = new double [n_nodes]; //!< sum distances to all joined elements */
317  double* sum_side_dist = new double [n_nodes]; //!< Sum distances to all joined sides */
318 
319  /** tmp variables, will be replaced by ini keys
320  * TODO include them into ini file*/
321  bool count_elements = true; //!< scalar is counted via elements*/
322  bool count_sides = true; //!< scalar is counted via sides */
323 
324 
326 
327  /** init arrays */
328  for (int i = 0; i < n_nodes; i++){
329  sum_elements[i] = 0;
330  sum_sides[i] = 0;
331  sum_ele_dist[i] = 0.0;
332  sum_side_dist[i] = 0.0;
333  scalars[i] = 0.0;
334  };
335 
336  /**first pass - calculate sums (weights)*/
337  if (count_elements){
338  FOR_ELEMENTS(mesh_, ele)
339  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
340  node = ele->node[li]; //!< get Node pointer from element */
341  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
342 
343  dist = sqrt(
344  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
345  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
346  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
347  );
348  sum_ele_dist[node_index] += dist;
349  sum_elements[node_index]++;
350  }
351  }
352  if (count_sides){
353  FOR_SIDES(mesh_, side) {
354  for (unsigned int li = 0; li < side->n_nodes(); li++) {
355  node = side->node(li);//!< get Node pointer from element */
356  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
357  dist = sqrt(
358  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
359  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
360  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
361  );
362 
363  sum_side_dist[node_index] += dist;
364  sum_sides[node_index]++;
365  }
366  }
367  }
368 
369  /**second pass - calculate scalar */
370  if (count_elements){
371  FOR_ELEMENTS(mesh_, ele)
372  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
373  node = ele->node[li];//!< get Node pointer from element */
374  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
375 
376  /**TODO - calculate it again or store it in prior pass*/
377  dist = sqrt(
378  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
379  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
380  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
381  );
382  scalars[node_index] += ele_pressure[ele.index()] *
383  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
384  (sum_elements[node_index] + sum_sides[node_index] - 1);
385  }
386  }
387  if (count_sides) {
388  FOR_SIDES(mesh_, side) {
389  for (unsigned int li = 0; li < side->n_nodes(); li++) {
390  node = side->node(li);//!< get Node pointer from element */
391  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
392 
393  /**TODO - calculate it again or store it in prior pass*/
394  dist = sqrt(
395  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
396  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
397  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
398  );
399 
400 
401  scalars[node_index] += dh.side_scalar( *side ) *
402  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
403  (sum_sides[node_index] + sum_elements[node_index] - 1);
404  }
405  }
406  }
407 
408  /** free memory */
409  delete [] sum_elements;
410  delete [] sum_sides;
411  delete [] sum_ele_dist;
412  delete [] sum_side_dist;
413 
414  make_corner_scalar(scalars);
415 }
416 
417 
418 /*
419  * Output of internal flow data.
420  */
422 {
423  START_TIMER("DarcyFlowMHOutput::output_internal_flow_data");
425 
426  if (! raw_output_file.is_open()) return;
427 
428  //char dbl_fmt[ 16 ]= "%.8g ";
429  // header
430  raw_output_file << "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n] ns_side_neighbors[n] neighbors[n*ns] n_vb_neighbors vb_neighbors[n_vb]\n";
431  raw_output_file << fmt::format("$FlowField\nT={}\n", darcy_flow->time().t());
432  raw_output_file << fmt::format("{}\n" , mesh_->n_elements() );
433 
434  ;
435  int cit = 0;
436 
437  // map from higher dim elements to its lower dim neighbors, using gmsh IDs: ele->id()
438  unsigned int undefined_ele_id = -1;
440  FOR_ELEMENTS( mesh_, ele ) {
441  if(ele->n_neighs_vb > 0){
442  for (unsigned int i = 0; i < ele->n_neighs_vb; i++){
443  ElementFullIter higher_ele = ele->neigh_vb[i]->side()->element();
444 
445  auto search = neigh_vb_map.find(higher_ele->id());
446  if(search != neigh_vb_map.end()){
447  // if found, add id to correct locel side idx
448  search->second[ele->neigh_vb[i]->side()->el_idx()] = ele->id();
449  }
450  else{
451  // if not found, create new vector, each side can have one vb neighbour
452  std::vector<unsigned int> higher_ele_side_ngh(higher_ele->n_sides(), undefined_ele_id);
453  higher_ele_side_ngh[ele->neigh_vb[i]->side()->el_idx()] = ele->id();
454  neigh_vb_map[higher_ele->id()] = higher_ele_side_ngh;
455  }
456  }
457  }
458  }
459 
460  FOR_ELEMENTS( mesh_, ele ) {
461  raw_output_file << fmt::format("{} {} ", ele.id(), ele_pressure[cit]);
462  for (unsigned int i = 0; i < 3; i++)
463  raw_output_file << ele_flux[3*cit+i] << " ";
464 
465  raw_output_file << ele->n_sides() << " ";
466 
467  for (unsigned int i = 0; i < ele->n_sides(); i++) {
468  raw_output_file << dh.side_scalar( *(ele->side(i) ) ) << " ";
469  }
470  for (unsigned int i = 0; i < ele->n_sides(); i++) {
471  raw_output_file << dh.side_flux( *(ele->side(i) ) ) << " ";
472  }
473 
474  auto search_neigh = neigh_vb_map.end();
475  for (unsigned int i = 0; i < ele->n_sides(); i++) {
476  unsigned int n_side_neighs = ele->side(i)->edge()->n_sides-1; //n_sides - the current one
477  // check vb neighbors (lower dimension)
478  if(n_side_neighs == 0){
479  //update search
480  if(search_neigh == neigh_vb_map.end())
481  search_neigh = neigh_vb_map.find(ele->id());
482 
483  if(search_neigh != neigh_vb_map.end())
484  if(search_neigh->second[i] != undefined_ele_id)
485  n_side_neighs = 1;
486  }
487  raw_output_file << n_side_neighs << " ";
488  }
489 
490  for (unsigned int i = 0; i < ele->n_sides(); i++) {
491  Edge* edge = ele->side(i)->edge();
492  if(ele->side(i)->edge()->n_sides > 1){
493  for (unsigned int j = 0; j < edge->n_sides; j++) {
494  if(edge->side(j) != ele->side(i))
495  raw_output_file << edge->side(j)->element()->id() << " ";
496  }
497  }
498  //check vb neighbour
499  else if(search_neigh != neigh_vb_map.end()
500  && search_neigh->second[i] != undefined_ele_id){
501  raw_output_file << search_neigh->second[i] << " ";
502  }
503  }
504 
505  // list higher dim neighbours
506  raw_output_file << ele->n_neighs_vb << " ";
507  for (unsigned int i = 0; i < ele->n_neighs_vb; i++)
508  raw_output_file << ele->neigh_vb[i]->side()->element()->id() << " ";
509 
510  raw_output_file << endl;
511  cit ++;
512  }
513  raw_output_file << "$EndFlowField\n" << endl;
514 }
515 
516 
518 #include "fem/fe_p.hh"
519 #include "fem/fe_values.hh"
520 #include "fem/mapping_p1.hh"
521 #include "fields/field_python.hh"
522 #include "fields/field_values.hh"
523 
525 
526 /*
527 * Calculate approximation of L2 norm for:
528  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
529  * 2) difference between RT velocities and analytical solution
530  * 3) difference of divergence
531  * */
532 
533 struct DiffData {
534  double pressure_error[2], velocity_error[2], div_error[2];
539 
540 
541  double * solution;
542  const MH_DofHandler * dh;
543 
544  //std::vector< std::vector<double> > *ele_flux;
548 };
549 
550 template <int dim>
552  FEValues<dim,3> &fe_values, FEValues<dim,3> &fv_rt,
553  ExactSolution &anal_sol, DiffData &result) {
554 
555  fv_rt.reinit(ele);
556  fe_values.reinit(ele);
557 
558  double conductivity = result.data_->conductivity.value(ele->centre(), ele->element_accessor() );
559  double cross = result.data_->cross_section.value(ele->centre(), ele->element_accessor() );
560 
561 
562  // get coefficients on the current element
563  vector<double> fluxes(dim+1);
564 // vector<double> pressure_traces(dim+1);
565 
566  for (unsigned int li = 0; li < ele->n_sides(); li++) {
567  fluxes[li] = result.dh->side_flux( *(ele->side( li ) ) );
568 // pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
569  }
570  double pressure_mean = result.dh->element_scalar(ele);
571 
572  arma::vec analytical(5);
573  arma::vec3 flux_in_q_point;
574  arma::vec3 anal_flux;
575 
576  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
577 
578  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
579  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
580  double mean_x_squared=0;
581  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
582  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
583  {
584  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
585  * arma::dot( ele->node[i_node]->point(), ele->node[j_node]->point());
586  }
587 
588  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
589  arma::vec3 q_point = fe_values.point(i_point);
590 
591 
592  analytical = anal_sol.value(q_point, ele->element_accessor() );
593  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
594 
595  // compute postprocesed pressure
596  diff = 0;
597  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
598  unsigned int oposite_node = RefElement<dim>::oposite_node(i_shape);
599 
600  diff += fluxes[ i_shape ] *
601  ( arma::dot( q_point, q_point )/ 2
602  - mean_x_squared / 2
603  - arma::dot( q_point, ele->node[oposite_node]->point() )
604  + arma::dot( ele->centre(), ele->node[oposite_node]->point() )
605  );
606  }
607 
608  diff = - (1.0 / conductivity) * diff / dim / ele->measure() / cross + pressure_mean ;
609  diff = ( diff - analytical[0]);
610  pressure_diff += diff * diff * fe_values.JxW(i_point);
611 
612 
613  // velocity difference
614  flux_in_q_point.zeros();
615  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
616  flux_in_q_point += fluxes[ i_shape ]
617  * fv_rt.shape_vector(i_shape, i_point)
618  / cross;
619  }
620 
621  flux_in_q_point -= anal_flux;
622  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
623 
624  // divergence diff
625  diff = 0;
626  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
627  diff = ( diff / ele->measure() / cross - analytical[4]);
628  divergence_diff += diff * diff * fe_values.JxW(i_point);
629 
630  }
631 
632 
633  result.velocity_diff[ele.index()] = velocity_diff;
634  result.velocity_error[dim-1] += velocity_diff;
635  if (dim == 2) {
636  result.mask_vel_error += (result.velocity_mask[ ele.index() ])? 0 : velocity_diff;
637  }
638 
639  result.pressure_diff[ele.index()] = pressure_diff;
640  result.pressure_error[dim-1] += pressure_diff;
641 
642  result.div_diff[ele.index()] = divergence_diff;
643  result.div_error[dim-1] += divergence_diff;
644 
645 }
646 
647 
648 
649 
650 
651 
653  DebugOut() << "l2 norm output\n";
654  ofstream os( FilePath("solution_error", FilePath::output_file) );
655 
656  const unsigned int order = 4; // order of Gauss quadrature
657 
658  // we create trivial Dofhandler , for P0 elements, to get access to, FEValues on individual elements
659  // this we use to integrate our own functions - difference of postprocessed pressure and analytical solution
660  FE_P_disc<0,1,3> fe_1d;
661  FE_P_disc<0,2,3> fe_2d;
662 
663  QGauss<1> quad_1d( order );
664  QGauss<2> quad_2d( order );
665 
666  MappingP1<1,3> mapp_1d;
667  MappingP1<2,3> mapp_2d;
668 
669  FEValues<1,3> fe_values_1d(mapp_1d, quad_1d, fe_1d, update_JxW_values | update_quadrature_points);
670  FEValues<2,3> fe_values_2d(mapp_2d, quad_2d, fe_2d, update_JxW_values | update_quadrature_points);
671 
672  // FEValues for velocity.
673  FE_RT0<1,3> fe_rt1d;
674  FE_RT0<2,3> fe_rt2d;
675  FEValues<1,3> fv_rt1d(mapp_1d,quad_1d, fe_rt1d, update_values | update_quadrature_points);
676  FEValues<2,3> fv_rt2d(mapp_2d,quad_2d, fe_rt2d, update_values | update_quadrature_points);
677 
678  FilePath source_file( "analytical_module.py", FilePath::input_file);
679  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
680  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
681 
682  ExactSolution anal_sol_2d(5);
683  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
684 
685 
686  DiffData result;
687 
688  // mask 2d elements crossing 1d
689  result.velocity_mask.resize(mesh_->n_elements(),0);
690  for(Intersection & isec : mesh_->intersections) {
691  result.velocity_mask[ mesh_->element.index( isec.slave_iter() ) ]++;
692  }
693 
694  result.pressure_diff.resize( mesh_->n_elements() );
695  result.velocity_diff.resize( mesh_->n_elements() );
696  result.div_diff.resize( mesh_->n_elements() );
697 
698  result.pressure_error[0] = 0;
699  result.velocity_error[0] = 0;
700  result.div_error[0] = 0;
701  result.pressure_error[1] = 0;
702  result.velocity_error[1] = 0;
703  result.div_error[1] = 0;
704  result.mask_vel_error=0;
705 
706  //result.ele_flux = &( ele_flux );
707 
709 
710  auto vel_diff_ptr = result.velocity_diff.create_field<3, FieldValue<3>::Scalar>(1);
712  auto pressure_diff_ptr = result.pressure_diff.create_field<3, FieldValue<3>::Scalar>(1);
713  output_fields.pressure_diff.set_field(mesh_->region_db().get_region_set("ALL"), pressure_diff_ptr, 0);
714  auto div_diff_ptr = result.div_diff.create_field<3, FieldValue<3>::Scalar>(1);
715  output_fields.div_diff.set_field(mesh_->region_db().get_region_set("ALL"), div_diff_ptr, 0);
716 
719 
720 
721  unsigned int solution_size;
722  darcy_flow->get_solution_vector(result.solution, solution_size);
723 
724  result.dh = &( darcy_flow->get_mh_dofhandler());
725  result.darcy = darcy_flow;
726  result.data_ = darcy_flow->data_.get();
727 
728  FOR_ELEMENTS( mesh_, ele) {
729 
730  switch (ele->dim()) {
731  case 1:
732 
733  l2_diff_local<1>( ele, fe_values_1d, fv_rt1d, anal_sol_1d, result);
734  break;
735  case 2:
736  l2_diff_local<2>( ele, fe_values_2d, fv_rt2d, anal_sol_2d, result);
737  break;
738  }
739  }
740 
741  os << "l2 norm output\n\n"
742  << "pressure error 1d: " << sqrt(result.pressure_error[0]) << endl
743  << "pressure error 2d: " << sqrt(result.pressure_error[1]) << endl
744  << "velocity error 1d: " << sqrt(result.velocity_error[0]) << endl
745  << "velocity error 2d: " << sqrt(result.velocity_error[1]) << endl
746  << "masked velocity error 2d: " << sqrt(result.mask_vel_error) <<endl
747  << "div error 1d: " << sqrt(result.div_error[0]) << endl
748  << "div error 2d: " << sqrt(result.div_error[1]);
749 }
750 
751 
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:187
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:61
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:426
Class for declaration of the input of type Bool.
Definition: type_base.hh:458
MappingP1< 1, 3 > map1
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:233
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:64
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:89
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:147
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
int n_sides
Definition: edges.h:36
Definition: edges.h:26
const RegionDB & region_db() const
Definition: mesh.h:155
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_
Edge * edge() const
Definition: side_impl.hh:65
unsigned int n_sides()
Definition: mesh.cc:175
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:303
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.
int id()
Returns id of the iterator. See VectorId documentation.
Definition: sys_vector.hh:169
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:124
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:141
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:292
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:362
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:490
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:54
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:196
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
ElementFullIter element() const
Definition: side_impl.hh:52
bool compute_errors_
Specific experimental error computing.
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:198
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:247
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:137
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:182
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:533
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
#define FOR_SIDES(_mesh_, it)
Definition: mesh.h:470
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:242
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)
SideIter side(const unsigned int i) const
Definition: edges.h:31
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:226
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:228
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