Flow123d  release_3.0.0-1117-gaafa50c
heat_model.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 heat_model.hh
15  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
16  * @author Jan Stebel
17  */
18 
19 #ifndef HEAT_MODEL_HH_
20 #define HEAT_MODEL_HH_
21 
22 #include <boost/exception/info.hpp> // for operator<<, error_inf...
23 #include <memory> // for shared_ptr
24 #include <string> // for string
25 #include <vector> // for vector
26 #include <armadillo>
29 #include "transport/substance.hh"
30 #include "fields/field_values.hh" // for FieldValue<>::Scalar
31 #include "fields/bc_field.hh"
32 #include "fields/bc_multi_field.hh"
33 #include "fields/field.hh"
34 #include "fields/multi_field.hh"
35 #include "fields/field_set.hh"
36 #include "input/type_base.hh" // for Array
37 #include "input/type_generic.hh" // for Instance
38 #include "input/type_record.hh" // for Record
39 #include "input/type_selection.hh" // for Selection
40 
41 class MH_DofHandler;
42 class Mesh;
43 class OutputTime;
44 namespace Input { class Record; }
45 template <int spacedim> class ElementAccessor;
46 
47 /*
48 class HeatProcessBase : public EquationBase
49 {
50 public:
51  typedef HeatProcessBase FactoryBaseType;
52 
53 
54  HeatProcessBase(Mesh &mesh, const Input::Record in_rec)
55  : EquationBase(mesh, in_rec)
56  {};
57 
58  /**
59  * This method takes sequential PETSc vector of side velocities and update
60  * transport matrix. The ordering is same as ordering of sides in the mesh.
61  * We just keep the pointer, but do not destroy the object.
62  *
63  * TODO: We should pass whole velocity field object (description of base functions and dof numbering) and vector.
64  *
65  virtual void set_velocity_field(const MH_DofHandler &dh) = 0;
66 
67  /// Common specification of the input record for secondary equations.
68  static Input::Type::Abstract & get_input_type() {
69  return Input::Type::Abstract("Heat",
70  "Equation for heat transfer.")
71  .close();
72  }
73 };*/
74 
75 /*
76 class HeatNothing : public HeatProcessBase {
77 public:
78  inline HeatNothing(Mesh &mesh_in)
79  : HeatProcessBase(mesh_in, Input::Record() )
80 
81  {
82  // make module solved for ever
83  time_=new TimeGovernor();
84  time_->next_time();
85  };
86 
87  inline virtual ~HeatNothing()
88  {}
89 
90  inline void set_velocity_field(const MH_DofHandler &dh) override {};
91 
92  inline virtual void output_data() override {};
93 
94 };
95 */
96 
97 class HeatTransferModel : public AdvectionDiffusionModel, public AdvectionProcessBase {
98 public:
99 
100  class ModelEqData : public FieldSet {
101  public:
102 
103  enum Heat_bc_types {
104  bc_inflow,
105  bc_dirichlet,
106  bc_total_flux,
107  bc_diffusive_flux
108  };
109 
110  /// Type of boundary condition (see also BC_Type)
111  BCField<3, FieldValue<3>::Enum > bc_type;
112  /// Dirichlet boundary condition for temperature.
113  BCMultiField<3, FieldValue<3>::Scalar> bc_dirichlet_value;
114  /// Flux value in total/diffusive flux b.c.
115  BCField<3, FieldValue<3>::Scalar > bc_flux;
116  /// Transition coefficient in total/diffusive flux b.c.
117  BCField<3, FieldValue<3>::Scalar > bc_robin_sigma;
118  /// Initial temperature.
119  Field<3, FieldValue<3>::Scalar> init_temperature;
120  /// Porosity of solid.
121  Field<3, FieldValue<3>::Scalar> porosity;
122  /// Water content passed from Darcy flow model
123  Field<3, FieldValue<3>::Scalar> water_content;
124  /// Density of fluid.
125  Field<3, FieldValue<3>::Scalar> fluid_density;
126  /// Heat capacity of fluid.
127  Field<3, FieldValue<3>::Scalar> fluid_heat_capacity;
128  /// Heat conductivity of fluid.
129  Field<3, FieldValue<3>::Scalar> fluid_heat_conductivity;
130  /// Density of solid.
131  Field<3, FieldValue<3>::Scalar> solid_density;
132  /// Heat capacity of solid.
133  Field<3, FieldValue<3>::Scalar> solid_heat_capacity;
134  /// Heat conductivity of solid.
135  Field<3, FieldValue<3>::Scalar> solid_heat_conductivity;
136  /// Longitudal heat dispersivity.
137  Field<3, FieldValue<3>::Scalar> disp_l;
138  /// Transversal heat dispersivity.
139  Field<3, FieldValue<3>::Scalar> disp_t;
140  /// Thermal source in fluid.
141  Field<3, FieldValue<3>::Scalar> fluid_thermal_source;
142  /// Thermal source in solid.
143  Field<3, FieldValue<3>::Scalar> solid_thermal_source;
144  /// Heat exchange rate in fluid.
145  Field<3, FieldValue<3>::Scalar> fluid_heat_exchange_rate;
146  /// Heat exchange rate in solid.
147  Field<3, FieldValue<3>::Scalar> solid_heat_exchange_rate;
148  /// Reference temperature in fluid.
149  Field<3, FieldValue<3>::Scalar> fluid_ref_temperature;
150  /// Reference temperature in solid.
151  Field<3, FieldValue<3>::Scalar> solid_ref_temperature;
152 
153  /// Pointer to DarcyFlow field cross_section
154  Field<3, FieldValue<3>::Scalar > cross_section;
155 
156 
157  MultiField<3, FieldValue<3>::Scalar> output_field;
158 
159 
160 
161  ModelEqData();
162 
163  static constexpr const char * name() { return "Heat_AdvectionDiffusion"; }
164 
165  static string default_output_field() { return "\"temperature\""; }
166 
167  static const Input::Type::Selection & get_bc_type_selection();
168 
169  static IT::Selection get_output_selection();
170  };
171 
172  typedef AdvectionProcessBase FactoryBaseType;
173 
174 
175  HeatTransferModel(Mesh &mesh, const Input::Record in_rec);
176 
177  void init_from_input(const Input::Record &in_rec) override {};
178 
179  void compute_mass_matrix_coefficient(const std::vector<arma::vec3 > &point_list,
180  const ElementAccessor<3> &ele_acc,
181  std::vector<double> &mm_coef) override;
182 
183  void compute_retardation_coefficient(const std::vector<arma::vec3 > &point_list,
184  const ElementAccessor<3> &ele_acc,
185  std::vector<std::vector<double> > &ret_coef) override {};
186 
187  void compute_advection_diffusion_coefficients(const std::vector<arma::vec3 > &point_list,
188  const std::vector<arma::vec3> &velocity,
189  const ElementAccessor<3> &ele_acc,
190  std::vector<std::vector<arma::vec3> > &ad_coef,
191  std::vector<std::vector<arma::mat33> > &dif_coef) override;
192 
193  void compute_init_cond(const std::vector<arma::vec3> &point_list,
194  const ElementAccessor<3> &ele_acc,
195  std::vector<std::vector<double> > &init_values) override;
196 
197  void get_bc_type(const ElementAccessor<3> &ele_acc,
198  arma::uvec &bc_types) override;
199 
200  void get_flux_bc_data(unsigned int index,
201  const std::vector<arma::vec3> &point_list,
202  const ElementAccessor<3> &ele_acc,
203  std::vector< double > &bc_flux,
204  std::vector< double > &bc_sigma,
205  std::vector< double > &bc_ref_value) override;
206 
207  void get_flux_bc_sigma(unsigned int index,
208  const std::vector<arma::vec3> &point_list,
209  const ElementAccessor<3> &ele_acc,
210  std::vector< double > &bc_sigma) override;
211 
212  void compute_source_coefficients(const std::vector<arma::vec3> &point_list,
213  const ElementAccessor<3> &ele_acc,
214  std::vector<std::vector<double> > &sources_conc,
215  std::vector<std::vector<double> > &sources_density,
216  std::vector<std::vector<double> > &sources_sigma) override;
217 
218  void compute_sources_sigma(const std::vector<arma::vec3> &point_list,
219  const ElementAccessor<3> &ele_acc,
220  std::vector<std::vector<double> > &sources_sigma) override;
221 
222  ~HeatTransferModel() override;
223 
224  /**
225  * @brief Updates the velocity field which determines some coefficients of the transport equation.
226  *
227  * @param dh mixed hybrid dof handler
228  *
229  * (So far it does not work since the flow module returns a vector of zeros.)
230  */
231  inline void set_velocity_field(const MH_DofHandler &dh) override
232  {
233  mh_dh = &dh;
234  flux_changed = true;
235  }
236 
237  /// Returns number of transported substances.
238  inline unsigned int n_substances()
239  { return 1; }
240 
241  /// Returns reference to the vector of substance names.
242  inline SubstanceList &substances()
243  { return substances_; }
244 
245 
246 protected:
247 
248  /// Derived class should implement getter for ModelEqData instance.
249  virtual ModelEqData &data() = 0;
250 
251  /**
252  * Create input type that can be passed to the derived class.
253  * @param implementation String characterizing the numerical method, e.g. DG, FEM, FVM.
254  * @param description Comment used to describe the record key.
255  * @return
256  */
257  static IT::Record get_input_type(const string &implementation, const string &description);
258 
259  void output_data() override;
260 
261  std::shared_ptr<OutputTime> &output_stream()
262  { return output_stream_; }
263 
264  virtual void calculate_cumulative_balance() = 0;
265 
266  /// Indicator of change in advection vector field.
267  bool flux_changed;
268 
269  /// Transported substances.
270  SubstanceList substances_;
271 
272  /**
273  * Temporary solution how to pass velocity field form the flow model.
274  * TODO: introduce FieldDiscrete -containing true DOFHandler and data vector and pass such object together with other
275  * data. Possibly make more general set_data method, allowing setting data given by name. needs support from EqDataBase.
276  */
277  const MH_DofHandler *mh_dh;
278 
279  /// List of indices used to call balance methods for a set of quantities.
280  vector<unsigned int> subst_idx;
281 
282  std::shared_ptr<OutputTime> output_stream_;
283 
284 
285 };
286 
287 
288 
289 
290 
291 #endif /* HEAT_MODEL_HH_ */
Abstract linear system class.
Definition: balance.hh:35
Definition: mesh.h:76
The class for outputting data during time.
Definition: output_time.hh:50
Discontinuous Galerkin method for equation of transport with dispersion.
Classes for storing substance data.