Flow123d  jenkins-Flow123d-windows32-release-multijob-51
heat_model.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
27  * @author Jan Stebel
28  */
29 
30 #include "input/input_type.hh"
31 #include "mesh/mesh.h"
32 #include "mesh/accessors.hh"
33 #include "flow/darcy_flow_mh.hh"
35 #include "heat_model.hh"
36 
37 
38 
39 using namespace std;
40 using namespace Input::Type;
41 
42 
43 
44 
45 
46 
47 
48 
49 
51 {
52 
53  *this+=bc_temperature
54  .name("bc_temperature")
55  .description("Boundary value of temperature.")
56  .units( UnitSI().K() )
57  .input_default("0.0")
58  .flags_add(in_rhs);
59 
60  *this+=init_temperature
61  .name("init_temperature")
62  .description("Initial temperature.")
63  .units( UnitSI().K() )
64  .input_default("0.0");
65 
66  *this+=porosity
67  .name("porosity")
68  .description("Porosity.")
69  .units( UnitSI::dimensionless() )
70  .input_default("1.0")
71  .flags_add(in_main_matrix & in_time_term);
72 
73  *this+=fluid_density
74  .name("fluid_density")
75  .description("Density of fluid.")
76  .units( UnitSI().kg().m(-3) )
77  .flags_add(in_main_matrix & in_time_term);
78 
79  *this+=fluid_heat_capacity
80  .name("fluid_heat_capacity")
81  .description("Heat capacity of fluid.")
82  .units( UnitSI::J() * UnitSI().kg(-1).K(-1) )
83  .flags_add(in_main_matrix & in_time_term);
84 
85  *this+=fluid_heat_conductivity
86  .name("fluid_heat_conductivity")
87  .description("Heat conductivity of fluid.")
88  .units( UnitSI::W() * UnitSI().m(-1).K(-1) )
89  .flags_add(in_main_matrix);
90 
91 
92  *this+=solid_density
93  .name("solid_density")
94  .description("Density of solid (rock).")
95  .units( UnitSI().kg().m(-3) )
96  .flags_add(in_time_term);
97 
98  *this+=solid_heat_capacity
99  .name("solid_heat_capacity")
100  .description("Heat capacity of solid (rock).")
101  .units( UnitSI::J() * UnitSI().kg(-1).K(-1) )
102  .flags_add(in_time_term);
103 
104  *this+=solid_heat_conductivity
105  .name("solid_heat_conductivity")
106  .description("Heat conductivity of solid (rock).")
107  .units( UnitSI::W() * UnitSI().m(-1).K(-1) )
108  .flags_add(in_main_matrix);
109 
110  *this+=disp_l
111  .name("disp_l")
112  .description("Longitudal heat dispersivity in fluid.")
113  .units( UnitSI().m() )
114  .input_default("0.0")
115  .flags_add(in_main_matrix);
116 
117  *this+=disp_t
118  .name("disp_t")
119  .description("Transversal heat dispersivity in fluid.")
120  .units( UnitSI().m() )
121  .input_default("0.0")
122  .flags_add(in_main_matrix);
123 
124  *this+=fluid_thermal_source
125  .name("fluid_thermal_source")
126  .description("Thermal source density in fluid.")
127  .units( UnitSI::W() * UnitSI().m(-3) )
128  .input_default("0.0")
129  .flags_add(in_rhs);
130 
131  *this+=solid_thermal_source
132  .name("solid_thermal_source")
133  .description("Thermal source density in solid.")
134  .units( UnitSI::W() * UnitSI().m(-3) )
135  .input_default("0.0")
136  .flags_add(in_rhs);
137 
138  *this+=fluid_heat_exchange_rate
139  .name("fluid_heat_exchange_rate")
140  .description("Heat exchange rate in fluid.")
141  .units( UnitSI().s(-1) )
142  .input_default("0.0")
143  .flags_add(in_rhs);
144 
145  *this+=solid_heat_exchange_rate
146  .name("solid_heat_exchange_rate")
147  .description("Heat exchange rate of source in solid.")
148  .units( UnitSI().s(-1) )
149  .input_default("0.0")
150  .flags_add(in_rhs);
151 
152  *this+=fluid_ref_temperature
153  .name("fluid_ref_temperature")
154  .description("Reference temperature of source in fluid.")
155  .units( UnitSI().K() )
156  .input_default("0.0")
157  .flags_add(in_rhs);
158 
159  *this+=solid_ref_temperature
160  .name("solid_ref_temperature")
161  .description("Reference temperature in solid.")
162  .units( UnitSI().K() )
163  .input_default("0.0")
164  .flags_add(in_rhs);
165 
166  *this+=cross_section
167  .name("cross_section")
168  .units( UnitSI().m(3).md() )
169  .flags(input_copy & in_time_term & in_main_matrix);
170 
171  *this+=output_field
172  .name("temperature")
173  .units( UnitSI().K() )
174  .flags(equation_result);
175 }
176 
177 
178 IT::Record &HeatTransferModel::get_input_type(const string &implementation, const string &description)
179 {
180  static IT::Record input_type = IT::Record(ModelEqData::name() + "_" + implementation, description + " for heat transfer.")
182 
183  return input_type;
184 
185 }
186 
187 
188 IT::Selection &HeatTransferModel::ModelEqData::get_output_selection_input_type(const string &implementation, const string &description)
189 {
190  static IT::Selection input_type = IT::Selection(ModelEqData::name() + "_" + implementation + "_Output", "Selection for output fields of " + description + " for heat transfer.");
191 
192  return input_type;
193 }
194 
195 
197  flux_changed(true)
198 {}
199 
200 
202 {
203  names.clear();
204  names.push_back("");
205 }
206 
207 
209  const ElementAccessor<3> &ele_acc,
210  std::vector<double> &mm_coef)
211 {
212  vector<double> elem_csec(point_list.size()),
213  por(point_list.size()),
214  f_rho(point_list.size()),
215  s_rho(point_list.size()),
216  f_c(point_list.size()),
217  s_c(point_list.size());
218 
219  data().cross_section.value_list(point_list, ele_acc, elem_csec);
220  data().porosity.value_list(point_list, ele_acc, por);
221  data().fluid_density.value_list(point_list, ele_acc, f_rho);
222  data().fluid_heat_capacity.value_list(point_list, ele_acc, f_c);
223  data().solid_density.value_list(point_list, ele_acc, s_rho);
224  data().solid_heat_capacity.value_list(point_list, ele_acc, s_c);
225 
226  for (unsigned int i=0; i<point_list.size(); i++)
227  mm_coef[i] = elem_csec[i]*(por[i]*f_rho[i]*f_c[i] + (1.-por[i])*s_rho[i]*s_c[i]);
228 }
229 
230 
232  const std::vector<arma::vec3> &velocity,
233  const ElementAccessor<3> &ele_acc,
236 {
237  const unsigned int qsize = point_list.size();
238  std::vector<double> f_rho(qsize), f_cap(qsize), f_cond(qsize),
239  s_cond(qsize), por(qsize), csection(qsize), disp_l(qsize), disp_t(qsize);
240 
241  data().fluid_density.value_list(point_list, ele_acc, f_rho);
242  data().fluid_heat_capacity.value_list(point_list, ele_acc, f_cap);
243  data().fluid_heat_conductivity.value_list(point_list, ele_acc, f_cond);
244  data().solid_heat_conductivity.value_list(point_list, ele_acc, s_cond);
245  data().disp_l.value_list(point_list, ele_acc, disp_l);
246  data().disp_t.value_list(point_list, ele_acc, disp_t);
247  data().porosity.value_list(point_list, ele_acc, por);
248  data().cross_section.value_list(point_list, ele_acc, csection);
249 
250  for (unsigned int k=0; k<qsize; k++) {
251  ad_coef[0][k] = velocity[k]*f_rho[k]*f_cap[k];
252 
253  // dispersive part of thermal diffusion
254  // Note that the velocity vector is in fact the Darcian flux,
255  // so to obtain |v| we have to divide vnorm by porosity and cross_section.
256  double vnorm = arma::norm(velocity[k], 2);
257  if (fabs(vnorm) > sqrt(numeric_limits<double>::epsilon()))
258  for (int i=0; i<3; i++)
259  for (int j=0; j<3; j++)
260  dif_coef[0][k](i,j) = (velocity[k][i]*velocity[k][j]/(vnorm*vnorm)*(disp_l[k]-disp_t[k]) + disp_t[k]*(i==j?1:0))
261  *vnorm*f_rho[k]*f_cond[k];
262  else
263  dif_coef[0][k].zeros();
264 
265  // conductive part of thermal diffusion
266  dif_coef[0][k] += csection[k]*(por[k]*f_cond[k] + (1.-por[k])*s_cond[k])*arma::eye(3,3);
267  }
268 }
269 
270 
272  const ElementAccessor<3> &ele_acc,
273  std::vector< arma::vec > &init_values)
274 {
275  vector<double> init_value(point_list.size());
276  data().init_temperature.value_list(point_list, ele_acc, init_value);
277  for (unsigned int i=0; i<point_list.size(); i++)
278  init_values[i] = init_value[i];
279 }
280 
281 
283  const ElementAccessor<3> &ele_acc,
284  std::vector< arma::vec > &bc_values)
285 {
286  vector<double> bc_value(point_list.size());
287  data().bc_temperature.value_list(point_list, ele_acc, bc_value);
288  for (unsigned int i=0; i<point_list.size(); i++)
289  bc_values[i] = bc_value[i];
290 }
291 
292 
294  const ElementAccessor<3> &ele_acc,
295  std::vector<arma::vec> &sources_value,
296  std::vector<arma::vec> &sources_density,
297  std::vector<arma::vec> &sources_sigma)
298 {
299  const unsigned int qsize = point_list.size();
300  std::vector<double> por(qsize), csection(qsize), f_rho(qsize), s_rho(qsize), f_cap(qsize), s_cap(qsize),
301  f_source(qsize), s_source(qsize), f_sigma(qsize), s_sigma(qsize), f_temp(qsize), s_temp(qsize);
302  data().porosity.value_list(point_list, ele_acc, por);
303  data().cross_section.value_list(point_list, ele_acc, csection);
304  data().fluid_density.value_list(point_list, ele_acc, f_rho);
305  data().solid_density.value_list(point_list, ele_acc, s_rho);
306  data().fluid_heat_capacity.value_list(point_list, ele_acc, f_cap);
307  data().solid_heat_capacity.value_list(point_list, ele_acc, s_cap);
308  data().fluid_thermal_source.value_list(point_list, ele_acc, f_source);
309  data().solid_thermal_source.value_list(point_list, ele_acc, s_source);
310  data().fluid_heat_exchange_rate.value_list(point_list, ele_acc, f_sigma);
311  data().solid_heat_exchange_rate.value_list(point_list, ele_acc, s_sigma);
312  data().fluid_ref_temperature.value_list(point_list, ele_acc, f_temp);
313  data().solid_ref_temperature.value_list(point_list, ele_acc, s_temp);
314 
315  for (unsigned int k=0; k<point_list.size(); k++)
316  {
317  sources_density[k].resize(1);
318  sources_sigma[k].resize(1);
319  sources_value[k].resize(1);
320 
321  sources_density[k][0] = csection[k]*(por[k]*f_source[k] + (1.-por[k])*s_source[k]);
322  sources_sigma[k][0] = csection[k]*(por[k]*f_rho[k]*f_cap[k]*f_sigma[k] + (1.-por[k])*s_rho[k]*s_cap[k]*s_sigma[k]);
323  if (fabs(sources_sigma[k][0]) > numeric_limits<double>::epsilon())
324  sources_value[k][0] = csection[k]*(por[k]*f_rho[k]*f_cap[k]*f_sigma[k]*f_temp[k]
325  + (1.-por[k])*s_rho[k]*s_cap[k]*s_sigma[k]*s_temp[k])/sources_sigma[k][0];
326  else
327  sources_value[k][0] = 0;
328  }
329 }
330 
331 
333  const ElementAccessor<3> &ele_acc,
334  std::vector<arma::vec> &sources_sigma)
335 {
336  const unsigned int qsize = point_list.size();
337  std::vector<double> por(qsize), csection(qsize), f_rho(qsize), s_rho(qsize), f_cap(qsize), s_cap(qsize),
338  f_source(qsize), s_source(qsize), f_sigma(qsize), s_sigma(qsize), f_temp(qsize), s_temp(qsize);
339  data().porosity.value_list(point_list, ele_acc, por);
340  data().cross_section.value_list(point_list, ele_acc, csection);
341  data().fluid_density.value_list(point_list, ele_acc, f_rho);
342  data().solid_density.value_list(point_list, ele_acc, s_rho);
343  data().fluid_heat_capacity.value_list(point_list, ele_acc, f_cap);
344  data().solid_heat_capacity.value_list(point_list, ele_acc, s_cap);
345  data().fluid_heat_exchange_rate.value_list(point_list, ele_acc, f_sigma);
346  data().solid_heat_exchange_rate.value_list(point_list, ele_acc, s_sigma);
347  for (unsigned int k=0; k<point_list.size(); k++)
348  {
349  sources_sigma[k].resize(1);
350  sources_sigma[k][0] = csection[k]*(por[k]*f_rho[k]*f_cap[k]*f_sigma[k] + (1.-por[k])*s_rho[k]*s_cap[k]*s_sigma[k]);
351  }
352 }
353 
354 
356 {}
357 
Field< 3, FieldValue< 3 >::Scalar > solid_heat_capacity
Heat capacity of solid.
Definition: heat_model.hh:59
void compute_sources_sigma(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &sources_sigma) override
Definition: heat_model.cc:332
Field< 3, FieldValue< 3 >::Scalar > fluid_density
Density of fluid.
Definition: heat_model.hh:51
BCField< 3, FieldValue< 3 >::Scalar > bc_temperature
Dirichlet boundary condition for temperature.
Definition: heat_model.hh:45
static Input::Type::AbstractRecord input_type
Common specification of the input record for secondary equations.
Field< 3, FieldValue< 3 >::Scalar > disp_t
Transversal heat dispersivity.
Definition: heat_model.hh:65
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Definition: heat_model.hh:80
???
Field< 3, FieldValue< 3 >::Scalar > fluid_heat_exchange_rate
Heat exchange rate in fluid.
Definition: heat_model.hh:71
Field< 3, FieldValue< 3 >::Scalar > solid_heat_exchange_rate
Heat exchange rate in solid.
Definition: heat_model.hh:73
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:169
void compute_init_cond(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &init_values) override
Definition: heat_model.cc:271
Field< 3, FieldValue< 3 >::Scalar > fluid_heat_capacity
Heat capacity of fluid.
Definition: heat_model.hh:53
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list) const
Definition: field.hh:413
Field< 3, FieldValue< 3 >::Scalar > porosity
Porosity of solid.
Definition: heat_model.hh:49
Field< 3, FieldValue< 3 >::Scalar > fluid_thermal_source
Thermal source in fluid.
Definition: heat_model.hh:67
Field< 3, FieldValue< 3 >::Scalar > fluid_heat_conductivity
Heat conductivity of fluid.
Definition: heat_model.hh:55
Field< 3, FieldValue< 3 >::Scalar > fluid_ref_temperature
Reference temperature in fluid.
Definition: heat_model.hh:75
void set_component_names(std::vector< string > &names, const Input::Record &in_rec) override
Read or set names of solution components.
Definition: heat_model.cc:201
Field< 3, FieldValue< 3 >::Scalar > init_temperature
Initial temperature.
Definition: heat_model.hh:47
Field< 3, FieldValue< 3 >::Scalar > solid_thermal_source
Thermal source in solid.
Definition: heat_model.hh:69
void compute_mass_matrix_coefficient(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &mm_coef) override
Definition: heat_model.cc:208
void compute_dirichlet_bc(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &bc_values) override
Definition: heat_model.cc:282
static UnitSI & W()
Returns Watt.
Definition: unit_si.cc:34
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
static UnitSI & J()
Returns Joule.
Definition: unit_si.cc:29
Field< 3, FieldValue< 3 >::Scalar > solid_heat_conductivity
Heat conductivity of solid.
Definition: heat_model.hh:61
const double epsilon
Definition: mathfce.h:6
Field< 3, FieldValue< 3 >::Scalar > solid_ref_temperature
Reference temperature in solid.
Definition: heat_model.hh:77
static IT::Selection & get_output_selection_input_type(const string &implementation, const string &description)
Definition: heat_model.cc:188
virtual ModelEqData & data()=0
Derived class should implement getter for ModelEqData instance.
void compute_advection_diffusion_coefficients(const std::vector< arma::vec3 > &point_list, const std::vector< arma::vec3 > &velocity, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< arma::vec3 > > &ad_coef, std::vector< std::vector< arma::mat33 > > &dif_coef) override
Definition: heat_model.cc:231
~HeatTransferModel() override
Definition: heat_model.cc:355
Field< 3, FieldValue< 3 >::Scalar > solid_density
Density of solid.
Definition: heat_model.hh:57
Discontinuous Galerkin method for equation of transport with dispersion.
mixed-hybrid model of linear Darcy flow, possibly unsteady.
Record type proxy class.
Definition: type_record.hh:161
void compute_source_coefficients(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &sources_conc, std::vector< arma::vec > &sources_density, std::vector< arma::vec > &sources_sigma) override
Definition: heat_model.cc:293
Class for representation SI units of Fields.
Definition: unit_si.hh:31
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:44
Template for classes storing finite set of named values.
static IT::Record & get_input_type(const string &implementation, const string &description)
Definition: heat_model.cc:178
Field< 3, FieldValue< 3 >::Scalar > disp_l
Longitudal heat dispersivity.
Definition: heat_model.hh:63