Flow123d  JS_before_hm-1626-gde32303
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.hh"
30 #include "flow/assembly_lmh.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 
80  return output_fields.make_output_type_from_record(rec,
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.data_ = nullptr;
194  if(DarcyMH* d = dynamic_cast<DarcyMH*>(darcy_flow))
195  {
196  diff_data.data_ = d->data_.get();
197  }
198  else if(DarcyLMH* d = dynamic_cast<DarcyLMH*>(darcy_flow))
199  {
200  diff_data.data_ = d->data_.get();
201  }
203 
204  { // init DOF handlers represents element DOFs
205  uint p_elem_component = 1;
206  diff_data.dh_ = std::make_shared<SubDOFHandlerMultiDim>(diff_data.data_->dh_, p_elem_component);
207  }
208 
209  // mask 2d elements crossing 1d
213  diff_data.velocity_mask[ isec.bulk_ele_idx() ]++;
214  }
215  }
216 
218 
219  diff_data.vel_diff_ptr = create_field_fe<3, FieldValue<3>::Scalar>(diff_data.dh_);
221  diff_data.pressure_diff_ptr = create_field_fe<3, FieldValue<3>::Scalar>(diff_data.dh_);
223  diff_data.div_diff_ptr = create_field_fe<3, FieldValue<3>::Scalar>(diff_data.dh_);
225 
228 }
229 
231 {}
232 
233 
234 
235 
236 
237 //=============================================================================
238 // CONVERT SOLUTION, CALCULATE BALANCES, ETC...
239 //=============================================================================
240 
241 
243 {
244  START_TIMER("Darcy fields output");
245 
246  {
247  START_TIMER("post-process output fields");
248 
249  {
251  }
252  }
253 
254  {
255  START_TIMER("evaluate output fields");
258  }
259 
260  if (compute_errors_)
261  {
262  START_TIMER("compute specific output fields");
264  }
265 
267  {
268  START_TIMER("evaluate output fields");
271  }
272 
273 
274 }
275 
276 
277 
278 /*
279  * Output of internal flow data.
280  */
282 {
283  START_TIMER("DarcyFlowMHOutput::output_internal_flow_data");
284 
285  if (! raw_output_file.is_open()) return;
286 
287  //char dbl_fmt[ 16 ]= "%.8g ";
288  // header
289  raw_output_file << "// fields:\n//ele_id ele_presure flux_in_barycenter[3] n_sides side_pressures[n] side_fluxes[n]\n";
290  raw_output_file << fmt::format("$FlowField\nT={}\n", darcy_flow->time().t());
291  raw_output_file << fmt::format("{}\n" , mesh_->n_elements() );
292 
293 
294  DarcyMH::EqData* data = nullptr;
295  if(DarcyMH* d = dynamic_cast<DarcyMH*>(darcy_flow))
296  {
297  data = d->data_.get();
298  }
299  else if(DarcyLMH* d = dynamic_cast<DarcyLMH*>(darcy_flow))
300  {
301  data = d->data_.get();
302  }
303  ASSERT_PTR(data);
304 
305  arma::vec3 flux_in_center;
306 
307  auto permutation_vec = data->dh_->mesh()->element_permutations();
308  for (unsigned int i_elem=0; i_elem<data->dh_->n_own_cells(); ++i_elem) {
309  ElementAccessor<3> ele(data->dh_->mesh(), permutation_vec[i_elem]);
310  DHCellAccessor dh_cell = data->dh_->cell_accessor_from_element( ele.idx() );
311  LocDofVec indices = dh_cell.get_loc_dof_indices();
312 
313  // pressure
314  raw_output_file << fmt::format("{} {} ", dh_cell.elm().index(), data->full_solution.get(indices[ele->n_sides()]));
315 
316  // velocity at element center
317  flux_in_center = data->field_ele_velocity.value(ele.centre(), ele);
318  for (unsigned int i = 0; i < 3; i++)
319  raw_output_file << flux_in_center[i] << " ";
320 
321  // number of sides
322  raw_output_file << ele->n_sides() << " ";
323 
324  // pressure on edges
325  unsigned int lid = ele->n_sides() + 1;
326  for (unsigned int i = 0; i < ele->n_sides(); i++, lid++) {
327  raw_output_file << data->full_solution.get(indices[lid]) << " ";
328  }
329  // fluxes on sides
330  for (unsigned int i = 0; i < ele->n_sides(); i++) {
331  raw_output_file << data->full_solution.get(indices[i]) << " ";
332  }
333 
334  raw_output_file << endl;
335  }
336 
337  raw_output_file << "$EndFlowField\n" << endl;
338 }
339 
340 
342 #include "fem/fe_p.hh"
343 #include "fem/fe_values.hh"
344 #include "fields/field_python.hh"
345 #include "fields/field_values.hh"
346 
348 
349 /*
350 * Calculate approximation of L2 norm for:
351  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
352  * 2) difference between RT velocities and analytical solution
353  * 3) difference of divergence
354  * */
355 
357  FEValues<3> &fe_values, FEValues<3> &fv_rt,
358  ExactSolution &anal_sol, DarcyFlowMHOutput::DiffData &result) {
359 
360  ASSERT_DBG( fe_values.dim() == fv_rt.dim());
361  unsigned int dim = fe_values.dim();
362 
363  ElementAccessor<3> ele = dh_cell.elm();
364  fv_rt.reinit(ele);
365  fe_values.reinit(ele);
366 
367  double conductivity = result.data_->conductivity.value(ele.centre(), ele );
368  double cross = result.data_->cross_section.value(ele.centre(), ele );
369 
370 
371  // get coefficients on the current element
372  vector<double> fluxes(dim+1);
373 // vector<double> pressure_traces(dim+1);
374 
375  for (unsigned int li = 0; li < ele->n_sides(); li++) {
376  fluxes[li] = diff_data.data_->full_solution.get( dh_cell.get_loc_dof_indices()[li] );
377 // pressure_traces[li] = result.dh->side_scalar( *(ele->side( li ) ) );
378  }
379  const uint ndofs = dh_cell.n_dofs();
380  // TODO: replace with DHCell getter when available for FESystem component
381  double pressure_mean = diff_data.data_->full_solution.get( dh_cell.get_loc_dof_indices()[ndofs/2] );
382 
383  arma::vec analytical(5);
384  arma::vec3 flux_in_q_point;
385  arma::vec3 anal_flux;
386 
387  double velocity_diff=0, divergence_diff=0, pressure_diff=0, diff;
388 
389  // 1d: mean_x_squared = 1/6 (v0^2 + v1^2 + v0*v1)
390  // 2d: mean_x_squared = 1/12 (v0^2 + v1^2 +v2^2 + v0*v1 + v0*v2 + v1*v2)
391  double mean_x_squared=0;
392  for(unsigned int i_node=0; i_node < ele->n_nodes(); i_node++ )
393  for(unsigned int j_node=0; j_node < ele->n_nodes(); j_node++ )
394  {
395  mean_x_squared += (i_node == j_node ? 2.0 : 1.0) / ( 6 * dim ) // multiply by 2 on diagonal
396  * arma::dot( *ele.node(i_node), *ele.node(j_node));
397  }
398 
399  for(unsigned int i_point=0; i_point < fe_values.n_points(); i_point++) {
400  arma::vec3 q_point = fe_values.point(i_point);
401 
402 
403  analytical = anal_sol.value(q_point, ele );
404  for(unsigned int i=0; i< 3; i++) anal_flux[i] = analytical[i+1];
405 
406  // compute postprocesed pressure
407  diff = 0;
408  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
409  unsigned int oposite_node = 0;
410  switch (dim) {
411  case 1: oposite_node = RefElement<1>::oposite_node(i_shape); break;
412  case 2: oposite_node = RefElement<2>::oposite_node(i_shape); break;
413  case 3: oposite_node = RefElement<3>::oposite_node(i_shape); break;
414  default: ASSERT(false)(dim).error("Unsupported FE dimension."); break;
415  }
416 
417  diff += fluxes[ i_shape ] *
418  ( arma::dot( q_point, q_point )/ 2
419  - mean_x_squared / 2
420  - arma::dot( q_point, *ele.node(oposite_node) )
421  + arma::dot( ele.centre(), *ele.node(oposite_node) )
422  );
423  }
424 
425  diff = - (1.0 / conductivity) * diff / dim / ele.measure() / cross + pressure_mean ;
426  diff = ( diff - analytical[0]);
427  pressure_diff += diff * diff * fe_values.JxW(i_point);
428 
429 
430  // velocity difference
431  flux_in_q_point.zeros();
432  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) {
433  flux_in_q_point += fluxes[ i_shape ]
434  * fv_rt.vector_view(0).value(i_shape, i_point)
435  / cross;
436  }
437 
438  flux_in_q_point -= anal_flux;
439  velocity_diff += dot(flux_in_q_point, flux_in_q_point) * fe_values.JxW(i_point);
440 
441  // divergence diff
442  diff = 0;
443  for(unsigned int i_shape=0; i_shape < ele->n_sides(); i_shape++) diff += fluxes[ i_shape ];
444  diff = ( diff / ele.measure() / cross - analytical[4]);
445  divergence_diff += diff * diff * fe_values.JxW(i_point);
446 
447  }
448 
449  // DHCell constructed with diff fields DH, get DOF indices of actual element
450  DHCellAccessor sub_dh_cell = dh_cell.cell_with_other_dh(result.dh_.get());
451  IntIdx idx = sub_dh_cell.get_loc_dof_indices()[0];
452 
453  auto velocity_data = result.vel_diff_ptr->vec();
454  velocity_data.set( idx, sqrt(velocity_diff) );
455  result.velocity_error[dim-1] += velocity_diff;
456  if (dim == 2 && result.velocity_mask.size() != 0 ) {
457  result.mask_vel_error += (result.velocity_mask[ ele.idx() ])? 0 : velocity_diff;
458  }
459 
460  auto pressure_data = result.pressure_diff_ptr->vec();
461  pressure_data.set( idx, sqrt(pressure_diff) );
462  result.pressure_error[dim-1] += pressure_diff;
463 
464  auto div_data = result.div_diff_ptr->vec();
465  div_data.set( idx, sqrt(divergence_diff) );
466  result.div_error[dim-1] += divergence_diff;
467 
468 }
469 
471 : order(4),
472  quad(QGauss::make_array(order)),
473  fe_p1(0), fe_p0(0),
474  fe_rt( )
475 {
477  fe_values = mixed_fe_values(quad, fe_p0, flags);
478  fv_rt = mixed_fe_values(quad, fe_rt, flags);
479 }
480 
481 
483  DebugOut() << "l2 norm output\n";
484  ofstream os( FilePath("solution_error", FilePath::output_file) );
485 
486  FilePath source_file( "analytical_module.py", FilePath::input_file);
487  ExactSolution anal_sol_1d(5); // components: pressure, flux vector 3d, divergence
488  anal_sol_1d.set_python_field_from_file( source_file, "all_values_1d");
489 
490  ExactSolution anal_sol_2d(5);
491  anal_sol_2d.set_python_field_from_file( source_file, "all_values_2d");
492 
493  ExactSolution anal_sol_3d(5);
494  anal_sol_3d.set_python_field_from_file( source_file, "all_values_3d");
495 
497  for(unsigned int j=0; j<3; j++){
498  diff_data.pressure_error[j] = 0;
499  diff_data.velocity_error[j] = 0;
500  diff_data.div_error[j] = 0;
501  }
502 
503  //diff_data.ele_flux = &( ele_flux );
504 
505  for (DHCellAccessor dh_cell : diff_data.data_->dh_->own_range()) {
506 
507  switch (dh_cell.dim()) {
508  case 1:
509  l2_diff_local( dh_cell, fe_data.fe_values[1], fe_data.fv_rt[1], anal_sol_1d, diff_data);
510  break;
511  case 2:
512  l2_diff_local( dh_cell, fe_data.fe_values[2], fe_data.fv_rt[2], anal_sol_2d, diff_data);
513  break;
514  case 3:
515  l2_diff_local( dh_cell, fe_data.fe_values[3], fe_data.fv_rt[3], anal_sol_3d, diff_data);
516  break;
517  }
518  }
519 
520  // square root for L2 norm
521  for(unsigned int j=0; j<3; j++){
524  diff_data.div_error[j] = sqrt(diff_data.div_error[j]);
525  }
527 
528  os << "l2 norm output\n\n"
529  << "pressure error 1d: " << diff_data.pressure_error[0] << endl
530  << "pressure error 2d: " << diff_data.pressure_error[1] << endl
531  << "pressure error 3d: " << diff_data.pressure_error[2] << endl
532  << "velocity error 1d: " << diff_data.velocity_error[0] << endl
533  << "velocity error 2d: " << diff_data.velocity_error[1] << endl
534  << "velocity error 3d: " << diff_data.velocity_error[2] << endl
535  << "masked velocity error 2d: " << diff_data.mask_vel_error <<endl
536  << "div error 1d: " << diff_data.div_error[0] << endl
537  << "div error 2d: " << diff_data.div_error[1] << endl
538  << "div error 3d: " << diff_data.div_error[2];
539 }
540 
541 
TimeGovernor & time()
Definition: equation.hh:149
void set_python_field_from_file(const FilePath &file_name, const string &func_name)
Shape function values.
Definition: update_flags.hh:87
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > vel_diff_ptr
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
FieldSet & eq_fieldset()
Definition: equation.hh:202
static auto subdomain(Mesh &mesh) -> IndexField
struct DarcyFlowMHOutput::DiffData diff_data
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:159
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
Output class for darcy_flow_mh model.
unsigned int n_nodes() const
Definition: elements.h:126
virtual void prepare_specific_output(Input::Record in_rec)
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:277
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
ArmaVec< double, N > vec
Definition: armor.hh:885
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
unsigned int uint
Classes with algorithms for computation of intersections of meshes.
std::shared_ptr< SubDOFHandlerMultiDim > dh_
virtual void prepare_output(Input::Record in_rec)
OutputSpecificFields output_specific_fields
Specific quantities for output in DarcyFlowMH - error estimates etc.
MixedMeshIntersections & mixed_intersections()
Definition: mesh.cc:904
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Class for declaration of the input of type Bool.
Definition: type_base.hh:452
int IntIdx
Definition: index_types.hh:25
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:267
std::shared_ptr< DOFHandlerMultiDim > dh_
full DOF handler represents DOFs of sides, elements and edges
void output(TimeStep step)
const Input::Type::Instance & make_output_type_from_record(Input::Type::Record &in_rec, const string &equation_name, const string &aditional_description="")
const Mesh * mesh() const
Returns pointer to mesh.
Definition: field_set.hh:370
virtual unsigned int n_elements() const
Returns count of boundary or bulk elements.
Definition: mesh.h:361
std::string format(CStringRef format_str, ArgList args)
Definition: format.h:3141
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:63
bool is_output_specific_fields
Output specific field stuff.
Standard quantities for output in DarcyFlowMH.
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
Fields computed from the mesh data.
Iterator< Ret > find(const string &key) const
Cell accessor allow iterate over DOF handler cells.
Helper class that stores data of generic types.
Definition: type_generic.hh:89
Class FEValues calculates finite element data on the actual cells such as shape function values...
FieldPython< 3, FieldValue< 3 >::Vector > ExactSolution
static const Input::Type::Instance & get_input_type(FieldSet &eq_data, const std::string &equation_name)
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
const TimeStep & step(int index=-1) const
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
double t() const
Field< 3, FieldValue< 3 >::Scalar > subdomain
unsigned int n_dofs() const
Return number of dofs on given cell.
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:541
std::shared_ptr< OutputTime > output_stream
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:200
Field< 3, FieldValue< 3 >::Scalar > pressure_diff
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
static unsigned int oposite_node(unsigned int sid)
Definition: ref_element.cc:327
Symmetric Gauss-Legendre quadrature formulae on simplices.
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Field< 3, FieldValue< 3 >::Scalar > velocity_diff
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
void open_stream(Stream &stream) const
Definition: file_path.cc:211
bool opt_val(const string &key, Ret &value) const
void set(FieldBasePtr field, double time, std::vector< std::string > region_set_names={"ALL"})
Definition: field.impl.hh:258
Transformed quadrature points.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Global macros to enhance readability and debugging, general constants.
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
static const Input::Type::Instance & get_input_type_specific()
unsigned int index() const
Definition: accessors.hh:191
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:291
const Ret val(const string &key) const
unsigned int n_sides() const
Definition: elements.h:132
static auto region_id(Mesh &mesh) -> IndexField
#define START_TIMER(tag)
Starts a timer with specified tag.
unsigned int dim() const
Return dimension of reference space.
Definition: fe_values.hh:308
DarcyFlowMHOutput(DarcyFlowInterface *flow, Input::Record in_rec)
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:452
Field< 3, FieldValue< 3 >::Scalar > div_diff
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:503
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > pressure_diff_ptr
ofstream raw_output_file
Raw data output file.
double measure() const
Computes the measure of the element.
VectorMPI full_solution
#define MPI_Comm_size
Definition: mpi.h:235
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
std::vector< FEValues< 3 > > fe_values
std::shared_ptr< TimeUnitConversion > get_unit_conversion() const
Getter for time unit conversion object.
#define MPI_Comm_rank
Definition: mpi.h:236
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
double get(unsigned int pos) const
Return value on given position.
Definition: vector_mpi.hh:105
FieldCommon & description(const string &description)
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
DarcyFlowInterface * darcy_flow
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
Field< 3, FieldValue< 3 >::Scalar > region_id
bool compute_errors_
Specific experimental error computing.
static Input::Type::Record & get_input_type()
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:164
#define ASSERT_DBG(expr)
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:270
FieldCommon & name(const string &name)
#define MPI_COMM_WORLD
Definition: mpi.h:123
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.
std::vector< FEValues< 3 > > fv_rt
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > div_diff_ptr
mixed-hybrid model of linear Darcy flow, possibly unsteady.
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:274
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...
Record type proxy class.
Definition: type_record.hh:182
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 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
Class for representation SI units of Fields.
Definition: unit_si.hh:40
static FileName output()
The factory function for declaring type FileName for input files.
Definition: type_base.cc:531
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
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:276
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:181
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
Definition: fe_values.cc:575
std::vector< IntersectionLocal< 1, 2 > > intersection_storage12_
Stores 1D-2D intersections.
unsigned int n_points() const
Returns the number of quadrature points.
Definition: fe_values.hh:296
Implementation of range helper class.
Definitions of Raviart-Thomas finite elements.
MortarMethod mortar_method_
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
void output()
Calculate values for output.
Transformed quadrature weights.
Mixed-hybrid model of linear Darcy flow, possibly unsteady.