Flow123d  JS_before_hm-1003-g4e68d2c
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 
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/fe_values.hh" // for FEValues
33 #include "quadrature/quadrature_lib.hh" // for QGauss
34 #include "fields/equation_output.hh" // for EquationOutput
35 #include "fields/field.hh" // for Field
36 #include "fields/field_set.hh" // for FieldSet
37 #include "fields/field_values.hh" // for FieldValue<>::Scalar, FieldV...
38 #include "fields/field_python.hh"
39 #include "la/vector_mpi.hh" // for VectorMPI
40 #include "input/type_base.hh" // for Array
41 #include "input/type_generic.hh" // for Instance
42 #include "petscvec.h" // for Vec, _p_Vec
43 #include "system/exceptions.hh" // for ExcAssertMsg::~ExcAssertMsg
44 
45 class DOFHandlerMultiDim;
46 class DarcyFlowInterface;
47 class Mesh;
48 class OutputTime;
49 namespace Input {
50  class Record;
51  namespace Type {
52  class Record;
53  }
54 }
55 
56 template<unsigned int spacedim> class FEValues;
57 template <int spacedim, class Value> class FieldFE;
58 
59 /**
60  * Actually this class only collect former code from postprocess.*
61  * This code depends on values stored in mesh and has to be changed to use fields or other data provided by
62  * interface to darcy_flow. We have to relay on the interface in order to allow different implementation of darcy_flow.
63  *
64  * Principal functionalities of current postprocess are:
65  * - move computed values into mesh - this can be removed after we get new output classes,
66  * other dependent code can use directly field interface
67  * - compute interpolation to nodes - this can be external and operate directly on fields
68  * - compute water balances over materials and boundary parts - is feasible through field interface
69  *
70  * Further functionality of this class should be:
71  * - prepare data for general output classes
72  *
73  */
75 public:
76 
77  /// Standard quantities for output in DarcyFlowMH.
78  class OutputFields : public EquationOutput {
79  public:
80 
81  OutputFields();
82 
85  };
86 
87  /// Specific quantities for output in DarcyFlowMH - error estimates etc.
89  public:
91 
95  };
96 
98  virtual ~DarcyFlowMHOutput();
99 
100  static const Input::Type::Instance & get_input_type(FieldSet& eq_data, const std::string &equation_name);
101  static const Input::Type::Instance & get_input_type_specific();
102 
103  /** \brief Calculate values for output. **/
104  void output();
105 
106  //const OutputFields &get_output_fields() { return output_fields; }
107 
108 
109 
110 protected:
112 
113  virtual void prepare_output(Input::Record in_rec);
114  virtual void prepare_specific_output(Input::Record in_rec);
115 
116  void output_internal_flow_data();
117 
118  /**
119  * Temporary hack.
120  * Calculate approximation of L2 norm for:
121  * 1) difference between regularized pressure and analytical solution (using FunctionPython)
122  * 2) difference between RT velocities and analytical solution
123  * 3) difference of divergence
124  *
125  * TODO:
126  * 1) implement field objects
127  * 2) implement DG_P2 finite elements
128  * 3) implement pressure postprocessing (result is DG_P2 field)
129  * 4) implement calculation of L2 norm for two field (compute the norm and values on individual elements as P0 field)
130  */
131  void compute_l2_difference();
132 
133 
136 
137  /// Specific experimental error computing.
139 
140 
141  /** Pressure head (in [m]) interpolated into nodes. Provides P1 approximation. Indexed by element-node numbering.*/
143  /** Pressure head (in [m]) in barycenters of elements (or equivalently mean pressure over every element). Indexed by element indexes in the mesh.*/
145  /** Piezo-metric head (in [m]) in barycenter of elements (or equivalently mean pressure over every element). Indexed by element indexes in the mesh.*/
147 
148  // integrals of squared differences on individual elements - error indicators, can be written out into VTK files
149  std::vector<double> l2_diff_pressure, l2_diff_velocity, l2_diff_divergence;
150 
151  std::shared_ptr<DOFHandlerMultiDim> dh_;
152  std::shared_ptr<DiscreteSpace> ds;
153 
156 
157  std::shared_ptr<OutputTime> output_stream;
158 
159  /// Raw data output file.
160  ofstream raw_output_file;
161 
162  /// Output specific field stuff
164  struct DiffData {
165  double pressure_error[3], velocity_error[3], div_error[3];
167 
168  // Temporary objects holding pointers to appropriate FieldFE
169  // TODO remove after final fix of equations
170  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar>> pressure_diff_ptr;
171  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar>> vel_diff_ptr;
172  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar>> div_diff_ptr;
173 
174  std::shared_ptr<SubDOFHandlerMultiDim> dh_;
175 
178  } diff_data;
179 
180  //MixedPtr<FE_P_disc> fe_p0;
181 
182  /// Struct containing all dim dependent FE classes needed for output
183  /// (and for computing solution error).
184  struct FEData{
185  FEData();
186 
187  const unsigned int order; // order of Gauss quadrature
190 
191  // following is used for calculation of postprocessed pressure difference
192  // and comparison to analytical solution
195 
196  // FEValues for velocity.
199  };
200 
202 
203  /// Computes L2 error on an element.
204  void l2_diff_local(DHCellAccessor dh_cell,
205  FEValues<3> &fe_values, FEValues<3> &fv_rt,
206  FieldPython<3, FieldValue<3>::Vector > &anal_sol, DiffData &result);
207 };
208 
209 
210 #endif /* DARCY_FLOW_MH_OUTPUT_HH_ */
211 
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > vel_diff_ptr
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:74
std::shared_ptr< SubDOFHandlerMultiDim > dh_
std::vector< double > l2_diff_velocity
OutputSpecificFields output_specific_fields
Specific quantities for output in DarcyFlowMH - error estimates etc.
Abstract linear system class.
Definition: balance.hh:40
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:92
bool is_output_specific_fields
Output specific field stuff.
Standard quantities for output in DarcyFlowMH.
Definition: mesh.h:78
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...
std::shared_ptr< DOFHandlerMultiDim > dh_
Field< 3, FieldValue< 3 >::Scalar > subdomain
std::shared_ptr< OutputTime > output_stream
Field< 3, FieldValue< 3 >::Scalar > pressure_diff
std::shared_ptr< DiscreteSpace > ds
Field< 3, FieldValue< 3 >::Scalar > velocity_diff
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:151
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
std::vector< FEValues< 3 > > fe_values
std::array< QGauss, 4 > array
DarcyFlowInterface * darcy_flow
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:59
std::vector< FEValues< 3 > > fv_rt
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > div_diff_ptr
const vector< unsigned int > & ElementSetRef
Definitions of Raviart-Thomas finite elements.
Definition: field.hh:60