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

Physics class for Heat Transfer. More...

#include <heat_transfer.h>

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

Public Member Functions

 HeatTransfer (const std::string &physics_name, const GetPot &input)
 
 ~HeatTransfer ()
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register postprocessing variables for HeatTransfer. 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 compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 Compute value of postprocessed quantities at libMesh::Point. More...
 
- Public Member Functions inherited from GRINS::HeatTransferBase< Conductivity >
 HeatTransferBase (const std::string &physics_name, const std::string &core_physics_name, const GetPot &input)
 
 ~HeatTransferBase ()
 
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 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 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 &)
 
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...
 

Private Member Functions

 HeatTransfer ()
 

Private Attributes

unsigned int _k_index
 Index from registering this postprocessed quantity. 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...
 
- Protected Attributes inherited from GRINS::HeatTransferBase< Conductivity >
unsigned int _dim
 Physical dimension of problem. More...
 
VelocityVariable_flow_vars
 
PressureFEVariable_press_var
 
PrimitiveTempFEVariables_temp_vars
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
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...
 
- 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::HeatTransfer< Conductivity >

Physics class for Heat Transfer.

Definition at line 40 of file heat_transfer.h.

Constructor & Destructor Documentation

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

Definition at line 45 of file heat_transfer.C.

References GRINS::Physics::_ic_handler.

46  : HeatTransferBase<K>(physics_name, PhysicsNaming::heat_transfer(), input),
47  _k_index(0)
48  {
49  this->_ic_handler = new GenericICHandler( physics_name, input );
50  }
static PhysicsName heat_transfer()
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
unsigned int _k_index
Index from registering this postprocessed quantity.
Definition: heat_transfer.h:74
template<class Conductivity >
GRINS::HeatTransfer< Conductivity >::~HeatTransfer ( )
inline

Definition at line 46 of file heat_transfer.h.

46 {};
template<class Conductivity >
GRINS::HeatTransfer< Conductivity >::HeatTransfer ( )
private

Member Function Documentation

template<class K >
void GRINS::HeatTransfer< K >::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtual

Compute value of postprocessed quantities at libMesh::Point.

Reimplemented from GRINS::Physics.

Definition at line 264 of file heat_transfer.C.

268  {
269  if( quantity_index == this->_k_index )
270  value = this->_k(point, context.get_time());
271  }
unsigned int _k_index
Index from registering this postprocessed quantity.
Definition: heat_transfer.h:74
Conductivity _k
Conductivity.
template<class K >
void GRINS::HeatTransfer< K >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 85 of file heat_transfer.C.

References GRINS::Physics::is_axisymmetric().

87  {
88  // The number of local degrees of freedom in each variable.
89  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
90  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
91 
92  //TODO: check n_T_dofs is same as n_u_dofs, n_v_dofs, n_w_dofs
93 
94  // We get some references to cell-specific data that
95  // will be used to assemble the linear system.
96 
97  // Element Jacobian * quadrature weights for interior integration.
98  const std::vector<libMesh::Real> &JxW =
99  context.get_element_fe(this->_temp_vars.T())->get_JxW();
100 
101  // The temperature shape functions at interior quadrature points.
102  const std::vector<std::vector<libMesh::Real> >& T_phi =
103  context.get_element_fe(this->_temp_vars.T())->get_phi();
104 
105  // The velocity shape functions at interior quadrature points.
106  const std::vector<std::vector<libMesh::Real> >& vel_phi =
107  context.get_element_fe(this->_flow_vars.u())->get_phi();
108 
109  // The temperature shape function gradients (in global coords.)
110  // at interior quadrature points.
111  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
112  context.get_element_fe(this->_temp_vars.T())->get_dphi();
113 
114  const std::vector<libMesh::Point>& u_qpoint =
115  context.get_element_fe(this->_flow_vars.u())->get_xyz();
116 
117  libMesh::DenseSubMatrix<libMesh::Number> &KTT = context.get_elem_jacobian(this->_temp_vars.T(), this->_temp_vars.T()); // R_{T},{T}
118 
119  libMesh::DenseSubMatrix<libMesh::Number> &KTu = context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.u()); // R_{T},{u}
120  libMesh::DenseSubMatrix<libMesh::Number> &KTv = context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.v()); // R_{T},{v}
121  libMesh::DenseSubMatrix<libMesh::Number>* KTw = NULL;
122 
123  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
124 
125  if( this->_flow_vars.dim() == 3 )
126  {
127  KTw = &context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.w()); // R_{T},{w}
128  }
129 
130  // Now we will build the element Jacobian and residual.
131  // Constructing the residual requires the solution and its
132  // gradient from the previous timestep. This must be
133  // calculated at each quadrature point by summing the
134  // solution degree-of-freedom values by the appropriate
135  // weight functions.
136  unsigned int n_qpoints = context.get_element_qrule().n_points();
137 
138  for (unsigned int qp=0; qp != n_qpoints; qp++)
139  {
140  // Compute the solution & its gradient at the old Newton iterate.
141  libMesh::Number u, v;
142  u = context.interior_value(this->_flow_vars.u(), qp);
143  v = context.interior_value(this->_flow_vars.v(), qp);
144 
145  libMesh::Gradient grad_T;
146  grad_T = context.interior_gradient(this->_temp_vars.T(), qp);
147 
148  libMesh::NumberVectorValue U (u,v);
149  if (this->_flow_vars.dim() == 3)
150  U(2) = context.interior_value(this->_flow_vars.w(), qp);
151 
152  const libMesh::Number r = u_qpoint[qp](0);
153 
154  libMesh::Real jac = JxW[qp];
155 
156  // Compute the conductivity at this qp
157  libMesh::Real _k_qp = this->_k(context, qp);
158 
160  {
161  jac *= r;
162  }
163 
164  // First, an i-loop over the degrees of freedom.
165  for (unsigned int i=0; i != n_T_dofs; i++)
166  {
167  FT(i) += jac *
168  (-this->_rho*this->_Cp*T_phi[i][qp]*(U*grad_T) // convection term
169  -_k_qp*(T_gradphi[i][qp]*grad_T) ); // diffusion term
170 
171  if (compute_jacobian)
172  {
173  for (unsigned int j=0; j != n_T_dofs; j++)
174  {
175  // TODO: precompute some terms like:
176  // this->_rho*this->_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])
177 
178  KTT(i,j) += jac * context.get_elem_solution_derivative() *
179  (-this->_rho*this->_Cp*T_phi[i][qp]*(U*T_gradphi[j][qp]) // convection term
180  -_k_qp*(T_gradphi[i][qp]*T_gradphi[j][qp])); // diffusion term
181  } // end of the inner dof (j) loop
182 
183  // Matrix contributions for the Tu, Tv and Tw couplings (n_T_dofs same as n_u_dofs, n_v_dofs and n_w_dofs)
184  for (unsigned int j=0; j != n_u_dofs; j++)
185  {
186  KTu(i,j) += jac * context.get_elem_solution_derivative()*(-this->_rho*this->_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(0)));
187  KTv(i,j) += jac * context.get_elem_solution_derivative()*(-this->_rho*this->_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(1)));
188  if (this->_flow_vars.dim() == 3)
189  (*KTw)(i,j) += jac * context.get_elem_solution_derivative()*(-this->_rho*this->_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(2)));
190  } // end of the inner dof (j) loop
191 
192  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
193 
194  } // end of the outer dof (i) loop
195  } // end of the quadrature point (qp) loop
196  }
static bool is_axisymmetric()
Definition: physics.h:132
VariableIndex T() const
VelocityVariable & _flow_vars
libMesh::Number _rho
Material parameters, read from input.
Conductivity _k
Conductivity.
unsigned int dim() const
Number of components.
PrimitiveTempFEVariables & _temp_vars
template<class K >
void GRINS::HeatTransfer< 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 200 of file heat_transfer.C.

References GRINS::Physics::is_axisymmetric().

201  {
202  // First we get some references to cell-specific data that
203  // will be used to assemble the linear system.
204 
205  // Element Jacobian * quadrature weights for interior integration
206  const std::vector<libMesh::Real> &JxW =
207  context.get_element_fe(this->_temp_vars.T())->get_JxW();
208 
209  // The shape functions at interior quadrature points.
210  const std::vector<std::vector<libMesh::Real> >& phi =
211  context.get_element_fe(this->_temp_vars.T())->get_phi();
212 
213  const std::vector<libMesh::Point>& u_qpoint =
214  context.get_element_fe(this->_flow_vars.u())->get_xyz();
215 
216  // The number of local degrees of freedom in each variable
217  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
218 
219  // The subvectors and submatrices we need to fill:
220  libMesh::DenseSubVector<libMesh::Real> &F = context.get_elem_residual(this->_temp_vars.T());
221 
222  libMesh::DenseSubMatrix<libMesh::Real> &M = context.get_elem_jacobian(this->_temp_vars.T(), this->_temp_vars.T());
223 
224  unsigned int n_qpoints = context.get_element_qrule().n_points();
225 
226  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
227  {
228  // For the mass residual, we need to be a little careful.
229  // The time integrator is handling the time-discretization
230  // for us so we need to supply M(u_fixed)*u' for the residual.
231  // u_fixed will be given by the fixed_interior_value function
232  // while u' will be given by the interior_rate function.
233  libMesh::Real T_dot;
234  context.interior_rate(this->_temp_vars.T(), qp, T_dot);
235 
236  const libMesh::Number r = u_qpoint[qp](0);
237 
238  libMesh::Real jac = JxW[qp];
239 
241  {
242  jac *= r;
243  }
244 
245  for (unsigned int i = 0; i != n_T_dofs; ++i)
246  {
247  F(i) -= this->_rho*this->_Cp*T_dot*phi[i][qp]*jac;
248 
249  if( compute_jacobian )
250  {
251  for (unsigned int j=0; j != n_T_dofs; j++)
252  {
253  // We're assuming rho, cp are constant w.r.t. T here.
254  M(i,j) -= this->_rho*this->_Cp*phi[j][qp]*phi[i][qp]*jac * context.get_elem_solution_rate_derivative();
255  }
256  }// End of check on Jacobian
257 
258  } // End of element dof loop
259 
260  } // End of the quadrature point loop
261  }
static bool is_axisymmetric()
Definition: physics.h:132
VariableIndex T() const
VelocityVariable & _flow_vars
libMesh::Number _rho
Material parameters, read from input.
PrimitiveTempFEVariables & _temp_vars
template<class K >
void GRINS::HeatTransfer< K >::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtual

Register postprocessing variables for HeatTransfer.

Reimplemented from GRINS::Physics.

Definition at line 53 of file heat_transfer.C.

References GRINS::PhysicsNaming::heat_transfer(), and GRINS::PostProcessedQuantities< NumericType >::register_quantity().

55  {
56  std::string section = "Physics/"+PhysicsNaming::heat_transfer()+"/output_vars";
57 
58  if( input.have_variable(section) )
59  {
60  unsigned int n_vars = input.vector_variable_size(section);
61 
62  for( unsigned int v = 0; v < n_vars; v++ )
63  {
64  std::string name = input(section,"DIE!",v);
65 
66  if( name == std::string("k") )
67  {
68  this->_k_index = postprocessing.register_quantity( name );
69  }
70  else
71  {
72  std::cerr << "Error: Invalid output_vars value for "+PhysicsNaming::heat_transfer() << std::endl
73  << " Found " << name << std::endl
74  << " Acceptable values are: k" << std::endl;
75  libmesh_error();
76  }
77  }
78  }
79 
80  return;
81  }
static PhysicsName heat_transfer()
unsigned int _k_index
Index from registering this postprocessed quantity.
Definition: heat_transfer.h:74

Member Data Documentation

template<class Conductivity >
unsigned int GRINS::HeatTransfer< Conductivity >::_k_index
private

Index from registering this postprocessed quantity.

Definition at line 74 of file 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