GRINS-0.7.0
List of all members | Public Member Functions | Private Member Functions
GRINS::HeatTransferSPGSMStabilization< Conductivity > Class Template Reference

Adds VMS-based stabilization to LowMachNavierStokes physics class. More...

#include <heat_transfer_spgsm_stab.h>

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

Public Member Functions

 HeatTransferSPGSMStabilization (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~HeatTransferSPGSMStabilization ()
 
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...
 
- Public Member Functions inherited from GRINS::HeatTransferStabilizationBase< Conductivity >
 HeatTransferStabilizationBase (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~HeatTransferStabilizationBase ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. 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 init_variables (libMesh::FEMSystem *system)
 Initialization Heat Transfer variables. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. 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...
 

Private Member Functions

 HeatTransferSPGSMStabilization ()
 

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)
 
- Protected Attributes inherited from GRINS::HeatTransferStabilizationBase< Conductivity >
HeatTransferStabilizationHelper _stab_helper
 
- Protected Attributes inherited from GRINS::HeatTransferBase< Conductivity >
unsigned int _dim
 Physical dimension of problem. More...
 
VelocityFEVariables _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::HeatTransferSPGSMStabilization< Conductivity >

Adds VMS-based stabilization to LowMachNavierStokes physics class.

Definition at line 36 of file heat_transfer_spgsm_stab.h.

Constructor & Destructor Documentation

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

Definition at line 38 of file heat_transfer_spgsm_stab.C.

40  : HeatTransferStabilizationBase<K>(physics_name,input)
41  {
42  return;
43  }

Definition at line 46 of file heat_transfer_spgsm_stab.C.

47  {
48  return;
49  }
template<class Conductivity >
GRINS::HeatTransferSPGSMStabilization< Conductivity >::HeatTransferSPGSMStabilization ( )
private

Member Function Documentation

template<class K >
void GRINS::HeatTransferSPGSMStabilization< 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 52 of file heat_transfer_spgsm_stab.C.

55  {
56 #ifdef GRINS_USE_GRVY_TIMERS
57  this->_timer->BeginTimer("HeatTransferSPGSMStabilization::element_time_derivative");
58 #endif
59 
60  // The number of local degrees of freedom in each variable.
61  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
62 
63  // Element Jacobian * quadrature weights for interior integration.
64  const std::vector<libMesh::Real> &JxW =
65  context.get_element_fe(this->_temp_vars.T())->get_JxW();
66 
67  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
68  context.get_element_fe(this->_temp_vars.T())->get_dphi();
69 
70  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
71 
72  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
73 
74  unsigned int n_qpoints = context.get_element_qrule().n_points();
75 
76  for (unsigned int qp=0; qp != n_qpoints; qp++)
77  {
78  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
79  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
80 
81  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
82  context.interior_value( this->_flow_vars.v(), qp ) );
83  if( this->_dim == 3 )
84  {
85  U(2) = context.interior_value( this->_flow_vars.w(), qp );
86  }
87 
88  // Compute Conductivity at this qp
89  libMesh::Real _k_qp = this->_k(context, qp);
90 
91  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, G, this->_rho, this->_Cp, _k_qp, U, this->_is_steady );
92 
93  libMesh::Real RE_s = this->_stab_helper.compute_res_energy_steady( context, qp, this->_rho, this->_Cp, _k_qp );
94 
95  for (unsigned int i=0; i != n_T_dofs; i++)
96  {
97  FT(i) += -tau_E*RE_s*this->_rho*this->_Cp*U*T_gradphi[i][qp]*JxW[qp];
98  }
99 
100  if( compute_jacobian )
101  {
102  libmesh_not_implemented();
103  }
104 
105  }
106 
107 #ifdef GRINS_USE_GRVY_TIMERS
108  this->_timer->EndTimer("HeatTransferSPGSMStabilization::element_time_derivative");
109 #endif
110  return;
111  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
libMesh::Number _rho
Material parameters, read from input.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
Conductivity _k
Conductivity.
HeatTransferStabilizationHelper _stab_helper
PrimitiveTempFEVariables _temp_vars
libMesh::Real compute_tau_energy(AssemblyContext &c, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Real cp, libMesh::Real k, libMesh::Gradient U, bool is_steady) const
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:277
VelocityFEVariables _flow_vars
unsigned int _dim
Physical dimension of problem.
libMesh::Real compute_res_energy_steady(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real Cp, const libMesh::Real k) const
template<class K >
void GRINS::HeatTransferSPGSMStabilization< 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 114 of file heat_transfer_spgsm_stab.C.

117  {
118 #ifdef GRINS_USE_GRVY_TIMERS
119  this->_timer->BeginTimer("HeatTransferSPGSMStabilization::mass_residual");
120 #endif
121 
122  // The number of local degrees of freedom in each variable.
123  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
124 
125  // Element Jacobian * quadrature weights for interior integration.
126  const std::vector<libMesh::Real> &JxW =
127  context.get_element_fe(this->_temp_vars.T())->get_JxW();
128 
129  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
130  context.get_element_fe(this->_temp_vars.T())->get_dphi();
131 
132  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
133 
134  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
135 
136  unsigned int n_qpoints = context.get_element_qrule().n_points();
137 
138  for (unsigned int qp=0; qp != n_qpoints; qp++)
139  {
140  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
141  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
142 
143  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
144  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
145  if( this->_dim == 3 )
146  {
147  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
148  }
149 
150  // Compute Conductivity at this qp
151  libMesh::Real _k_qp = this->_k(context, qp);
152 
153  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, G, this->_rho, this->_Cp, _k_qp, U, false );
154 
155  libMesh::Real RE_t = this->_stab_helper.compute_res_energy_transient( context, qp, this->_rho, this->_Cp );
156 
157  for (unsigned int i=0; i != n_T_dofs; i++)
158  {
159  FT(i) -= tau_E*RE_t*this->_rho*this->_Cp*U*T_gradphi[i][qp]*JxW[qp];
160  }
161 
162  }
163 
164 #ifdef GRINS_USE_GRVY_TIMERS
165  this->_timer->EndTimer("HeatTransferSPGSMStabilization::mass_residual");
166 #endif
167  return;
168  }
libMesh::Real compute_res_energy_transient(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real Cp) const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
libMesh::Number _rho
Material parameters, read from input.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
Conductivity _k
Conductivity.
HeatTransferStabilizationHelper _stab_helper
PrimitiveTempFEVariables _temp_vars
libMesh::Real compute_tau_energy(AssemblyContext &c, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Real cp, libMesh::Real k, libMesh::Gradient U, bool is_steady) const
VelocityFEVariables _flow_vars
unsigned int _dim
Physical dimension of problem.

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