GRINS-0.8.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-2017 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(GRINSPrivate::VariableWarehouse::get_variable_subclass<PrimitiveTempFEVariables>(VariablesParsing::temp_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
50  _rho(0.0),
51  _Cp(0.0),
52  _k(input,MaterialsParsing::material_name(input,PhysicsNaming::heat_conduction()))
53  {
55 
57 
58  // This is deleted in the base class
59  this->_ic_handler = new GenericICHandler( physics_name, input );
60 
62  }
63 
64  template<class K>
65  void HeatConduction<K>::set_time_evolving_vars( libMesh::FEMSystem* system )
66  {
67  // Tell the system to march temperature forward in time
68  system->time_evolving(_temp_vars.T(), 1);
69  }
70 
71  template<class K>
73  {
74  // We should prerequest all the data
75  // we will need to build the linear system
76  // or evaluate a quantity of interest.
77  context.get_element_fe(_temp_vars.T())->get_JxW();
78  context.get_element_fe(_temp_vars.T())->get_phi();
79  context.get_element_fe(_temp_vars.T())->get_dphi();
80  context.get_element_fe(_temp_vars.T())->get_xyz();
81 
82  context.get_side_fe(_temp_vars.T())->get_JxW();
83  context.get_side_fe(_temp_vars.T())->get_phi();
84  context.get_side_fe(_temp_vars.T())->get_dphi();
85  context.get_side_fe(_temp_vars.T())->get_xyz();
86  }
87 
88  template<class K>
90  ( bool compute_jacobian,
91  AssemblyContext & context )
92  {
93  // The number of local degrees of freedom in each variable.
94  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
95 
96  // We get some references to cell-specific data that
97  // will be used to assemble the linear system.
98 
99  // Element Jacobian * quadrature weights for interior integration.
100  const std::vector<libMesh::Real> &JxW =
101  context.get_element_fe(_temp_vars.T())->get_JxW();
102 
103  // The temperature shape function gradients (in global coords.)
104  // at interior quadrature points.
105  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
106  context.get_element_fe(_temp_vars.T())->get_dphi();
107 
108  // The subvectors and submatrices we need to fill:
109  //
110  // K_{\alpha \beta} = R_{\alpha},{\beta} = \partial{ R_{\alpha} } / \partial{ {\beta} } (where R denotes residual)
111  // e.g., for \alpha = T and \beta = v we get: K_{Tu} = R_{T},{u}
112  //
113 
114  libMesh::DenseSubMatrix<libMesh::Number> &KTT = context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T()); // R_{T},{T}
115 
116  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(_temp_vars.T()); // R_{T}
117 
118  // Now we will build the element Jacobian and residual.
119  // Constructing the residual requires the solution and its
120  // gradient from the previous timestep. This must be
121  // calculated at each quadrature point by summing the
122  // solution degree-of-freedom values by the appropriate
123  // weight functions.
124  unsigned int n_qpoints = context.get_element_qrule().n_points();
125 
126  for (unsigned int qp=0; qp != n_qpoints; qp++)
127  {
128  // Compute the solution & its gradient at the old Newton iterate.
129  libMesh::Gradient grad_T;
130  grad_T = context.interior_gradient(_temp_vars.T(), qp);
131 
132  // Compute the conductivity at this qp
133  libMesh::Real _k_qp = this->_k(context, qp);
134 
135  // First, an i-loop over the degrees of freedom.
136  for (unsigned int i=0; i != n_T_dofs; i++)
137  {
138  FT(i) += JxW[qp]*(-_k_qp*(T_gradphi[i][qp]*grad_T));
139 
140  if (compute_jacobian)
141  {
142  for (unsigned int j=0; j != n_T_dofs; j++)
143  {
144  // TODO: precompute some terms like:
145  // _rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])
146 
147  KTT(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
148  ( -_k_qp*(T_gradphi[i][qp]*T_gradphi[j][qp]) ); // diffusion term
149  } // end of the inner dof (j) loop
150 
151  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
152 
153  } // end of the outer dof (i) loop
154  } // end of the quadrature point (qp) loop
155 
156  return;
157  }
158 
159  template<class K>
161  ( bool compute_jacobian, AssemblyContext & context )
162  {
163  // First we get some references to cell-specific data that
164  // will be used to assemble the linear system.
165 
166  // Element Jacobian * quadrature weights for interior integration
167  const std::vector<libMesh::Real> &JxW =
168  context.get_element_fe(_temp_vars.T())->get_JxW();
169 
170  // The shape functions at interior quadrature points.
171  const std::vector<std::vector<libMesh::Real> >& phi =
172  context.get_element_fe(_temp_vars.T())->get_phi();
173 
174  // The number of local degrees of freedom in each variable
175  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
176 
177  // The subvectors and submatrices we need to fill:
178  libMesh::DenseSubVector<libMesh::Real> &F =
179  context.get_elem_residual(_temp_vars.T());
180 
181  libMesh::DenseSubMatrix<libMesh::Real> &M =
182  context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T());
183 
184  unsigned int n_qpoints = context.get_element_qrule().n_points();
185 
186  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
187  {
188  // For the mass residual, we need to be a little careful.
189  // The time integrator is handling the time-discretization
190  // for us so we need to supply M(u_fixed)*u' for the residual.
191  // u_fixed will be given by the fixed_interior_value function
192  // while u' will be given by the interior_rate function.
193  libMesh::Real T_dot;
194  context.interior_rate(_temp_vars.T(), qp, T_dot);
195 
196  for (unsigned int i = 0; i != n_T_dofs; ++i)
197  {
198  F(i) -= JxW[qp]*(_rho*_Cp*T_dot*phi[i][qp] );
199 
200  if( compute_jacobian )
201  {
202  for (unsigned int j=0; j != n_T_dofs; j++)
203  {
204  // We're assuming rho, cp are constant w.r.t. T here.
205  M(i,j) -=
206  context.get_elem_solution_rate_derivative()
207  * JxW[qp]*_rho*_Cp*phi[j][qp]*phi[i][qp] ;
208  }
209  }// End of check on Jacobian
210 
211  } // End of element dof loop
212 
213  } // End of the quadrature point loop
214 
215  return;
216  }
217 
218  template<class K>
220  ( const std::string & param_name,
222  const
223  {
224  ParameterUser::register_parameter(param_name, param_pointer);
225  _k.register_parameter(param_name, param_pointer);
226  }
227 
228 
229 } // namespace GRINS
230 
231 // Instantiate
232 INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatConduction);
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
libMesh::Number _Cp
PrimitiveTempFEVariables & _temp_vars
virtual void mass_residual(bool compute_jacobian, AssemblyContext &context)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
libMesh::Number _rho
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
void check_var_subdomain_consistency(const FEVariablesBase &var) const
Check that var is enabled on at least the subdomains this Physics is.
Definition: physics.C:174
INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatConduction)
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.

Generated on Tue Dec 19 2017 12:47:28 for GRINS-0.8.0 by  doxygen 1.8.9.1