GRINS-0.8.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 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)
 Time dependent part(s) of physics for element interiors. More...
 
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. 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 void init_variables (libMesh::FEMSystem *)
 Initialize variables for this physics. More...
 
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 preassembly (MultiphysicsSystem &)
 Perform any necessary setup before element assembly begins. More...
 
virtual void reinit (MultiphysicsSystem &)
 Any reinitialization that needs to be done. More...
 
virtual void side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 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 (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
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

const VelocityVariable_flow_vars
 
const PressureFEVariable_press_var
 
const 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...
 

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)
 
void check_var_subdomain_consistency (const FEVariablesBase &var) const
 Check that var is enabled on at least the subdomains this Physics is. More...
 
- 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 48 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, and GRINS::AxisymmetricHeatTransfer< Conductivity >::read_input_options().

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))),
56  {
57  this->read_input_options(input);
58 
59  this->_ic_handler = new GenericICHandler( physics_name, input );
60  }
void read_input_options(const GetPot &input)
Read options from GetPot input file.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
static std::string velocity_variable_name(const GetPot &input, const std::string &subsection_name, const SECTION_TYPE section_type)
const PrimitiveTempFEVariables & _temp_vars
static std::string press_variable_name(const GetPot &input, const std::string &subsection_name, const SECTION_TYPE section_type)
const VelocityVariable & _flow_vars
static std::string temp_variable_name(const GetPot &input, const std::string &subsection_name, const SECTION_TYPE section_type)
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.
const PressureFEVariable & _press_var
template<class Conductivity >
GRINS::AxisymmetricHeatTransfer< Conductivity >::~AxisymmetricHeatTransfer ( )
inline

Definition at line 54 of file axisym_heat_transfer.h.

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

Member Function Documentation

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

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

Reimplemented from GRINS::Physics.

Definition at line 96 of file axisym_heat_transfer.C.

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  }
VariableIndex T() const
const PrimitiveTempFEVariables & _temp_vars
const VelocityVariable & _flow_vars
libMesh::Number _rho
Material parameters, read from input.
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 78 of file axisym_heat_transfer.C.

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  }
VariableIndex T() const
const PrimitiveTempFEVariables & _temp_vars
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::mass_residual ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 212 of file axisym_heat_transfer.C.

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

Read options from GetPot input file.

Definition at line 63 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().

64  {
66 
68  }
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 273 of file axisym_heat_transfer.C.

References GRINS::ParameterUser::register_parameter().

276  {
277  ParameterUser::register_parameter(param_name, param_pointer);
278  _k.register_parameter(param_name, param_pointer);
279  }
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 >::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtual

Sets velocity variables to be time-evolving.

Reimplemented from GRINS::Physics.

Definition at line 71 of file axisym_heat_transfer.C.

72  {
73  // Tell the system to march temperature forward in time
74  system->time_evolving(this->_temp_vars.T(), 1);
75  }
VariableIndex T() const
const PrimitiveTempFEVariables & _temp_vars

Member Data Documentation

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

Definition at line 91 of file axisym_heat_transfer.h.

template<class Conductivity >
const VelocityVariable& GRINS::AxisymmetricHeatTransfer< Conductivity >::_flow_vars
protected

Definition at line 82 of file axisym_heat_transfer.h.

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

Definition at line 93 of file axisym_heat_transfer.h.

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

Definition at line 83 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 91 of file axisym_heat_transfer.h.

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

Definition at line 84 of file axisym_heat_transfer.h.


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

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