Flow123d
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  ADD_FIELD(bc_temperature, "Boundary value of temperature.", "0.0");
53 
54  ADD_FIELD(init_temperature, "Initial temperature.", "0.0");
55  ADD_FIELD(porosity, "Porosity.", "1.0" );
56  ADD_FIELD(fluid_density, "Density of fluid.");
57  ADD_FIELD(fluid_heat_capacity, "Heat capacity of fluid.");
58  ADD_FIELD(fluid_heat_conductivity, "Heat conductivity of fluid.");
59  ADD_FIELD(solid_density, "Density of solid (rock).");
60  ADD_FIELD(solid_heat_capacity, "Heat capacity of solid (rock).");
61  ADD_FIELD(solid_heat_conductivity, "Heat conductivity of solid (rock).");
62  ADD_FIELD(disp_l, "Longitudal heat dispersivity in fluid.", "0.0" );
63  ADD_FIELD(disp_t, "Transversal heat dispersivity in fluid.", "0.0" );
64  ADD_FIELD(fluid_thermal_source, "Thermal source density in fluid.", "0.0");
65  ADD_FIELD(solid_thermal_source, "Thermal source density in solid.", "0.0");
66  ADD_FIELD(fluid_heat_exchange_rate, "Heat exchange rate in fluid.", "0.0");
67  ADD_FIELD(solid_heat_exchange_rate, "Heat exchange rate in solid.", "0.0");
68  ADD_FIELD(fluid_ref_temperature, "Reference temperature in fluid.", "0.0");
69  ADD_FIELD(solid_ref_temperature, "Reference temperature in solid.", "0.0");
70 
71  *this += cross_section.name("cross_section").just_copy();
72  output_fields += output_field.name("temperature").units("M/L^3");
73 }
74 
75 
76 IT::Record &HeatTransferModel::get_input_type(const string &implementation, const string &description)
77 {
78  static IT::Record input_type = IT::Record(ModelEqData::name() + "_" + implementation, description + " for heat transfer.")
80 
81  return input_type;
82 
83 }
84 
85 
86 IT::Selection &HeatTransferModel::ModelEqData::get_output_selection_input_type(const string &implementation, const string &description)
87 {
88  static IT::Selection input_type = IT::Selection(ModelEqData::name() + "_" + implementation + "_Output", "Selection for output fields of " + description + " for heat transfer.");
89 
90  return input_type;
91 }
92 
93 
95  flux_changed(true)
96 {}
97 
98 
99 void HeatTransferModel::init_data(unsigned int n_subst_)
100 {}
101 
102 /*
103 void HeatTransferModel::set_cross_section_field(const Field< 3, FieldValue<3>::Scalar > &cross_section)
104 {
105  data().cross_section.copy_from(cross_section);
106 }*/
107 
108 
110 {
111  names.clear();
112  names.push_back("");
113 }
114 
115 
117 {
118  return (data().porosity.changed() ||
119  data().fluid_density.changed() ||
121  data().solid_density.changed() ||
124 }
125 
126 
128 {
129  return (flux_changed ||
130  data().fluid_density.changed() ||
132  data().porosity.changed() ||
135  data().disp_l.changed() ||
136  data().disp_t.changed() ||
138 }
139 
140 
142 {
143  return (flux_changed ||
144  data().fluid_thermal_source.changed() ||
151 }
152 
154  const ElementAccessor<3> &ele_acc,
155  std::vector<double> &mm_coef)
156 {
157  vector<double> elem_csec(point_list.size()),
158  por(point_list.size()),
159  f_rho(point_list.size()),
160  s_rho(point_list.size()),
161  f_c(point_list.size()),
162  s_c(point_list.size());
163 
164  data().cross_section.value_list(point_list, ele_acc, elem_csec);
165  data().porosity.value_list(point_list, ele_acc, por);
166  data().fluid_density.value_list(point_list, ele_acc, f_rho);
167  data().fluid_heat_capacity.value_list(point_list, ele_acc, f_c);
168  data().solid_density.value_list(point_list, ele_acc, s_rho);
169  data().solid_heat_capacity.value_list(point_list, ele_acc, s_c);
170 
171  for (unsigned int i=0; i<point_list.size(); i++)
172  mm_coef[i] = elem_csec[i]*(por[i]*f_rho[i]*f_c[i] + (1.-por[i])*s_rho[i]*s_c[i]);
173 }
174 
175 
177  const std::vector<arma::vec3> &velocity,
178  const ElementAccessor<3> &ele_acc,
181 {
182  const unsigned int qsize = point_list.size();
183  std::vector<double> f_rho(qsize), f_cap(qsize), f_cond(qsize),
184  s_cond(qsize), por(qsize), csection(qsize), disp_l(qsize), disp_t(qsize);
185 
186  data().fluid_density.value_list(point_list, ele_acc, f_rho);
187  data().fluid_heat_capacity.value_list(point_list, ele_acc, f_cap);
188  data().fluid_heat_conductivity.value_list(point_list, ele_acc, f_cond);
189  data().solid_heat_conductivity.value_list(point_list, ele_acc, s_cond);
190  data().disp_l.value_list(point_list, ele_acc, disp_l);
191  data().disp_t.value_list(point_list, ele_acc, disp_t);
192  data().porosity.value_list(point_list, ele_acc, por);
193  data().cross_section.value_list(point_list, ele_acc, csection);
194 
195  for (unsigned int k=0; k<qsize; k++) {
196  ad_coef[0][k] = velocity[k]*f_rho[k]*f_cap[k];
197 
198  // dispersive part of thermal diffusion
199  // Note that the velocity vector is in fact the Darcian flux,
200  // so to obtain |v| we have to divide vnorm by porosity and cross_section.
201  double vnorm = arma::norm(velocity[k], 2);
202  if (fabs(vnorm) > sqrt(numeric_limits<double>::epsilon()))
203  for (int i=0; i<3; i++)
204  for (int j=0; j<3; j++)
205  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))
206  *vnorm*f_rho[k]*f_cond[k];
207  else
208  dif_coef[0][k].zeros();
209 
210  // conductive part of thermal diffusion
211  dif_coef[0][k] += csection[k]*(por[k]*f_cond[k] + (1.-por[k])*s_cond[k])*arma::eye(3,3);
212  }
213 }
214 
215 
217  const ElementAccessor<3> &ele_acc,
218  std::vector< arma::vec > &init_values)
219 {
220  vector<double> init_value(point_list.size());
221  data().init_temperature.value_list(point_list, ele_acc, init_value);
222  for (int i=0; i<point_list.size(); i++)
223  init_values[i] = init_value[i];
224 }
225 
226 
228  const ElementAccessor<3> &ele_acc,
229  std::vector< arma::vec > &bc_values)
230 {
231  vector<double> bc_value(point_list.size());
232  data().bc_temperature.value_list(point_list, ele_acc, bc_value);
233  for (int i=0; i<point_list.size(); i++)
234  bc_values[i] = bc_value[i];
235 }
236 
237 
239  const ElementAccessor<3> &ele_acc,
240  std::vector<arma::vec> &sources_value,
241  std::vector<arma::vec> &sources_density,
242  std::vector<arma::vec> &sources_sigma)
243 {
244  const unsigned int qsize = point_list.size();
245  std::vector<double> por(qsize), csection(qsize), f_rho(qsize), s_rho(qsize), f_cap(qsize), s_cap(qsize),
246  f_source(qsize), s_source(qsize), f_sigma(qsize), s_sigma(qsize), f_temp(qsize), s_temp(qsize);
247  data().porosity.value_list(point_list, ele_acc, por);
248  data().cross_section.value_list(point_list, ele_acc, csection);
249  data().fluid_density.value_list(point_list, ele_acc, f_rho);
250  data().solid_density.value_list(point_list, ele_acc, s_rho);
251  data().fluid_heat_capacity.value_list(point_list, ele_acc, f_cap);
252  data().solid_heat_capacity.value_list(point_list, ele_acc, s_cap);
253  data().fluid_thermal_source.value_list(point_list, ele_acc, f_source);
254  data().solid_thermal_source.value_list(point_list, ele_acc, s_source);
255  data().fluid_heat_exchange_rate.value_list(point_list, ele_acc, f_sigma);
256  data().solid_heat_exchange_rate.value_list(point_list, ele_acc, s_sigma);
257  data().fluid_ref_temperature.value_list(point_list, ele_acc, f_temp);
258  data().solid_ref_temperature.value_list(point_list, ele_acc, s_temp);
259 
260  for (int k=0; k<point_list.size(); k++)
261  {
262  sources_density[k].resize(1);
263  sources_sigma[k].resize(1);
264  sources_value[k].resize(1);
265 
266  sources_density[k][0] = csection[k]*(por[k]*f_source[k] + (1.-por[k])*s_source[k]);
267  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]);
268  if (fabs(sources_sigma[k][0]) > numeric_limits<double>::epsilon())
269  sources_value[k][0] = csection[k]*(por[k]*f_rho[k]*f_cap[k]*f_sigma[k]*f_temp[k]
270  + (1.-por[k])*s_rho[k]*s_cap[k]*s_sigma[k]*s_temp[k])/sources_sigma[k][0];
271  else
272  sources_value[k][0] = 0;
273  }
274 }
275 
276 
278  const ElementAccessor<3> &ele_acc,
279  std::vector<arma::vec> &sources_sigma)
280 {
281  const unsigned int qsize = point_list.size();
282  std::vector<double> por(qsize), csection(qsize), f_rho(qsize), s_rho(qsize), f_cap(qsize), s_cap(qsize),
283  f_source(qsize), s_source(qsize), f_sigma(qsize), s_sigma(qsize), f_temp(qsize), s_temp(qsize);
284  data().porosity.value_list(point_list, ele_acc, por);
285  data().cross_section.value_list(point_list, ele_acc, csection);
286  data().fluid_density.value_list(point_list, ele_acc, f_rho);
287  data().solid_density.value_list(point_list, ele_acc, s_rho);
288  data().fluid_heat_capacity.value_list(point_list, ele_acc, f_cap);
289  data().solid_heat_capacity.value_list(point_list, ele_acc, s_cap);
290  data().fluid_heat_exchange_rate.value_list(point_list, ele_acc, f_sigma);
291  data().solid_heat_exchange_rate.value_list(point_list, ele_acc, s_sigma);
292  for (int k=0; k<point_list.size(); k++)
293  {
294  sources_sigma[k].resize(1);
295  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]);
296  }
297 }
298 
299 
301 {}
302