GRINS-0.7.0
List of all members | Public Member Functions | Protected Attributes | Private Member Functions
GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity > Class Template Reference

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

#include <low_mach_navier_stokes_stab_base.h>

Inheritance diagram for GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity >:
Collaboration graph
[legend]

Public Member Functions

 LowMachNavierStokesStabilizationBase (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~LowMachNavierStokesStabilizationBase ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
libMesh::Real compute_res_continuity_steady (AssemblyContext &context, unsigned int qp) const
 
libMesh::Real compute_res_continuity_transient (AssemblyContext &context, unsigned int qp) const
 
libMesh::RealGradient compute_res_momentum_steady (AssemblyContext &context, unsigned int qp) const
 
libMesh::RealGradient compute_res_momentum_transient (AssemblyContext &context, unsigned int qp) const
 
libMesh::Real compute_res_energy_steady (AssemblyContext &context, unsigned int qp) const
 
libMesh::Real compute_res_energy_transient (AssemblyContext &context, unsigned int qp) const
 
- Public Member Functions inherited from GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >
 LowMachNavierStokesBase (const PhysicsName &physics_name, const std::string &core_physics_name, const GetPot &input)
 
 ~LowMachNavierStokesBase ()
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialize variables for this physics. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
libMesh::Real T (const libMesh::Point &p, const AssemblyContext &c) const
 
libMesh::Real rho (libMesh::Real T, libMesh::Real p0) const
 
libMesh::Real d_rho_dT (libMesh::Real T, libMesh::Real p0) const
 
libMesh::Real get_p0_steady (const AssemblyContext &c, unsigned int qp) const
 
libMesh::Real get_p0_steady_side (const AssemblyContext &c, unsigned int qp) const
 
libMesh::Real get_p0_steady (const AssemblyContext &c, const libMesh::Point &p) const
 
libMesh::Real get_p0_transient (AssemblyContext &c, unsigned int qp) const
 
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 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 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 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 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...
 

Protected Attributes

LowMachNavierStokesStabilizationHelper _stab_helper
 
- Protected Attributes inherited from GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >
libMesh::Number _p0_over_R
 Thermodynamic pressure divided by gas constant. More...
 
libMesh::Number _p0
 
libMesh::Number _R
 
libMesh::Number _T0
 
unsigned int _dim
 Physical dimension of problem. More...
 
VelocityFEVariables _flow_vars
 
PressureFEVariable _press_var
 
PrimitiveTempFEVariables _temp_vars
 
libMesh::UniquePtr< ThermoPressureFEVariable_p0_var
 
Viscosity _mu
 Viscosity object. More...
 
SpecificHeat _cp
 Specific heat object. More...
 
ThermalConductivity _k
 Thermal conductivity object. More...
 
libMesh::Point _g
 Gravity vector. More...
 
bool _enable_thermo_press_calc
 Flag to enable thermodynamic pressure calculation. 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

 LowMachNavierStokesStabilizationBase ()
 

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)
 
- 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 SpecificHeat, class ThermalConductivity>
class GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity >

Adds VMS-based stabilization to LowMachNavierStokes physics class.

Definition at line 37 of file low_mach_navier_stokes_stab_base.h.

Constructor & Destructor Documentation

template<class Mu , class SH , class TC >
GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::LowMachNavierStokesStabilizationBase ( const GRINS::PhysicsName physics_name,
const GetPot &  input 
)

Definition at line 39 of file low_mach_navier_stokes_stab_base.C.

41  : LowMachNavierStokesBase<Mu,SH,TC>(physics_name,PhysicsNaming::low_mach_navier_stokes(),input),
42  _stab_helper( physics_name+"StabHelper", input )
43  {
44  return;
45  }
static PhysicsName low_mach_navier_stokes()
LowMachNavierStokesStabilizationHelper _stab_helper
template<class Mu , class SH , class TC >
GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::~LowMachNavierStokesStabilizationBase ( )
virtual

Definition at line 48 of file low_mach_navier_stokes_stab_base.C.

49  {
50  return;
51  }
template<class Viscosity , class SpecificHeat , class ThermalConductivity >
GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity >::LowMachNavierStokesStabilizationBase ( )
private

Member Function Documentation

template<class Mu , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::compute_res_continuity_steady ( AssemblyContext context,
unsigned int  qp 
) const

Definition at line 70 of file low_mach_navier_stokes_stab_base.C.

72  {
73  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
74  libMesh::RealGradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
75 
76  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
77  context.fixed_interior_value(this->_flow_vars.v(), qp) );
78 
79  libMesh::RealGradient grad_u, grad_v;
80 
81  grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
82  grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
83 
84  libMesh::Real divU = grad_u(0) + grad_v(1);
85 
86  if( this->_dim == 3 )
87  {
88  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
89  divU += (context.fixed_interior_gradient(this->_flow_vars.w(), qp))(2);
90  }
91 
92  return divU - (U*grad_T)/T;
93  }
unsigned int _dim
Physical dimension of problem.
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
template<class Mu , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::compute_res_continuity_transient ( AssemblyContext context,
unsigned int  qp 
) const

Definition at line 96 of file low_mach_navier_stokes_stab_base.C.

98  {
99  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
100  libMesh::Real T_dot = context.interior_value(this->_temp_vars.T(), qp);
101 
102  libMesh::Real RC_t = -T_dot/T;
103 
104  if( this->_enable_thermo_press_calc )
105  {
106  libMesh::Real p0 = context.fixed_interior_value(this->_p0_var->p0(), qp);
107  libMesh::Real p0_dot = context.interior_value(this->_p0_var->p0(), qp);
108 
109  RC_t += p0_dot/p0;
110  }
111 
112  return RC_t;
113  }
libMesh::UniquePtr< ThermoPressureFEVariable > _p0_var
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
template<class Mu , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::compute_res_energy_steady ( AssemblyContext context,
unsigned int  qp 
) const

Definition at line 220 of file low_mach_navier_stokes_stab_base.C.

222  {
223  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
224  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
225  libMesh::Tensor hess_T = context.fixed_interior_hessian(this->_temp_vars.T(), qp);
226 
227  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
228  libMesh::Real rho_cp = rho*this->_cp(T);
229 
230  libMesh::RealGradient rhocpU( rho_cp*context.fixed_interior_value(this->_flow_vars.u(), qp),
231  rho_cp*context.fixed_interior_value(this->_flow_vars.v(), qp) );
232  if(this->_dim == 3)
233  rhocpU(2) = rho_cp*context.fixed_interior_value(this->_flow_vars.w(), qp);
234 
235  libMesh::Real hess_term = hess_T(0,0) + hess_T(1,1);
236 #if LIBMESH_DIM > 2
237  hess_term += hess_T(2,2);
238 #endif
239 
240  return rhocpU*grad_T - this->_k.deriv(T)*(grad_T*grad_T) - this->_k(T)*(hess_term);
241  }
unsigned int _dim
Physical dimension of problem.
SpecificHeat _cp
Specific heat object.
ThermalConductivity _k
Thermal conductivity object.
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
template<class Mu , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::compute_res_energy_transient ( AssemblyContext context,
unsigned int  qp 
) const

Definition at line 245 of file low_mach_navier_stokes_stab_base.C.

247  {
248  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
249  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
250  libMesh::Real rho_cp = rho*this->_cp(T);
251  libMesh::Real T_dot = context.interior_value(this->_temp_vars.T(), qp);
252 
253  libMesh::Real RE_t = rho_cp*T_dot;
254 
255  if( this->_enable_thermo_press_calc )
256  {
257  RE_t -= context.interior_value(this->_p0_var->p0(), qp);
258  }
259 
260  return RE_t;
261  }
libMesh::UniquePtr< ThermoPressureFEVariable > _p0_var
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
SpecificHeat _cp
Specific heat object.
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
template<class Mu , class SH , class TC >
libMesh::RealGradient GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::compute_res_momentum_steady ( AssemblyContext context,
unsigned int  qp 
) const

Definition at line 116 of file low_mach_navier_stokes_stab_base.C.

118  {
119  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
120 
121  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
122 
123  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
124  context.fixed_interior_value(this->_flow_vars.v(), qp) );
125  if(this->_dim == 3)
126  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
127 
128  libMesh::RealGradient grad_p = context.fixed_interior_gradient(this->_press_var.p(), qp);
129 
130  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
131  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
132 
133  libMesh::RealTensor hess_u = context.fixed_interior_hessian(this->_flow_vars.u(), qp);
134  libMesh::RealTensor hess_v = context.fixed_interior_hessian(this->_flow_vars.v(), qp);
135 
136  libMesh::RealGradient rhoUdotGradU;
137  libMesh::RealGradient divGradU;
138  libMesh::RealGradient divGradUT;
139  libMesh::RealGradient divdivU;
140 
141  if( this->_dim < 3 )
142  {
143  rhoUdotGradU = rho*_stab_helper.UdotGradU( U, grad_u, grad_v );
144  divGradU = _stab_helper.div_GradU( hess_u, hess_v );
145  divGradUT = _stab_helper.div_GradU_T( hess_u, hess_v );
146  divdivU = _stab_helper.div_divU_I( hess_u, hess_v );
147  }
148  else
149  {
150  libMesh::RealGradient grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
151  libMesh::RealTensor hess_w = context.fixed_interior_hessian(this->_flow_vars.w(), qp);
152 
153  rhoUdotGradU = rho*_stab_helper.UdotGradU( U, grad_u, grad_v, grad_w );
154 
155  divGradU = _stab_helper.div_GradU( hess_u, hess_v, hess_w );
156  divGradUT = _stab_helper.div_GradU_T( hess_u, hess_v, hess_w );
157  divdivU = _stab_helper.div_divU_I( hess_u, hess_v, hess_w );
158  }
159 
160  libMesh::RealGradient divT = this->_mu(T)*(divGradU + divGradUT - 2.0/3.0*divdivU);
161 
162  if( this->_mu.deriv(T) != 0.0 )
163  {
164  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
165 
166  libMesh::Gradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
167  libMesh::Gradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
168 
169  libMesh::Gradient gradTgradu( grad_T*grad_u, grad_T*grad_v );
170 
171  libMesh::Gradient gradTgraduT( grad_T(0)*grad_u(0) + grad_T(1)*grad_u(1),
172  grad_T(0)*grad_v(0) + grad_T(1)*grad_v(1) );
173 
174  libMesh::Real divU = grad_u(0) + grad_v(1);
175 
176  libMesh::Gradient gradTdivU( grad_T(0)*divU, grad_T(1)*divU );
177 
178  if(this->_dim == 3)
179  {
180  libMesh::Gradient grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
181 
182  gradTgradu(2) = grad_T*grad_w;
183 
184  gradTgraduT(0) += grad_T(2)*grad_u(2);
185  gradTgraduT(1) += grad_T(2)*grad_v(2);
186  gradTgraduT(2) = grad_T(0)*grad_w(0) + grad_T(1)*grad_w(1) + grad_T(2)*grad_w(2);
187 
188  divU += grad_w(2);
189  gradTdivU(0) += grad_T(0)*grad_w(2);
190  gradTdivU(1) += grad_T(1)*grad_w(2);
191  gradTdivU(2) += grad_T(2)*divU;
192  }
193 
194  divT += this->_mu.deriv(T)*( gradTgradu + gradTgraduT - 2.0/3.0*gradTdivU );
195  }
196 
197  libMesh::RealGradient rhog( rho*this->_g(0), rho*this->_g(1) );
198  if(this->_dim == 3)
199  rhog(2) = rho*this->_g(2);
200 
201  return rhoUdotGradU + grad_p - divT - rhog;
202  }
libMesh::RealGradient div_divU_I(libMesh::RealTensor &hess_u, libMesh::RealTensor &hess_v) const
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient div_GradU(libMesh::RealTensor &hess_u, libMesh::RealTensor &hess_v) const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
VariableIndex p() const
libMesh::Point _g
Gravity vector.
libMesh::RealGradient UdotGradU(libMesh::Gradient &U, libMesh::Gradient &grad_u, libMesh::Gradient &grad_v) const
libMesh::RealGradient div_GradU_T(libMesh::RealTensor &hess_u, libMesh::RealTensor &hess_v) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
template<class Mu , class SH , class TC >
libMesh::RealGradient GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::compute_res_momentum_transient ( AssemblyContext context,
unsigned int  qp 
) const

Definition at line 205 of file low_mach_navier_stokes_stab_base.C.

207  {
208  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
209  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
210 
211  libMesh::RealGradient u_dot( context.interior_value(this->_flow_vars.u(), qp), context.interior_value(this->_flow_vars.v(), qp) );
212 
213  if(this->_dim == 3)
214  u_dot(2) = context.interior_value(this->_flow_vars.w(), qp);
215 
216  return rho*u_dot;
217  }
unsigned int _dim
Physical dimension of problem.
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >.

Definition at line 54 of file low_mach_navier_stokes_stab_base.C.

References GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::init_context().

55  {
56  // First call base class
58 
59  // We need pressure derivatives
60  context.get_element_fe(this->_press_var.p())->get_dphi();
61 
62  // We also need second derivatives, so initialize those.
63  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
64  context.get_element_fe(this->_temp_vars.T())->get_d2phi();
65 
66  return;
67  }
VariableIndex p() const
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.

Member Data Documentation

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
LowMachNavierStokesStabilizationHelper GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity >::_stab_helper
protected

Definition at line 69 of file low_mach_navier_stokes_stab_base.h.


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