GRINS-0.7.0
heat_conduction.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_conduction.h"
28 
29 // GRINS
30 #include "grins/common.h"
31 #include "grins/assembly_context.h"
34 #include "grins/physics_naming.h"
38 
39 // libMesh
40 #include "libmesh/quadrature.h"
41 #include "libmesh/fem_system.h"
42 
43 namespace GRINS
44 {
45 
46  template<class K>
47  HeatConduction<K>::HeatConduction( const GRINS::PhysicsName& physics_name, const GetPot& input )
48  : Physics(physics_name,input),
49  _temp_vars(input,PhysicsNaming::heat_conduction()),
50  _rho(0.0),
51  _Cp(0.0),
52  _k(input,MaterialsParsing::material_name(input,PhysicsNaming::heat_conduction()))
53  {
55 
57 
58  this->register_variables();
59 
60  // This is deleted in the base class
61  this->_ic_handler = new GenericICHandler( physics_name, input );
62  }
63 
64  template<class K>
66  {
67  return;
68  }
69 
70  template<class K>
72  {
74  this->_temp_vars);
75  }
76 
77  template<class K>
78  void HeatConduction<K>::init_variables( libMesh::FEMSystem* system )
79  {
80  // Get libMesh to assign an index for each variable
81  this->_dim = system->get_mesh().mesh_dimension();
82 
83  _temp_vars.init(system);
84 
85  return;
86  }
87 
88  template<class K>
89  void HeatConduction<K>::set_time_evolving_vars( libMesh::FEMSystem* system )
90  {
91  // Tell the system to march temperature forward in time
92  system->time_evolving(_temp_vars.T());
93 
94  return;
95  }
96 
97  template<class K>
99  {
100  // We should prerequest all the data
101  // we will need to build the linear system
102  // or evaluate a quantity of interest.
103  context.get_element_fe(_temp_vars.T())->get_JxW();
104  context.get_element_fe(_temp_vars.T())->get_phi();
105  context.get_element_fe(_temp_vars.T())->get_dphi();
106  context.get_element_fe(_temp_vars.T())->get_xyz();
107 
108  context.get_side_fe(_temp_vars.T())->get_JxW();
109  context.get_side_fe(_temp_vars.T())->get_phi();
110  context.get_side_fe(_temp_vars.T())->get_dphi();
111  context.get_side_fe(_temp_vars.T())->get_xyz();
112 
113  return;
114  }
115 
116  template<class K>
117  void HeatConduction<K>::element_time_derivative( bool compute_jacobian,
118  AssemblyContext& context,
119  CachedValues& /*cache*/ )
120  {
121  // The number of local degrees of freedom in each variable.
122  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
123 
124  // We get some references to cell-specific data that
125  // will be used to assemble the linear system.
126 
127  // Element Jacobian * quadrature weights for interior integration.
128  const std::vector<libMesh::Real> &JxW =
129  context.get_element_fe(_temp_vars.T())->get_JxW();
130 
131  // The temperature shape function gradients (in global coords.)
132  // at interior quadrature points.
133  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
134  context.get_element_fe(_temp_vars.T())->get_dphi();
135 
136  // The subvectors and submatrices we need to fill:
137  //
138  // K_{\alpha \beta} = R_{\alpha},{\beta} = \partial{ R_{\alpha} } / \partial{ {\beta} } (where R denotes residual)
139  // e.g., for \alpha = T and \beta = v we get: K_{Tu} = R_{T},{u}
140  //
141 
142  libMesh::DenseSubMatrix<libMesh::Number> &KTT = context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T()); // R_{T},{T}
143 
144  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(_temp_vars.T()); // R_{T}
145 
146  // Now we will build the element Jacobian and residual.
147  // Constructing the residual requires the solution and its
148  // gradient from the previous timestep. This must be
149  // calculated at each quadrature point by summing the
150  // solution degree-of-freedom values by the appropriate
151  // weight functions.
152  unsigned int n_qpoints = context.get_element_qrule().n_points();
153 
154  for (unsigned int qp=0; qp != n_qpoints; qp++)
155  {
156  // Compute the solution & its gradient at the old Newton iterate.
157  libMesh::Gradient grad_T;
158  grad_T = context.interior_gradient(_temp_vars.T(), qp);
159 
160  // Compute the conductivity at this qp
161  libMesh::Real _k_qp = this->_k(context, qp);
162 
163  // First, an i-loop over the degrees of freedom.
164  for (unsigned int i=0; i != n_T_dofs; i++)
165  {
166  FT(i) += JxW[qp]*(-_k_qp*(T_gradphi[i][qp]*grad_T));
167 
168  if (compute_jacobian)
169  {
170  for (unsigned int j=0; j != n_T_dofs; j++)
171  {
172  // TODO: precompute some terms like:
173  // _rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])
174 
175  KTT(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
176  ( -_k_qp*(T_gradphi[i][qp]*T_gradphi[j][qp]) ); // diffusion term
177  } // end of the inner dof (j) loop
178 
179  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
180 
181  } // end of the outer dof (i) loop
182  } // end of the quadrature point (qp) loop
183 
184  return;
185  }
186 
187  template<class K>
188  void HeatConduction<K>::mass_residual( bool compute_jacobian,
189  AssemblyContext& context,
190  CachedValues& /*cache*/ )
191  {
192  // First we get some references to cell-specific data that
193  // will be used to assemble the linear system.
194 
195  // Element Jacobian * quadrature weights for interior integration
196  const std::vector<libMesh::Real> &JxW =
197  context.get_element_fe(_temp_vars.T())->get_JxW();
198 
199  // The shape functions at interior quadrature points.
200  const std::vector<std::vector<libMesh::Real> >& phi =
201  context.get_element_fe(_temp_vars.T())->get_phi();
202 
203  // The number of local degrees of freedom in each variable
204  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
205 
206  // The subvectors and submatrices we need to fill:
207  libMesh::DenseSubVector<libMesh::Real> &F =
208  context.get_elem_residual(_temp_vars.T());
209 
210  libMesh::DenseSubMatrix<libMesh::Real> &M =
211  context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T());
212 
213  unsigned int n_qpoints = context.get_element_qrule().n_points();
214 
215  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
216  {
217  // For the mass residual, we need to be a little careful.
218  // The time integrator is handling the time-discretization
219  // for us so we need to supply M(u_fixed)*u' for the residual.
220  // u_fixed will be given by the fixed_interior_value function
221  // while u' will be given by the interior_rate function.
222  libMesh::Real T_dot;
223  context.interior_rate(_temp_vars.T(), qp, T_dot);
224 
225  for (unsigned int i = 0; i != n_T_dofs; ++i)
226  {
227  F(i) -= JxW[qp]*(_rho*_Cp*T_dot*phi[i][qp] );
228 
229  if( compute_jacobian )
230  {
231  for (unsigned int j=0; j != n_T_dofs; j++)
232  {
233  // We're assuming rho, cp are constant w.r.t. T here.
234  M(i,j) -=
235  context.get_elem_solution_rate_derivative()
236  * JxW[qp]*_rho*_Cp*phi[j][qp]*phi[i][qp] ;
237  }
238  }// End of check on Jacobian
239 
240  } // End of element dof loop
241 
242  } // End of the quadrature point loop
243 
244  return;
245  }
246 
247  template<class K>
249  ( const std::string & param_name,
251  const
252  {
253  ParameterUser::register_parameter(param_name, param_pointer);
254  _k.register_parameter(param_name, param_pointer);
255  }
256 
257 
258 } // namespace GRINS
259 
260 // Instantiate
261 INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatConduction);
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...
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:269
libMesh::Number _Cp
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:107
libMesh::Number _rho
static std::string temperature_section()
INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatConduction)
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
static void check_and_register_variable(const std::string &var_name, const FEVariablesBase &variable)
First check if var_name is registered and then register.
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
Base class for reading and handling initial conditions for physics classes.
static void read_density(const std::string &core_physics_name, const GetPot &input, ParameterUser &params, libMesh::Real &rho)
Helper function to reading density from input.
GRINS namespace.
static void read_specific_heat(const std::string &core_physics_name, const GetPot &input, ParameterUser &params, libMesh::Real &cp)
Helper function to reading scalar specific heat from input.
Helper functions for parsing material properties.
std::string PhysicsName
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Set which variables are time evolving.
static PhysicsName heat_conduction()
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
virtual void init_variables(libMesh::FEMSystem *system)
Initialize variables for this physics.

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