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

#include <heat_conduction.h>

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

Public Member Functions

 HeatConduction (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
 ~HeatConduction ()
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Set which variables are 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

PrimitiveTempFEVariables_temp_vars
 
libMesh::Number _rho
 
libMesh::Number _Cp
 
Conductivity _k
 Conductivity. More...
 
- 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

 HeatConduction ()
 

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::HeatConduction< Conductivity >

Definition at line 36 of file heat_conduction.h.

Constructor & Destructor Documentation

template<class K >
GRINS::HeatConduction< K >::HeatConduction ( const GRINS::PhysicsName physics_name,
const GetPot &  input 
)

Definition at line 47 of file heat_conduction.C.

References GRINS::HeatConduction< Conductivity >::_Cp, GRINS::Physics::_ic_handler, GRINS::HeatConduction< Conductivity >::_rho, GRINS::HeatConduction< Conductivity >::_temp_vars, GRINS::Physics::check_var_subdomain_consistency(), GRINS::PhysicsNaming::heat_conduction(), GRINS::MaterialsParsing::read_density(), and GRINS::MaterialsParsing::read_specific_heat().

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),
53  {
55 
57 
58  // This is deleted in the base class
59  this->_ic_handler = new GenericICHandler( physics_name, input );
60 
62  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
libMesh::Number _Cp
PrimitiveTempFEVariables & _temp_vars
libMesh::Number _rho
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
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.
static std::string temp_variable_name(const GetPot &input, const std::string &subsection_name, const SECTION_TYPE section_type)
static std::string material_name(const GetPot &input, const std::string &physics)
Get the name of the material in the Physics/physics section.
static PhysicsName heat_conduction()
Conductivity _k
Conductivity.
template<class Conductivity >
GRINS::HeatConduction< Conductivity >::~HeatConduction ( )
inline

Definition at line 42 of file heat_conduction.h.

42 {};
template<class Conductivity >
GRINS::HeatConduction< Conductivity >::HeatConduction ( )
private

Member Function Documentation

template<class K >
void GRINS::HeatConduction< K >::element_time_derivative ( bool  compute_jacobian,
AssemblyContext context 
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 90 of file heat_conduction.C.

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  }
PrimitiveTempFEVariables & _temp_vars
VariableIndex T() const
Conductivity _k
Conductivity.
template<class K >
void GRINS::HeatConduction< K >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::Physics.

Definition at line 72 of file heat_conduction.C.

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  }
PrimitiveTempFEVariables & _temp_vars
VariableIndex T() const
template<class K >
void GRINS::HeatConduction< K >::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 161 of file heat_conduction.C.

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  }
libMesh::Number _Cp
PrimitiveTempFEVariables & _temp_vars
VariableIndex T() const
libMesh::Number _rho
template<class K >
void GRINS::HeatConduction< K >::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 220 of file heat_conduction.C.

References GRINS::ParameterUser::register_parameter().

223  {
224  ParameterUser::register_parameter(param_name, param_pointer);
225  _k.register_parameter(param_name, param_pointer);
226  }
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.
Conductivity _k
Conductivity.
template<class K >
void GRINS::HeatConduction< K >::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtual

Set which variables are time evolving.

Set those variables which evolve in time (as opposed to variables that behave like constraints). This is done separately from init_variables() because the MultiphysicsSystem must initialize its base class before time_evolving variables can be set. Default implementation is no time evolving variables.

Reimplemented from GRINS::Physics.

Definition at line 65 of file heat_conduction.C.

66  {
67  // Tell the system to march temperature forward in time
68  system->time_evolving(_temp_vars.T(), 1);
69  }
PrimitiveTempFEVariables & _temp_vars
VariableIndex T() const

Member Data Documentation

template<class Conductivity >
libMesh::Number GRINS::HeatConduction< Conductivity >::_Cp
protected
template<class Conductivity >
Conductivity GRINS::HeatConduction< Conductivity >::_k
protected

Conductivity.

Definition at line 73 of file heat_conduction.h.

template<class Conductivity >
libMesh::Number GRINS::HeatConduction< Conductivity >::_rho
protected
template<class Conductivity >
PrimitiveTempFEVariables& GRINS::HeatConduction< Conductivity >::_temp_vars
protected

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