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