GRINS-0.8.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 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 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 element_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for element interiors. 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 mass_residual (bool, AssemblyContext &)
 Mass 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

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
 
VelocityVariable_flow_vars
 
PressureFEVariable_press_var
 
PrimitiveTempFEVariables_temp_vars
 
ThermoPressureVariable_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)
 
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 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  {}
static PhysicsName low_mach_navier_stokes()
LowMachNavierStokesStabilizationHelper _stab_helper
template<class Viscosity , class SpecificHeat , class ThermalConductivity >
virtual GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity >::~LowMachNavierStokesStabilizationBase ( )
inlinevirtual

Definition at line 44 of file low_mach_navier_stokes_stab_base.h.

44 {};
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 60 of file low_mach_navier_stokes_stab_base.C.

62  {
63  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
64  libMesh::RealGradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
65 
66  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
67  context.fixed_interior_value(this->_flow_vars.v(), qp) );
68 
69  libMesh::RealGradient grad_u, grad_v;
70 
71  grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
72  grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
73 
74  libMesh::Real divU = grad_u(0) + grad_v(1);
75 
76 
77  if( this->_flow_vars.dim() == 3 )
78  {
79  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
80  divU += (context.fixed_interior_gradient(this->_flow_vars.w(), qp))(2);
81  }
82 
83  return divU - (U*grad_T)/T;
84  }
VariableIndex T() const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
unsigned int dim() const
Number of components.
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 87 of file low_mach_navier_stokes_stab_base.C.

89  {
90  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
91  libMesh::Real T_dot = context.interior_value(this->_temp_vars.T(), qp);
92 
93  libMesh::Real RC_t = -T_dot/T;
94 
95  if( this->_enable_thermo_press_calc )
96  {
97  libMesh::Real p0 = context.fixed_interior_value(this->_p0_var->p0(), qp);
98  libMesh::Real p0_dot = context.interior_value(this->_p0_var->p0(), qp);
99 
100  RC_t += p0_dot/p0;
101  }
102 
103  return RC_t;
104  }
VariableIndex T() const
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
VariableIndex p0() 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 212 of file low_mach_navier_stokes_stab_base.C.

214  {
215  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
216  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
217  libMesh::Tensor hess_T = context.fixed_interior_hessian(this->_temp_vars.T(), qp);
218 
219  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
220  libMesh::Real rho_cp = rho*this->_cp(T);
221 
222  libMesh::RealGradient rhocpU( rho_cp*context.fixed_interior_value(this->_flow_vars.u(), qp),
223  rho_cp*context.fixed_interior_value(this->_flow_vars.v(), qp) );
224 
225  if(this->_flow_vars.dim() == 3)
226  rhocpU(2) = rho_cp*context.fixed_interior_value(this->_flow_vars.w(), qp);
227 
228  libMesh::Real hess_term = hess_T(0,0) + hess_T(1,1);
229 #if LIBMESH_DIM > 2
230  hess_term += hess_T(2,2);
231 #endif
232 
233  return rhocpU*grad_T - this->_k.deriv(T)*(grad_T*grad_T) - this->_k(T)*(hess_term);
234  }
VariableIndex T() const
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
unsigned int dim() const
Number of components.
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 238 of file low_mach_navier_stokes_stab_base.C.

240  {
241  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
242  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
243  libMesh::Real rho_cp = rho*this->_cp(T);
244  libMesh::Real T_dot = context.interior_value(this->_temp_vars.T(), qp);
245 
246  libMesh::Real RE_t = rho_cp*T_dot;
247 
248  if( this->_enable_thermo_press_calc )
249  {
250  RE_t -= context.interior_value(this->_p0_var->p0(), qp);
251  }
252 
253  return RE_t;
254  }
VariableIndex T() const
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
VariableIndex p0() 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 107 of file low_mach_navier_stokes_stab_base.C.

109  {
110  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
111 
112  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
113 
114  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
115  context.fixed_interior_value(this->_flow_vars.v(), qp) );
116 
117  if(this->_flow_vars.dim() == 3)
118  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
119 
120  libMesh::RealGradient grad_p = context.fixed_interior_gradient(this->_press_var.p(), qp);
121 
122  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
123  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
124 
125  libMesh::RealTensor hess_u = context.fixed_interior_hessian(this->_flow_vars.u(), qp);
126  libMesh::RealTensor hess_v = context.fixed_interior_hessian(this->_flow_vars.v(), qp);
127 
128  libMesh::RealGradient rhoUdotGradU;
129  libMesh::RealGradient divGradU;
130  libMesh::RealGradient divGradUT;
131  libMesh::RealGradient divdivU;
132 
133  if( this->_flow_vars.dim() < 3 )
134  {
135  rhoUdotGradU = rho*_stab_helper.UdotGradU( U, grad_u, grad_v );
136  divGradU = _stab_helper.div_GradU( hess_u, hess_v );
137  divGradUT = _stab_helper.div_GradU_T( hess_u, hess_v );
138  divdivU = _stab_helper.div_divU_I( hess_u, hess_v );
139  }
140  else
141  {
142  libMesh::RealGradient grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
143  libMesh::RealTensor hess_w = context.fixed_interior_hessian(this->_flow_vars.w(), qp);
144 
145  rhoUdotGradU = rho*_stab_helper.UdotGradU( U, grad_u, grad_v, grad_w );
146 
147  divGradU = _stab_helper.div_GradU( hess_u, hess_v, hess_w );
148  divGradUT = _stab_helper.div_GradU_T( hess_u, hess_v, hess_w );
149  divdivU = _stab_helper.div_divU_I( hess_u, hess_v, hess_w );
150  }
151 
152  libMesh::RealGradient divT = this->_mu(T)*(divGradU + divGradUT - 2.0/3.0*divdivU);
153 
154  if( this->_mu.deriv(T) != 0.0 )
155  {
156  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
157 
158  libMesh::Gradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
159  libMesh::Gradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
160 
161  libMesh::Gradient gradTgradu( grad_T*grad_u, grad_T*grad_v );
162 
163  libMesh::Gradient gradTgraduT( grad_T(0)*grad_u(0) + grad_T(1)*grad_u(1),
164  grad_T(0)*grad_v(0) + grad_T(1)*grad_v(1) );
165 
166  libMesh::Real divU = grad_u(0) + grad_v(1);
167 
168  libMesh::Gradient gradTdivU( grad_T(0)*divU, grad_T(1)*divU );
169 
170  if(this->_flow_vars.dim() == 3)
171  {
172  libMesh::Gradient grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
173 
174  gradTgradu(2) = grad_T*grad_w;
175 
176  gradTgraduT(0) += grad_T(2)*grad_u(2);
177  gradTgraduT(1) += grad_T(2)*grad_v(2);
178  gradTgraduT(2) = grad_T(0)*grad_w(0) + grad_T(1)*grad_w(1) + grad_T(2)*grad_w(2);
179 
180  divU += grad_w(2);
181  gradTdivU(0) += grad_T(0)*grad_w(2);
182  gradTdivU(1) += grad_T(1)*grad_w(2);
183  gradTdivU(2) += grad_T(2)*divU;
184  }
185 
186  divT += this->_mu.deriv(T)*( gradTgradu + gradTgraduT - 2.0/3.0*gradTdivU );
187  }
188 
189  libMesh::RealGradient rhog( rho*this->_g(0), rho*this->_g(1) );
190  if(this->_flow_vars.dim() == 3)
191  rhog(2) = rho*this->_g(2);
192 
193  return rhoUdotGradU + grad_p - divT - rhog;
194  }
libMesh::RealGradient div_GradU_T(libMesh::RealTensor &hess_u) const
VariableIndex T() const
libMesh::RealGradient UdotGradU(libMesh::Gradient &U, libMesh::Gradient &grad_u) const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
VariableIndex p() const
libMesh::Point _g
Gravity vector.
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
libMesh::RealGradient div_GradU(libMesh::RealTensor &hess_u) const
unsigned int dim() const
Number of components.
libMesh::RealGradient div_divU_I(libMesh::RealTensor &hess_u) 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 197 of file low_mach_navier_stokes_stab_base.C.

199  {
200  libMesh::Real T = context.fixed_interior_value(this->_temp_vars.T(), qp);
201  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
202 
203  libMesh::RealGradient u_dot( context.interior_value(this->_flow_vars.u(), qp), context.interior_value(this->_flow_vars.v(), qp) );
204 
205  if(this->_flow_vars.dim() == 3)
206  u_dot(2) = context.interior_value(this->_flow_vars.w(), qp);
207 
208  return rho*u_dot;
209  }
VariableIndex T() const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
unsigned int dim() const
Number of components.
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 46 of file low_mach_navier_stokes_stab_base.C.

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

47  {
48  // First call base class
50 
51  // We need pressure derivatives
52  context.get_element_fe(this->_press_var.p())->get_dphi();
53 
54  // We also need second derivatives, so initialize those.
55  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
56  context.get_element_fe(this->_temp_vars.T())->get_d2phi();
57  }
VariableIndex T() const
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 Tue Dec 19 2017 12:47:31 for GRINS-0.8.0 by  doxygen 1.8.9.1