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

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

#include <low_mach_navier_stokes_vms_stab.h>

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

Public Member Functions

 LowMachNavierStokesVMSStabilization (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~LowMachNavierStokesVMSStabilization ()
 
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...
 
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
 
virtual void read_input_options (const GetPot &input)
 Read options from GetPot input file. More...
 
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 (const AssemblyContext &c, const libMesh::Point &p) const
 
libMesh::Real get_p0_steady_side (const AssemblyContext &c, unsigned int qp) const
 
libMesh::Real get_p0_transient (AssemblyContext &c, unsigned int qp) const
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. 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 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 nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
void init_bcs (libMesh::FEMSystem *system)
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
void attach_neumann_bound_func (GRINS::NBCContainer &neumann_bcs)
 
void attach_dirichlet_bound_func (const GRINS::DBCContainer &dirichlet_bc)
 
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_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)
 
BCHandlingBaseget_bc_handler ()
 
ICHandlingBaseget_ic_handler ()
 
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...
 

Protected Member Functions

void assemble_continuity_time_deriv (bool compute_jacobian, AssemblyContext &context)
 
void assemble_momentum_time_deriv (bool compute_jacobian, AssemblyContext &context)
 
void assemble_energy_time_deriv (bool compute_jacobian, AssemblyContext &context)
 
void assemble_continuity_mass_residual (bool compute_jacobian, AssemblyContext &context)
 
void assemble_momentum_mass_residual (bool compute_jacobian, AssemblyContext &context)
 
void assemble_energy_mass_residual (bool compute_jacobian, AssemblyContext &context)
 

Protected Attributes

LowMachNavierStokesStabilizationHelper _stab_helper
 
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...
 
VariableIndex _u_var
 Indices for each (owned) variable;. More...
 
VariableIndex _v_var
 
VariableIndex _w_var
 
VariableIndex _p_var
 
VariableIndex _T_var
 
VariableIndex _p0_var
 
std::string _u_var_name
 Names of each (owned) variable in the system. More...
 
std::string _v_var_name
 
std::string _w_var_name
 
std::string _p_var_name
 
std::string _T_var_name
 
std::string _p0_var_name
 
GRINSEnums::FEFamily _V_FE_family
 Element type, read from input. More...
 
GRINSEnums::FEFamily _P_FE_family
 
GRINSEnums::FEFamily _T_FE_family
 
GRINSEnums::Order _V_order
 Element orders, read from input. More...
 
GRINSEnums::Order _P_order
 
GRINSEnums::Order _T_order
 
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...
 
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::BCHandlingBase_bc_handler
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 
bool _is_axisymmetric
 

Static Protected Attributes

static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 

Private Member Functions

 LowMachNavierStokesVMSStabilization ()
 

Detailed Description

template<class Viscosity, class SpecificHeat, class ThermalConductivity>
class GRINS::LowMachNavierStokesVMSStabilization< Viscosity, SpecificHeat, ThermalConductivity >

Adds VMS-based stabilization to LowMachNavierStokes physics class.

Definition at line 36 of file low_mach_navier_stokes_vms_stab.h.

Constructor & Destructor Documentation

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

Definition at line 42 of file low_mach_navier_stokes_vms_stab.C.

44  : LowMachNavierStokesStabilizationBase<Mu,SH,TC>(physics_name,input)
45  {
46  return;
47  }
template<class Mu , class SH , class TC >
GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::~LowMachNavierStokesVMSStabilization ( )
virtual

Definition at line 50 of file low_mach_navier_stokes_vms_stab.C.

51  {
52  return;
53  }
template<class Viscosity , class SpecificHeat , class ThermalConductivity >
GRINS::LowMachNavierStokesVMSStabilization< Viscosity, SpecificHeat, ThermalConductivity >::LowMachNavierStokesVMSStabilization ( )
private

Member Function Documentation

template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::assemble_continuity_mass_residual ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 299 of file low_mach_navier_stokes_vms_stab.C.

301  {
302  // The number of local degrees of freedom in each variable.
303  const unsigned int n_p_dofs = context.get_dof_indices(this->_p_var).size();
304 
305  // Element Jacobian * quadrature weights for interior integration.
306  const std::vector<libMesh::Real> &JxW =
307  context.get_element_fe(this->_u_var)->get_JxW();
308 
309  // The pressure shape functions at interior quadrature points.
310  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
311  context.get_element_fe(this->_p_var)->get_dphi();
312 
313  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_p_var); // R_{p}
314 
315  unsigned int n_qpoints = context.get_element_qrule().n_points();
316 
317  for (unsigned int qp=0; qp != n_qpoints; qp++)
318  {
319  libMesh::FEBase* fe = context.get_element_fe(this->_u_var);
320 
321  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
322  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
323 
324  libMesh::Real T = context.fixed_interior_value( this->_T_var, qp );
325  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
326 
327  libMesh::Real mu = this->_mu(T);
328 
329  libMesh::RealGradient U( context.fixed_interior_value( this->_u_var, qp ),
330  context.fixed_interior_value( this->_v_var, qp ) );
331  if( this->_dim == 3 )
332  U(2) = context.fixed_interior_value( this->_w_var, qp );
333 
334  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
335  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
336 
337  // Now a loop over the pressure degrees of freedom. This
338  // computes the contributions of the continuity equation.
339  for (unsigned int i=0; i != n_p_dofs; i++)
340  {
341  Fp(i) += tau_M*RM_t*p_dphi[i][qp]*JxW[qp];
342  }
343  }
344 
345  return;
346  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
VariableIndex _u_var
Indices for each (owned) variable;.
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient compute_res_momentum_transient(AssemblyContext &context, unsigned int qp) const
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
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
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::assemble_continuity_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 94 of file low_mach_navier_stokes_vms_stab.C.

96  {
97  // The number of local degrees of freedom in each variable.
98  const unsigned int n_p_dofs = context.get_dof_indices(this->_p_var).size();
99 
100  // Element Jacobian * quadrature weights for interior integration.
101  const std::vector<libMesh::Real> &JxW =
102  context.get_element_fe(this->_u_var)->get_JxW();
103 
104  // The pressure shape functions at interior quadrature points.
105  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
106  context.get_element_fe(this->_p_var)->get_dphi();
107 
108  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_p_var); // R_{p}
109 
110  unsigned int n_qpoints = context.get_element_qrule().n_points();
111 
112  for (unsigned int qp=0; qp != n_qpoints; qp++)
113  {
114  libMesh::FEBase* fe = context.get_element_fe(this->_u_var);
115 
116  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
117  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
118 
119  libMesh::Real T = context.interior_value( this->_T_var, qp );
120 
121  libMesh::Real mu = this->_mu(T);
122 
123  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
124 
125  libMesh::RealGradient U( context.interior_value( this->_u_var, qp ),
126  context.interior_value( this->_v_var, qp ) );
127  if( this->_dim == 3 )
128  U(2) = context.interior_value( this->_w_var, qp );
129 
130  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
131 
132  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
133 
134  // Now a loop over the pressure degrees of freedom. This
135  // computes the contributions of the continuity equation.
136  for (unsigned int i=0; i != n_p_dofs; i++)
137  {
138  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*JxW[qp];
139  }
140 
141  }
142 
143  return;
144  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
VariableIndex _u_var
Indices for each (owned) variable;.
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) 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_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
libMesh::RealGradient compute_res_momentum_steady(AssemblyContext &context, unsigned int qp) const
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::assemble_energy_mass_residual ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 438 of file low_mach_navier_stokes_vms_stab.C.

440  {
441  // The number of local degrees of freedom in each variable.
442  const unsigned int n_T_dofs = context.get_dof_indices(this->_T_var).size();
443 
444  // Element Jacobian * quadrature weights for interior integration.
445  const std::vector<libMesh::Real> &JxW =
446  context.get_element_fe(this->_T_var)->get_JxW();
447 
448  // The temperature shape functions at interior quadrature points.
449  const std::vector<std::vector<libMesh::Real> >& T_phi =
450  context.get_element_fe(this->_T_var)->get_phi();
451 
452  // The temperature shape functions gradients at interior quadrature points.
453  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
454  context.get_element_fe(this->_T_var)->get_dphi();
455 
456  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_T_var); // R_{T}
457 
458  unsigned int n_qpoints = context.get_element_qrule().n_points();
459 
460  for (unsigned int qp=0; qp != n_qpoints; qp++)
461  {
462  libMesh::Number u, v;
463  u = context.fixed_interior_value(this->_u_var, qp);
464  v = context.fixed_interior_value(this->_v_var, qp);
465 
466  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_T_var, qp);
467 
468  libMesh::NumberVectorValue U(u,v);
469  if (this->_dim == 3)
470  U(2) = context.fixed_interior_value(this->_w_var, qp); // w
471 
472  libMesh::Real T = context.fixed_interior_value( this->_T_var, qp );
473  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
474 
475  libMesh::Real mu = this->_mu(T);
476  libMesh::Real k = this->_k(T);
477  libMesh::Real cp = this->_cp(T);
478 
479  libMesh::Number rho_cp = rho*this->_cp(T);
480 
481  libMesh::FEBase* fe = context.get_element_fe(this->_u_var);
482 
483  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
484  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
485 
486  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
487  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
488 
489  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
490  libMesh::Real RE_t = this->compute_res_energy_transient( context, qp );
491 
492  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
493  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
494 
495  for (unsigned int i=0; i != n_T_dofs; i++)
496  {
497  FT(i) -= ( -rho_cp*tau_M*RM_t*grad_T*T_phi[i][qp]
498  +rho_cp*tau_E*RE_t*U*T_gradphi[i][qp]
499  - rho_cp*tau_E*(RE_s+RE_t)*tau_M*RM_t*T_gradphi[i][qp]
500  - rho_cp*tau_E*RE_t*tau_M*RM_s*T_gradphi[i][qp] )*JxW[qp];
501  }
502 
503  }
504 
505  return;
506  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
VariableIndex _u_var
Indices for each (owned) variable;.
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient compute_res_momentum_transient(AssemblyContext &context, unsigned int qp) const
SpecificHeat _cp
Specific heat object.
ThermalConductivity _k
Thermal conductivity object.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) 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
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
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
libMesh::RealGradient compute_res_momentum_steady(AssemblyContext &context, unsigned int qp) const
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
libMesh::Real compute_tau_energy(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real k, libMesh::Real cp, bool is_steady) const
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::assemble_energy_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 232 of file low_mach_navier_stokes_vms_stab.C.

234  {
235  // The number of local degrees of freedom in each variable.
236  const unsigned int n_T_dofs = context.get_dof_indices(this->_T_var).size();
237 
238  // Element Jacobian * quadrature weights for interior integration.
239  const std::vector<libMesh::Real> &JxW =
240  context.get_element_fe(this->_T_var)->get_JxW();
241 
242  // The temperature shape functions at interior quadrature points.
243  const std::vector<std::vector<libMesh::Real> >& T_phi =
244  context.get_element_fe(this->_T_var)->get_phi();
245 
246  // The temperature shape functions gradients at interior quadrature points.
247  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
248  context.get_element_fe(this->_T_var)->get_dphi();
249 
250  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_T_var); // R_{T}
251 
252  unsigned int n_qpoints = context.get_element_qrule().n_points();
253 
254  for (unsigned int qp=0; qp != n_qpoints; qp++)
255  {
256  libMesh::Number u, v;
257  u = context.interior_value(this->_u_var, qp);
258  v = context.interior_value(this->_v_var, qp);
259 
260  libMesh::Gradient grad_T = context.interior_gradient(this->_T_var, qp);
261 
262  libMesh::NumberVectorValue U(u,v);
263  if (this->_dim == 3)
264  U(2) = context.interior_value(this->_w_var, qp); // w
265 
266  libMesh::Real T = context.interior_value( this->_T_var, qp );
267  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
268 
269  libMesh::Real mu = this->_mu(T);
270  libMesh::Real k = this->_k(T);
271  libMesh::Real cp = this->_cp(T);
272 
273  libMesh::Number rho_cp = rho*this->_cp(T);
274 
275  libMesh::FEBase* fe = context.get_element_fe(this->_u_var);
276 
277  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
278  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
279 
280  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
281  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, this->_is_steady );
282 
283  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
284  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
285 
286  for (unsigned int i=0; i != n_T_dofs; i++)
287  {
288  FT(i) += ( rho_cp*tau_M*RM_s*grad_T*T_phi[i][qp]
289  - rho_cp*tau_E*RE_s*U*T_gradphi[i][qp]
290  + rho_cp*tau_E*RE_s*tau_M*RM_s*T_gradphi[i][qp] )*JxW[qp];
291  }
292 
293  }
294 
295  return;
296  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
VariableIndex _u_var
Indices for each (owned) variable;.
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int _dim
Physical dimension of problem.
SpecificHeat _cp
Specific heat object.
ThermalConductivity _k
Thermal conductivity object.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real compute_res_energy_steady(AssemblyContext &context, unsigned int qp) 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_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
libMesh::RealGradient compute_res_momentum_steady(AssemblyContext &context, unsigned int qp) const
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
libMesh::Real compute_tau_energy(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real k, libMesh::Real cp, bool is_steady) const
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::assemble_momentum_mass_residual ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 349 of file low_mach_navier_stokes_vms_stab.C.

351  {
352  // The number of local degrees of freedom in each variable.
353  const unsigned int n_u_dofs = context.get_dof_indices(this->_u_var).size();
354 
355  // Check number of dofs is same for _u_var, v_var and w_var.
356  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_v_var).size());
357  if (this->_dim == 3)
358  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_w_var).size());
359 
360  // Element Jacobian * quadrature weights for interior integration.
361  const std::vector<libMesh::Real> &JxW =
362  context.get_element_fe(this->_u_var)->get_JxW();
363 
364  // The pressure shape functions at interior quadrature points.
365  const std::vector<std::vector<libMesh::Real> >& u_phi =
366  context.get_element_fe(this->_u_var)->get_phi();
367 
368  // The velocity shape function gradients at interior quadrature points.
369  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
370  context.get_element_fe(this->_u_var)->get_dphi();
371 
372  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_u_var); // R_{u}
373  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_v_var); // R_{v}
374  libMesh::DenseSubVector<libMesh::Number> &Fw = context.get_elem_residual(this->_w_var); // R_{w}
375 
376  unsigned int n_qpoints = context.get_element_qrule().n_points();
377  for (unsigned int qp=0; qp != n_qpoints; qp++)
378  {
379  libMesh::Real T = context.fixed_interior_value( this->_T_var, qp );
380  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
381 
382  libMesh::Real mu = this->_mu(T);
383 
384  libMesh::RealGradient U( context.fixed_interior_value(this->_u_var, qp),
385  context.fixed_interior_value(this->_v_var, qp) );
386 
387  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_u_var, qp);
388  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_v_var, qp);
389  libMesh::RealGradient grad_w;
390 
391  if( this->_dim == 3 )
392  {
393  U(2) = context.fixed_interior_value(this->_w_var, qp);
394  grad_w = context.fixed_interior_gradient(this->_w_var, qp);
395  }
396 
397  libMesh::FEBase* fe = context.get_element_fe(this->_u_var);
398 
399  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
400  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
401 
402  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
403  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
404 
405  libMesh::Real RC_t = this->compute_res_continuity_transient( context, qp );
406  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
407  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
408 
409  for (unsigned int i=0; i != n_u_dofs; i++)
410  {
411  Fu(i) -= ( tau_C*RC_t*u_gradphi[i][qp](0)
412  + tau_M*RM_t(0)*rho*U*u_gradphi[i][qp]
413  - rho*tau_M*RM_t*grad_u*u_phi[i][qp]
414  - tau_M*(RM_s(0)+RM_t(0))*rho*tau_M*RM_t*u_gradphi[i][qp]
415  - tau_M*RM_t(0)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
416 
417  Fv(i) -= ( tau_C*RC_t*u_gradphi[i][qp](1)
418  - rho*tau_M*RM_t*grad_v*u_phi[i][qp]
419  + tau_M*RM_t(1)*rho*U*u_gradphi[i][qp]
420  - tau_M*(RM_s(1)+RM_t(1))*rho*tau_M*RM_t*u_gradphi[i][qp]
421  - tau_M*RM_t(1)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
422 
423  if( this->_dim == 3 )
424  {
425  Fw(i) -= ( tau_C*RC_t*u_gradphi[i][qp](2)
426  - rho*tau_M*RM_t*grad_w*u_phi[i][qp]
427  + tau_M*RM_t(2)*rho*U*u_gradphi[i][qp]
428  - tau_M*(RM_s(2)+RM_t(2))*rho*tau_M*RM_t*u_gradphi[i][qp]
429  - tau_M*RM_t(2)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
430  }
431  }
432 
433  }
434  return;
435  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
libMesh::Real compute_tau_continuity(libMesh::Real tau_C, libMesh::RealGradient &g) const
VariableIndex _u_var
Indices for each (owned) variable;.
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient compute_res_momentum_transient(AssemblyContext &context, unsigned int qp) const
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
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
libMesh::RealGradient compute_res_momentum_steady(AssemblyContext &context, unsigned int qp) const
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real compute_res_continuity_transient(AssemblyContext &context, unsigned int qp) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::assemble_momentum_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 147 of file low_mach_navier_stokes_vms_stab.C.

149  {
150  // The number of local degrees of freedom in each variable.
151  const unsigned int n_u_dofs = context.get_dof_indices(this->_u_var).size();
152 
153  // Check number of dofs is same for _u_var, v_var and w_var.
154  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_v_var).size());
155  if (this->_dim == 3)
156  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_w_var).size());
157 
158  // Element Jacobian * quadrature weights for interior integration.
159  const std::vector<libMesh::Real> &JxW =
160  context.get_element_fe(this->_u_var)->get_JxW();
161 
162  // The pressure shape functions at interior quadrature points.
163  const std::vector<std::vector<libMesh::Real> >& u_phi =
164  context.get_element_fe(this->_u_var)->get_phi();
165 
166  // The velocity shape function gradients at interior quadrature points.
167  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
168  context.get_element_fe(this->_u_var)->get_dphi();
169 
170  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_u_var); // R_{u}
171  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_v_var); // R_{v}
172  libMesh::DenseSubVector<libMesh::Number> &Fw = context.get_elem_residual(this->_w_var); // R_{w}
173 
174  unsigned int n_qpoints = context.get_element_qrule().n_points();
175 
176  for (unsigned int qp=0; qp != n_qpoints; qp++)
177  {
178  libMesh::Real T = context.interior_value( this->_T_var, qp );
179  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
180 
181  libMesh::RealGradient U( context.interior_value(this->_u_var, qp),
182  context.interior_value(this->_v_var, qp) );
183 
184  libMesh::RealGradient grad_u = context.interior_gradient(this->_u_var, qp);
185  libMesh::RealGradient grad_v = context.interior_gradient(this->_v_var, qp);
186  libMesh::RealGradient grad_w;
187 
188  if( this->_dim == 3 )
189  {
190  U(2) = context.interior_value(this->_w_var, qp);
191  grad_w = context.interior_gradient(this->_w_var, qp);
192  }
193 
194  libMesh::FEBase* fe = context.get_element_fe(this->_u_var);
195 
196  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
197  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
198  libMesh::Real mu = this->_mu(T);
199 
200  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
201  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
202 
203  libMesh::Real RC_s = this->compute_res_continuity_steady( context, qp );
204  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
205 
206  for (unsigned int i=0; i != n_u_dofs; i++)
207  {
208  Fu(i) += ( -tau_C*RC_s*u_gradphi[i][qp](0)
209  - tau_M*RM_s(0)*rho*U*u_gradphi[i][qp]
210  + rho*tau_M*RM_s*grad_u*u_phi[i][qp]
211  + tau_M*RM_s(0)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
212 
213  Fv(i) += ( -tau_C*RC_s*u_gradphi[i][qp](1)
214  - tau_M*RM_s(1)*rho*U*u_gradphi[i][qp]
215  + rho*tau_M*RM_s*grad_v*u_phi[i][qp]
216  + tau_M*RM_s(1)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
217 
218  if( this->_dim == 3 )
219  {
220  Fw(i) += ( -tau_C*RC_s*u_gradphi[i][qp](2)
221  - tau_M*RM_s(2)*rho*U*u_gradphi[i][qp]
222  + rho*tau_M*RM_s*grad_w*u_phi[i][qp]
223  + tau_M*RM_s(2)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
224  }
225  }
226 
227  }
228  return;
229  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
libMesh::Real compute_tau_continuity(libMesh::Real tau_C, libMesh::RealGradient &g) const
VariableIndex _u_var
Indices for each (owned) variable;.
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) 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_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
libMesh::RealGradient compute_res_momentum_steady(AssemblyContext &context, unsigned int qp) const
libMesh::Real compute_res_continuity_steady(AssemblyContext &context, unsigned int qp) const
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
void GRINS::Physics::attach_dirichlet_bound_func ( const GRINS::DBCContainer dirichlet_bc)
inherited

Definition at line 150 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_dirichlet_bound_func().

151  {
152  _bc_handler->attach_dirichlet_bound_func( dirichlet_bc );
153  return;
154  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_dirichlet_bound_func(const GRINS::DBCContainer &dirichlet_bc)
void GRINS::Physics::attach_neumann_bound_func ( GRINS::NBCContainer neumann_bcs)
inherited

Definition at line 144 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_neumann_bound_func().

145  {
146  _bc_handler->attach_neumann_bound_func( neumann_bcs );
147  return;
148  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_neumann_bound_func(GRINS::NBCContainer &neumann_bcs)
void GRINS::Physics::auxiliary_init ( MultiphysicsSystem system)
virtualinherited

Any auxillary initialization a Physics class may need.

This is called after all variables are added, so this method can safely query the MultiphysicsSystem about variable information.

Definition at line 113 of file physics.C.

114  {
115  return;
116  }
void GRINS::Physics::compute_element_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 185 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::element_constraint().

187  {
188  return;
189  }
void GRINS::Physics::compute_element_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited
void GRINS::Physics::compute_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 203 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::mass_residual().

205  {
206  return;
207  }
void GRINS::Physics::compute_nonlocal_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 197 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

199  {
200  return;
201  }
void GRINS::Physics::compute_nonlocal_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 209 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

211  {
212  return;
213  }
void GRINS::Physics::compute_nonlocal_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 179 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

181  {
182  return;
183  }
void GRINS::Physics::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtualinherited
template<class Mu , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::compute_res_continuity_steady ( AssemblyContext context,
unsigned int  qp 
) const
inherited

Definition at line 70 of file low_mach_navier_stokes_stab_base.C.

72  {
73  libMesh::Real T = context.fixed_interior_value(this->_T_var, qp);
74  libMesh::RealGradient grad_T = context.fixed_interior_gradient(this->_T_var, qp);
75 
76  libMesh::RealGradient U( context.fixed_interior_value(this->_u_var, qp),
77  context.fixed_interior_value(this->_v_var, qp) );
78 
79  libMesh::RealGradient grad_u, grad_v;
80 
81  grad_u = context.fixed_interior_gradient(this->_u_var, qp);
82  grad_v = context.fixed_interior_gradient(this->_v_var, 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->_w_var, qp);
89  divU += (context.fixed_interior_gradient(this->_w_var, qp))(2);
90  }
91 
92  return divU - (U*grad_T)/T;
93  }
VariableIndex _u_var
Indices for each (owned) variable;.
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
inherited

Definition at line 96 of file low_mach_navier_stokes_stab_base.C.

98  {
99  libMesh::Real T = context.fixed_interior_value(this->_T_var, qp);
100  libMesh::Real T_dot = context.interior_value(this->_T_var, 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, qp);
107  libMesh::Real p0_dot = context.interior_value(this->_p0_var, qp);
108 
109  RC_t += p0_dot/p0;
110  }
111 
112  return RC_t;
113  }
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
inherited

Definition at line 220 of file low_mach_navier_stokes_stab_base.C.

222  {
223  libMesh::Real T = context.fixed_interior_value(this->_T_var, qp);
224  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_T_var, qp);
225  libMesh::Tensor hess_T = context.fixed_interior_hessian(this->_T_var, 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->_u_var, qp),
231  rho_cp*context.fixed_interior_value(this->_v_var, qp) );
232  if(this->_dim == 3)
233  rhocpU(2) = rho_cp*context.fixed_interior_value(this->_w_var, 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
inherited

Definition at line 245 of file low_mach_navier_stokes_stab_base.C.

247  {
248  libMesh::Real T = context.fixed_interior_value(this->_T_var, 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->_T_var, 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, qp);
258  }
259 
260  return RE_t;
261  }
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
inherited

Definition at line 116 of file low_mach_navier_stokes_stab_base.C.

118  {
119  libMesh::Real T = context.fixed_interior_value(this->_T_var, 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->_u_var, qp),
124  context.fixed_interior_value(this->_v_var, qp) );
125  if(this->_dim == 3)
126  U(2) = context.fixed_interior_value(this->_w_var, qp);
127 
128  libMesh::RealGradient grad_p = context.fixed_interior_gradient(this->_p_var, qp);
129 
130  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_u_var, qp);
131  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_v_var, qp);
132 
133  libMesh::RealTensor hess_u = context.fixed_interior_hessian(this->_u_var, qp);
134  libMesh::RealTensor hess_v = context.fixed_interior_hessian(this->_v_var, 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->_w_var, qp);
151  libMesh::RealTensor hess_w = context.fixed_interior_hessian(this->_w_var, 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->_T_var, qp);
165 
166  libMesh::Gradient grad_u = context.fixed_interior_gradient(this->_u_var, qp);
167  libMesh::Gradient grad_v = context.fixed_interior_gradient(this->_v_var, 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->_w_var, 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
VariableIndex _u_var
Indices for each (owned) variable;.
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
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
inherited

Definition at line 205 of file low_mach_navier_stokes_stab_base.C.

207  {
208  libMesh::Real T = context.fixed_interior_value(this->_T_var, qp);
209  libMesh::Real rho = this->rho(T, this->get_p0_transient(context,qp) );
210 
211  libMesh::RealGradient u_dot( context.interior_value(this->_u_var, qp), context.interior_value(this->_v_var, qp) );
212 
213  if(this->_dim == 3)
214  u_dot(2) = context.interior_value(this->_w_var, 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
void GRINS::Physics::compute_side_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 191 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

193  {
194  return;
195  }
void GRINS::Physics::compute_side_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Reimplemented in GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 173 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_time_derivative().

175  {
176  return;
177  }
template<class V , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesBase< V, SH, TC >::d_rho_dT ( libMesh::Real  T,
libMesh::Real  p0 
) const
inlineinherited

Definition at line 151 of file low_mach_navier_stokes_base.h.

152  {
153  return -p0/(this->_R*(T*T));
154  }
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
void GRINS::Physics::element_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::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 56 of file low_mach_navier_stokes_vms_stab.C.

59  {
60 #ifdef GRINS_USE_GRVY_TIMERS
61  this->_timer->BeginTimer("LowMachNavierStokesVMSStabilization::element_time_derivative");
62 #endif
63 
64  this->assemble_continuity_time_deriv( compute_jacobian, context );
65  this->assemble_momentum_time_deriv( compute_jacobian, context );
66  this->assemble_energy_time_deriv( compute_jacobian, context );
67 
68 #ifdef GRINS_USE_GRVY_TIMERS
69  this->_timer->EndTimer("LowMachNavierStokesVMSStabilization::element_time_derivative");
70 #endif
71  return;
72  }
void assemble_energy_time_deriv(bool compute_jacobian, AssemblyContext &context)
void assemble_momentum_time_deriv(bool compute_jacobian, AssemblyContext &context)
void assemble_continuity_time_deriv(bool compute_jacobian, AssemblyContext &context)
bool GRINS::Physics::enabled_on_elem ( const libMesh::Elem *  elem)
virtualinherited

Find if current physics is active on supplied element.

Definition at line 83 of file physics.C.

References GRINS::Physics::_enabled_subdomains.

84  {
85  // Check if enabled_subdomains flag has been set and if we're
86  // looking at a real element (rather than a nonlocal evaluation)
87  if( !elem || _enabled_subdomains.empty() )
88  return true;
89 
90  // Check if current physics is enabled on elem
91  if( _enabled_subdomains.find( elem->subdomain_id() ) == _enabled_subdomains.end() )
92  return false;
93 
94  return true;
95  }
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
BCHandlingBase * GRINS::Physics::get_bc_handler ( )
inlineinherited

Definition at line 282 of file physics.h.

References GRINS::Physics::_bc_handler.

283  {
284  return _bc_handler;
285  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
ICHandlingBase * GRINS::Physics::get_ic_handler ( )
inlineinherited

Definition at line 288 of file physics.h.

References GRINS::Physics::_ic_handler.

289  {
290  return _ic_handler;
291  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
template<class V , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesBase< V, SH, TC >::get_p0_steady ( const AssemblyContext c,
unsigned int  qp 
) const
inlineinherited

Definition at line 158 of file low_mach_navier_stokes_base.h.

160  {
161  libMesh::Real p0;
162  if( this->_enable_thermo_press_calc )
163  {
164  p0 = c.interior_value( _p0_var, qp );
165  }
166  else
167  {
168  p0 = _p0;
169  }
170  return p0;
171  }
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
template<class V , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesBase< V, SH, TC >::get_p0_steady ( const AssemblyContext c,
const libMesh::Point &  p 
) const
inlineinherited

Definition at line 192 of file low_mach_navier_stokes_base.h.

194  {
195  libMesh::Real p0;
196  if( this->_enable_thermo_press_calc )
197  {
198  p0 = c.point_value( _p0_var, p );
199  }
200  else
201  {
202  p0 = _p0;
203  }
204  return p0;
205  }
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
template<class V , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesBase< V, SH, TC >::get_p0_steady_side ( const AssemblyContext c,
unsigned int  qp 
) const
inlineinherited

Definition at line 175 of file low_mach_navier_stokes_base.h.

177  {
178  libMesh::Real p0;
179  if( this->_enable_thermo_press_calc )
180  {
181  p0 = c.side_value( _p0_var, qp );
182  }
183  else
184  {
185  p0 = _p0;
186  }
187  return p0;
188  }
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
template<class V , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesBase< V, SH, TC >::get_p0_transient ( AssemblyContext c,
unsigned int  qp 
) const
inlineinherited

Definition at line 209 of file low_mach_navier_stokes_base.h.

210  {
211  libMesh::Real p0;
212  if( this->_enable_thermo_press_calc )
213  {
214  p0 = c.fixed_interior_value( _p0_var, qp );
215  }
216  else
217  {
218  p0 = _p0;
219  }
220  return p0;
221  }
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
void GRINS::Physics::init_bcs ( libMesh::FEMSystem *  system)
inherited

Definition at line 118 of file physics.C.

References GRINS::Physics::_bc_handler, GRINS::BCHandlingBase::init_bc_data(), GRINS::BCHandlingBase::init_dirichlet_bc_func_objs(), GRINS::BCHandlingBase::init_dirichlet_bcs(), and GRINS::BCHandlingBase::init_periodic_bcs().

119  {
120  // Only need to init BC's if the physics actually created a handler
121  if( _bc_handler )
122  {
123  _bc_handler->init_bc_data( *system );
124  _bc_handler->init_dirichlet_bcs( system );
126  _bc_handler->init_periodic_bcs( system );
127  }
128 
129  return;
130  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
virtual void init_dirichlet_bcs(libMesh::FEMSystem *system) const
virtual void init_periodic_bcs(libMesh::FEMSystem *system) const
virtual void init_dirichlet_bc_func_objs(libMesh::FEMSystem *system) const
virtual void init_bc_data(const libMesh::FEMSystem &system)
Override this method to initialize any system-dependent data.
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesStabilizationBase< Mu, SH, TC >::init_context ( AssemblyContext context)
virtualinherited

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->_p_var)->get_dphi();
61 
62  // We also need second derivatives, so initialize those.
63  context.get_element_fe(this->_u_var)->get_d2phi();
64  context.get_element_fe(this->_T_var)->get_d2phi();
65 
66  return;
67  }
VariableIndex _u_var
Indices for each (owned) variable;.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
void GRINS::Physics::init_ics ( libMesh::FEMSystem *  system,
libMesh::CompositeFunction< libMesh::Number > &  all_ics 
)
inherited

Definition at line 133 of file physics.C.

References GRINS::Physics::_ic_handler, and GRINS::ICHandlingBase::init_ic_data().

135  {
136  if( _ic_handler )
137  {
138  _ic_handler->init_ic_data( *system, all_ics );
139  }
140 
141  return;
142  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
virtual void init_ic_data(const libMesh::FEMSystem &system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
Override this method to initialize any system-dependent data.
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesBase< Mu, SH, TC >::init_variables ( libMesh::FEMSystem *  system)
virtualinherited

Initialize variables for this physics.

Implements GRINS::Physics.

Definition at line 130 of file low_mach_navier_stokes_base.C.

131  {
132  // Get libMesh to assign an index for each variable
133  this->_dim = system->get_mesh().mesh_dimension();
134 
135  _u_var = system->add_variable( _u_var_name, this->_V_order, _V_FE_family);
136  _v_var = system->add_variable( _v_var_name, this->_V_order, _V_FE_family);
137 
138  if (_dim == 3)
139  _w_var = system->add_variable( _w_var_name, this->_V_order, _V_FE_family);
140  else
141  _w_var = _u_var;
142 
143  _p_var = system->add_variable( _p_var_name, this->_P_order, _P_FE_family);
144  _T_var = system->add_variable( _T_var_name, this->_T_order, _T_FE_family);
145 
146  /* If we need to compute the thermodynamic pressure, we force this to be a first
147  order scalar variable. */
149  _p0_var = system->add_variable( _p0_var_name, libMesh::FIRST, libMesh::SCALAR);
150 
151  return;
152  }
GRINSEnums::FEFamily _V_FE_family
Element type, read from input.
VariableIndex _u_var
Indices for each (owned) variable;.
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
unsigned int _dim
Physical dimension of problem.
std::string _u_var_name
Names of each (owned) variable in the system.
GRINSEnums::Order _V_order
Element orders, read from input.
bool GRINS::Physics::is_steady ( ) const
inherited

Returns whether or not this physics is being solved with a steady solver.

Definition at line 103 of file physics.C.

References GRINS::Physics::_is_steady.

Referenced by GRINS::Physics::set_is_steady().

104  {
105  return _is_steady;
106  }
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::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 75 of file low_mach_navier_stokes_vms_stab.C.

78  {
79 #ifdef GRINS_USE_GRVY_TIMERS
80  this->_timer->BeginTimer("LowMachNavierStokesVMSStabilization::mass_residual");
81 #endif
82 
83  this->assemble_continuity_mass_residual( compute_jacobian, context );
84  this->assemble_momentum_mass_residual( compute_jacobian, context );
85  this->assemble_energy_mass_residual( compute_jacobian, context );
86 
87 #ifdef GRINS_USE_GRVY_TIMERS
88  this->_timer->EndTimer("LowMachNavierStokesVMSStabilization::mass_residual");
89 #endif
90  return;
91  }
void assemble_momentum_mass_residual(bool compute_jacobian, AssemblyContext &context)
void assemble_continuity_mass_residual(bool compute_jacobian, AssemblyContext &context)
void assemble_energy_mass_residual(bool compute_jacobian, AssemblyContext &context)
void GRINS::Physics::nonlocal_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for scalar variables.

Reimplemented in GRINS::ScalarODE.

Definition at line 250 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

253  {
254  return;
255  }
void GRINS::Physics::nonlocal_mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Mass matrix part(s) for scalar variables.

Reimplemented in GRINS::ScalarODE, and GRINS::AveragedTurbine< Viscosity >.

Definition at line 264 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

267  {
268  return;
269  }
void GRINS::Physics::nonlocal_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Time dependent part(s) of physics for scalar variables.

Reimplemented in GRINS::AveragedTurbine< Viscosity >, and GRINS::ScalarODE.

Definition at line 229 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

232  {
233  return;
234  }
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesBase< Mu, SH, TC >::read_input_options ( const GetPot &  input)
virtualinherited

Read options from GetPot input file.

Reimplemented from GRINS::Physics.

Reimplemented in GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >.

Definition at line 63 of file low_mach_navier_stokes_base.C.

References GRINS::low_mach_navier_stokes, GRINS::p_var_name_default, GRINS::T_var_name_default, GRINS::u_var_name_default, GRINS::v_var_name_default, and GRINS::w_var_name_default.

Referenced by GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::LowMachNavierStokesBase().

64  {
65  // Read FE info
66  this->_V_FE_family =
67  libMesh::Utility::string_to_enum<GRINSEnums::FEFamily>( input("Physics/"+low_mach_navier_stokes+"/V_FE_family", "LAGRANGE") );
68 
69  this->_P_FE_family =
70  libMesh::Utility::string_to_enum<GRINSEnums::FEFamily>( input("Physics/"+low_mach_navier_stokes+"/P_FE_family", "LAGRANGE") );
71 
72  this->_T_FE_family =
73  libMesh::Utility::string_to_enum<GRINSEnums::FEFamily>( input("Physics/"+low_mach_navier_stokes+"/T_FE_family", "LAGRANGE") );
74 
75  this->_V_order =
76  libMesh::Utility::string_to_enum<GRINSEnums::Order>( input("Physics/"+low_mach_navier_stokes+"/V_order", "SECOND") );
77 
78  this->_P_order =
79  libMesh::Utility::string_to_enum<GRINSEnums::Order>( input("Physics/"+low_mach_navier_stokes+"/P_order", "FIRST") );
80 
81  this->_T_order =
82  libMesh::Utility::string_to_enum<GRINSEnums::Order>( input("Physics/"+low_mach_navier_stokes+"/T_order", "SECOND") );
83 
84  // Read variable naming info
85  this->_u_var_name = input("Physics/VariableNames/u_velocity", u_var_name_default );
86  this->_v_var_name = input("Physics/VariableNames/v_velocity", v_var_name_default );
87  this->_w_var_name = input("Physics/VariableNames/w_velocity", w_var_name_default );
88  this->_p_var_name = input("Physics/VariableNames/pressure", p_var_name_default );
89  this->_T_var_name = input("Physics/VariableNames/temperature", T_var_name_default );
90 
91  // Read thermodynamic state info
92  this->set_parameter
93  (_p0, input, "Physics/"+low_mach_navier_stokes+"/p0", 0.0 ); /* thermodynamic pressure */
94  this->set_parameter
95  (_T0, input, "Physics/"+low_mach_navier_stokes+"/T0", 0.0 ); /* Reference temperature */
96  this->set_parameter
97  (_R, input, "Physics/"+low_mach_navier_stokes+"/R", 0.0 ); /* gas constant */
98 
99  if( _R <= 0.0 )
100  {
101  std::cerr << "=========================================" << std::endl
102  << " Error: Gas constant R must be positive. " << std::endl
103  << " Detected value R = " << _R << std::endl
104  << "=========================================" << std::endl;
105  libmesh_error();
106  }
107 
108  _p0_over_R = _p0/_R;
109 
110  _enable_thermo_press_calc = input("Physics/"+low_mach_navier_stokes+"/enable_thermo_press_calc", false );
111 
113  {
114  _p0_var_name = input("Physics/VariableNames/thermo_presure", "p0" );
115  }
116 
117  // Read gravity vector
118  unsigned int g_dim = input.vector_variable_size("Physics/"+low_mach_navier_stokes+"/g");
119 
120  _g(0) = input("Physics/"+low_mach_navier_stokes+"/g", 0.0, 0 );
121  _g(1) = input("Physics/"+low_mach_navier_stokes+"/g", 0.0, 1 );
122 
123  if( g_dim == 3)
124  _g(2) = input("Physics/"+low_mach_navier_stokes+"/g", 0.0, 2 );
125 
126  return;
127  }
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.
GRINSEnums::FEFamily _V_FE_family
Element type, read from input.
libMesh::Number _p0_over_R
Thermodynamic pressure divided by gas constant.
const std::string p_var_name_default
pressure
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
const std::string v_var_name_default
y-velocity
const std::string T_var_name_default
temperature
libMesh::Point _g
Gravity vector.
std::string _u_var_name
Names of each (owned) variable in the system.
GRINSEnums::Order _V_order
Element orders, read from input.
const std::string u_var_name_default
Default physics variable names.
const PhysicsName low_mach_navier_stokes
const std::string w_var_name_default
z-velocity
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesBase< Mu, SH, TC >::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiPointer< libMesh::Number > &  param_pointer 
) const
virtualinherited

Each subclass will register its copy of an independent.

Reimplemented from GRINS::ParameterUser.

Definition at line 198 of file low_mach_navier_stokes_base.C.

References GRINS::ParameterUser::register_parameter().

201  {
202  ParameterUser::register_parameter(param_name, param_pointer);
203  _mu.register_parameter(param_name, param_pointer);
204  _cp.register_parameter(param_name, param_pointer);
205  _k.register_parameter(param_name, param_pointer);
206  }
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
SpecificHeat _cp
Specific heat object.
ThermalConductivity _k
Thermal conductivity object.
void GRINS::Physics::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtualinherited

Register name of postprocessed quantity with PostProcessedQuantities.

Each Physics class will need to cache an unsigned int corresponding to each postprocessed quantity. This will be used in computing the values and putting them in the CachedVariables object.

Reimplemented in GRINS::ParsedVelocitySource< Viscosity >, GRINS::VelocityPenalty< Viscosity >, GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::HeatTransfer< Conductivity >, GRINS::ElasticCable< StressStrainLaw >, GRINS::ElasticMembrane< StressStrainLaw >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 161 of file physics.C.

163  {
164  return;
165  }
template<class V , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesBase< V, SH, TC >::rho ( libMesh::Real  T,
libMesh::Real  p0 
) const
inlineinherited

Definition at line 144 of file low_mach_navier_stokes_base.h.

145  {
146  return p0/(this->_R*T);
147  }
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
void GRINS::Physics::set_is_steady ( bool  is_steady)
inherited

Sets whether this physics is to be solved with a steady solver or not.

Since the member variable is static, only needs to be called on a single physics.

Definition at line 97 of file physics.C.

References GRINS::Physics::_is_steady, and GRINS::Physics::is_steady().

98  {
100  return;
101  }
bool is_steady() const
Returns whether or not this physics is being solved with a steady solver.
Definition: physics.C:103
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
void GRINS::ParameterUser::set_parameter ( libMesh::Number &  param_variable,
const GetPot &  input,
const std::string &  param_name,
libMesh::Number  param_default 
)
virtualinherited

Each subclass can simultaneously read a parameter value from.

Definition at line 35 of file parameter_user.C.

References GRINS::ParameterUser::_my_name, and GRINS::ParameterUser::_my_parameters.

Referenced by GRINS::AveragedFanAdjointStabilization< Viscosity >::AveragedFanAdjointStabilization(), GRINS::AveragedTurbineAdjointStabilization< Viscosity >::AveragedTurbineAdjointStabilization(), GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::BoussinesqBuoyancyAdjointStabilization(), GRINS::BoussinesqBuoyancyBase::BoussinesqBuoyancyBase(), GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::BoussinesqBuoyancySPGSMStabilization(), GRINS::ConstantConductivity::ConstantConductivity(), GRINS::ConstantPrandtlConductivity::ConstantPrandtlConductivity(), GRINS::ConstantSourceFunction::ConstantSourceFunction(), GRINS::ConstantSourceTerm::ConstantSourceTerm(), GRINS::ConstantSpecificHeat::ConstantSpecificHeat(), GRINS::ConstantViscosity::ConstantViscosity(), GRINS::ElasticCable< StressStrainLaw >::ElasticCable(), GRINS::ElasticCableConstantGravity::ElasticCableConstantGravity(), GRINS::ElasticMembrane< StressStrainLaw >::ElasticMembrane(), GRINS::ElasticMembraneConstantPressure::ElasticMembraneConstantPressure(), GRINS::HeatConduction< Conductivity >::HeatConduction(), GRINS::HeatTransferBase< Conductivity >::HeatTransferBase(), GRINS::IncompressibleNavierStokesBase< Viscosity >::IncompressibleNavierStokesBase(), GRINS::AverageNusseltNumber::init(), GRINS::MooneyRivlin::MooneyRivlin(), GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >::ReactingLowMachNavierStokesBase(), GRINS::HookesLaw1D::read_input_options(), GRINS::HookesLaw::read_input_options(), GRINS::AxisymmetricBoussinesqBuoyancy::read_input_options(), and GRINS::VelocityDragAdjointStabilization< Viscosity >::VelocityDragAdjointStabilization().

39  {
40  param_variable = input(param_name, param_default);
41 
42  libmesh_assert_msg(!_my_parameters.count(param_name),
43  "ERROR: " << _my_name << " double-registered parameter " <<
44  param_name);
45 
46  _my_parameters[param_name] = &param_variable;
47  }
std::map< std::string, libMesh::Number * > _my_parameters
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesBase< Mu, SH, TC >::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtualinherited

Sets velocity variables to be time-evolving.

Reimplemented from GRINS::Physics.

Definition at line 155 of file low_mach_navier_stokes_base.C.

156  {
157  const unsigned int dim = system->get_mesh().mesh_dimension();
158 
159  system->time_evolving(_u_var);
160  system->time_evolving(_v_var);
161 
162  if (dim == 3)
163  system->time_evolving(_w_var);
164 
165  system->time_evolving(_T_var);
166  system->time_evolving(_p_var);
167 
169  system->time_evolving(_p0_var);
170 
171  return;
172  }
VariableIndex _u_var
Indices for each (owned) variable;.
bool _enable_thermo_press_calc
Flag to enable thermodynamic pressure calculation.
void GRINS::Physics::side_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for boundaries of elements on the domain boundary.

Reimplemented in GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 243 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

246  {
247  return;
248  }
void GRINS::Physics::side_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited
template<class V , class SH , class TC >
libMesh::Real GRINS::LowMachNavierStokesBase< V, SH, TC >::T ( const libMesh::Point &  p,
const AssemblyContext c 
) const
inlineinherited

Definition at line 137 of file low_mach_navier_stokes_base.h.

138  {
139  return c.point_value(_T_var,p);
140  }

Member Data Documentation

GRINS::BCHandlingBase* GRINS::Physics::_bc_handler
protectedinherited
template<class Viscosity , class SpecificHeat , class ThermalConductivity >
SpecificHeat GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_cp
protectedinherited

Specific heat object.

Definition at line 118 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
unsigned int GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_dim
protectedinherited

Physical dimension of problem.

Definition at line 95 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
bool GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_enable_thermo_press_calc
protectedinherited

Flag to enable thermodynamic pressure calculation.

Definition at line 127 of file low_mach_navier_stokes_base.h.

std::set<libMesh::subdomain_id_type> GRINS::Physics::_enabled_subdomains
protectedinherited

Subdomains on which the current Physics class is enabled.

Definition at line 261 of file physics.h.

Referenced by GRINS::Physics::enabled_on_elem(), and GRINS::Physics::read_input_options().

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
libMesh::Point GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_g
protectedinherited

Gravity vector.

Definition at line 124 of file low_mach_navier_stokes_base.h.

GRINS::ICHandlingBase* GRINS::Physics::_ic_handler
protectedinherited
bool GRINS::Physics::_is_axisymmetric
protectedinherited
bool GRINS::Physics::_is_steady = false
staticprotectedinherited

Caches whether or not the solver that's being used is steady or not.

This is need, for example, in flow stabilization as the tau terms change depending on whether the solver is steady or unsteady.

Definition at line 266 of file physics.h.

Referenced by GRINS::Physics::is_steady(), and GRINS::Physics::set_is_steady().

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
ThermalConductivity GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_k
protectedinherited

Thermal conductivity object.

Definition at line 121 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
Viscosity GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_mu
protectedinherited

Viscosity object.

Definition at line 115 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
libMesh::Number GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_p0
protectedinherited

Definition at line 92 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
libMesh::Number GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_p0_over_R
protectedinherited

Thermodynamic pressure divided by gas constant.

Definition at line 90 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
VariableIndex GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_p0_var
protectedinherited

Definition at line 103 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
std::string GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_p0_var_name
protectedinherited

Definition at line 106 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
GRINSEnums::FEFamily GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_P_FE_family
protectedinherited

Definition at line 109 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
GRINSEnums::Order GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_P_order
protectedinherited

Definition at line 112 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
VariableIndex GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_p_var
protectedinherited

Definition at line 101 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
std::string GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_p_var_name
protectedinherited

Definition at line 106 of file low_mach_navier_stokes_base.h.

const PhysicsName GRINS::Physics::_physics_name
protectedinherited

Name of the physics object. Used for reading physics specific inputs.

We use a reference because the physics names are const global objects in GRINS namespace

No, we use a copy, because otherwise as soon as the memory in std::set<std::string> requested_physics gets overwritten we get in trouble.

Definition at line 254 of file physics.h.

Referenced by GRINS::SourceTermBase::parse_var_info(), and GRINS::Physics::read_input_options().

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
libMesh::Number GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_R
protectedinherited

Definition at line 92 of file low_mach_navier_stokes_base.h.

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

Definition at line 69 of file low_mach_navier_stokes_stab_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
libMesh::Number GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_T0
protectedinherited

Definition at line 92 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
GRINSEnums::FEFamily GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_T_FE_family
protectedinherited

Definition at line 109 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
GRINSEnums::Order GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_T_order
protectedinherited

Definition at line 112 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
VariableIndex GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_T_var
protectedinherited

Definition at line 102 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
std::string GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_T_var_name
protectedinherited

Definition at line 106 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
VariableIndex GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_u_var
protectedinherited

Indices for each (owned) variable;.

Definition at line 98 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
std::string GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_u_var_name
protectedinherited

Names of each (owned) variable in the system.

Definition at line 106 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
GRINSEnums::FEFamily GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_V_FE_family
protectedinherited

Element type, read from input.

Definition at line 109 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
GRINSEnums::Order GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_V_order
protectedinherited

Element orders, read from input.

Definition at line 112 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
VariableIndex GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_v_var
protectedinherited

Definition at line 99 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
std::string GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_v_var_name
protectedinherited

Definition at line 106 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
VariableIndex GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_w_var
protectedinherited

Definition at line 100 of file low_mach_navier_stokes_base.h.

template<class Viscosity , class SpecificHeat , class ThermalConductivity >
std::string GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::_w_var_name
protectedinherited

Definition at line 106 of file low_mach_navier_stokes_base.h.


The documentation for this class was generated from the following files:

Generated on Mon Jun 22 2015 21:32:23 for GRINS-0.6.0 by  doxygen 1.8.9.1