Flow123d  release_3.0.0-885-g06276d1
darcy_flow_mh_output.hh
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.hh
15  * @brief Output class for darcy_flow_mh model.
16  * @author Jan Brezina
17  */
18 
19 #ifndef DARCY_FLOW_MH_OUTPUT_HH_
20 #define DARCY_FLOW_MH_OUTPUT_HH_
21 
22 #include <stdio.h> // for sprintf
23 #include <string.h> // for memcpy
24 #include <boost/exception/info.hpp> // for operator<<, error_info::erro...
25 #include <cmath> // for pow
26 #include <iosfwd> // for ofstream
27 #include <memory> // for shared_ptr
28 #include <vector> // for vector
29 #include <armadillo>
30 #include "fem/fe_p.hh" // for FE_P_disc
31 #include "fem/fe_rt.hh" // for FE_RT0
32 #include "fem/mapping_p1.hh" // for MappingP1
33 #include "fem/fe_values.hh" // for FEValues
34 #include "quadrature/quadrature_lib.hh" // for QGauss
35 #include "fields/equation_output.hh" // for EquationOutput
36 #include "fields/field.hh" // for Field
37 #include "fields/field_set.hh" // for FieldSet
38 #include "fields/field_values.hh" // for FieldValue<>::Scalar, FieldV...
39 #include "fields/field_python.hh"
40 #include "la/vector_mpi.hh" // for VectorMPI
41 #include "input/type_base.hh" // for Array
42 #include "input/type_generic.hh" // for Instance
43 #include "petscvec.h" // for Vec, _p_Vec
44 #include "system/exceptions.hh" // for ExcAssertMsg::~ExcAssertMsg
45 
46 class DOFHandlerMultiDim;
47 class DarcyMH;
48 class Mesh;
49 class OutputTime;
50 namespace Input {
51  class Record;
52  namespace Type {
53  class Record;
54  }
55 }
56 
57 template<unsigned int dim, unsigned int spacedim> class FEValues;
58 template <int spacedim, class Value> class FieldFE;
59 
60 /**
61  * Actually this class only collect former code from postprocess.*
62  * This code depends on values stored in mesh and has to be changed to use fields or other data provided by
63  * interface to darcy_flow. We have to relay on the interface in order to allow different implementation of darcy_flow.
64  *
65  * Principal functionalities of current postprocess are:
66  * - move computed values into mesh - this can be removed after we get new output classes,
67  * other dependent code can use directly field interface
68  * - compute interpolation to nodes - this can be external and operate directly on fields
69  * - compute water balances over materials and boundary parts - is feasible through field interface
70  *
71  * Further functionality of this class should be:
72  * - prepare data for general output classes
73  *
74  */
76 public:
77 
78  /// Standard quantities for output in DarcyFlowMH.
79  class OutputFields : public EquationOutput {
80  public:
81 
82  OutputFields();
83 
90  };
91 
92  /// Specific quantities for output in DarcyFlowMH - error estimates etc.
94  public:
96 
100  };
101 
103  virtual ~DarcyFlowMHOutput();
104 
105  static const Input::Type::Instance & get_input_type();
106  static const Input::Type::Instance & get_input_type_specific();
107 
108  /** \brief Calculate values for output. **/
109  void output();
110 
111  //const OutputFields &get_output_fields() { return output_fields; }
112 
113 
114 
115 protected:
117 
118  virtual void prepare_output(Input::Record in_rec);
119  virtual void prepare_specific_output(Input::Record in_rec);
120 
121  void make_side_flux();
122  void make_element_scalar(ElementSetRef element_indices);
123 
124  /** Computes fluxes at the barycenters of elements.
125  * TODO:
126  * We use FEValues to get fluxes at the barycenters of elements,
127  * but we still use MHDofHandler. Once we are able to make output routines
128  * parallel, we can use simply FieldFE for velocity here.
129  */
130  void make_element_vector(ElementSetRef element_indices);
131 
132  //void make_sides_scalar();
133  /**
134  * \brief Calculate nodes scalar,
135  * store it in double* node_scalars instead of node->scalar
136  * */
137  void make_node_scalar_param(ElementSetRef element_indices);
138  //void make_node_scalar();
139  void make_corner_scalar(vector<double> &node_scalar);
140  //void make_neighbour_flux();
141  void output_internal_flow_data();
142 
143  /**
144  * Temporary hack.
145  * Calculate approximation of L2 norm for:
146  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
147  * 2) difference between RT velocities and analytical solution
148  * 3) difference of divergence
149  *
150  * TODO:
151  * 1) implement field objects
152  * 2) implement DG_P2 finite elements
153  * 3) implement pressure postprocessing (result is DG_P2 field)
154  * 4) implement calculation of L2 norm for two field (compute the norm and values on individual elements as P0 field)
155  */
156  void compute_l2_difference();
157 
158 
161 
162  /// Specific experimental error computing.
164 
165 
166  /** Pressure head (in [m]) interpolated into nodes. Provides P1 approximation. Indexed by element-node numbering.*/
168  /** Pressure head (in [m]) in barycenters of elements (or equivalently mean pressure over every element). Indexed by element indexes in the mesh.*/
170  /** Piezo-metric head (in [m]) in barycenter of elements (or equivalently mean pressure over every element). Indexed by element indexes in the mesh.*/
172 
173  /** Average flux in barycenter of every element. Indexed as elements in the mesh. */
174  // TODO: Definitely we need more general (templated) implementation of Output that accept arbitrary containers. So
175  // that we can pass there directly vector< arma:: vec3 >
177 
178  // Temporary objects holding pointers to appropriate FieldFE
179  // TODO remove after final fix of equations
180  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar>> ele_pressure_ptr; ///< Field of pressure head in barycenters of elements.
181  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar>> ele_piezo_head_ptr; ///< Field of piezo-metric head in barycenters of elements.
182  std::shared_ptr<FieldFE<3, FieldValue<3>::VectorFixed>> ele_flux_ptr; ///< Field of flux in barycenter of every element.
183 
184  // A vector of all element indexes
186 
187  // integrals of squared differences on individual elements - error indicators, can be written out into VTK files
188  std::vector<double> l2_diff_pressure, l2_diff_velocity, l2_diff_divergence;
189 
190  std::shared_ptr<DOFHandlerMultiDim> dh_;
191  FE_P_disc<0> fe0; //TODO temporary solution - add support of FEData<0>
192  std::shared_ptr<DiscreteSpace> ds;
193 
196 
197  std::shared_ptr<OutputTime> output_stream;
198 
199  /// Raw data output file.
200  ofstream raw_output_file;
201 
202  /// Output specific field stuff
204  struct DiffData {
205  double pressure_error[3], velocity_error[3], div_error[3];
210 
211  // Temporary objects holding pointers to appropriate FieldFE
212  // TODO remove after final fix of equations
213  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar>> pressure_diff_ptr;
214  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar>> vel_diff_ptr;
215  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar>> div_diff_ptr;
216 
217  double * solution;
218  const MH_DofHandler * dh;
219 
223  } diff_data;
224 
225 
226  /// Struct containing all dim dependent FE classes needed for output
227  /// (and for computing solution error).
228  template<int dim> struct FEData{
229  FEData();
230 
231  // we create trivial Dofhandler , for P0 elements, to get access to, FEValues on individual elements
232  // this we use to integrate our own functions - difference of postprocessed pressure and analytical solution
235 
236  const unsigned int order; // order of Gauss quadrature
238 
240 
242 
243  // FEValues for velocity.
246  };
247 
251 
252  /// Computes L2 error on an element.
253  template <int dim>
254  void l2_diff_local(ElementAccessor<3> &ele,
255  FEValues<dim,3> &fe_values, FEValues<dim,3> &fv_rt,
256  FieldPython<3, FieldValue<3>::Vector > &anal_sol, DiffData &result);
257 };
258 
259 
260 #endif /* DARCY_FLOW_MH_OUTPUT_HH_ */
261 
Field< 3, FieldValue< 3 >::VectorFixed > field_ele_flux
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > vel_diff_ptr
std::vector< double > l2_diff_velocity
OutputSpecificFields output_specific_fields
Specific quantities for output in DarcyFlowMH - error estimates etc.
Field< 3, FieldValue< 3 >::Scalar > field_node_pressure
Abstract linear system class.
Definition: balance.hh:35
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:83
bool is_output_specific_fields
Output specific field stuff.
Standard quantities for output in DarcyFlowMH.
Definition: mesh.h:80
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...
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > ele_pressure_ptr
Field of pressure head in barycenters of elements.
std::shared_ptr< DOFHandlerMultiDim > dh_
Field< 3, FieldValue< 3 >::Scalar > subdomain
std::vector< unsigned int > all_element_idx_
std::shared_ptr< OutputTime > output_stream
Field< 3, FieldValue< 3 >::Scalar > pressure_diff
std::shared_ptr< DiscreteSpace > ds
Symmetric Gauss-Legendre quadrature formulae on simplices.
Field< 3, FieldValue< 3 >::Scalar > velocity_diff
Field< 3, FieldValue< 3 >::Scalar > field_ele_pressure
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:60
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:152
Field< 3, FieldValue< 3 >::Scalar > div_diff
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > pressure_diff_ptr
ofstream raw_output_file
Raw data output file.
The class for outputting data during time.
Definition: output_time.hh:50
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > ele_piezo_head_ptr
Field of piezo-metric head in barycenters of elements.
Field< 3, FieldValue< 3 >::Scalar > region_id
bool compute_errors_
Specific experimental error computing.
Definitions of particular quadrature rules on simplices.
Calculates finite element data on the actual cell.
Definition: fe_values.hh:532
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > div_diff_ptr
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > ele_flux_ptr
Field of flux in barycenter of every element.
const vector< unsigned int > & ElementSetRef
Definitions of Raviart-Thomas finite elements.
Definition: field.hh:56
Mixed-hybrid model of linear Darcy flow, possibly unsteady.