GRINS-0.7.0
List of all members | Public Member Functions | Protected Attributes | Private Member Functions
GRINS::AxisymmetricHeatTransfer< Conductivity > Class Template Reference

Physics class for Axisymmetric Heat Transfer. More...

#include <axisym_heat_transfer.h>

Inheritance diagram for GRINS::AxisymmetricHeatTransfer< Conductivity >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::AxisymmetricHeatTransfer< Conductivity >:
Collaboration graph
[legend]

Public Member Functions

 AxisymmetricHeatTransfer (const std::string &physics_name, const GetPot &input)
 
 ~AxisymmetricHeatTransfer ()
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialization AxisymmetricHeatTransfer variables. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for element interiors. More...
 
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. More...
 
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. More...
 
- Public Member Functions inherited from GRINS::Physics
 Physics (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~Physics ()
 
virtual bool enabled_on_elem (const libMesh::Elem *elem)
 Find if current physics is active on supplied element. More...
 
void set_is_steady (bool is_steady)
 Sets whether this physics is to be solved with a steady solver or not. More...
 
bool is_steady () const
 Returns whether or not this physics is being solved with a steady solver. More...
 
virtual void auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. More...
 
virtual void side_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_element_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_damping_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

Protected Attributes

unsigned int _dim
 Physical dimension of problem. More...
 
VelocityFEVariables _flow_vars
 
PressureFEVariable _press_var
 
PrimitiveTempFEVariables _temp_vars
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
libMesh::Number _Cp
 
Conductivity _k
 
- Protected Attributes inherited from GRINS::Physics
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 

Private Member Functions

 AxisymmetricHeatTransfer ()
 
void read_input_options (const GetPot &input)
 Read options from GetPot input file. More...
 
void register_variables ()
 

Additional Inherited Members

- Static Public Member Functions inherited from GRINS::Physics
static void set_is_axisymmetric (bool is_axisymmetric)
 Set whether we should treat the problem as axisymmetric. More...
 
static bool is_axisymmetric ()
 
- Static Public Attributes inherited from GRINS::ParameterUser
static std::string zero_vector_function = std::string("{0}")
 A parseable function string with LIBMESH_DIM components, all 0. More...
 
- Protected Member Functions inherited from GRINS::Physics
libMesh::UniquePtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe (const libMesh::Elem *elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
 
void parse_enabled_subdomains (const GetPot &input, const std::string &physics_name)
 
- Static Protected Attributes inherited from GRINS::Physics
static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 
static bool _is_axisymmetric = false
 Caches whether we are solving an axisymmetric problem or not. More...
 

Detailed Description

template<class Conductivity>
class GRINS::AxisymmetricHeatTransfer< Conductivity >

Physics class for Axisymmetric Heat Transfer.

Definition at line 49 of file axisym_heat_transfer.h.

Constructor & Destructor Documentation

template<class Conductivity >
GRINS::AxisymmetricHeatTransfer< Conductivity >::AxisymmetricHeatTransfer ( const std::string &  physics_name,
const GetPot &  input 
)

Definition at line 49 of file axisym_heat_transfer.C.

References GRINS::Physics::_ic_handler, GRINS::AxisymmetricHeatTransfer< Conductivity >::read_input_options(), and GRINS::AxisymmetricHeatTransfer< Conductivity >::register_variables().

51  : Physics(physics_name, input),
53  _press_var(input,PhysicsNaming::incompressible_navier_stokes(), true /*is_constraint_var*/),
56  {
57  this->read_input_options(input);
58 
59  this->_ic_handler = new GenericICHandler( physics_name, input );
60 
61  this->register_variables();
62  }
void read_input_options(const GetPot &input)
Read options from GetPot input file.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:269
static PhysicsName incompressible_navier_stokes()
PrimitiveTempFEVariables _temp_vars
static PhysicsName axisymmetric_heat_transfer()
static std::string material_name(const GetPot &input, const std::string &physics)
Get the name of the material in the Physics/physics section.
template<class Conductivity >
GRINS::AxisymmetricHeatTransfer< Conductivity >::~AxisymmetricHeatTransfer ( )
inline

Definition at line 55 of file axisym_heat_transfer.h.

55 {};
template<class Conductivity >
GRINS::AxisymmetricHeatTransfer< Conductivity >::AxisymmetricHeatTransfer ( )
private

Member Function Documentation

template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::element_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Time dependent part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 125 of file axisym_heat_transfer.C.

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  }
libMesh::Number _rho
Material parameters, read from input.
PrimitiveTempFEVariables _temp_vars
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::Physics.

Definition at line 103 of file axisym_heat_transfer.C.

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  }
PrimitiveTempFEVariables _temp_vars
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::init_variables ( libMesh::FEMSystem *  system)
virtual

Initialization AxisymmetricHeatTransfer variables.

Add velocity and pressure variables to system.

Implements GRINS::Physics.

Definition at line 84 of file axisym_heat_transfer.C.

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  }
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
unsigned int _dim
Physical dimension of problem.
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
PrimitiveTempFEVariables _temp_vars
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part.

Reimplemented from GRINS::Physics.

Definition at line 252 of file axisym_heat_transfer.C.

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  }
libMesh::Number _rho
Material parameters, read from input.
PrimitiveTempFEVariables _temp_vars
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::read_input_options ( const GetPot &  input)
private

Read options from GetPot input file.

Definition at line 65 of file axisym_heat_transfer.C.

References GRINS::PhysicsNaming::axisymmetric_heat_transfer(), GRINS::MaterialsParsing::read_density(), and GRINS::MaterialsParsing::read_specific_heat().

Referenced by GRINS::AxisymmetricHeatTransfer< Conductivity >::AxisymmetricHeatTransfer().

66  {
68 
70  }
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.
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.
libMesh::Number _rho
Material parameters, read from input.
static PhysicsName axisymmetric_heat_transfer()
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiAccessor< libMesh::Number > &  param_pointer 
) const
virtual

Each subclass will register its copy of an independent.

Reimplemented from GRINS::ParameterUser.

Definition at line 325 of file axisym_heat_transfer.C.

References GRINS::ParameterUser::register_parameter().

328  {
329  ParameterUser::register_parameter(param_name, param_pointer);
330  _k.register_parameter(param_name, param_pointer);
331  }
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.
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::register_variables ( )
private

Definition at line 73 of file axisym_heat_transfer.C.

References GRINS::GRINSPrivate::VariableWarehouse::check_and_register_variable(), GRINS::VariablesParsing::pressure_section(), GRINS::VariablesParsing::temperature_section(), and GRINS::VariablesParsing::velocity_section().

Referenced by GRINS::AxisymmetricHeatTransfer< Conductivity >::AxisymmetricHeatTransfer().

74  {
76  this->_press_var);
78  this->_flow_vars);
80  this->_temp_vars);
81  }
static std::string temperature_section()
static void check_and_register_variable(const std::string &var_name, const FEVariablesBase &variable)
First check if var_name is registered and then register.
static std::string velocity_section()
static std::string pressure_section()
PrimitiveTempFEVariables _temp_vars
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtual

Sets velocity variables to be time-evolving.

Reimplemented from GRINS::Physics.

Definition at line 95 of file axisym_heat_transfer.C.

96  {
97  // Tell the system to march temperature forward in time
98  system->time_evolving(this->_temp_vars.T());
99  return;
100  }
PrimitiveTempFEVariables _temp_vars

Member Data Documentation

template<class Conductivity >
libMesh::Number GRINS::AxisymmetricHeatTransfer< Conductivity >::_Cp
protected

Definition at line 104 of file axisym_heat_transfer.h.

template<class Conductivity >
unsigned int GRINS::AxisymmetricHeatTransfer< Conductivity >::_dim
protected

Physical dimension of problem.

Todo:
Make this static member of base class?

Definition at line 93 of file axisym_heat_transfer.h.

template<class Conductivity >
VelocityFEVariables GRINS::AxisymmetricHeatTransfer< Conductivity >::_flow_vars
protected

Definition at line 95 of file axisym_heat_transfer.h.

template<class Conductivity >
Conductivity GRINS::AxisymmetricHeatTransfer< Conductivity >::_k
protected

Definition at line 106 of file axisym_heat_transfer.h.

template<class Conductivity >
PressureFEVariable GRINS::AxisymmetricHeatTransfer< Conductivity >::_press_var
protected

Definition at line 96 of file axisym_heat_transfer.h.

template<class Conductivity >
libMesh::Number GRINS::AxisymmetricHeatTransfer< Conductivity >::_rho
protected

Material parameters, read from input.

Todo:
Need to generalize material parameters. Right now they are assumed constant
Todo:
Shouldn't this rho be the same as the one in the flow? Need to figure out how to have those shared

Definition at line 104 of file axisym_heat_transfer.h.

template<class Conductivity >
PrimitiveTempFEVariables GRINS::AxisymmetricHeatTransfer< Conductivity >::_temp_vars
protected

Definition at line 97 of file axisym_heat_transfer.h.


The documentation for this class was generated from the following files:

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