GRINS-0.7.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 init_variables (libMesh::FEMSystem *system)
 Initialize variables for this physics. More...
 
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, 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
 
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 ()
 
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::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::PhysicsNaming::heat_conduction(), GRINS::MaterialsParsing::read_density(), GRINS::MaterialsParsing::read_specific_heat(), and GRINS::HeatConduction< Conductivity >::register_variables().

48  : Physics(physics_name,input),
50  _rho(0.0),
51  _Cp(0.0),
53  {
55 
57 
58  this->register_variables();
59 
60  // This is deleted in the base class
61  this->_ic_handler = new GenericICHandler( physics_name, input );
62  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:269
libMesh::Number _Cp
libMesh::Number _rho
PrimitiveTempFEVariables _temp_vars
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 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 K >
GRINS::HeatConduction< K >::~HeatConduction ( )

Definition at line 65 of file heat_conduction.C.

66  {
67  return;
68  }
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,
CachedValues cache 
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 117 of file heat_conduction.C.

120  {
121  // The number of local degrees of freedom in each variable.
122  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
123 
124  // We get some references to cell-specific data that
125  // will be used to assemble the linear system.
126 
127  // Element Jacobian * quadrature weights for interior integration.
128  const std::vector<libMesh::Real> &JxW =
129  context.get_element_fe(_temp_vars.T())->get_JxW();
130 
131  // The temperature shape function gradients (in global coords.)
132  // at interior quadrature points.
133  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
134  context.get_element_fe(_temp_vars.T())->get_dphi();
135 
136  // The subvectors and submatrices we need to fill:
137  //
138  // K_{\alpha \beta} = R_{\alpha},{\beta} = \partial{ R_{\alpha} } / \partial{ {\beta} } (where R denotes residual)
139  // e.g., for \alpha = T and \beta = v we get: K_{Tu} = R_{T},{u}
140  //
141 
142  libMesh::DenseSubMatrix<libMesh::Number> &KTT = context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T()); // R_{T},{T}
143 
144  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(_temp_vars.T()); // R_{T}
145 
146  // Now we will build the element Jacobian and residual.
147  // Constructing the residual requires the solution and its
148  // gradient from the previous timestep. This must be
149  // calculated at each quadrature point by summing the
150  // solution degree-of-freedom values by the appropriate
151  // weight functions.
152  unsigned int n_qpoints = context.get_element_qrule().n_points();
153 
154  for (unsigned int qp=0; qp != n_qpoints; qp++)
155  {
156  // Compute the solution & its gradient at the old Newton iterate.
157  libMesh::Gradient grad_T;
158  grad_T = context.interior_gradient(_temp_vars.T(), qp);
159 
160  // Compute the conductivity at this qp
161  libMesh::Real _k_qp = this->_k(context, qp);
162 
163  // First, an i-loop over the degrees of freedom.
164  for (unsigned int i=0; i != n_T_dofs; i++)
165  {
166  FT(i) += JxW[qp]*(-_k_qp*(T_gradphi[i][qp]*grad_T));
167 
168  if (compute_jacobian)
169  {
170  for (unsigned int j=0; j != n_T_dofs; j++)
171  {
172  // TODO: precompute some terms like:
173  // _rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])
174 
175  KTT(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
176  ( -_k_qp*(T_gradphi[i][qp]*T_gradphi[j][qp]) ); // diffusion term
177  } // end of the inner dof (j) loop
178 
179  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
180 
181  } // end of the outer dof (i) loop
182  } // end of the quadrature point (qp) loop
183 
184  return;
185  }
PrimitiveTempFEVariables _temp_vars
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 98 of file heat_conduction.C.

99  {
100  // We should prerequest all the data
101  // we will need to build the linear system
102  // or evaluate a quantity of interest.
103  context.get_element_fe(_temp_vars.T())->get_JxW();
104  context.get_element_fe(_temp_vars.T())->get_phi();
105  context.get_element_fe(_temp_vars.T())->get_dphi();
106  context.get_element_fe(_temp_vars.T())->get_xyz();
107 
108  context.get_side_fe(_temp_vars.T())->get_JxW();
109  context.get_side_fe(_temp_vars.T())->get_phi();
110  context.get_side_fe(_temp_vars.T())->get_dphi();
111  context.get_side_fe(_temp_vars.T())->get_xyz();
112 
113  return;
114  }
PrimitiveTempFEVariables _temp_vars
template<class K >
void GRINS::HeatConduction< K >::init_variables ( libMesh::FEMSystem *  system)
virtual

Initialize variables for this physics.

Implements GRINS::Physics.

Definition at line 78 of file heat_conduction.C.

79  {
80  // Get libMesh to assign an index for each variable
81  this->_dim = system->get_mesh().mesh_dimension();
82 
83  _temp_vars.init(system);
84 
85  return;
86  }
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
PrimitiveTempFEVariables _temp_vars
template<class K >
void GRINS::HeatConduction< K >::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 188 of file heat_conduction.C.

191  {
192  // First we get some references to cell-specific data that
193  // will be used to assemble the linear system.
194 
195  // Element Jacobian * quadrature weights for interior integration
196  const std::vector<libMesh::Real> &JxW =
197  context.get_element_fe(_temp_vars.T())->get_JxW();
198 
199  // The shape functions at interior quadrature points.
200  const std::vector<std::vector<libMesh::Real> >& phi =
201  context.get_element_fe(_temp_vars.T())->get_phi();
202 
203  // The number of local degrees of freedom in each variable
204  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
205 
206  // The subvectors and submatrices we need to fill:
207  libMesh::DenseSubVector<libMesh::Real> &F =
208  context.get_elem_residual(_temp_vars.T());
209 
210  libMesh::DenseSubMatrix<libMesh::Real> &M =
211  context.get_elem_jacobian(_temp_vars.T(), _temp_vars.T());
212 
213  unsigned int n_qpoints = context.get_element_qrule().n_points();
214 
215  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
216  {
217  // For the mass residual, we need to be a little careful.
218  // The time integrator is handling the time-discretization
219  // for us so we need to supply M(u_fixed)*u' for the residual.
220  // u_fixed will be given by the fixed_interior_value function
221  // while u' will be given by the interior_rate function.
222  libMesh::Real T_dot;
223  context.interior_rate(_temp_vars.T(), qp, T_dot);
224 
225  for (unsigned int i = 0; i != n_T_dofs; ++i)
226  {
227  F(i) -= JxW[qp]*(_rho*_Cp*T_dot*phi[i][qp] );
228 
229  if( compute_jacobian )
230  {
231  for (unsigned int j=0; j != n_T_dofs; j++)
232  {
233  // We're assuming rho, cp are constant w.r.t. T here.
234  M(i,j) -=
235  context.get_elem_solution_rate_derivative()
236  * JxW[qp]*_rho*_Cp*phi[j][qp]*phi[i][qp] ;
237  }
238  }// End of check on Jacobian
239 
240  } // End of element dof loop
241 
242  } // End of the quadrature point loop
243 
244  return;
245  }
libMesh::Number _Cp
libMesh::Number _rho
PrimitiveTempFEVariables _temp_vars
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 249 of file heat_conduction.C.

References GRINS::ParameterUser::register_parameter().

252  {
253  ParameterUser::register_parameter(param_name, param_pointer);
254  _k.register_parameter(param_name, param_pointer);
255  }
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 >::register_variables ( )
private

Definition at line 71 of file heat_conduction.C.

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

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

72  {
74  this->_temp_vars);
75  }
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.
PrimitiveTempFEVariables _temp_vars
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 89 of file heat_conduction.C.

90  {
91  // Tell the system to march temperature forward in time
92  system->time_evolving(_temp_vars.T());
93 
94  return;
95  }
PrimitiveTempFEVariables _temp_vars

Member Data Documentation

template<class Conductivity >
libMesh::Number GRINS::HeatConduction< Conductivity >::_Cp
protected
template<class Conductivity >
unsigned int GRINS::HeatConduction< Conductivity >::_dim
protected

Definition at line 73 of file heat_conduction.h.

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

Conductivity.

Definition at line 80 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

Definition at line 75 of file heat_conduction.h.


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

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