GRINS-0.8.0
axisym_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-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
28 
29 // GRINS
30 #include "grins/assembly_context.h"
32 #include "grins/grins_enums.h"
37 
38 // libMesh
39 #include "libmesh/utility.h"
40 #include "libmesh/string_to_enum.h"
41 #include "libmesh/getpot.h"
42 #include "libmesh/fem_system.h"
43 #include "libmesh/quadrature.h"
44 
45 namespace GRINS
46 {
47 
48  template< class Conductivity>
50  const GetPot& input)
51  : Physics(physics_name, input),
52  _flow_vars(GRINSPrivate::VariableWarehouse::get_variable_subclass<VelocityVariable>(VariablesParsing::velocity_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
53  _press_var(GRINSPrivate::VariableWarehouse::get_variable_subclass<PressureFEVariable>(VariablesParsing::press_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
54  _temp_vars(GRINSPrivate::VariableWarehouse::get_variable_subclass<PrimitiveTempFEVariables>(VariablesParsing::temp_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
55  _k(input,MaterialsParsing::material_name(input,PhysicsNaming::axisymmetric_heat_transfer()))
56  {
57  this->read_input_options(input);
58 
59  this->_ic_handler = new GenericICHandler( physics_name, input );
60  }
61 
62  template< class Conductivity>
64  {
66 
68  }
69 
70  template< class Conductivity>
72  {
73  // Tell the system to march temperature forward in time
74  system->time_evolving(this->_temp_vars.T(), 1);
75  }
76 
77  template< class Conductivity>
79  {
80  // We should prerequest all the data
81  // we will need to build the linear system
82  // or evaluate a quantity of interest.
83  context.get_element_fe(_temp_vars.T())->get_JxW();
84  context.get_element_fe(_temp_vars.T())->get_phi();
85  context.get_element_fe(_temp_vars.T())->get_dphi();
86  context.get_element_fe(_temp_vars.T())->get_xyz();
87 
88  context.get_side_fe(_temp_vars.T())->get_JxW();
89  context.get_side_fe(_temp_vars.T())->get_phi();
90  context.get_side_fe(_temp_vars.T())->get_dphi();
91  context.get_side_fe(_temp_vars.T())->get_xyz();
92  }
93 
94  template< class Conductivity>
96  ( bool compute_jacobian, AssemblyContext & context )
97  {
98  // The number of local degrees of freedom in each variable.
99  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
100  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
101 
102  //TODO: check n_T_dofs is same as n_u_dofs, n_v_dofs, n_w_dofs
103 
104  // We get some references to cell-specific data that
105  // will be used to assemble the linear system.
106 
107  // Element Jacobian * quadrature weights for interior integration.
108  const std::vector<libMesh::Real> &JxW =
109  context.get_element_fe(_temp_vars.T())->get_JxW();
110 
111  // The temperature shape functions at interior quadrature points.
112  const std::vector<std::vector<libMesh::Real> >& T_phi =
113  context.get_element_fe(_temp_vars.T())->get_phi();
114 
115  // The velocity shape functions at interior quadrature points.
116  const std::vector<std::vector<libMesh::Real> >& vel_phi =
117  context.get_element_fe(_flow_vars.u())->get_phi();
118 
119  // The temperature shape function gradients (in global coords.)
120  // at interior quadrature points.
121  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
122  context.get_element_fe(_temp_vars.T())->get_dphi();
123 
124  // Physical location of the quadrature points
125  const std::vector<libMesh::Point>& u_qpoint =
126  context.get_element_fe(_flow_vars.u())->get_xyz();
127 
128  // The subvectors and submatrices we need to fill:
129  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(_temp_vars.T()); // R_{T}
130 
131  libMesh::DenseSubMatrix<libMesh::Number> &KTT = context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T()); // R_{T},{T}
132 
133  libMesh::DenseSubMatrix<libMesh::Number> &KTr = context.get_elem_jacobian(_temp_vars.T(), _flow_vars.u()); // R_{T},{r}
134  libMesh::DenseSubMatrix<libMesh::Number> &KTz = context.get_elem_jacobian(_temp_vars.T(), _flow_vars.v()); // R_{T},{z}
135 
136 
137  // Now we will build the element Jacobian and residual.
138  // Constructing the residual requires the solution and its
139  // gradient from the previous timestep. This must be
140  // calculated at each quadrature point by summing the
141  // solution degree-of-freedom values by the appropriate
142  // weight functions.
143  unsigned int n_qpoints = context.get_element_qrule().n_points();
144 
145  for (unsigned int qp=0; qp != n_qpoints; qp++)
146  {
147  const libMesh::Number r = u_qpoint[qp](0);
148 
149  // Compute the solution & its gradient at the old Newton iterate.
150  libMesh::Number u_r, u_z;
151  u_r = context.interior_value(_flow_vars.u(), qp);
152  u_z = context.interior_value(_flow_vars.v(), qp);
153 
154  libMesh::Gradient grad_T;
155  grad_T = context.interior_gradient(_temp_vars.T(), qp);
156 
157  libMesh::NumberVectorValue U (u_r,u_z);
158 
159  libMesh::Number k = this->_k( context, qp );
160 
161  // FIXME - once we have T-dependent k, we'll need its
162  // derivatives in Jacobians
163  // libMesh::Number dk_dT = this->_k.deriv( T );
164 
165  // First, an i-loop over the degrees of freedom.
166  for (unsigned int i=0; i != n_T_dofs; i++)
167  {
168  FT(i) += JxW[qp]*r*
169  (-_rho*_Cp*T_phi[i][qp]*(U*grad_T) // convection term
170  -k*(T_gradphi[i][qp]*grad_T) ); // diffusion term
171 
172  if (compute_jacobian)
173  {
174  libmesh_assert (context.get_elem_solution_derivative() == 1.0);
175 
176  for (unsigned int j=0; j != n_T_dofs; j++)
177  {
178  // TODO: precompute some terms like:
179  // _rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])
180 
181  KTT(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*
182  (-_rho*_Cp*T_phi[i][qp]*(U*T_gradphi[j][qp]) // convection term
183  -k*(T_gradphi[i][qp]*T_gradphi[j][qp])); // diffusion term
184  } // end of the inner dof (j) loop
185 
186 #if 0
187  if( dk_dT != 0.0 )
188  {
189  for (unsigned int j=0; j != n_T_dofs; j++)
190  {
191  // TODO: precompute some terms like:
192  KTT(i,j) -= JxW[qp] * context.get_elem_solution_derivative() *r*( dk_dT*T_phi[j][qp]*T_gradphi[i][qp]*grad_T );
193  }
194  }
195 #endif
196 
197  // Matrix contributions for the Tu, Tv and Tw couplings (n_T_dofs same as n_u_dofs, n_v_dofs and n_w_dofs)
198  for (unsigned int j=0; j != n_u_dofs; j++)
199  {
200  KTr(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(0)));
201  KTz(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(1)));
202  } // end of the inner dof (j) loop
203 
204  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
205 
206  } // end of the outer dof (i) loop
207  } // end of the quadrature point (qp) loop
208  }
209 
210  template< class Conductivity>
212  ( bool compute_jacobian, AssemblyContext & context )
213  {
214  // First we get some references to cell-specific data that
215  // will be used to assemble the linear system.
216 
217  // Element Jacobian * quadrature weights for interior integration
218  const std::vector<libMesh::Real> &JxW =
219  context.get_element_fe(_temp_vars.T())->get_JxW();
220 
221  // The shape functions at interior quadrature points.
222  const std::vector<std::vector<libMesh::Real> >& phi =
223  context.get_element_fe(_temp_vars.T())->get_phi();
224 
225  // The number of local degrees of freedom in each variable
226  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
227 
228  // Physical location of the quadrature points
229  const std::vector<libMesh::Point>& u_qpoint =
230  context.get_element_fe(_flow_vars.u())->get_xyz();
231 
232  // The subvectors and submatrices we need to fill:
233  libMesh::DenseSubVector<libMesh::Real> &F =
234  context.get_elem_residual(_temp_vars.T());
235 
236  libMesh::DenseSubMatrix<libMesh::Real> &M =
237  context.get_elem_jacobian(_temp_vars.T(), _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  const libMesh::Number r = u_qpoint[qp](0);
244 
245  // For the mass residual, we need to be a little careful.
246  // The time integrator is handling the time-discretization
247  // for us so we need to supply M(u_fixed)*u' for the residual.
248  // u_fixed will be given by the fixed_interior_value function
249  // while u' will be given by the interior_rate function.
250  libMesh::Real T_dot;
251  context.interior_rate(_temp_vars.T(), qp, T_dot);
252 
253  for (unsigned int i = 0; i != n_T_dofs; ++i)
254  {
255  F(i) -= JxW[qp]*r*(_rho*_Cp*T_dot*phi[i][qp] );
256 
257  if( compute_jacobian )
258  {
259  for (unsigned int j=0; j != n_T_dofs; j++)
260  {
261  // We're assuming rho, cp are constant w.r.t. T here.
262  M(i,j) -= JxW[qp] * context.get_elem_solution_rate_derivative() *r*_rho*_Cp*phi[j][qp]*phi[i][qp] ;
263  }
264  }// End of check on Jacobian
265 
266  } // End of element dof loop
267 
268  } // End of the quadrature point loop
269  }
270 
271  template< class Conductivity>
273  ( const std::string & param_name,
275  const
276  {
277  ParameterUser::register_parameter(param_name, param_pointer);
278  _k.register_parameter(param_name, param_pointer);
279  }
280 
281 
282 } // namespace GRINS
283 
284 // Instantiate
285 INSTANTIATE_HEAT_TRANSFER_SUBCLASS(AxisymmetricHeatTransfer);
void read_input_options(const GetPot &input)
Read options from GetPot input file.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
INSTANTIATE_HEAT_TRANSFER_SUBCLASS(AxisymmetricHeatTransfer)
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
Base class for reading and handling initial conditions for physics classes.
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Sets velocity variables to be time-evolving.
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.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
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.
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 mass_residual(bool compute_jacobian, AssemblyContext &context)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
Helper functions for parsing material properties.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
static PhysicsName axisymmetric_heat_transfer()
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