GRINS-0.6.0
Public Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | List of all members
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 read_input_options (const GetPot &input)
 Read options from GetPot input file. More...
 
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 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 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::ParameterMultiPointer< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. 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 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 nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
void init_bcs (libMesh::FEMSystem *system)
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
void attach_neumann_bound_func (GRINS::NBCContainer &neumann_bcs)
 
void attach_dirichlet_bound_func (const GRINS::DBCContainer &dirichlet_bc)
 
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_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)
 
BCHandlingBaseget_bc_handler ()
 
ICHandlingBaseget_ic_handler ()
 
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...
 

Protected Attributes

unsigned int _dim
 Physical dimension of problem. More...
 
VariableIndex _T_var
 Index for temperature field. More...
 
VariableIndex _u_r_var
 Index for r-velocity field. More...
 
VariableIndex _u_z_var
 Index for z-velocity field. More...
 
std::string _T_var_name
 Name for temperature variable. More...
 
std::string _u_r_var_name
 Name of r-velocity. More...
 
std::string _u_z_var_name
 Name of z-velocity. More...
 
GRINSEnums::FEFamily _T_FE_family
 Element type, read from input. More...
 
GRINSEnums::FEFamily _V_FE_family
 
GRINSEnums::Order _T_order
 Temperature element order, read from input. More...
 
GRINSEnums::Order _V_order
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
libMesh::Number _Cp
 
Conductivity _k
 
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::BCHandlingBase_bc_handler
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 
bool _is_axisymmetric
 

Static Protected Attributes

static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 

Private Member Functions

 AxisymmetricHeatTransfer ()
 

Detailed Description

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

Physics class for Axisymmetric Heat Transfer.

Definition at line 46 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 47 of file axisym_heat_transfer.C.

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

49  : Physics(physics_name, input),
50  _k(input)
51  {
52  this->read_input_options(input);
53 
54  // This is deleted in the base class
55  this->_bc_handler = new AxisymmetricHeatTransferBCHandling( physics_name, input );
56  this->_ic_handler = new GenericICHandler( physics_name, input );
57 
58  return;
59  }
virtual void read_input_options(const GetPot &input)
Read options from GetPot input file.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
template<class Conductivity >
GRINS::AxisymmetricHeatTransfer< Conductivity >::~AxisymmetricHeatTransfer ( )

Definition at line 62 of file axisym_heat_transfer.C.

63  {
64  return;
65  }
template<class Conductivity >
GRINS::AxisymmetricHeatTransfer< Conductivity >::AxisymmetricHeatTransfer ( )
private

Member Function Documentation

void GRINS::Physics::attach_dirichlet_bound_func ( const GRINS::DBCContainer dirichlet_bc)
inherited

Definition at line 150 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_dirichlet_bound_func().

151  {
152  _bc_handler->attach_dirichlet_bound_func( dirichlet_bc );
153  return;
154  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_dirichlet_bound_func(const GRINS::DBCContainer &dirichlet_bc)
void GRINS::Physics::attach_neumann_bound_func ( GRINS::NBCContainer neumann_bcs)
inherited

Definition at line 144 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_neumann_bound_func().

145  {
146  _bc_handler->attach_neumann_bound_func( neumann_bcs );
147  return;
148  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_neumann_bound_func(GRINS::NBCContainer &neumann_bcs)
void GRINS::Physics::auxiliary_init ( MultiphysicsSystem system)
virtualinherited

Any auxillary initialization a Physics class may need.

This is called after all variables are added, so this method can safely query the MultiphysicsSystem about variable information.

Definition at line 113 of file physics.C.

114  {
115  return;
116  }
void GRINS::Physics::compute_element_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 185 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::element_constraint().

187  {
188  return;
189  }
void GRINS::Physics::compute_element_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited
void GRINS::Physics::compute_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 203 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::mass_residual().

205  {
206  return;
207  }
void GRINS::Physics::compute_nonlocal_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 197 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

199  {
200  return;
201  }
void GRINS::Physics::compute_nonlocal_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 209 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

211  {
212  return;
213  }
void GRINS::Physics::compute_nonlocal_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 179 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

181  {
182  return;
183  }
void GRINS::Physics::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtualinherited
void GRINS::Physics::compute_side_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 191 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

193  {
194  return;
195  }
void GRINS::Physics::compute_side_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Reimplemented in GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 173 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_time_derivative().

175  {
176  return;
177  }
void GRINS::Physics::element_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited
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 146 of file axisym_heat_transfer.C.

149  {
150 #ifdef GRINS_USE_GRVY_TIMERS
151  this->_timer->BeginTimer("AxisymmetricHeatTransfer::element_time_derivative");
152 #endif
153 
154  // The number of local degrees of freedom in each variable.
155  const unsigned int n_T_dofs = context.get_dof_indices(_T_var).size();
156  const unsigned int n_u_dofs = context.get_dof_indices(_u_r_var).size();
157 
158  //TODO: check n_T_dofs is same as n_u_dofs, n_v_dofs, n_w_dofs
159 
160  // We get some references to cell-specific data that
161  // will be used to assemble the linear system.
162 
163  // Element Jacobian * quadrature weights for interior integration.
164  const std::vector<libMesh::Real> &JxW =
165  context.get_element_fe(_T_var)->get_JxW();
166 
167  // The temperature shape functions at interior quadrature points.
168  const std::vector<std::vector<libMesh::Real> >& T_phi =
169  context.get_element_fe(_T_var)->get_phi();
170 
171  // The velocity shape functions at interior quadrature points.
172  const std::vector<std::vector<libMesh::Real> >& vel_phi =
173  context.get_element_fe(_u_r_var)->get_phi();
174 
175  // The temperature shape function gradients (in global coords.)
176  // at interior quadrature points.
177  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
178  context.get_element_fe(_T_var)->get_dphi();
179 
180  // Physical location of the quadrature points
181  const std::vector<libMesh::Point>& u_qpoint =
182  context.get_element_fe(_u_r_var)->get_xyz();
183 
184  // The subvectors and submatrices we need to fill:
185  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(_T_var); // R_{T}
186 
187  libMesh::DenseSubMatrix<libMesh::Number> &KTT = context.get_elem_jacobian(_T_var, _T_var); // R_{T},{T}
188 
189  libMesh::DenseSubMatrix<libMesh::Number> &KTr = context.get_elem_jacobian(_T_var, _u_r_var); // R_{T},{r}
190  libMesh::DenseSubMatrix<libMesh::Number> &KTz = context.get_elem_jacobian(_T_var, _u_z_var); // R_{T},{z}
191 
192 
193  // Now we will build the element Jacobian and residual.
194  // Constructing the residual requires the solution and its
195  // gradient from the previous timestep. This must be
196  // calculated at each quadrature point by summing the
197  // solution degree-of-freedom values by the appropriate
198  // weight functions.
199  unsigned int n_qpoints = context.get_element_qrule().n_points();
200 
201  for (unsigned int qp=0; qp != n_qpoints; qp++)
202  {
203  const libMesh::Number r = u_qpoint[qp](0);
204 
205  // Compute the solution & its gradient at the old Newton iterate.
206  libMesh::Number u_r, u_z;
207  u_r = context.interior_value(_u_r_var, qp);
208  u_z = context.interior_value(_u_z_var, qp);
209 
210  libMesh::Gradient grad_T;
211  grad_T = context.interior_gradient(_T_var, qp);
212 
213  libMesh::NumberVectorValue U (u_r,u_z);
214 
215  libMesh::Number k = this->_k( context, qp );
216 
217  // FIXME - once we have T-dependent k, we'll need its
218  // derivatives in Jacobians
219  // libMesh::Number dk_dT = this->_k.deriv( T );
220 
221  // First, an i-loop over the degrees of freedom.
222  for (unsigned int i=0; i != n_T_dofs; i++)
223  {
224  FT(i) += JxW[qp]*r*
225  (-_rho*_Cp*T_phi[i][qp]*(U*grad_T) // convection term
226  -k*(T_gradphi[i][qp]*grad_T) ); // diffusion term
227 
228  if (compute_jacobian)
229  {
230  libmesh_assert (context.get_elem_solution_derivative() == 1.0);
231 
232  for (unsigned int j=0; j != n_T_dofs; j++)
233  {
234  // TODO: precompute some terms like:
235  // _rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])
236 
237  KTT(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*
238  (-_rho*_Cp*T_phi[i][qp]*(U*T_gradphi[j][qp]) // convection term
239  -k*(T_gradphi[i][qp]*T_gradphi[j][qp])); // diffusion term
240  } // end of the inner dof (j) loop
241 
242 #if 0
243  if( dk_dT != 0.0 )
244  {
245  for (unsigned int j=0; j != n_T_dofs; j++)
246  {
247  // TODO: precompute some terms like:
248  KTT(i,j) -= JxW[qp] * context.get_elem_solution_derivative() *r*( dk_dT*T_phi[j][qp]*T_gradphi[i][qp]*grad_T );
249  }
250  }
251 #endif
252 
253  // Matrix contributions for the Tu, Tv and Tw couplings (n_T_dofs same as n_u_dofs, n_v_dofs and n_w_dofs)
254  for (unsigned int j=0; j != n_u_dofs; j++)
255  {
256  KTr(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(0)));
257  KTz(i,j) += JxW[qp] * context.get_elem_solution_derivative() *r*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(1)));
258  } // end of the inner dof (j) loop
259 
260  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
261 
262  } // end of the outer dof (i) loop
263  } // end of the quadrature point (qp) loop
264 
265 #ifdef GRINS_USE_GRVY_TIMERS
266  this->_timer->EndTimer("AxisymmetricHeatTransfer::element_time_derivative");
267 #endif
268 
269  return;
270  }
VariableIndex _u_r_var
Index for r-velocity field.
libMesh::Number _rho
Material parameters, read from input.
VariableIndex _u_z_var
Index for z-velocity field.
VariableIndex _T_var
Index for temperature field.
bool GRINS::Physics::enabled_on_elem ( const libMesh::Elem *  elem)
virtualinherited

Find if current physics is active on supplied element.

Definition at line 83 of file physics.C.

References GRINS::Physics::_enabled_subdomains.

84  {
85  // Check if enabled_subdomains flag has been set and if we're
86  // looking at a real element (rather than a nonlocal evaluation)
87  if( !elem || _enabled_subdomains.empty() )
88  return true;
89 
90  // Check if current physics is enabled on elem
91  if( _enabled_subdomains.find( elem->subdomain_id() ) == _enabled_subdomains.end() )
92  return false;
93 
94  return true;
95  }
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
BCHandlingBase * GRINS::Physics::get_bc_handler ( )
inlineinherited

Definition at line 282 of file physics.h.

References GRINS::Physics::_bc_handler.

283  {
284  return _bc_handler;
285  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
ICHandlingBase * GRINS::Physics::get_ic_handler ( )
inlineinherited

Definition at line 288 of file physics.h.

References GRINS::Physics::_ic_handler.

289  {
290  return _ic_handler;
291  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
void GRINS::Physics::init_bcs ( libMesh::FEMSystem *  system)
inherited

Definition at line 118 of file physics.C.

References GRINS::Physics::_bc_handler, GRINS::BCHandlingBase::init_bc_data(), GRINS::BCHandlingBase::init_dirichlet_bc_func_objs(), GRINS::BCHandlingBase::init_dirichlet_bcs(), and GRINS::BCHandlingBase::init_periodic_bcs().

119  {
120  // Only need to init BC's if the physics actually created a handler
121  if( _bc_handler )
122  {
123  _bc_handler->init_bc_data( *system );
124  _bc_handler->init_dirichlet_bcs( system );
126  _bc_handler->init_periodic_bcs( system );
127  }
128 
129  return;
130  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
virtual void init_dirichlet_bcs(libMesh::FEMSystem *system) const
virtual void init_periodic_bcs(libMesh::FEMSystem *system) const
virtual void init_dirichlet_bc_func_objs(libMesh::FEMSystem *system) const
virtual void init_bc_data(const libMesh::FEMSystem &system)
Override this method to initialize any system-dependent data.
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 124 of file axisym_heat_transfer.C.

125  {
126  // We should prerequest all the data
127  // we will need to build the linear system
128  // or evaluate a quantity of interest.
129  context.get_element_fe(_T_var)->get_JxW();
130  context.get_element_fe(_T_var)->get_phi();
131  context.get_element_fe(_T_var)->get_dphi();
132  context.get_element_fe(_T_var)->get_xyz();
133 
134  context.get_side_fe(_T_var)->get_JxW();
135  context.get_side_fe(_T_var)->get_phi();
136  context.get_side_fe(_T_var)->get_dphi();
137  context.get_side_fe(_T_var)->get_xyz();
138 
139  // _u_var is registered so can we assume things related to _u_var
140  // are available in FEMContext
141 
142  return;
143  }
VariableIndex _T_var
Index for temperature field.
void GRINS::Physics::init_ics ( libMesh::FEMSystem *  system,
libMesh::CompositeFunction< libMesh::Number > &  all_ics 
)
inherited

Definition at line 133 of file physics.C.

References GRINS::Physics::_ic_handler, and GRINS::ICHandlingBase::init_ic_data().

135  {
136  if( _ic_handler )
137  {
138  _ic_handler->init_ic_data( *system, all_ics );
139  }
140 
141  return;
142  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
virtual void init_ic_data(const libMesh::FEMSystem &system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
Override this method to initialize any system-dependent data.
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 101 of file axisym_heat_transfer.C.

102  {
103  // Get libMesh to assign an index for each variable
104  this->_dim = system->get_mesh().mesh_dimension();
105 
106  _T_var = system->add_variable( _T_var_name, _T_order, _T_FE_family);
107 
108  // If these are already added, then we just get the index.
109  _u_r_var = system->add_variable(_u_r_var_name, _V_order, _V_FE_family);
110  _u_z_var = system->add_variable(_u_z_var_name, _V_order, _V_FE_family);
111 
112  return;
113  }
GRINSEnums::Order _T_order
Temperature element order, read from input.
std::string _u_z_var_name
Name of z-velocity.
unsigned int _dim
Physical dimension of problem.
std::string _u_r_var_name
Name of r-velocity.
VariableIndex _u_r_var
Index for r-velocity field.
VariableIndex _u_z_var
Index for z-velocity field.
GRINSEnums::FEFamily _T_FE_family
Element type, read from input.
VariableIndex _T_var
Index for temperature field.
std::string _T_var_name
Name for temperature variable.
bool GRINS::Physics::is_steady ( ) const
inherited

Returns whether or not this physics is being solved with a steady solver.

Definition at line 103 of file physics.C.

References GRINS::Physics::_is_steady.

Referenced by GRINS::Physics::set_is_steady().

104  {
105  return _is_steady;
106  }
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
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 299 of file axisym_heat_transfer.C.

302  {
303 #ifdef GRINS_USE_GRVY_TIMERS
304  this->_timer->BeginTimer("AxisymmetricHeatTransfer::mass_residual");
305 #endif
306 
307  // First we get some references to cell-specific data that
308  // will be used to assemble the linear system.
309 
310  // Element Jacobian * quadrature weights for interior integration
311  const std::vector<libMesh::Real> &JxW =
312  context.get_element_fe(_T_var)->get_JxW();
313 
314  // The shape functions at interior quadrature points.
315  const std::vector<std::vector<libMesh::Real> >& phi =
316  context.get_element_fe(_T_var)->get_phi();
317 
318  // The number of local degrees of freedom in each variable
319  const unsigned int n_T_dofs = context.get_dof_indices(_T_var).size();
320 
321  // Physical location of the quadrature points
322  const std::vector<libMesh::Point>& u_qpoint =
323  context.get_element_fe(_u_r_var)->get_xyz();
324 
325  // The subvectors and submatrices we need to fill:
326  libMesh::DenseSubVector<libMesh::Real> &F =
327  context.get_elem_residual(_T_var);
328 
329  libMesh::DenseSubMatrix<libMesh::Real> &M =
330  context.get_elem_jacobian(_T_var, _T_var);
331 
332  unsigned int n_qpoints = context.get_element_qrule().n_points();
333 
334  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
335  {
336  const libMesh::Number r = u_qpoint[qp](0);
337 
338  // For the mass residual, we need to be a little careful.
339  // The time integrator is handling the time-discretization
340  // for us so we need to supply M(u_fixed)*u' for the residual.
341  // u_fixed will be given by the fixed_interior_value function
342  // while u' will be given by the interior_rate function.
343  libMesh::Real T_dot;
344  context.interior_rate(_T_var, qp, T_dot);
345 
346  for (unsigned int i = 0; i != n_T_dofs; ++i)
347  {
348  F(i) -= JxW[qp]*r*(_rho*_Cp*T_dot*phi[i][qp] );
349 
350  if( compute_jacobian )
351  {
352  for (unsigned int j=0; j != n_T_dofs; j++)
353  {
354  // We're assuming rho, cp are constant w.r.t. T here.
355  M(i,j) -= JxW[qp] * context.get_elem_solution_rate_derivative() *r*_rho*_Cp*phi[j][qp]*phi[i][qp] ;
356  }
357  }// End of check on Jacobian
358 
359  } // End of element dof loop
360 
361  } // End of the quadrature point loop
362 
363 #ifdef GRINS_USE_GRVY_TIMERS
364  this->_timer->EndTimer("AxisymmetricHeatTransfer::mass_residual");
365 #endif
366 
367  return;
368  }
VariableIndex _u_r_var
Index for r-velocity field.
libMesh::Number _rho
Material parameters, read from input.
VariableIndex _T_var
Index for temperature field.
void GRINS::Physics::nonlocal_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for scalar variables.

Reimplemented in GRINS::ScalarODE.

Definition at line 250 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

253  {
254  return;
255  }
void GRINS::Physics::nonlocal_mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Mass matrix part(s) for scalar variables.

Reimplemented in GRINS::ScalarODE, and GRINS::AveragedTurbine< Viscosity >.

Definition at line 264 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

267  {
268  return;
269  }
void GRINS::Physics::nonlocal_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Time dependent part(s) of physics for scalar variables.

Reimplemented in GRINS::AveragedTurbine< Viscosity >, and GRINS::ScalarODE.

Definition at line 229 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

232  {
233  return;
234  }
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::read_input_options ( const GetPot &  input)
virtual

Read options from GetPot input file.

Reimplemented from GRINS::Physics.

Definition at line 68 of file axisym_heat_transfer.C.

References GRINS::axisymmetric_heat_transfer, GRINS::incompressible_navier_stokes, GRINS::T_var_name_default, GRINS::u_r_var_name_default, and GRINS::u_z_var_name_default.

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

69  {
70  this->_T_FE_family =
71  libMesh::Utility::string_to_enum<GRINSEnums::FEFamily>( input("Physics/"+axisymmetric_heat_transfer+"/FE_family", "LAGRANGE") );
72 
73  this->_T_order =
74  libMesh::Utility::string_to_enum<GRINSEnums::Order>( input("Physics/"+axisymmetric_heat_transfer+"/T_order", "SECOND") );
75 
76  this->_V_FE_family =
77  libMesh::Utility::string_to_enum<GRINSEnums::FEFamily>( input("Physics/"+incompressible_navier_stokes+"/FE_family", "LAGRANGE") );
78 
79  this->_V_order =
80  libMesh::Utility::string_to_enum<GRINSEnums::Order>( input("Physics/"+incompressible_navier_stokes+"/V_order", "SECOND") );
81 
82  //TODO: same as Incompressible NS
83  this->set_parameter
84  (this->_rho, input,
85  "Physics/"+axisymmetric_heat_transfer+"/rho", 1.0);
86 
87  this->set_parameter
88  (this->_Cp, input,
89  "Physics/"+axisymmetric_heat_transfer+"/Cp", 1.0);
90 
91  this->_T_var_name = input("Physics/VariableNames/Temperature", T_var_name_default );
92 
93  // registered/non-owned variable names
94  this->_u_r_var_name = input("Physics/VariableNames/r_velocity", u_r_var_name_default );
95  this->_u_z_var_name = input("Physics/VariableNames/z_velocity", u_z_var_name_default );
96 
97  return;
98  }
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.
const PhysicsName incompressible_navier_stokes
GRINSEnums::Order _T_order
Temperature element order, read from input.
const PhysicsName axisymmetric_heat_transfer
std::string _u_z_var_name
Name of z-velocity.
std::string _u_r_var_name
Name of r-velocity.
const std::string T_var_name_default
temperature
libMesh::Number _rho
Material parameters, read from input.
GRINSEnums::FEFamily _T_FE_family
Element type, read from input.
const std::string u_z_var_name_default
z-velocity (axisymmetric case)
std::string _T_var_name
Name for temperature variable.
const std::string u_r_var_name_default
r-velocity (axisymmetric case)
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiPointer< libMesh::Number > &  param_pointer 
) const
virtual

Each subclass will register its copy of an independent.

Reimplemented from GRINS::ParameterUser.

Definition at line 372 of file axisym_heat_transfer.C.

References GRINS::ParameterUser::register_parameter().

375  {
376  ParameterUser::register_parameter(param_name, param_pointer);
377  _k.register_parameter(param_name, param_pointer);
378  }
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
void GRINS::Physics::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtualinherited

Register name of postprocessed quantity with PostProcessedQuantities.

Each Physics class will need to cache an unsigned int corresponding to each postprocessed quantity. This will be used in computing the values and putting them in the CachedVariables object.

Reimplemented in GRINS::ParsedVelocitySource< Viscosity >, GRINS::VelocityPenalty< Viscosity >, GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::HeatTransfer< Conductivity >, GRINS::ElasticCable< StressStrainLaw >, GRINS::ElasticMembrane< StressStrainLaw >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 161 of file physics.C.

163  {
164  return;
165  }
void GRINS::Physics::set_is_steady ( bool  is_steady)
inherited

Sets whether this physics is to be solved with a steady solver or not.

Since the member variable is static, only needs to be called on a single physics.

Definition at line 97 of file physics.C.

References GRINS::Physics::_is_steady, and GRINS::Physics::is_steady().

98  {
100  return;
101  }
bool is_steady() const
Returns whether or not this physics is being solved with a steady solver.
Definition: physics.C:103
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
void GRINS::ParameterUser::set_parameter ( libMesh::Number &  param_variable,
const GetPot &  input,
const std::string &  param_name,
libMesh::Number  param_default 
)
virtualinherited

Each subclass can simultaneously read a parameter value from.

Definition at line 35 of file parameter_user.C.

References GRINS::ParameterUser::_my_name, and GRINS::ParameterUser::_my_parameters.

Referenced by GRINS::AveragedFanAdjointStabilization< Viscosity >::AveragedFanAdjointStabilization(), GRINS::AveragedTurbineAdjointStabilization< Viscosity >::AveragedTurbineAdjointStabilization(), GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::BoussinesqBuoyancyAdjointStabilization(), GRINS::BoussinesqBuoyancyBase::BoussinesqBuoyancyBase(), GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::BoussinesqBuoyancySPGSMStabilization(), GRINS::ConstantConductivity::ConstantConductivity(), GRINS::ConstantPrandtlConductivity::ConstantPrandtlConductivity(), GRINS::ConstantSourceFunction::ConstantSourceFunction(), GRINS::ConstantSourceTerm::ConstantSourceTerm(), GRINS::ConstantSpecificHeat::ConstantSpecificHeat(), GRINS::ConstantViscosity::ConstantViscosity(), GRINS::ElasticCable< StressStrainLaw >::ElasticCable(), GRINS::ElasticCableConstantGravity::ElasticCableConstantGravity(), GRINS::ElasticMembrane< StressStrainLaw >::ElasticMembrane(), GRINS::ElasticMembraneConstantPressure::ElasticMembraneConstantPressure(), GRINS::HeatConduction< Conductivity >::HeatConduction(), GRINS::HeatTransferBase< Conductivity >::HeatTransferBase(), GRINS::IncompressibleNavierStokesBase< Viscosity >::IncompressibleNavierStokesBase(), GRINS::AverageNusseltNumber::init(), GRINS::MooneyRivlin::MooneyRivlin(), GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >::ReactingLowMachNavierStokesBase(), GRINS::HookesLaw1D::read_input_options(), GRINS::HookesLaw::read_input_options(), GRINS::AxisymmetricBoussinesqBuoyancy::read_input_options(), and GRINS::VelocityDragAdjointStabilization< Viscosity >::VelocityDragAdjointStabilization().

39  {
40  param_variable = input(param_name, param_default);
41 
42  libmesh_assert_msg(!_my_parameters.count(param_name),
43  "ERROR: " << _my_name << " double-registered parameter " <<
44  param_name);
45 
46  _my_parameters[param_name] = &param_variable;
47  }
std::map< std::string, libMesh::Number * > _my_parameters
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 116 of file axisym_heat_transfer.C.

117  {
118  // Tell the system to march temperature forward in time
119  system->time_evolving(_T_var);
120  return;
121  }
VariableIndex _T_var
Index for temperature field.
void GRINS::Physics::side_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for boundaries of elements on the domain boundary.

Reimplemented in GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 243 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

246  {
247  return;
248  }
template<class Conductivity >
void GRINS::AxisymmetricHeatTransfer< Conductivity >::side_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Time dependent part(s) of physics for boundaries of elements on the domain boundary.

Reimplemented from GRINS::Physics.

Definition at line 273 of file axisym_heat_transfer.C.

276  {
277 #ifdef GRINS_USE_GRVY_TIMERS
278  this->_timer->BeginTimer("AxisymmetricHeatTransfer::side_time_derivative");
279 #endif
280 
281  std::vector<BoundaryID> ids = context.side_boundary_ids();
282 
283  for( std::vector<BoundaryID>::const_iterator it = ids.begin();
284  it != ids.end(); it++ )
285  {
286  libmesh_assert (*it != libMesh::BoundaryInfo::invalid_id);
287 
288  _bc_handler->apply_neumann_bcs( context, cache, compute_jacobian, *it );
289  }
290 
291 #ifdef GRINS_USE_GRVY_TIMERS
292  this->_timer->EndTimer("AxisymmetricHeatTransfer::side_time_derivative");
293 #endif
294 
295  return;
296  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
virtual void apply_neumann_bcs(AssemblyContext &context, const GRINS::CachedValues &cache, const bool request_jacobian, const GRINS::BoundaryID bc_id) const

Member Data Documentation

GRINS::BCHandlingBase* GRINS::Physics::_bc_handler
protectedinherited
template<class Conductivity >
libMesh::Number GRINS::AxisymmetricHeatTransfer< Conductivity >::_Cp
protected

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

std::set<libMesh::subdomain_id_type> GRINS::Physics::_enabled_subdomains
protectedinherited

Subdomains on which the current Physics class is enabled.

Definition at line 261 of file physics.h.

Referenced by GRINS::Physics::enabled_on_elem(), and GRINS::Physics::read_input_options().

GRINS::ICHandlingBase* GRINS::Physics::_ic_handler
protectedinherited
bool GRINS::Physics::_is_axisymmetric
protectedinherited
bool GRINS::Physics::_is_steady = false
staticprotectedinherited

Caches whether or not the solver that's being used is steady or not.

This is need, for example, in flow stabilization as the tau terms change depending on whether the solver is steady or unsteady.

Definition at line 266 of file physics.h.

Referenced by GRINS::Physics::is_steady(), and GRINS::Physics::set_is_steady().

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

Definition at line 132 of file axisym_heat_transfer.h.

const PhysicsName GRINS::Physics::_physics_name
protectedinherited

Name of the physics object. Used for reading physics specific inputs.

We use a reference because the physics names are const global objects in GRINS namespace

No, we use a copy, because otherwise as soon as the memory in std::set<std::string> requested_physics gets overwritten we get in trouble.

Definition at line 254 of file physics.h.

Referenced by GRINS::SourceTermBase::parse_var_info(), and GRINS::Physics::read_input_options().

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 130 of file axisym_heat_transfer.h.

template<class Conductivity >
GRINSEnums::FEFamily GRINS::AxisymmetricHeatTransfer< Conductivity >::_T_FE_family
protected

Element type, read from input.

Definition at line 120 of file axisym_heat_transfer.h.

template<class Conductivity >
GRINSEnums::Order GRINS::AxisymmetricHeatTransfer< Conductivity >::_T_order
protected

Temperature element order, read from input.

Definition at line 123 of file axisym_heat_transfer.h.

template<class Conductivity >
VariableIndex GRINS::AxisymmetricHeatTransfer< Conductivity >::_T_var
protected

Index for temperature field.

Definition at line 101 of file axisym_heat_transfer.h.

template<class Conductivity >
std::string GRINS::AxisymmetricHeatTransfer< Conductivity >::_T_var_name
protected

Name for temperature variable.

Definition at line 111 of file axisym_heat_transfer.h.

template<class Conductivity >
VariableIndex GRINS::AxisymmetricHeatTransfer< Conductivity >::_u_r_var
protected

Index for r-velocity field.

Definition at line 104 of file axisym_heat_transfer.h.

template<class Conductivity >
std::string GRINS::AxisymmetricHeatTransfer< Conductivity >::_u_r_var_name
protected

Name of r-velocity.

Definition at line 114 of file axisym_heat_transfer.h.

template<class Conductivity >
VariableIndex GRINS::AxisymmetricHeatTransfer< Conductivity >::_u_z_var
protected

Index for z-velocity field.

Definition at line 107 of file axisym_heat_transfer.h.

template<class Conductivity >
std::string GRINS::AxisymmetricHeatTransfer< Conductivity >::_u_z_var_name
protected

Name of z-velocity.

Definition at line 117 of file axisym_heat_transfer.h.

template<class Conductivity >
GRINSEnums::FEFamily GRINS::AxisymmetricHeatTransfer< Conductivity >::_V_FE_family
protected

Definition at line 120 of file axisym_heat_transfer.h.

template<class Conductivity >
GRINSEnums::Order GRINS::AxisymmetricHeatTransfer< Conductivity >::_V_order
protected

Definition at line 123 of file axisym_heat_transfer.h.


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

Generated on Mon Jun 22 2015 21:32:22 for GRINS-0.6.0 by  doxygen 1.8.9.1