GRINS-0.7.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-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
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(input, PhysicsNaming::incompressible_navier_stokes()),
53  _press_var(input,PhysicsNaming::incompressible_navier_stokes(), true /*is_constraint_var*/),
54  _temp_vars(input, PhysicsNaming::axisymmetric_heat_transfer()),
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  this->register_variables();
62  }
63 
64  template< class Conductivity>
66  {
68 
70  }
71 
72  template< class Conductivity>
74  {
76  this->_press_var);
78  this->_flow_vars);
80  this->_temp_vars);
81  }
82 
83  template< class Conductivity>
84  void AxisymmetricHeatTransfer<Conductivity>::init_variables( libMesh::FEMSystem* system )
85  {
86  // Get libMesh to assign an index for each variable
87  this->_dim = system->get_mesh().mesh_dimension();
88 
89  this->_temp_vars.init(system);
90  this->_flow_vars.init(system);
91  this->_press_var.init(system);
92  }
93 
94  template< class Conductivity>
96  {
97  // Tell the system to march temperature forward in time
98  system->time_evolving(this->_temp_vars.T());
99  return;
100  }
101 
102  template< class Conductivity>
104  {
105  // We should prerequest all the data
106  // we will need to build the linear system
107  // or evaluate a quantity of interest.
108  context.get_element_fe(_temp_vars.T())->get_JxW();
109  context.get_element_fe(_temp_vars.T())->get_phi();
110  context.get_element_fe(_temp_vars.T())->get_dphi();
111  context.get_element_fe(_temp_vars.T())->get_xyz();
112 
113  context.get_side_fe(_temp_vars.T())->get_JxW();
114  context.get_side_fe(_temp_vars.T())->get_phi();
115  context.get_side_fe(_temp_vars.T())->get_dphi();
116  context.get_side_fe(_temp_vars.T())->get_xyz();
117 
118  // _u_var is registered so can we assume things related to _u_var
119  // are available in FEMContext
120 
121  return;
122  }
123 
124  template< class Conductivity>
126  AssemblyContext& context,
127  CachedValues& /*cache*/ )
128  {
129 #ifdef GRINS_USE_GRVY_TIMERS
130  this->_timer->BeginTimer("AxisymmetricHeatTransfer::element_time_derivative");
131 #endif
132 
133  // The number of local degrees of freedom in each variable.
134  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
135  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
136 
137  //TODO: check n_T_dofs is same as n_u_dofs, n_v_dofs, n_w_dofs
138 
139  // We get some references to cell-specific data that
140  // will be used to assemble the linear system.
141 
142  // Element Jacobian * quadrature weights for interior integration.
143  const std::vector<libMesh::Real> &JxW =
144  context.get_element_fe(_temp_vars.T())->get_JxW();
145 
146  // The temperature shape functions at interior quadrature points.
147  const std::vector<std::vector<libMesh::Real> >& T_phi =
148  context.get_element_fe(_temp_vars.T())->get_phi();
149 
150  // The velocity shape functions at interior quadrature points.
151  const std::vector<std::vector<libMesh::Real> >& vel_phi =
152  context.get_element_fe(_flow_vars.u())->get_phi();
153 
154  // The temperature shape function gradients (in global coords.)
155  // at interior quadrature points.
156  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
157  context.get_element_fe(_temp_vars.T())->get_dphi();
158 
159  // Physical location of the quadrature points
160  const std::vector<libMesh::Point>& u_qpoint =
161  context.get_element_fe(_flow_vars.u())->get_xyz();
162 
163  // The subvectors and submatrices we need to fill:
164  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(_temp_vars.T()); // R_{T}
165 
166  libMesh::DenseSubMatrix<libMesh::Number> &KTT = context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T()); // R_{T},{T}
167 
168  libMesh::DenseSubMatrix<libMesh::Number> &KTr = context.get_elem_jacobian(_temp_vars.T(), _flow_vars.u()); // R_{T},{r}
169  libMesh::DenseSubMatrix<libMesh::Number> &KTz = context.get_elem_jacobian(_temp_vars.T(), _flow_vars.v()); // R_{T},{z}
170 
171 
172  // Now we will build the element Jacobian and residual.
173  // Constructing the residual requires the solution and its
174  // gradient from the previous timestep. This must be
175  // calculated at each quadrature point by summing the
176  // solution degree-of-freedom values by the appropriate
177  // weight functions.
178  unsigned int n_qpoints = context.get_element_qrule().n_points();
179 
180  for (unsigned int qp=0; qp != n_qpoints; qp++)
181  {
182  const libMesh::Number r = u_qpoint[qp](0);
183 
184  // Compute the solution & its gradient at the old Newton iterate.
185  libMesh::Number u_r, u_z;
186  u_r = context.interior_value(_flow_vars.u(), qp);
187  u_z = context.interior_value(_flow_vars.v(), qp);
188 
189  libMesh::Gradient grad_T;
190  grad_T = context.interior_gradient(_temp_vars.T(), qp);
191 
192  libMesh::NumberVectorValue U (u_r,u_z);
193 
194  libMesh::Number k = this->_k( context, qp );
195 
196  // FIXME - once we have T-dependent k, we'll need its
197  // derivatives in Jacobians
198  // libMesh::Number dk_dT = this->_k.deriv( T );
199 
200  // First, an i-loop over the degrees of freedom.
201  for (unsigned int i=0; i != n_T_dofs; i++)
202  {
203  FT(i) += JxW[qp]*r*
204  (-_rho*_Cp*T_phi[i][qp]*(U*grad_T) // convection term
205  -k*(T_gradphi[i][qp]*grad_T) ); // diffusion term
206 
207  if (compute_jacobian)
208  {
209  libmesh_assert (context.get_elem_solution_derivative() == 1.0);
210 
211  for (unsigned int j=0; j != n_T_dofs; j++)
212  {
213  // TODO: precompute some terms like:
214  // _rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])
215 
216  KTT(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*
217  (-_rho*_Cp*T_phi[i][qp]*(U*T_gradphi[j][qp]) // convection term
218  -k*(T_gradphi[i][qp]*T_gradphi[j][qp])); // diffusion term
219  } // end of the inner dof (j) loop
220 
221 #if 0
222  if( dk_dT != 0.0 )
223  {
224  for (unsigned int j=0; j != n_T_dofs; j++)
225  {
226  // TODO: precompute some terms like:
227  KTT(i,j) -= JxW[qp] * context.get_elem_solution_derivative() *r*( dk_dT*T_phi[j][qp]*T_gradphi[i][qp]*grad_T );
228  }
229  }
230 #endif
231 
232  // Matrix contributions for the Tu, Tv and Tw couplings (n_T_dofs same as n_u_dofs, n_v_dofs and n_w_dofs)
233  for (unsigned int j=0; j != n_u_dofs; j++)
234  {
235  KTr(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(0)));
236  KTz(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(1)));
237  } // end of the inner dof (j) loop
238 
239  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
240 
241  } // end of the outer dof (i) loop
242  } // end of the quadrature point (qp) loop
243 
244 #ifdef GRINS_USE_GRVY_TIMERS
245  this->_timer->EndTimer("AxisymmetricHeatTransfer::element_time_derivative");
246 #endif
247 
248  return;
249  }
250 
251  template< class Conductivity>
253  AssemblyContext& context,
254  CachedValues& /*cache*/ )
255  {
256 #ifdef GRINS_USE_GRVY_TIMERS
257  this->_timer->BeginTimer("AxisymmetricHeatTransfer::mass_residual");
258 #endif
259 
260  // First we get some references to cell-specific data that
261  // will be used to assemble the linear system.
262 
263  // Element Jacobian * quadrature weights for interior integration
264  const std::vector<libMesh::Real> &JxW =
265  context.get_element_fe(_temp_vars.T())->get_JxW();
266 
267  // The shape functions at interior quadrature points.
268  const std::vector<std::vector<libMesh::Real> >& phi =
269  context.get_element_fe(_temp_vars.T())->get_phi();
270 
271  // The number of local degrees of freedom in each variable
272  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
273 
274  // Physical location of the quadrature points
275  const std::vector<libMesh::Point>& u_qpoint =
276  context.get_element_fe(_flow_vars.u())->get_xyz();
277 
278  // The subvectors and submatrices we need to fill:
279  libMesh::DenseSubVector<libMesh::Real> &F =
280  context.get_elem_residual(_temp_vars.T());
281 
282  libMesh::DenseSubMatrix<libMesh::Real> &M =
283  context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T());
284 
285  unsigned int n_qpoints = context.get_element_qrule().n_points();
286 
287  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
288  {
289  const libMesh::Number r = u_qpoint[qp](0);
290 
291  // For the mass residual, we need to be a little careful.
292  // The time integrator is handling the time-discretization
293  // for us so we need to supply M(u_fixed)*u' for the residual.
294  // u_fixed will be given by the fixed_interior_value function
295  // while u' will be given by the interior_rate function.
296  libMesh::Real T_dot;
297  context.interior_rate(_temp_vars.T(), qp, T_dot);
298 
299  for (unsigned int i = 0; i != n_T_dofs; ++i)
300  {
301  F(i) -= JxW[qp]*r*(_rho*_Cp*T_dot*phi[i][qp] );
302 
303  if( compute_jacobian )
304  {
305  for (unsigned int j=0; j != n_T_dofs; j++)
306  {
307  // We're assuming rho, cp are constant w.r.t. T here.
308  M(i,j) -= JxW[qp] * context.get_elem_solution_rate_derivative() *r*_rho*_Cp*phi[j][qp]*phi[i][qp] ;
309  }
310  }// End of check on Jacobian
311 
312  } // End of element dof loop
313 
314  } // End of the quadrature point loop
315 
316 #ifdef GRINS_USE_GRVY_TIMERS
317  this->_timer->EndTimer("AxisymmetricHeatTransfer::mass_residual");
318 #endif
319 
320  return;
321  }
322 
323  template< class Conductivity>
325  ( const std::string & param_name,
327  const
328  {
329  ParameterUser::register_parameter(param_name, param_pointer);
330  _k.register_parameter(param_name, param_pointer);
331  }
332 
333 
334 } // namespace GRINS
335 
336 // Instantiate
337 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:269
INSTANTIATE_HEAT_TRANSFER_SUBCLASS(AxisymmetricHeatTransfer)
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:107
static std::string temperature_section()
virtual void init_variables(libMesh::FEMSystem *system)
Initialization AxisymmetricHeatTransfer variables.
static void check_and_register_variable(const std::string &var_name, const FEVariablesBase &variable)
First check if var_name is registered and then register.
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.
static std::string velocity_section()
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.
Helper functions for parsing material properties.
static std::string pressure_section()
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...
static PhysicsName axisymmetric_heat_transfer()
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
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 Thu Jun 2 2016 21:52:27 for GRINS-0.7.0 by  doxygen 1.8.10