GRINS-0.7.0
heat_transfer.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2016 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 // This class
27 #include "grins/heat_transfer.h"
28 
29 // GRINS
30 #include "grins_config.h"
31 #include "grins/assembly_context.h"
35 
36 // libMesh
37 #include "libmesh/getpot.h"
38 #include "libmesh/quadrature.h"
39 #include "libmesh/boundary_info.h"
40 
41 namespace GRINS
42 {
43 
44  template<class K>
45  HeatTransfer<K>::HeatTransfer( const std::string& physics_name, const GetPot& input )
46  : HeatTransferBase<K>(physics_name, PhysicsNaming::heat_transfer(), input),
47  _k_index(0)
48  {
49  this->_ic_handler = new GenericICHandler( physics_name, input );
50  }
51 
52  template<class K>
55  {
56  std::string section = "Physics/"+PhysicsNaming::heat_transfer()+"/output_vars";
57 
58  if( input.have_variable(section) )
59  {
60  unsigned int n_vars = input.vector_variable_size(section);
61 
62  for( unsigned int v = 0; v < n_vars; v++ )
63  {
64  std::string name = input(section,"DIE!",v);
65 
66  if( name == std::string("k") )
67  {
68  this->_k_index = postprocessing.register_quantity( name );
69  }
70  else
71  {
72  std::cerr << "Error: Invalid output_vars value for "+PhysicsNaming::heat_transfer() << std::endl
73  << " Found " << name << std::endl
74  << " Acceptable values are: k" << std::endl;
75  libmesh_error();
76  }
77  }
78  }
79 
80  return;
81  }
82 
83  template<class K>
84  void HeatTransfer<K>::element_time_derivative( bool compute_jacobian,
85  AssemblyContext& context,
86  CachedValues& /*cache*/ )
87  {
88 #ifdef GRINS_USE_GRVY_TIMERS
89  this->_timer->BeginTimer("HeatTransfer::element_time_derivative");
90 #endif
91 
92  // The number of local degrees of freedom in each variable.
93  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
94  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
95 
96  //TODO: check n_T_dofs is same as n_u_dofs, n_v_dofs, n_w_dofs
97 
98  // We get some references to cell-specific data that
99  // will be used to assemble the linear system.
100 
101  // Element Jacobian * quadrature weights for interior integration.
102  const std::vector<libMesh::Real> &JxW =
103  context.get_element_fe(this->_temp_vars.T())->get_JxW();
104 
105  // The temperature shape functions at interior quadrature points.
106  const std::vector<std::vector<libMesh::Real> >& T_phi =
107  context.get_element_fe(this->_temp_vars.T())->get_phi();
108 
109  // The velocity shape functions at interior quadrature points.
110  const std::vector<std::vector<libMesh::Real> >& vel_phi =
111  context.get_element_fe(this->_flow_vars.u())->get_phi();
112 
113  // The temperature shape function gradients (in global coords.)
114  // at interior quadrature points.
115  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
116  context.get_element_fe(this->_temp_vars.T())->get_dphi();
117 
118  const std::vector<libMesh::Point>& u_qpoint =
119  context.get_element_fe(this->_flow_vars.u())->get_xyz();
120 
121  libMesh::DenseSubMatrix<libMesh::Number> &KTT = context.get_elem_jacobian(this->_temp_vars.T(), this->_temp_vars.T()); // R_{T},{T}
122 
123  libMesh::DenseSubMatrix<libMesh::Number> &KTu = context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.u()); // R_{T},{u}
124  libMesh::DenseSubMatrix<libMesh::Number> &KTv = context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.v()); // R_{T},{v}
125  libMesh::DenseSubMatrix<libMesh::Number>* KTw = NULL;
126 
127  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
128 
129  if( this->_dim == 3 )
130  {
131  KTw = &context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.w()); // R_{T},{w}
132  }
133 
134  // Now we will build the element Jacobian and residual.
135  // Constructing the residual requires the solution and its
136  // gradient from the previous timestep. This must be
137  // calculated at each quadrature point by summing the
138  // solution degree-of-freedom values by the appropriate
139  // weight functions.
140  unsigned int n_qpoints = context.get_element_qrule().n_points();
141 
142  for (unsigned int qp=0; qp != n_qpoints; qp++)
143  {
144  // Compute the solution & its gradient at the old Newton iterate.
145  libMesh::Number u, v;
146  u = context.interior_value(this->_flow_vars.u(), qp);
147  v = context.interior_value(this->_flow_vars.v(), qp);
148 
149  libMesh::Gradient grad_T;
150  grad_T = context.interior_gradient(this->_temp_vars.T(), qp);
151 
152  libMesh::NumberVectorValue U (u,v);
153  if (this->_dim == 3)
154  U(2) = context.interior_value(this->_flow_vars.w(), qp);
155 
156  const libMesh::Number r = u_qpoint[qp](0);
157 
158  libMesh::Real jac = JxW[qp];
159 
160  // Compute the conductivity at this qp
161  libMesh::Real _k_qp = this->_k(context, qp);
162 
164  {
165  jac *= r;
166  }
167 
168  // First, an i-loop over the degrees of freedom.
169  for (unsigned int i=0; i != n_T_dofs; i++)
170  {
171  FT(i) += jac *
172  (-this->_rho*this->_Cp*T_phi[i][qp]*(U*grad_T) // convection term
173  -_k_qp*(T_gradphi[i][qp]*grad_T) ); // diffusion term
174 
175  if (compute_jacobian)
176  {
177  for (unsigned int j=0; j != n_T_dofs; j++)
178  {
179  // TODO: precompute some terms like:
180  // this->_rho*this->_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])
181 
182  KTT(i,j) += jac * context.get_elem_solution_derivative() *
183  (-this->_rho*this->_Cp*T_phi[i][qp]*(U*T_gradphi[j][qp]) // convection term
184  -_k_qp*(T_gradphi[i][qp]*T_gradphi[j][qp])); // diffusion term
185  } // end of the inner dof (j) loop
186 
187  // Matrix contributions for the Tu, Tv and Tw couplings (n_T_dofs same as n_u_dofs, n_v_dofs and n_w_dofs)
188  for (unsigned int j=0; j != n_u_dofs; j++)
189  {
190  KTu(i,j) += jac * context.get_elem_solution_derivative()*(-this->_rho*this->_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(0)));
191  KTv(i,j) += jac * context.get_elem_solution_derivative()*(-this->_rho*this->_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(1)));
192  if (this->_dim == 3)
193  (*KTw)(i,j) += jac * context.get_elem_solution_derivative()*(-this->_rho*this->_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(2)));
194  } // end of the inner dof (j) loop
195 
196  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
197 
198  } // end of the outer dof (i) loop
199  } // end of the quadrature point (qp) loop
200 
201 #ifdef GRINS_USE_GRVY_TIMERS
202  this->_timer->EndTimer("HeatTransfer::element_time_derivative");
203 #endif
204 
205  return;
206  }
207 
208  template<class K>
209  void HeatTransfer<K>::mass_residual( bool compute_jacobian,
210  AssemblyContext& context,
211  CachedValues& /*cache*/ )
212  {
213 #ifdef GRINS_USE_GRVY_TIMERS
214  this->_timer->BeginTimer("HeatTransfer::mass_residual");
215 #endif
216 
217  // First we get some references to cell-specific data that
218  // will be used to assemble the linear system.
219 
220  // Element Jacobian * quadrature weights for interior integration
221  const std::vector<libMesh::Real> &JxW =
222  context.get_element_fe(this->_temp_vars.T())->get_JxW();
223 
224  // The shape functions at interior quadrature points.
225  const std::vector<std::vector<libMesh::Real> >& phi =
226  context.get_element_fe(this->_temp_vars.T())->get_phi();
227 
228  const std::vector<libMesh::Point>& u_qpoint =
229  context.get_element_fe(this->_flow_vars.u())->get_xyz();
230 
231  // The number of local degrees of freedom in each variable
232  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
233 
234  // The subvectors and submatrices we need to fill:
235  libMesh::DenseSubVector<libMesh::Real> &F = context.get_elem_residual(this->_temp_vars.T());
236 
237  libMesh::DenseSubMatrix<libMesh::Real> &M = context.get_elem_jacobian(this->_temp_vars.T(), this->_temp_vars.T());
238 
239  unsigned int n_qpoints = context.get_element_qrule().n_points();
240 
241  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
242  {
243  // For the mass residual, we need to be a little careful.
244  // The time integrator is handling the time-discretization
245  // for us so we need to supply M(u_fixed)*u' for the residual.
246  // u_fixed will be given by the fixed_interior_value function
247  // while u' will be given by the interior_rate function.
248  libMesh::Real T_dot;
249  context.interior_rate(this->_temp_vars.T(), qp, T_dot);
250 
251  const libMesh::Number r = u_qpoint[qp](0);
252 
253  libMesh::Real jac = JxW[qp];
254 
256  {
257  jac *= r;
258  }
259 
260  for (unsigned int i = 0; i != n_T_dofs; ++i)
261  {
262  F(i) -= this->_rho*this->_Cp*T_dot*phi[i][qp]*jac;
263 
264  if( compute_jacobian )
265  {
266  for (unsigned int j=0; j != n_T_dofs; j++)
267  {
268  // We're assuming rho, cp are constant w.r.t. T here.
269  M(i,j) -= this->_rho*this->_Cp*phi[j][qp]*phi[i][qp]*jac * context.get_elem_solution_rate_derivative();
270  }
271  }// End of check on Jacobian
272 
273  } // End of element dof loop
274 
275  } // End of the quadrature point loop
276 
277 #ifdef GRINS_USE_GRVY_TIMERS
278  this->_timer->EndTimer("HeatTransfer::mass_residual");
279 #endif
280 
281  return;
282  }
283 
284  template<class K>
285  void HeatTransfer<K>::compute_postprocessed_quantity( unsigned int quantity_index,
286  const AssemblyContext& context,
287  const libMesh::Point& point,
288  libMesh::Real& value )
289  {
290  if( quantity_index == this->_k_index )
291  {
292  value = this->_k(point, context.get_time());
293  }
294 
295  return;
296  }
297 
298 } // namespace GRINS
299 
300 // Instantiate
static bool is_axisymmetric()
Definition: physics.h:133
static PhysicsName heat_transfer()
unsigned int register_quantity(std::string name)
Register quantity to be postprocessed.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:269
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
Definition: heat_transfer.C:84
virtual void mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
Physics class for Heat Transfer.
Base class for reading and handling initial conditions for physics classes.
GRINS namespace.
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
Compute value of postprocessed quantities at libMesh::Point.
INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatTransfer)
virtual void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Register postprocessing variables for HeatTransfer.
Definition: heat_transfer.C:53

Generated on Thu Jun 2 2016 21:52:28 for GRINS-0.7.0 by  doxygen 1.8.10