GRINS-0.8.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)
 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...
 
- 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 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 void init_variables (libMesh::FEMSystem *)
 Initialize variables for this physics. More...
 
virtual bool enabled_on_elem (const libMesh::Elem *elem)
 Find if current physics is active on supplied element. More...
 
void set_is_steady (bool is_steady)
 Sets whether this physics is to be solved with a steady solver or not. More...
 
bool is_steady () const
 Returns whether or not this physics is being solved with a steady solver. More...
 
virtual void auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. More...
 
virtual void preassembly (MultiphysicsSystem &)
 Perform any necessary setup before element assembly begins. More...
 
virtual void reinit (MultiphysicsSystem &)
 Any reinitialization that needs to be done. More...
 
virtual void side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

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)
 
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::HeatTransferStabilizationBase< Conductivity >
HeatTransferStabilizationHelper _stab_helper
 
- 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::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  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 53 of file heat_transfer_spgsm_stab.C.

54  {
55  // The number of local degrees of freedom in each variable.
56  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
57 
58  // Element Jacobian * quadrature weights for interior integration.
59  const std::vector<libMesh::Real> &JxW =
60  context.get_element_fe(this->_temp_vars.T())->get_JxW();
61 
62  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
63  context.get_element_fe(this->_temp_vars.T())->get_dphi();
64 
65  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
66 
67  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
68 
69  unsigned int n_qpoints = context.get_element_qrule().n_points();
70 
71  for (unsigned int qp=0; qp != n_qpoints; qp++)
72  {
73  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
74  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
75 
76  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
77  context.interior_value( this->_flow_vars.v(), qp ) );
78  if( this->_flow_vars.dim() == 3 )
79  {
80  U(2) = context.interior_value( this->_flow_vars.w(), qp );
81  }
82 
83  // Compute Conductivity at this qp
84  libMesh::Real _k_qp = this->_k(context, qp);
85 
86  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, G, this->_rho, this->_Cp, _k_qp, U, this->_is_steady );
87 
88  libMesh::Real RE_s = this->_stab_helper.compute_res_energy_steady( context, qp, this->_rho, this->_Cp, _k_qp );
89 
90  for (unsigned int i=0; i != n_T_dofs; i++)
91  {
92  FT(i) += -tau_E*RE_s*this->_rho*this->_Cp*U*T_gradphi[i][qp]*JxW[qp];
93  }
94 
95  if( compute_jacobian )
96  {
97  libmesh_not_implemented();
98  }
99 
100  }
101  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
VelocityVariable & _flow_vars
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
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:266
libMesh::Real compute_res_energy_steady(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real Cp, const libMesh::Real k) const
unsigned int dim() const
Number of components.
PrimitiveTempFEVariables & _temp_vars
template<class K >
void GRINS::HeatTransferSPGSMStabilization< 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 104 of file heat_transfer_spgsm_stab.C.

106  {
107  if( compute_jacobian )
108  libmesh_not_implemented();
109 
110  // The number of local degrees of freedom in each variable.
111  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
112 
113  // Element Jacobian * quadrature weights for interior integration.
114  const std::vector<libMesh::Real> &JxW =
115  context.get_element_fe(this->_temp_vars.T())->get_JxW();
116 
117  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
118  context.get_element_fe(this->_temp_vars.T())->get_dphi();
119 
120  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
121 
122  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
123 
124  unsigned int n_qpoints = context.get_element_qrule().n_points();
125 
126  for (unsigned int qp=0; qp != n_qpoints; qp++)
127  {
128  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
129  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
130 
131  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
132  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
133  if( this->_flow_vars.dim() == 3 )
134  {
135  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
136  }
137 
138  // Compute Conductivity at this qp
139  libMesh::Real _k_qp = this->_k(context, qp);
140 
141  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, G, this->_rho, this->_Cp, _k_qp, U, false );
142 
143  libMesh::Real RE_t = this->_stab_helper.compute_res_energy_transient( context, qp, this->_rho, this->_Cp );
144 
145  for (unsigned int i=0; i != n_T_dofs; i++)
146  {
147  FT(i) -= tau_E*RE_t*this->_rho*this->_Cp*U*T_gradphi[i][qp]*JxW[qp];
148  }
149 
150  }
151  }
libMesh::Real compute_res_energy_transient(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real Cp) const
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
VelocityVariable & _flow_vars
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
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
unsigned int dim() const
Number of components.
PrimitiveTempFEVariables & _temp_vars

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