Flow123d  last_with_con_2.0.0-4-g42e6930
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  raw_output_file_path.open_stream(raw_output_file);
151  }
152 
153  }
154  }
155 }
156 
157 
158 
160 {
161  chkerr(VecDestroy(&vec_corner_pressure));
162 
163  delete dh;
164 };
165 
166 
167 
168 
169 
170 //=============================================================================
171 // CONVERT SOLUTION, CALCULATE BALANCES, ETC...
172 //=============================================================================
173 
174 
176 {
177  START_TIMER("Darcy fields output");
178 
179  ElementSetRef observed_elements = output_stream->observe()->observed_elements();
180  {
181  START_TIMER("post-process output fields");
182 
184 
188  else
189  make_element_scalar(observed_elements);
190 
193  else
194  make_element_vector(observed_elements);
195 
198  //else
199  // make_node_scalar_param(observed_elements);
200 
201  // Internal output only if both ele_pressure and ele_flux are output.
205 
207  }
208 
209  {
210  START_TIMER("evaluate output fields");
212  }
213 
214  {
215  START_TIMER("write time frame");
216  output_stream->write_time_frame();
217  }
218 
219 
220 }
221 
222 
223 //=============================================================================
224 // FILL TH "SCALAR" FIELD FOR ALL ELEMENTS IN THE MESH
225 //=============================================================================
226 
228 {
229  START_TIMER("DarcyFlowMHOutput::make_element_scalar");
230  unsigned int sol_size;
231  double *sol;
232 
233  darcy_flow->get_solution_vector(sol, sol_size);
234  unsigned int soi = mesh_->n_sides();
235  for(unsigned int i_ele : element_indices) {
236  ElementFullIter ele = mesh_->element(i_ele);
237  ele_pressure[i_ele] = sol[ soi + i_ele];
238  ele_piezo_head[i_ele] = sol[soi + i_ele ] + ele->centre()[Mesh::z_coord];
239  }
240 }
241 
242 
243 
244 /****
245  * compute Darcian velocity in centre of elements
246  *
247  */
249  START_TIMER("DarcyFlowMHOutput::make_element_vector");
250  // need to call this to create mh solution vector
252 
253  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(darcy_flow->data_);
254  arma::vec3 flux_in_center;
255  for(unsigned int i_ele : element_indices) {
256  ElementFullIter ele = mesh_->element(i_ele);
257 
258  flux_in_center = multidim_assembler[ele->dim() -1]->make_element_vector(ele);
259 
260  // place it in the sequential vector
261  for(unsigned int j=0; j<3; j++) ele_flux[3*i_ele + j]=flux_in_center[j];
262  }
263 }
264 
265 
267 {
268  START_TIMER("DarcyFlowMHOutput::make_corner_scalar");
269  unsigned int ndofs = max(dh->fe<1>()->n_dofs(), max(dh->fe<2>()->n_dofs(), dh->fe<3>()->n_dofs()));
270  unsigned int indices[ndofs];
271  unsigned int i_node;
272  FOR_ELEMENTS(mesh_, ele)
273  {
274  dh->get_dof_indices(ele, indices);
275  FOR_ELEMENT_NODES(ele, i_node)
276  {
277  corner_pressure[indices[i_node]] = node_scalar[mesh_->node_vector.index(ele->node[i_node])];
278  }
279  }
280 }
281 
282 
283 //=============================================================================
284 // pressure interpolation
285 //
286 // some sort of weighted average over elements and sides/edges of an node
287 //
288 // TODO:
289 // questions:
290 // - explain details of averaging and motivation for this type
291 // - edge scalars are stronger since thay are computed twice by each side
292 // - why division by (nod.aux-1)
293 // - ?implement without TED, TSD global arrays with single loop over elements, sides and one loop over nodes for normalization
294 //
295 //=============================================================================
296 
298  START_TIMER("DarcyFlowMHOutput::make_node_scalar_param");
299 
300  vector<double> scalars(mesh_->n_nodes());
301 
302  double dist; //!< tmp variable for storing particular distance node --> element, node --> side*/
303 
304  /** Iterators */
305  const Node * node;
306 
307  int n_nodes = mesh_->node_vector.size(); //!< number of nodes in the mesh */
308  int node_index = 0; //!< index of each node */
309 
310  int* sum_elements = new int [n_nodes]; //!< sum elements joined to node */
311  int* sum_sides = new int [n_nodes]; //!< sum sides joined to node */
312  double* sum_ele_dist = new double [n_nodes]; //!< sum distances to all joined elements */
313  double* sum_side_dist = new double [n_nodes]; //!< Sum distances to all joined sides */
314 
315  /** tmp variables, will be replaced by ini keys
316  * TODO include them into ini file*/
317  bool count_elements = true; //!< scalar is counted via elements*/
318  bool count_sides = true; //!< scalar is counted via sides */
319 
320 
322 
323  /** init arrays */
324  for (int i = 0; i < n_nodes; i++){
325  sum_elements[i] = 0;
326  sum_sides[i] = 0;
327  sum_ele_dist[i] = 0.0;
328  sum_side_dist[i] = 0.0;
329  scalars[i] = 0.0;
330  };
331 
332  /**first pass - calculate sums (weights)*/
333  if (count_elements){
334  FOR_ELEMENTS(mesh_, ele)
335  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
336  node = ele->node[li]; //!< get Node pointer from element */
337  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
338 
339  dist = sqrt(
340  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
341  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
342  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
343  );
344  sum_ele_dist[node_index] += dist;
345  sum_elements[node_index]++;
346  }
347  }
348  if (count_sides){
349  FOR_SIDES(mesh_, side) {
350  for (unsigned int li = 0; li < side->n_nodes(); li++) {
351  node = side->node(li);//!< get Node pointer from element */
352  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
353  dist = sqrt(
354  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
355  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
356  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
357  );
358 
359  sum_side_dist[node_index] += dist;
360  sum_sides[node_index]++;
361  }
362  }
363  }
364 
365  /**second pass - calculate scalar */
366  if (count_elements){
367  FOR_ELEMENTS(mesh_, ele)
368  for (unsigned int li = 0; li < ele->n_nodes(); li++) {
369  node = ele->node[li];//!< get Node pointer from element */
370  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
371 
372  /**TODO - calculate it again or store it in prior pass*/
373  dist = sqrt(
374  ((node->getX() - ele->centre()[ 0 ])*(node->getX() - ele->centre()[ 0 ])) +
375  ((node->getY() - ele->centre()[ 1 ])*(node->getY() - ele->centre()[ 1 ])) +
376  ((node->getZ() - ele->centre()[ 2 ])*(node->getZ() - ele->centre()[ 2 ]))
377  );
378  scalars[node_index] += ele_pressure[ele.index()] *
379  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
380  (sum_elements[node_index] + sum_sides[node_index] - 1);
381  }
382  }
383  if (count_sides) {
384  FOR_SIDES(mesh_, side) {
385  for (unsigned int li = 0; li < side->n_nodes(); li++) {
386  node = side->node(li);//!< get Node pointer from element */
387  node_index = mesh_->node_vector.index(node); //!< get nod index from mesh */
388 
389  /**TODO - calculate it again or store it in prior pass*/
390  dist = sqrt(
391  ((node->getX() - side->centre()[ 0 ])*(node->getX() - side->centre()[ 0 ])) +
392  ((node->getY() - side->centre()[ 1 ])*(node->getY() - side->centre()[ 1 ])) +
393  ((node->getZ() - side->centre()[ 2 ])*(node->getZ() - side->centre()[ 2 ]))
394  );
395 
396 
397  scalars[node_index] += dh.side_scalar( *side ) *
398  (1 - dist / (sum_ele_dist[node_index] + sum_side_dist[node_index])) /
399  (sum_sides[node_index] + sum_elements[node_index] - 1);
400  }
401  }
402  }
403 
404  /** free memory */
405  delete [] sum_elements;
406  delete [] sum_sides;
407  delete [] sum_ele_dist;
408  delete [] sum_side_dist;
409 
410  make_corner_scalar(scalars);
411 }
412 
413 
414 /*
415  * Output of internal flow data.
416  */
418 {
419  START_TIMER("DarcyFlowMHOutput::output_internal_flow_data");
421 
422  if (! raw_output_file.is_open()) return;
423 
424  //char dbl_fmt[ 16 ]= "%.8g ";
425  // header
426  raw_output_file << "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n]\n";
427  raw_output_file << fmt::format("$FlowField\nT={}\n", darcy_flow->time().t());
428  raw_output_file << fmt::format("{}\n" , mesh_->n_elements() );
429 
430  ;
431  int cit = 0;
432  FOR_ELEMENTS( mesh_, ele ) {
433  raw_output_file << fmt::format("{} {} ", ele.id(), ele_pressure[cit]);
434  for (unsigned int i = 0; i < 3; i++)
435  raw_output_file << ele_flux[3*cit+i] << " ";
436 
437  raw_output_file << ele->n_sides() << " ";
438  std::vector< std::vector<unsigned int > > old_to_new_side =
439  { {0, 1},
440  {0, 1, 2},
441  {0, 1, 2, 3}
442  };
443  for (unsigned int i = 0; i < ele->n_sides(); i++) {
444  unsigned int i_new_side = old_to_new_side[ele->dim()-1][i];
445  raw_output_file << dh.side_scalar( *(ele->side(i_new_side) ) ) << " ";
446  }
447  for (unsigned int i = 0; i < ele->n_sides(); i++) {
448  unsigned int i_new_side = old_to_new_side[ele->dim()-1][i];
449  raw_output_file << dh.side_flux( *(ele->side(i_new_side) ) ) << " ";
450  }
451 
452  raw_output_file << endl;
453  cit ++;
454  }
455  raw_output_file << "$EndFlowField\n" << endl;
456 }
457 
458 
460 #include "fem/fe_p.hh"
461 #include "fem/fe_values.hh"
462 #include "fem/mapping_p1.hh"
463 #include "fields/field_python.hh"
464 #include "fields/field_values.hh"
465 
467 
468 /*
469 * Calculate approximation of L2 norm for:
470  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
471  * 2) difference between RT velocities and analytical solution
472  * 3) difference of divergence
473  * */
474 
475 struct DiffData {
476  double pressure_error[2], velocity_error[2], div_error[2];
481 
482 
483  double * solution;
484  const MH_DofHandler * dh;
485 
486  //std::vector< std::vector<double> > *ele_flux;
490 };
491 
492 template <int dim>
494  FEValues<dim,3> &fe_values, FEValues<dim,3> &fv_rt,
495  ExactSolution &anal_sol, DiffData &result) {
496 
497  fv_rt.reinit(ele);
498  fe_values.reinit(ele);
499 
500  double conductivity = result.data_->conductivity.value(ele->centre(), ele->element_accessor() );
501  double cross = result.data_->cross_section.value(ele->centre(), ele->element_accessor() );
502 
503 
504  // get coefficients on the current element
505  vector<double> fluxes(dim+1);
506 // vector<double> pressure_traces(dim+1);
507 
508  for (unsigned int li = 0; li < ele->n_sides(); li++) {
509  fluxes[li] = result.dh->side_flux( *(ele->side( li ) ) );
510 // pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
511  }
512  double pressure_mean = result.dh->element_scalar(ele);
513 
514  arma::vec analytical(5);
515  arma::vec3 flux_in_q_point;
516  arma::vec3 anal_flux;
517 
518  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
519 
520  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
521  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
522  double mean_x_squared=0;
523  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
524  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
525  {
526  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
527  * arma::dot( ele->node[i_node]->point(), ele->node[j_node]->point());
528  }
529 
530  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
531  arma::vec3 q_point = fe_values.point(i_point);
532 
533 
534  analytical = anal_sol.value(q_point, ele->element_accessor() );
535  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
536 
537  // compute postprocesed pressure
538  diff = 0;
539  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
540  unsigned int oposite_node = RefElement<dim>::oposite_node(i_shape);
541 
542  diff += fluxes[ i_shape ] *
543  ( arma::dot( q_point, q_point )/ 2
544  - mean_x_squared / 2
545  - arma::dot( q_point, ele->node[oposite_node]->point() )
546  + arma::dot( ele->centre(), ele->node[oposite_node]->point() )
547  );
548  }
549 
550  diff = - (1.0 / conductivity) * diff / dim / ele->measure() / cross + pressure_mean ;
551  diff = ( diff - analytical[0]);
552  pressure_diff += diff * diff * fe_values.JxW(i_point);
553 
554 
555  // velocity difference
556  flux_in_q_point.zeros();
557  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
558  flux_in_q_point += fluxes[ i_shape ]
559  * fv_rt.shape_vector(i_shape, i_point)
560  / cross;
561  }
562 
563  flux_in_q_point -= anal_flux;
564  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
565 
566  // divergence diff
567  diff = 0;
568  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
569  diff = ( diff / ele->measure() / cross - analytical[4]);
570  divergence_diff += diff * diff * fe_values.JxW(i_point);
571 
572  }
573 
574 
575  result.velocity_diff[ele.index()] = velocity_diff;
576  result.velocity_error[dim-1] += velocity_diff;
577  if (dim == 2) {
578  result.mask_vel_error += (result.velocity_mask[ ele.index() ])? 0 : velocity_diff;
579  }
580 
581  result.pressure_diff[ele.index()] = pressure_diff;
582  result.pressure_error[dim-1] += pressure_diff;
583 
584  result.div_diff[ele.index()] = divergence_diff;
585  result.div_error[dim-1] += divergence_diff;
586 
587 }
588 
589 
590 
591 
592 
593 
595  DebugOut() << "l2 norm output\n";
596  ofstream os( FilePath("solution_error", FilePath::output_file) );
597 
598  const unsigned int order = 4; // order of Gauss quadrature
599 
600  // we create trivial Dofhandler , for P0 elements, to get access to, FEValues on individual elements
601  // this we use to integrate our own functions - difference of postprocessed pressure and analytical solution
602  FE_P_disc<0,1,3> fe_1d;
603  FE_P_disc<0,2,3> fe_2d;
604 
605  QGauss<1> quad_1d( order );
606  QGauss<2> quad_2d( order );
607 
608  MappingP1<1,3> mapp_1d;
609  MappingP1<2,3> mapp_2d;
610 
611  FEValues<1,3> fe_values_1d(mapp_1d, quad_1d, fe_1d, update_JxW_values | update_quadrature_points);
612  FEValues<2,3> fe_values_2d(mapp_2d, quad_2d, fe_2d, update_JxW_values | update_quadrature_points);
613 
614  // FEValues for velocity.
615  FE_RT0<1,3> fe_rt1d;
616  FE_RT0<2,3> fe_rt2d;
617  FEValues<1,3> fv_rt1d(mapp_1d,quad_1d, fe_rt1d, update_values | update_quadrature_points);
618  FEValues<2,3> fv_rt2d(mapp_2d,quad_2d, fe_rt2d, update_values | update_quadrature_points);
619 
620  FilePath source_file( "analytical_module.py", FilePath::input_file);
621  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
622  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
623 
624  ExactSolution anal_sol_2d(5);
625  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
626 
627 
628  DiffData result;
629 
630  // mask 2d elements crossing 1d
631  result.velocity_mask.resize(mesh_->n_elements(),0);
632  for(Intersection & isec : mesh_->intersections) {
633  result.velocity_mask[ mesh_->element.index( isec.slave_iter() ) ]++;
634  }
635 
636  result.pressure_diff.resize( mesh_->n_elements() );
637  result.velocity_diff.resize( mesh_->n_elements() );
638  result.div_diff.resize( mesh_->n_elements() );
639 
640  result.pressure_error[0] = 0;
641  result.velocity_error[0] = 0;
642  result.div_error[0] = 0;
643  result.pressure_error[1] = 0;
644  result.velocity_error[1] = 0;
645  result.div_error[1] = 0;
646  result.mask_vel_error=0;
647 
648  //result.ele_flux = &( ele_flux );
649 
651 
652  auto vel_diff_ptr = result.velocity_diff.create_field<3, FieldValue<3>::Scalar>(1);
654  auto pressure_diff_ptr = result.pressure_diff.create_field<3, FieldValue<3>::Scalar>(1);
655  output_fields.pressure_diff.set_field(mesh_->region_db().get_region_set("ALL"), pressure_diff_ptr, 0);
656  auto div_diff_ptr = result.div_diff.create_field<3, FieldValue<3>::Scalar>(1);
657  output_fields.div_diff.set_field(mesh_->region_db().get_region_set("ALL"), div_diff_ptr, 0);
658 
661 
662 
663  unsigned int solution_size;
664  darcy_flow->get_solution_vector(result.solution, solution_size);
665 
666  result.dh = &( darcy_flow->get_mh_dofhandler());
667  result.darcy = darcy_flow;
668  result.data_ = darcy_flow->data_.get();
669 
670  FOR_ELEMENTS( mesh_, ele) {
671 
672  switch (ele->dim()) {
673  case 1:
674 
675  l2_diff_local<1>( ele, fe_values_1d, fv_rt1d, anal_sol_1d, result);
676  break;
677  case 2:
678  l2_diff_local<2>( ele, fe_values_2d, fv_rt2d, anal_sol_2d, result);
679  break;
680  }
681  }
682 
683  os << "l2 norm output\n\n"
684  << "pressure error 1d: " << sqrt(result.pressure_error[0]) << endl
685  << "pressure error 2d: " << sqrt(result.pressure_error[1]) << endl
686  << "velocity error 1d: " << sqrt(result.velocity_error[0]) << endl
687  << "velocity error 2d: " << sqrt(result.velocity_error[1]) << endl
688  << "masked velocity error 2d: " << sqrt(result.mask_vel_error) <<endl
689  << "div error 1d: " << sqrt(result.div_error[0]) << endl
690  << "div error 2d: " << sqrt(result.div_error[1]);
691 }
692 
TimeGovernor & time()
Definition: equation.hh:146
FieldSet & data()
Definition: equation.hh:191
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.
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
#define FOR_ELEMENT_NODES(i, j)
Definition: elements.h:188
double getY() const
Definition: nodes.hh:59
Definition: nodes.hh:32
Field< 3, FieldValue< 3 >::Scalar > field_node_pressure
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:50
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:410
Class for declaration of the input of type Bool.
Definition: type_base.hh:434
MappingP1< 1, 3 > map1
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:231
VectorSeqDouble velocity_diff
void set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.hh:195
void output(TimeStep step)
void make_element_vector(ElementSetRef element_indices)
std::string format(CStringRef format_str, ArgList args)
Definition: format.h:3141
DOFHandlerMultiDim * dh
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
void distribute_dofs(FiniteElement< 1, 3 > &fe1d, FiniteElement< 2, 3 > &fe2d, FiniteElement< 3, 3 > &fe3d, const unsigned int offset=0)
Distributes degrees of freedom on the mesh needed for the given finite elements.
Definition: dofhandler.cc:245
Fields computed from the mesh data.
Iterator< Ret > find(const string &key) const
void make_corner_scalar(vector< double > &node_scalar)
Helper class that stores data of generic types.
Definition: type_generic.hh:88
Class FEValues calculates finite element data on the actual cells such as shape function values...
int index() const
Definition: sys_vector.hh:78
VectorSeqDouble pressure_diff
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:142
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
const RegionDB & region_db() const
Definition: mesh.h:152
Field< 3, FieldValue< 3 >::Integer > subdomain
const TimeStep & step(int index=-1) const
double t() const
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
std::vector< unsigned int > all_element_idx_
unsigned int n_sides()
Definition: mesh.cc:172
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:286
double div_error[2]
FiniteElement< dim, 3 > * fe() const
Returns finite element object for given space dimension.
void make_node_scalar_param(ElementSetRef element_indices)
Calculate nodes scalar, store it in double* node_scalars instead of node->scalar. ...
static unsigned int oposite_node(unsigned int sid)
Definition: ref_element.hh:251
Symmetric Gauss-Legendre quadrature formulae on simplices.
double getZ() const
Definition: nodes.hh:61
unsigned int size() const
Returns size of the container. This is independent of the allocated space.
Definition: sys_vector.hh:391
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:113
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:138
FE_P_disc< 1, 2, 3 > fe2
const MH_DofHandler & get_mh_dofhandler()
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:277
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:342
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:468
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:187
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
static FileName output()
The factory function for declaring type FileName for input files.
Definition: type_base.hh:612
bool compute_errors_
Specific experimental error computing.
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:195
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:157
std::vector< int > velocity_mask
vector< Intersection > intersections
Definition: mesh.h:240
std::shared_ptr< FieldElementwise< spacedim, Value > > create_field(unsigned int n_comp)
Create and return shared pointer to FieldElementwise object.
Definitions of particular quadrature rules on simplices.
FieldCommon & name(const string &name)
Definition: field_common.hh:86
#define MPI_COMM_WORLD
Definition: mpi.h:123
Field< 3, FieldValue< 3 >::Scalar > div_diff
FE_P_disc< 1, 3, 3 > fe3
unsigned int n_nodes() const
Definition: mesh.h:134
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:171
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 UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:53
#define FOR_SIDES(_mesh_, it)
Definition: mesh.h:454
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:240
void make_element_scalar(ElementSetRef element_indices)
void l2_diff_local(ElementFullIter &ele, FEValues< dim, 3 > &fe_values, FEValues< dim, 3 > &fv_rt, ExactSolution &anal_sol, DiffData &result)
Definitions of Raviart-Thomas finite elements.
const unsigned int n_global_dofs() const
Getter for the number of all mesh dofs required by the given finite element.
Definition: dofhandler.hh:61
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:219
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:221
void output()
Calculate values for output.
Transformed quadrature weights.
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