Flow123d  JS_before_hm-2139-g05f84414c
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"
28 #include "flow/darcy_flow_lmh.hh"
29 #include "flow/assembly_mh_old.hh"
30 #include "flow/assembly_lmh_old.hh"
32 
33 #include "io/output_time.hh"
34 #include "io/observe.hh"
35 #include "system/system.hh"
36 #include "system/sys_profiler.hh"
37 #include "system/index_types.hh"
38 
39 #include "fields/field_set.hh"
40 #include "fem/dofhandler.hh"
41 #include "fem/fe_values.hh"
42 #include "fem/fe_rt.hh"
43 #include "fem/fe_values_views.hh"
45 #include "fields/field_fe.hh"
47 #include "fields/generic_field.hh"
48 
49 #include "mesh/mesh.h"
50 #include "mesh/partitioning.hh"
51 #include "mesh/accessors.hh"
52 #include "mesh/node_accessor.hh"
53 #include "mesh/range_wrapper.hh"
54 
55 // #include "coupling/balance.hh"
58 
59 namespace it = Input::Type;
60 
61 
62 const it::Instance & DarcyFlowMHOutput::get_input_type(FieldSet& eq_data, const std::string &equation_name) {
64  output_fields += eq_data;
65  return output_fields.make_output_type(equation_name, "");
66 }
67 
68 
70 
71  static it::Record& rec = it::Record("Output_DarcyMHSpecific", "Specific Darcy flow MH output.")
73  .declare_key("compute_errors", it::Bool(), it::Default("false"),
74  "SPECIAL PURPOSE. Computes error norms of the solution, particulary suited for non-compatible coupling models.")
76  "Output file with raw data from MH module.")
77  .close();
78 
81  "Flow_Darcy_MH_specific",
82  "");
83 }
84 
85 
88 {
89 
90 // *this += field_ele_pressure.name("pressure_p0_old").units(UnitSI().m()) // TODO remove: obsolete field
91 // .flags(FieldFlag::equation_result)
92 // .description("Pressure solution - P0 interpolation.");
93 // *this += field_node_pressure.name("pressure_p1").units(UnitSI().m())
94 // .flags(FieldFlag::equation_result)
95 // .description("Pressure solution - P1 interpolation.");
96 // *this += field_ele_piezo_head.name("piezo_head_p0_old").units(UnitSI().m()) // TODO remove: obsolete field
97 // .flags(FieldFlag::equation_result)
98 // .description("Piezo head solution - P0 interpolation.");
99 // *this += field_ele_flux.name("velocity_p0_old").units(UnitSI().m()) // TODO remove: obsolete field
100 // .flags(FieldFlag::equation_result)
101 // .description("Velocity solution - P0 interpolation.");
102  *this += subdomain.name("subdomain")
105  .description("Subdomain ids of the domain decomposition.");
106  *this += region_id.name("region_id")
109  .description("Region ids.");
110 }
111 
112 
114 : EquationOutput()
115 {
116  *this += pressure_diff.name("pressure_diff").units(UnitSI().m())
118  .description("Error norm of the pressure solution. [Experimental]");
119  *this += velocity_diff.name("velocity_diff").units(UnitSI().m().s(-1))
121  .description("Error norm of the velocity solution. [Experimental]");
122  *this += div_diff.name("div_diff").units(UnitSI().s(-1))
124  .description("Error norm of the divergence of the velocity solution. [Experimental]");
125 }
126 
128 : darcy_flow(flow),
129  mesh_(&darcy_flow->mesh()),
130  compute_errors_(false),
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 
149  // optionally open raw output file
150  FilePath raw_output_file_path;
151  if (in_rec_specific->opt_val("raw_flow_output", raw_output_file_path))
152  {
153  int mpi_size;
154  MPI_Comm_size(MPI_COMM_WORLD, &mpi_size);
155  if(mpi_size > 1)
156  {
157  WarningOut() << "Raw output is not available in parallel computation. MPI size: " << mpi_size << "\n";
158  }
159  else
160  {
161  MessageOut() << "Opening raw flow output: " << raw_output_file_path << "\n";
162  try {
163  raw_output_file_path.open_stream(raw_output_file);
164  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, (*in_rec_specific))
165  }
166  }
167  }
168 
169  auto fields_array = in_rec_specific->val<Input::Array>("fields");
170  if(fields_array.size() > 0){
172  prepare_specific_output(*in_rec_specific);
173  }
174  }
175 }
176 
178 {
179  // we need to add data from the flow equation at this point, not in constructor of OutputFields
182 
185 
186  //output_stream->add_admissible_field_names(in_rec_output.val<Input::Array>("fields"));
187  //output_stream->mark_output_times(darcy_flow->time());
189 }
190 
192 {
193  diff_data.eq_fields_ = nullptr;
194  diff_data.eq_data_ = nullptr;
195  if(DarcyMH* d = dynamic_cast<DarcyMH*>(darcy_flow))
196  {
197  diff_data.eq_fields_ = d->eq_fields_.get();
198  diff_data.eq_data_ = d->eq_data_.get();
199  }
200  else if(DarcyLMH* d = dynamic_cast<DarcyLMH*>(darcy_flow))
201  {
202  diff_data.eq_fields_ = d->eq_fields_.get();
203  diff_data.eq_data_ = d->eq_data_.get();
204  }
206 
207  { // init DOF handlers represents element DOFs
208  uint p_elem_component = 1;
209  diff_data.dh_ = std::make_shared<SubDOFHandlerMultiDim>(diff_data.eq_data_->dh_, p_elem_component);
210  }
211 
212  // mask 2d elements crossing 1d
216  diff_data.velocity_mask[ isec.bulk_ele_idx() ]++;
217  }
218  }
219 
221 
222  diff_data.vel_diff_ptr = create_field_fe<3, FieldValue<3>::Scalar>(diff_data.dh_);
224  diff_data.pressure_diff_ptr = create_field_fe<3, FieldValue<3>::Scalar>(diff_data.dh_);
226  diff_data.div_diff_ptr = create_field_fe<3, FieldValue<3>::Scalar>(diff_data.dh_);
228 
229  darcy_flow->time().step().use_fparser_ = true;
232 }
233 
235 {}
236 
237 
238 
239 
240 
241 //=============================================================================
242 // CONVERT SOLUTION, CALCULATE BALANCES, ETC...
243 //=============================================================================
244 
245 
247 {
248  START_TIMER("Darcy fields output");
249 
250  {
251  START_TIMER("post-process output fields");
252 
253  {
255  }
256  }
257 
258  {
259  START_TIMER("evaluate output fields");
260  darcy_flow->time().step().use_fparser_ = true;
263  }
264 
265  if (compute_errors_)
266  {
267  START_TIMER("compute specific output fields");
269  }
270 
272  {
273  START_TIMER("evaluate output fields");
274  darcy_flow->time().step().use_fparser_ = true;
277  }
278 
279 
280 }
281 
282 
283 
284 /*
285  * Output of internal flow data.
286  */
288 {
289  START_TIMER("DarcyFlowMHOutput::output_internal_flow_data");
290 
291  if (! raw_output_file.is_open()) return;
292 
293  //char dbl_fmt[ 16 ]= "%.8g ";
294  // header
295  raw_output_file << "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n]\n";
296  raw_output_file << fmt::format("$FlowField\nT={}\n", darcy_flow->time().t());
297  raw_output_file << fmt::format("{}\n" , mesh_->n_elements() );
298 
299 
300  DarcyMH::EqFields* eq_fields = nullptr;
301  DarcyMH::EqData* eq_data = nullptr;
302  if(DarcyMH* d = dynamic_cast<DarcyMH*>(darcy_flow))
303  {
304  eq_fields = d->eq_fields_.get();
305  eq_data = d->eq_data_.get();
306  }
307  else if(DarcyLMH* d = dynamic_cast<DarcyLMH*>(darcy_flow))
308  {
309  eq_fields = d->eq_fields_.get();
310  eq_data = d->eq_data_.get();
311  }
312  ASSERT_PTR(eq_data);
313 
314  arma::vec3 flux_in_center;
315 
316  auto permutation_vec = eq_data->dh_->mesh()->element_permutations();
317  for (unsigned int i_elem=0; i_elem<eq_data->dh_->n_own_cells(); ++i_elem) {
318  ElementAccessor<3> ele(eq_data->dh_->mesh(), permutation_vec[i_elem]);
319  DHCellAccessor dh_cell = eq_data->dh_->cell_accessor_from_element( ele.idx() );
320  LocDofVec indices = dh_cell.get_loc_dof_indices();
321 
322  std::stringstream ss;
323  // pressure
324  ss << fmt::format("{} {} ", dh_cell.elm().input_id(), eq_data->full_solution.get(indices[ele->n_sides()]));
325 
326  // velocity at element center
327  flux_in_center = eq_fields->field_ele_velocity.value(ele.centre(), ele);
328  for (unsigned int i = 0; i < 3; i++)
329  ss << flux_in_center[i] << " ";
330 
331  // number of sides
332  ss << ele->n_sides() << " ";
333 
334  // use node permutation to permute sides
335  auto &new_to_old_node = ele.orig_nodes_order();
336  std::vector<uint> old_to_new_side(ele->n_sides());
337  for (unsigned int i = 0; i < ele->n_sides(); i++) {
338  // According to RefElement<dim>::opposite_node()
339  uint new_opp_node = ele->n_sides() - i - 1;
340  uint old_opp_node = new_to_old_node[new_opp_node];
341  uint old_iside = ele->n_sides() - old_opp_node - 1;
342  old_to_new_side[old_iside] = i;
343  }
344 
345  // pressure on edges
346  // unsigned int lid = ele->n_sides() + 1;
347  for (unsigned int i = 0; i < ele->n_sides(); i++) {
348  uint new_lid = ele->n_sides() + 1 + old_to_new_side[i];
349  ss << eq_data->full_solution.get(indices[new_lid]) << " ";
350  }
351  // fluxes on sides
352  for (unsigned int i = 0; i < ele->n_sides(); i++) {
353  uint new_iside = old_to_new_side[i];
354  ss << eq_data->full_solution.get(indices[new_iside]) << " ";
355  }
356 
357  // remove last white space
358  string line = ss.str();
359  raw_output_file << line.substr(0, line.size()-1) << endl;
360  }
361 
362  raw_output_file << "$EndFlowField\n" << endl;
363 }
364 
365 
367 #include "fem/fe_p.hh"
368 #include "fem/fe_values.hh"
369 #include "fields/field_python.hh"
370 #include "fields/field_values.hh"
371 
373 
374 /*
375 * Calculate approximation of L2 norm for:
376  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
377  * 2) difference between RT velocities and analytical solution
378  * 3) difference of divergence
379  * */
380 
382  FEValues<3> &fe_values, FEValues<3> &fv_rt,
383  ExactSolution &anal_sol, DarcyFlowMHOutput::DiffData &result) {
384 
385  ASSERT_DBG( fe_values.dim() == fv_rt.dim());
386  unsigned int dim = fe_values.dim();
387 
388  ElementAccessor<3> ele = dh_cell.elm();
389  fv_rt.reinit(ele);
390  fe_values.reinit(ele);
391 
392  double conductivity = result.eq_fields_->conductivity.value(ele.centre(), ele );
393  double cross = result.eq_fields_->cross_section.value(ele.centre(), ele );
394 
395 
396  // get coefficients on the current element
397  vector<double> fluxes(dim+1);
398 // vector<double> pressure_traces(dim+1);
399 
400  for (unsigned int li = 0; li < ele->n_sides(); li++) {
401  fluxes[li] = diff_data.eq_data_->full_solution.get( dh_cell.get_loc_dof_indices()[li] );
402 // pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
403  }
404  const uint ndofs = dh_cell.n_dofs();
405  // TODO: replace with DHCell getter when available for FESystem component
406  double pressure_mean = diff_data.eq_data_->full_solution.get( dh_cell.get_loc_dof_indices()[ndofs/2] );
407 
408  arma::vec analytical(5);
409  arma::vec3 flux_in_q_point;
410  arma::vec3 anal_flux;
411 
412  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
413 
414  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
415  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
416  double mean_x_squared=0;
417  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
418  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
419  {
420  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
421  * arma::dot( *ele.node(i_node), *ele.node(j_node));
422  }
423 
424  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
425  arma::vec3 q_point = fe_values.point(i_point);
426 
427 
428  analytical = anal_sol.value(q_point, ele );
429  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
430 
431  // compute postprocesed pressure
432  diff = 0;
433  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
434  unsigned int oposite_node = 0;
435  switch (dim) {
436  case 1: oposite_node = RefElement<1>::oposite_node(i_shape); break;
437  case 2: oposite_node = RefElement<2>::oposite_node(i_shape); break;
438  case 3: oposite_node = RefElement<3>::oposite_node(i_shape); break;
439  default: ASSERT(false)(dim).error("Unsupported FE dimension."); break;
440  }
441 
442  diff += fluxes[ i_shape ] *
443  ( arma::dot( q_point, q_point )/ 2
444  - mean_x_squared / 2
445  - arma::dot( q_point, *ele.node(oposite_node) )
446  + arma::dot( ele.centre(), *ele.node(oposite_node) )
447  );
448  }
449 
450  diff = - (1.0 / conductivity) * diff / dim / ele.measure() / cross + pressure_mean ;
451  diff = ( diff - analytical[0]);
452  pressure_diff += diff * diff * fe_values.JxW(i_point);
453 
454 
455  // velocity difference
456  flux_in_q_point.zeros();
457  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
458  flux_in_q_point += fluxes[ i_shape ]
459  * fv_rt.vector_view(0).value(i_shape, i_point)
460  / cross;
461  }
462 
463  flux_in_q_point -= anal_flux;
464  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
465 
466  // divergence diff
467  diff = 0;
468  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
469  diff = ( diff / ele.measure() / cross - analytical[4]);
470  divergence_diff += diff * diff * fe_values.JxW(i_point);
471 
472  }
473 
474  // DHCell constructed with diff fields DH, get DOF indices of actual element
475  DHCellAccessor sub_dh_cell = dh_cell.cell_with_other_dh(result.dh_.get());
476  IntIdx idx = sub_dh_cell.get_loc_dof_indices()[0];
477 
478  auto velocity_data = result.vel_diff_ptr->vec();
479  velocity_data.set( idx, sqrt(velocity_diff) );
480  result.velocity_error[dim-1] += velocity_diff;
481  if (dim == 2 && result.velocity_mask.size() != 0 ) {
482  result.mask_vel_error += (result.velocity_mask[ ele.idx() ])? 0 : velocity_diff;
483  }
484 
485  auto pressure_data = result.pressure_diff_ptr->vec();
486  pressure_data.set( idx, sqrt(pressure_diff) );
487  result.pressure_error[dim-1] += pressure_diff;
488 
489  auto div_data = result.div_diff_ptr->vec();
490  div_data.set( idx, sqrt(divergence_diff) );
491  result.div_error[dim-1] += divergence_diff;
492 
493 }
494 
496 : order(4),
497  quad(QGauss::make_array(order)),
498  fe_p1(0), fe_p0(0),
499  fe_rt( )
500 {
502  fe_values = mixed_fe_values(quad, fe_p0, flags);
503  fv_rt = mixed_fe_values(quad, fe_rt, flags);
504 }
505 
506 
508  DebugOut() << "l2 norm output\n";
509  ofstream os( FilePath("solution_error", FilePath::output_file) );
510 
511  FilePath source_file( "analytical_module.py", FilePath::input_file);
512  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
513  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
514 
515  ExactSolution anal_sol_2d(5);
516  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
517 
518  ExactSolution anal_sol_3d(5);
519  anal_sol_3d.set_python_field_from_file( source_file, "all_values_3d");
520 
522  for(unsigned int j=0; j<3; j++){
523  diff_data.pressure_error[j] = 0;
524  diff_data.velocity_error[j] = 0;
525  diff_data.div_error[j] = 0;
526  }
527 
528  //diff_data.ele_flux = &( ele_flux );
529 
530  for (DHCellAccessor dh_cell : diff_data.eq_data_->dh_->own_range()) {
531 
532  switch (dh_cell.dim()) {
533  case 1:
534  l2_diff_local( dh_cell, fe_data.fe_values[1], fe_data.fv_rt[1], anal_sol_1d, diff_data);
535  break;
536  case 2:
537  l2_diff_local( dh_cell, fe_data.fe_values[2], fe_data.fv_rt[2], anal_sol_2d, diff_data);
538  break;
539  case 3:
540  l2_diff_local( dh_cell, fe_data.fe_values[3], fe_data.fv_rt[3], anal_sol_3d, diff_data);
541  break;
542  }
543  }
544 
545  // square root for L2 norm
546  for(unsigned int j=0; j<3; j++){
549  diff_data.div_error[j] = sqrt(diff_data.div_error[j]);
550  }
552 
553  os << "l2 norm output\n\n"
554  << "pressure error 1d: " << diff_data.pressure_error[0] << endl
555  << "pressure error 2d: " << diff_data.pressure_error[1] << endl
556  << "pressure error 3d: " << diff_data.pressure_error[2] << endl
557  << "velocity error 1d: " << diff_data.velocity_error[0] << endl
558  << "velocity error 2d: " << diff_data.velocity_error[1] << endl
559  << "velocity error 3d: " << diff_data.velocity_error[2] << endl
560  << "masked velocity error 2d: " << diff_data.mask_vel_error <<endl
561  << "div error 1d: " << diff_data.div_error[0] << endl
562  << "div error 2d: " << diff_data.div_error[1] << endl
563  << "div error 3d: " << diff_data.div_error[2];
564 }
565 
566 
FieldCommon::units
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
Definition: field_common.hh:152
LimitSide::right
@ right
DarcyFlowMHOutput::OutputSpecificFields::OutputSpecificFields
OutputSpecificFields()
Definition: darcy_flow_mh_output.cc:113
OutputTime::create_output_stream
static std::shared_ptr< OutputTime > create_output_stream(const std::string &equation_name, const Input::Record &in_rec, const std::shared_ptr< TimeUnitConversion > &time_unit_conv)
This method delete all object instances of class OutputTime stored in output_streams vector.
Definition: output_time.cc:187
LocDofVec
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
FEValues::point
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:236
Input::Type::Bool
Class for declaration of the input of type Bool.
Definition: type_base.hh:452
DarcyFlowMHOutput::prepare_specific_output
virtual void prepare_specific_output(Input::Record in_rec)
Definition: darcy_flow_mh_output.cc:191
Armor::vec
ArmaVec< double, N > vec
Definition: armor.hh:885
fe_rt.hh
Definitions of Raviart-Thomas finite elements.
DarcyFlowMHOutput::DiffData::eq_data_
DarcyMH::EqData * eq_data_
Definition: darcy_flow_mh_output.hh:178
DarcyFlowMHOutput::output
void output()
Calculate values for output.
Definition: darcy_flow_mh_output.cc:246
FEValues::vector_view
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:277
UnitSI::dimensionless
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
FEValues::n_points
unsigned int n_points() const
Returns the number of quadrature points.
Definition: fe_values.hh:296
FieldSet::set_time
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:245
fe_values_views.hh
ExactSolution
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
Definition: darcy_flow_mh_output.cc:372
EquationOutput::initialize
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
Definition: equation_output.cc:126
field_python.hh
TimeGovernor::get_unit_conversion
std::shared_ptr< TimeUnitConversion > get_unit_conversion() const
Getter for time unit conversion object.
Definition: time_governor.cc:810
DarcyFlowMHOutput::is_output_specific_fields
bool is_output_specific_fields
Output specific field stuff.
Definition: darcy_flow_mh_output.hh:163
DarcyFlowMHOutput::output_specific_fields
OutputSpecificFields output_specific_fields
Definition: darcy_flow_mh_output.hh:155
ASSERT
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:347
DarcyFlowMHOutput::DiffData::eq_fields_
DarcyMH::EqFields * eq_fields_
Definition: darcy_flow_mh_output.hh:177
DarcyFlowMHOutput::compute_l2_difference
void compute_l2_difference()
Definition: darcy_flow_mh_output.cc:507
FEValues::JxW
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:223
DarcyFlowMHOutput::FEData::FEData
FEData()
Definition: darcy_flow_mh_output.cc:495
Input::Record::val
const Ret val(const string &key) const
Definition: accessors_impl.hh:31
IntIdx
int IntIdx
Definition: index_types.hh:25
EquationBase::time
TimeGovernor & time()
Definition: equation.hh:149
DarcyFlowMHOutput::prepare_output
virtual void prepare_output(Input::Record in_rec)
Definition: darcy_flow_mh_output.cc:177
fe_values.hh
Class FEValues calculates finite element data on the actual cells such as shape function values,...
EquationBase::eq_fieldset
FieldSet & eq_fieldset()
Definition: equation.hh:204
FilePath
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
ASSERT_DBG
#define ASSERT_DBG(expr)
Definition: include_fadbad.hh:28
FEValues::dim
unsigned int dim() const
Return dimension of reference space.
Definition: fe_values.hh:308
field_set.hh
update_values
@ update_values
Shape function values.
Definition: update_flags.hh:87
ElementAccessor::orig_nodes_order
auto & orig_nodes_order() const
Definition: accessors.hh:246
std::vector< uint >
ElementAccessor< 3 >
Input::Type::FileName::output
static FileName output()
The factory function for declaring type FileName for input files.
Definition: type_base.cc:531
system.hh
update_quadrature_points
@ update_quadrature_points
Transformed quadrature points.
Definition: update_flags.hh:102
DarcyFlowMHOutput::mesh_
Mesh * mesh_
Definition: darcy_flow_mh_output.hh:135
arma::vec3
Definition: doxy_dummy_defs.hh:17
FieldCommon::flags
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Definition: field_common.hh:191
field_fe.hh
dofhandler.hh
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
uint
unsigned int uint
Definition: mh_dofhandler.hh:101
fmt::format
std::string format(CStringRef format_str, ArgList args)
Definition: format.h:3141
TimeStep::use_fparser_
bool use_fparser_
Definition: time_governor.hh:235
darcy_flow_lmh.hh
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
FilePath::open_stream
void open_stream(Stream &stream) const
Definition: file_path.cc:211
DarcyFlowMHOutput::output_stream
std::shared_ptr< OutputTime > output_stream
Definition: darcy_flow_mh_output.hh:157
flow
Definition: output_msh.hh:24
index_types.hh
DarcyMH::EqData::dh_
std::shared_ptr< DOFHandlerMultiDim > dh_
full DOF handler represents DOFs of sides, elements and edges
Definition: darcy_flow_mh.hh:218
DarcyFlowMHOutput::OutputSpecificFields::pressure_diff
Field< 3, FieldValue< 3 >::Scalar > pressure_diff
Definition: darcy_flow_mh_output.hh:93
FieldFlag::equation_external_output
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:58
DarcyFlowMHOutput::OutputFields
Standard quantities for output in DarcyFlowMH.
Definition: darcy_flow_mh_output.hh:78
DarcyFlowMHOutput::OutputFields::region_id
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: darcy_flow_mh_output.hh:84
MPI_Comm_rank
#define MPI_Comm_rank
Definition: mpi.h:236
FEValues< 3 >
EquationOutput::make_output_type
const Input::Type::Instance & make_output_type(const string &equation_name, const string &aditional_description="")
Definition: equation_output.cc:109
fe_p.hh
Definitions of basic Lagrangean finite elements with polynomial shape functions.
DarcyMH::EqData
Definition: darcy_flow_mh.hh:204
DarcyFlowMHOutput::get_input_type
static const Input::Type::Instance & get_input_type(FieldSet &eq_data, const std::string &equation_name)
Definition: darcy_flow_mh_output.cc:62
mixed_fe_values
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
Definition: fe_values.cc:569
Input::Type::Default
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
MPI_Comm_size
#define MPI_Comm_size
Definition: mpi.h:235
DarcyMH::EqFields::cross_section
Field< 3, FieldValue< 3 >::Scalar > cross_section
Definition: darcy_flow_mh.hh:171
DarcyFlowInterface::NoMortar
@ NoMortar
Definition: darcy_flow_interface.hh:31
DarcyFlowMHOutput::DiffData::div_diff_ptr
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > div_diff_ptr
Definition: darcy_flow_mh_output.hh:172
accessors.hh
DarcyFlowMHOutput::DiffData::velocity_error
double velocity_error[3]
Definition: darcy_flow_mh_output.hh:165
DarcyFlowMHOutput::output_fields
OutputFields output_fields
Definition: darcy_flow_mh_output.hh:154
DarcyFlowMHOutput::DiffData::mask_vel_error
double mask_vel_error
Definition: darcy_flow_mh_output.hh:166
DHCellAccessor::cell_with_other_dh
DHCellAccessor cell_with_other_dh(const DOFHandlerMultiDim *dh) const
Create new accessor with same local idx and given DOF handler. Actual and given DOF handler must be c...
Definition: dh_cell_accessor.hh:135
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
FieldFlag::equation_result
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
sys_profiler.hh
DarcyMH::EqFields::conductivity
Field< 3, FieldValue< 3 >::Scalar > conductivity
Definition: darcy_flow_mh.hh:170
darcy_flow_mh_output.hh
Output class for darcy_flow_mh model.
DarcyFlowMHOutput::OutputSpecificFields
Specific quantities for output in DarcyFlowMH - error estimates etc.
Definition: darcy_flow_mh_output.hh:88
quadrature_lib.hh
Definitions of particular quadrature rules on simplices.
DarcyFlowMHOutput::DiffData::vel_diff_ptr
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > vel_diff_ptr
Definition: darcy_flow_mh_output.hh:171
mixed_mesh_intersections.hh
DarcyFlowMHOutput::FEData::quad
QGauss::array quad
Definition: darcy_flow_mh_output.hh:189
EquationOutput
Definition: equation_output.hh:44
Field::set
void set(FieldBasePtr field, double time, std::vector< std::string > region_set_names={"ALL"})
Definition: field.impl.hh:244
DarcyMH::EqFields
Definition: darcy_flow_mh.hh:143
INPUT_CATCH
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:63
DarcyFlowMHOutput::OutputSpecificFields::velocity_diff
Field< 3, FieldValue< 3 >::Scalar > velocity_diff
Definition: darcy_flow_mh_output.hh:92
DarcyMH::EqFields::field_ele_velocity
Field< 3, FieldValue< 3 >::VectorFixed > field_ele_velocity
Definition: darcy_flow_mh.hh:188
TimeGovernor::step
const TimeStep & step(int index=-1) const
Definition: time_governor.cc:756
output_time.hh
DarcyFlowMHOutput::OutputSpecificFields::div_diff
Field< 3, FieldValue< 3 >::Scalar > div_diff
Definition: darcy_flow_mh_output.hh:94
DarcyFlowMHOutput::OutputFields::OutputFields
OutputFields()
Definition: darcy_flow_mh_output.cc:86
DHCellAccessor::elm
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Definition: dh_cell_accessor.hh:71
field_values.hh
GenericField::subdomain
static auto subdomain(Mesh &mesh) -> IndexField
Definition: generic_field.impl.hh:58
UnitSI
Class for representation SI units of Fields.
Definition: unit_si.hh:40
DarcyFlowMHOutput::FEData::fv_rt
std::vector< FEValues< 3 > > fv_rt
Definition: darcy_flow_mh_output.hh:199
DarcyFlowMHOutput::DiffData::pressure_diff_ptr
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > pressure_diff_ptr
Definition: darcy_flow_mh_output.hh:170
RefElement::oposite_node
static unsigned int oposite_node(unsigned int sid)
Definition: ref_element.cc:209
EquationOutput::get_input_type
static Input::Type::Record & get_input_type()
Definition: equation_output.cc:26
Input::Record::opt_val
bool opt_val(const string &key, Ret &value) const
Definition: accessors_impl.hh:107
DarcyFlowMHOutput::DarcyFlowMHOutput
DarcyFlowMHOutput(DarcyFlowInterface *flow, Input::Record in_rec)
Definition: darcy_flow_mh_output.cc:127
Input::Type::Instance
Helper class that stores data of generic types.
Definition: type_generic.hh:89
Input::Type::Record::declare_key
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
DarcyFlowMHOutput::fe_data
FEData fe_data
Definition: darcy_flow_mh_output.hh:202
DarcyMH::EqData::full_solution
VectorMPI full_solution
Definition: darcy_flow_mh.hh:237
DarcyFlowMHOutput::DiffData::div_error
double div_error[3]
Definition: darcy_flow_mh_output.hh:165
DarcyFlowMHOutput::get_input_type_specific
static const Input::Type::Instance & get_input_type_specific()
Definition: darcy_flow_mh_output.cc:69
FieldSet
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:159
mesh.h
DHCellAccessor::n_dofs
unsigned int n_dofs() const
Return number of dofs on given cell.
Definition: dh_cell_accessor.hh:420
DarcyFlowMHOutput::output_internal_flow_data
void output_internal_flow_data()
Definition: darcy_flow_mh_output.cc:287
Input::Record::find
Iterator< Ret > find(const string &key) const
Definition: accessors_impl.hh:91
ElementAccessor::node
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:233
Input::Type::Record::close
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
generic_field.hh
DarcyFlowMHOutput::DiffData
Definition: darcy_flow_mh_output.hh:164
ElementAccessor::input_id
unsigned int input_id() const
Return the element ID in the input mesh. Should be only used for special output.
Definition: accessors.hh:223
Input::Type
Definition: balance.hh:41
partitioning.hh
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
Field::value
virtual const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:449
QGauss
Symmetric Gauss-Legendre quadrature formulae on simplices.
Definition: quadrature_lib.hh:34
observe.hh
GenericField::region_id
static auto region_id(Mesh &mesh) -> IndexField
Definition: generic_field.impl.hh:39
DarcyFlowMHOutput::raw_output_file
ofstream raw_output_file
Raw data output file.
Definition: darcy_flow_mh_output.hh:160
EquationOutput::output
void output(TimeStep step)
Definition: equation_output.cc:245
FieldPython
Definition: field_python.hh:46
MixedMeshIntersections::intersection_storage12_
std::vector< IntersectionLocal< 1, 2 > > intersection_storage12_
Stores 1D-2D intersections.
Definition: mixed_mesh_intersections.hh:74
DarcyFlowMHOutput::diff_data
struct DarcyFlowMHOutput::DiffData diff_data
Element::n_sides
unsigned int n_sides() const
Definition: elements.h:131
Input::Type::Record::copy_keys
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
DHCellAccessor
Cell accessor allow iterate over DOF handler cells.
Definition: dh_cell_accessor.hh:43
DarcyFlowMHOutput::DiffData::velocity_mask
std::vector< int > velocity_mask
Definition: darcy_flow_mh_output.hh:176
Input::Array
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
global_defs.h
Global macros to enhance readability and debugging, general constants.
DarcyFlowMHOutput::FEData::fe_values
std::vector< FEValues< 3 > > fe_values
Definition: darcy_flow_mh_output.hh:195
FieldSet::set_mesh
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:281
WarningOut
#define WarningOut()
Macro defining 'warning' record of log.
Definition: logger.hh:278
DarcyFlowMHOutput::darcy_flow
DarcyFlowInterface * darcy_flow
Definition: darcy_flow_mh_output.hh:134
FilePath::input_file
@ input_file
Definition: file_path.hh:68
DarcyFlowMHOutput::FEData::fe_p0
MixedPtr< FE_P_disc > fe_p0
Definition: darcy_flow_mh_output.hh:194
fe_value_handler.hh
DarcyFlowMHOutput::compute_errors_
bool compute_errors_
Specific experimental error computing.
Definition: darcy_flow_mh_output.hh:138
DarcyLMH
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
Definition: darcy_flow_lmh.hh:128
IntersectionLocal< 1, 2 >
darcy_flow_mh.hh
mixed-hybrid model of linear Darcy flow, possibly unsteady.
update_JxW_values
@ update_JxW_values
Transformed quadrature weights.
Definition: update_flags.hh:114
EquationOutput::make_output_type_from_record
const Input::Type::Instance & make_output_type_from_record(Input::Type::Record &in_rec, const string &equation_name, const string &aditional_description="")
Definition: equation_output.cc:114
MPI_COMM_WORLD
#define MPI_COMM_WORLD
Definition: mpi.h:123
DHCellAccessor::get_loc_dof_indices
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
Definition: dh_cell_accessor.hh:88
assembly_lmh_old.hh
DarcyMH::EqData::mortar_method_
MortarMethod mortar_method_
Definition: darcy_flow_mh.hh:225
FieldPython::value
virtual const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_python.impl.hh:157
ElementAccessor::idx
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
Definition: accessors.hh:218
FilePath::output_file
@ output_file
Definition: file_path.hh:69
UpdateFlags
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:67
DarcyFlowMHOutput::DiffData::dh_
std::shared_ptr< SubDOFHandlerMultiDim > dh_
Definition: darcy_flow_mh_output.hh:174
VectorMPI::get
double get(unsigned int pos) const
Return value on given position.
Definition: vector_mpi.hh:105
FEValues::reinit
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:536
DarcyFlowMHOutput::OutputFields::subdomain
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: darcy_flow_mh_output.hh:83
FieldPython::set_python_field_from_file
void set_python_field_from_file(const FilePath &file_name, const string &func_name)
Definition: field_python.impl.hh:109
DarcyFlowMHOutput::~DarcyFlowMHOutput
virtual ~DarcyFlowMHOutput()
Definition: darcy_flow_mh_output.cc:234
MeshBase::n_elements
unsigned int n_elements() const
Definition: mesh.h:111
FieldCommon::description
FieldCommon & description(const string &description)
Definition: field_common.hh:127
DarcyFlowMHOutput::DiffData::pressure_error
double pressure_error[3]
Definition: darcy_flow_mh_output.hh:165
Mesh::mixed_intersections
MixedMeshIntersections & mixed_intersections()
Definition: mesh.cc:849
DebugOut
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
ASSERT_PTR
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
DarcyFlowMHOutput::l2_diff_local
void l2_diff_local(DHCellAccessor dh_cell, FEValues< 3 > &fe_values, FEValues< 3 > &fv_rt, FieldPython< 3, FieldValue< 3 >::Vector > &anal_sol, DiffData &result)
Computes L2 error on an element.
Definition: darcy_flow_mh_output.cc:381
Input::Type::Default::optional
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
DarcyMH
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
Definition: darcy_flow_mh.hh:125
START_TIMER
#define START_TIMER(tag)
Starts a timer with specified tag.
Definition: sys_profiler.hh:115
DarcyFlowMHOutput::FEData::fe_rt
MixedPtr< FE_RT0 > fe_rt
Definition: darcy_flow_mh_output.hh:198
DarcyFlowInterface
Definition: darcy_flow_interface.hh:15
assembly_mh_old.hh
ElementAccessor::measure
double measure() const
Computes the measure of the element.
Definition: accessors_impl.hh:67
Element::n_nodes
unsigned int n_nodes() const
Definition: elements.h:125
FieldCommon::name
FieldCommon & name(const string &name)
Definition: field_common.hh:120
range_wrapper.hh
Implementation of range helper class.
ElementAccessor::centre
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Definition: accessors_impl.hh:89
MessageOut
#define MessageOut()
Macro defining 'message' record of log.
Definition: logger.hh:275
intersection_local.hh
Classes with algorithms for computation of intersections of meshes.
node_accessor.hh
TimeGovernor::t
double t() const
Definition: time_governor.hh:542