GRINS-0.8.0
List of all members | Public Member Functions | Protected Attributes | Private Member Functions
GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity > Class Template Reference

Adds Boussinesq bouyancy SPGSM stabilization source term. More...

#include <boussinesq_buoyancy_spgsm_stab.h>

Inheritance diagram for GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >:
Collaboration graph
[legend]

Public Member Functions

 BoussinesqBuoyancySPGSMStabilization (const std::string &physics_name, const GetPot &input)
 
 ~BoussinesqBuoyancySPGSMStabilization ()
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context)
 Constraint 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 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::BoussinesqBuoyancyBase
 BoussinesqBuoyancyBase (const std::string &physics_name, const GetPot &input)
 
 ~BoussinesqBuoyancyBase ()
 
- 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 set_time_evolving_vars (libMesh::FEMSystem *system)
 Set which variables are time evolving. 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 init_context (AssemblyContext &context)
 Initialize context for added physics variables. 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 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...
 

Protected Attributes

IncompressibleNavierStokesStabilizationHelper _flow_stab_helper
 
HeatTransferStabilizationHelper _temp_stab_helper
 
libMesh::Number _Cp
 
libMesh::Number _k
 
Viscosity _mu
 Viscosity object. More...
 
- Protected Attributes inherited from GRINS::BoussinesqBuoyancyBase
const VelocityVariable_flow_vars
 
const PressureFEVariable_press_var
 
const PrimitiveTempFEVariables_temp_vars
 
libMesh::Number _rho
 $ \rho = $ density More...
 
libMesh::Number _T_ref
 $ T_0 = $ reference temperature More...
 
libMesh::Number _beta_T
 $ \beta_T = $ coefficient of thermal expansion More...
 
libMesh::Point _g
 Gravitational vector. 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

 BoussinesqBuoyancySPGSMStabilization ()
 

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::BoussinesqBuoyancyBase
void read_property (const GetPot &input, const std::string &old_option, const std::string &property, libMesh::Real &value)
 Helper function for parsing/maintaing backward compatibility. More...
 
void no_input_warning (const GetPot &input, const std::string &old_option, const std::string &material, const std::string &property)
 Helper function for parsing/maintaing backward compatibility. 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...
 
- 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 Viscosity>
class GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >

Adds Boussinesq bouyancy SPGSM stabilization source term.

This class implements the SPGSM stabilization term for the BoussinesqBuoyancy Physics. Intended to be used with IncompressibleNavierStokesSPGSMStabilization and HeatTransferSPGSMStabilization. This, essentially, is just adding SUPG terms to the momentum equation and PSPG terms to the pressure equation.

Definition at line 45 of file boussinesq_buoyancy_spgsm_stab.h.

Constructor & Destructor Documentation

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

Definition at line 42 of file boussinesq_buoyancy_spgsm_stab.C.

References GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::_Cp, GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::_k, GRINS::PhysicsNaming::boussinesq_buoyancy(), GRINS::PhysicsNaming::heat_transfer(), GRINS::MaterialsParsing::read_specific_heat(), and GRINS::ParameterUser::set_parameter().

43  : BoussinesqBuoyancyBase(physics_name,input),
44  _flow_stab_helper(physics_name+"FlowStabHelper", input),
45  _temp_stab_helper(physics_name+"TempStabHelper", input),
46  _Cp(0.0),
47  _k(1.0),
49  {
51 
52  this->set_parameter
53  (_k, input, "Physics/"+PhysicsNaming::heat_transfer()+"/k", _k);
54  }
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.
static PhysicsName heat_transfer()
IncompressibleNavierStokesStabilizationHelper _flow_stab_helper
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 PhysicsName boussinesq_buoyancy()
static std::string material_name(const GetPot &input, const std::string &physics)
Get the name of the material in the Physics/physics section.

Definition at line 57 of file boussinesq_buoyancy_spgsm_stab.C.

58  {
59  return;
60  }
template<class Viscosity >
GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::BoussinesqBuoyancySPGSMStabilization ( )
private

Member Function Documentation

template<class Mu >
void GRINS::BoussinesqBuoyancySPGSMStabilization< Mu >::element_constraint ( bool  ,
AssemblyContext  
)
virtual

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 153 of file boussinesq_buoyancy_spgsm_stab.C.

155  {
156  // The number of local degrees of freedom in each variable.
157  const unsigned int n_p_dofs = context.get_dof_indices(_press_var.p()).size();
158 
159  // Element Jacobian * quadrature weights for interior integration.
160  const std::vector<libMesh::Real> &JxW =
161  context.get_element_fe(_flow_vars.u())->get_JxW();
162 
163  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
164  context.get_element_fe(this->_press_var.p())->get_dphi();
165 
166  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
167 
168  // Now we will build the element Jacobian and residual.
169  // Constructing the residual requires the solution and its
170  // gradient from the previous timestep. This must be
171  // calculated at each quadrature point by summing the
172  // solution degree-of-freedom values by the appropriate
173  // weight functions.
174  unsigned int n_qpoints = context.get_element_qrule().n_points();
175 
176  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
177 
178  for (unsigned int qp=0; qp != n_qpoints; qp++)
179  {
180  libMesh::RealGradient g = this->_flow_stab_helper.compute_g( fe, context, qp );
181  libMesh::RealTensor G = this->_flow_stab_helper.compute_G( fe, context, qp );
182 
183  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
184  context.interior_value( this->_flow_vars.v(), qp ) );
185  if( this->_flow_vars.dim() == 3 )
186  U(2) = context.interior_value( this->_flow_vars.w(), qp );
187 
188  // Compute the viscosity at this qp
189  libMesh::Real mu_qp = this->_mu(context, qp);
190 
191  libMesh::Real tau_M = this->_flow_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, mu_qp, this->_is_steady );
192 
193  // Compute the solution & its gradient at the old Newton iterate.
194  libMesh::Number T;
195  T = context.interior_value(_temp_vars.T(), qp);
196 
197  libMesh::RealGradient residual = _rho*_beta_T*(T-_T_ref)*_g;
198 
199  // First, an i-loop over the velocity degrees of freedom.
200  // We know that n_u_dofs == n_v_dofs so we can compute contributions
201  // for both at the same time.
202  for (unsigned int i=0; i != n_p_dofs; i++)
203  {
204  Fp(i) += -tau_M*residual*p_dphi[i][qp]*JxW[qp];
205 
206  if (compute_jacobian)
207  {
208  libmesh_not_implemented();
209  } // End compute_jacobian check
210 
211  } // End i dof loop
212 
213  } // End quadrature loop
214  }
libMesh::Point _g
Gravitational vector.
const PressureFEVariable & _press_var
VariableIndex T() const
const VelocityVariable & _flow_vars
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
VariableIndex p() const
IncompressibleNavierStokesStabilizationHelper _flow_stab_helper
libMesh::Number _T_ref
reference temperature
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
const PrimitiveTempFEVariables & _temp_vars
libMesh::Number _beta_T
coefficient of thermal expansion
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
template<class Mu >
void GRINS::BoussinesqBuoyancySPGSMStabilization< Mu >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 64 of file boussinesq_buoyancy_spgsm_stab.C.

65  {
66  // The number of local degrees of freedom in each variable.
67  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
68 
69  // Element Jacobian * quadrature weights for interior integration.
70  const std::vector<libMesh::Real> &JxW =
71  context.get_element_fe(_flow_vars.u())->get_JxW();
72 
73  /*
74  const std::vector<std::vector<libMesh::Real> >& u_phi =
75  context.get_element_fe(this->_flow_vars.u())->get_phi();
76  */
77 
78  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
79  context.get_element_fe(this->_flow_vars.u())->get_dphi();
80 
81  // Get residuals
82  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u()); // R_{u}
83  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v()); // R_{v}
84  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
85 
86  if(this->_flow_vars.dim() == 3)
87  {
88  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
89  }
90 
91  // Now we will build the element Jacobian and residual.
92  // Constructing the residual requires the solution and its
93  // gradient from the previous timestep. This must be
94  // calculated at each quadrature point by summing the
95  // solution degree-of-freedom values by the appropriate
96  // weight functions.
97  unsigned int n_qpoints = context.get_element_qrule().n_points();
98 
99  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
100 
101  for (unsigned int qp=0; qp != n_qpoints; qp++)
102  {
103  libMesh::RealGradient g = this->_flow_stab_helper.compute_g( fe, context, qp );
104  libMesh::RealTensor G = this->_flow_stab_helper.compute_G( fe, context, qp );
105 
106  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
107  context.interior_value( this->_flow_vars.v(), qp ) );
108  if( this->_flow_vars.dim() == 3 )
109  {
110  U(2) = context.interior_value( this->_flow_vars.w(), qp );
111  }
112 
113  // Compute the viscosity at this qp
114  libMesh::Real mu_qp = this->_mu(context, qp);
115 
116  libMesh::Real tau_M = this->_flow_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, mu_qp, this->_is_steady );
117 
118  //libMesh::Real tau_E = this->_temp_stab_helper.compute_tau_energy( context, G, _rho, _Cp, _k, U, this->_is_steady );
119 
120  //libMesh::Real RE = this->_temp_stab_helper.compute_res_energy_steady( context, qp, _rho, _Cp, _k );
121 
122  // Compute the solution & its gradient at the old Newton iterate.
123  libMesh::Number T;
124  T = context.interior_value(_temp_vars.T(), qp);
125 
126  libMesh::RealGradient residual = _rho*_beta_T*(T-_T_ref)*_g;
127 
128  for (unsigned int i=0; i != n_u_dofs; i++)
129  {
130  Fu(i) += ( -tau_M*residual(0)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
131  // + _rho*_beta_T*tau_E*RE*_g(0)*u_phi[i][qp] )*JxW[qp];
132 
133  Fv(i) += ( -tau_M*residual(1)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
134  // + _rho*_beta_T*tau_E*RE*_g(1)*u_phi[i][qp] )*JxW[qp];
135 
136  if (this->_flow_vars.dim() == 3)
137  {
138  (*Fw)(i) += ( -tau_M*residual(2)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
139  // + _rho*_beta_T*tau_E*RE*_g(2)*u_phi[i][qp] )*JxW[qp];
140  }
141 
142  if (compute_jacobian)
143  {
144  libmesh_not_implemented();
145  } // End compute_jacobian check
146 
147  } // End i dof loop
148  } // End quadrature loop
149  }
libMesh::Point _g
Gravitational vector.
VariableIndex T() const
const VelocityVariable & _flow_vars
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
IncompressibleNavierStokesStabilizationHelper _flow_stab_helper
libMesh::Number _T_ref
reference temperature
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
const PrimitiveTempFEVariables & _temp_vars
libMesh::Number _beta_T
coefficient of thermal expansion
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
template<class Mu >
void GRINS::BoussinesqBuoyancySPGSMStabilization< Mu >::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 217 of file boussinesq_buoyancy_spgsm_stab.C.

219  {
220  /*
221  // The number of local degrees of freedom in each variable.
222  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
223 
224  // Element Jacobian * quadrature weights for interior integration.
225  const std::vector<libMesh::Real> &JxW =
226  context.get_element_fe(_flow_vars.u())->get_JxW();
227 
228  const std::vector<std::vector<libMesh::Real> >& u_phi =
229  context.get_element_fe(this->_flow_vars.u())->get_phi();
230 
231  // Get residuals
232  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u()); // R_{u}
233  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v()); // R_{v}
234  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
235  if(this->_flow_vars.dim() == 3)
236  {
237  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
238  }
239 
240  // Now we will build the element Jacobian and residual.
241  // Constructing the residual requires the solution and its
242  // gradient from the previous timestep. This must be
243  // calculated at each quadrature point by summing the
244  // solution degree-of-freedom values by the appropriate
245  // weight functions.
246  unsigned int n_qpoints = context.get_element_qrule().n_points();
247 
248  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
249 
250  for (unsigned int qp=0; qp != n_qpoints; qp++)
251  {
252  libMesh::RealGradient g = this->_flow_stab_helper.compute_g( fe, context, qp );
253  libMesh::RealTensor G = this->_flow_stab_helper.compute_G( fe, context, qp );
254 
255  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
256  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
257  if( this->_flow_vars.dim() == 3 )
258  {
259  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
260  }
261 
262  libMesh::Real tau_E = this->_temp_stab_helper.compute_tau_energy( context, G, _rho, _Cp, _k, U, false );
263 
264  libMesh::Real RE = this->_temp_stab_helper.compute_res_energy_transient( context, qp, _rho, _Cp );
265 
266 
267  for (unsigned int i=0; i != n_u_dofs; i++)
268  {
269  Fu(i) += -_rho*_beta_T*tau_E*RE*_g(0)*u_phi[i][qp]*JxW[qp];
270 
271  Fv(i) += -_rho*_beta_T*tau_E*RE*_g(1)*u_phi[i][qp]*JxW[qp];
272 
273  if (this->_flow_vars.dim() == 3)
274  {
275  (*Fw)(i) += -_rho*_beta_T*tau_E*RE*_g(2)*u_phi[i][qp]*JxW[qp];
276  }
277 
278  if (compute_jacobian)
279  {
280  libmesh_not_implemented();
281  } // End compute_jacobian check
282 
283  } // End i dof loop
284  } // End quadrature loop
285  */
286  }
template<class Mu >
void GRINS::BoussinesqBuoyancySPGSMStabilization< Mu >::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 290 of file boussinesq_buoyancy_spgsm_stab.C.

References GRINS::ParameterUser::register_parameter().

293  {
294  ParameterUser::register_parameter(param_name, param_pointer);
295  _mu.register_parameter(param_name, param_pointer);
296  }
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.

Member Data Documentation

template<class Viscosity >
libMesh::Number GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::_Cp
protected
template<class Viscosity >
IncompressibleNavierStokesStabilizationHelper GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::_flow_stab_helper
protected

Definition at line 73 of file boussinesq_buoyancy_spgsm_stab.h.

template<class Viscosity >
libMesh::Number GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::_k
protected
template<class Viscosity >
Viscosity GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::_mu
protected

Viscosity object.

Definition at line 79 of file boussinesq_buoyancy_spgsm_stab.h.

template<class Viscosity >
HeatTransferStabilizationHelper GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::_temp_stab_helper
protected

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