GRINS-0.8.0
List of all members | Public Member Functions | Protected Member Functions | Private Member Functions
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)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
- Public Member Functions inherited from GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity >
 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 side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

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 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...
 

Private Member Functions

 LowMachNavierStokesVMSStabilization ()
 

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 Attributes inherited from GRINS::LowMachNavierStokesStabilizationBase< Viscosity, SpecificHeat, ThermalConductivity >
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...
 
- 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::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 285 of file low_mach_navier_stokes_vms_stab.C.

287  {
288  // The number of local degrees of freedom in each variable.
289  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
290 
291  // Element Jacobian * quadrature weights for interior integration.
292  const std::vector<libMesh::Real> &JxW =
293  context.get_element_fe(this->_flow_vars.u())->get_JxW();
294 
295  // The pressure shape functions at interior quadrature points.
296  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
297  context.get_element_fe(this->_press_var.p())->get_dphi();
298 
299  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
300 
301  unsigned int n_qpoints = context.get_element_qrule().n_points();
302 
303  for (unsigned int qp=0; qp != n_qpoints; qp++)
304  {
305  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
306 
307  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
308  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
309 
310  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
311  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
312 
313  libMesh::Real mu = this->_mu(T);
314 
315  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
316  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
317  if( this->_flow_vars.dim() == 3 )
318  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
319 
320  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
321  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
322 
323  // Now a loop over the pressure degrees of freedom. This
324  // computes the contributions of the continuity equation.
325  for (unsigned int i=0; i != n_p_dofs; i++)
326  {
327  Fp(i) += tau_M*RM_t*p_dphi[i][qp]*JxW[qp];
328  }
329  }
330 
331  return;
332  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
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
VariableIndex p() const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
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 75 of file low_mach_navier_stokes_vms_stab.C.

77  {
78  // The number of local degrees of freedom in each variable.
79  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
80 
81  // Element Jacobian * quadrature weights for interior integration.
82  const std::vector<libMesh::Real> &JxW =
83  context.get_element_fe(this->_flow_vars.u())->get_JxW();
84 
85  // The pressure shape functions at interior quadrature points.
86  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
87  context.get_element_fe(this->_press_var.p())->get_dphi();
88 
89  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
90 
91  unsigned int n_qpoints = context.get_element_qrule().n_points();
92 
93  for (unsigned int qp=0; qp != n_qpoints; qp++)
94  {
95  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
96 
97  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
98  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
99 
100  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
101 
102  libMesh::Real mu = this->_mu(T);
103 
104  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
105 
106  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
107  context.interior_value( this->_flow_vars.v(), qp ) );
108  if( this->_flow_vars.dim() == 3 )
109  U(2) = context.interior_value( this->_flow_vars.w(), qp );
110 
111  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
112 
113  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
114 
115  // Now a loop over the pressure degrees of freedom. This
116  // computes the contributions of the continuity equation.
117  for (unsigned int i=0; i != n_p_dofs; i++)
118  {
119  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*JxW[qp];
120  }
121 
122  }
123 
124  return;
125  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::Real get_p0_steady(const AssemblyContext &c, 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
VariableIndex p() const
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
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 429 of file low_mach_navier_stokes_vms_stab.C.

431  {
432  // The number of local degrees of freedom in each variable.
433  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
434 
435  // Element Jacobian * quadrature weights for interior integration.
436  const std::vector<libMesh::Real> &JxW =
437  context.get_element_fe(this->_temp_vars.T())->get_JxW();
438 
439  // The temperature shape functions at interior quadrature points.
440  const std::vector<std::vector<libMesh::Real> >& T_phi =
441  context.get_element_fe(this->_temp_vars.T())->get_phi();
442 
443  // The temperature shape functions gradients at interior quadrature points.
444  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
445  context.get_element_fe(this->_temp_vars.T())->get_dphi();
446 
447  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
448 
449  unsigned int n_qpoints = context.get_element_qrule().n_points();
450 
451  for (unsigned int qp=0; qp != n_qpoints; qp++)
452  {
453  libMesh::Number u, v;
454  u = context.fixed_interior_value(this->_flow_vars.u(), qp);
455  v = context.fixed_interior_value(this->_flow_vars.v(), qp);
456 
457  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
458 
459  libMesh::NumberVectorValue U(u,v);
460  if (this->_flow_vars.dim() == 3)
461  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp); // w
462 
463  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
464  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
465 
466  libMesh::Real mu = this->_mu(T);
467  libMesh::Real k = this->_k(T);
468  libMesh::Real cp = this->_cp(T);
469 
470  libMesh::Number rho_cp = rho*this->_cp(T);
471 
472  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
473 
474  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
475  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
476 
477  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
478  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
479 
480  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
481  libMesh::Real RE_t = this->compute_res_energy_transient( context, qp );
482 
483  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
484  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
485 
486  for (unsigned int i=0; i != n_T_dofs; i++)
487  {
488  FT(i) -= ( -rho_cp*tau_M*RM_t*grad_T*T_phi[i][qp]
489  +rho_cp*tau_E*RE_t*U*T_gradphi[i][qp]
490  - rho_cp*tau_E*(RE_s+RE_t)*tau_M*RM_t*T_gradphi[i][qp]
491  - rho_cp*tau_E*RE_t*tau_M*RM_s*T_gradphi[i][qp] )*JxW[qp];
492  }
493 
494  }
495 
496  return;
497  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
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
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
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 218 of file low_mach_navier_stokes_vms_stab.C.

220  {
221  // The number of local degrees of freedom in each variable.
222  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
223 
224  // Element Jacobian * quadrature weights for interior integration.
225  const std::vector<libMesh::Real> &JxW =
226  context.get_element_fe(this->_temp_vars.T())->get_JxW();
227 
228  // The temperature shape functions at interior quadrature points.
229  const std::vector<std::vector<libMesh::Real> >& T_phi =
230  context.get_element_fe(this->_temp_vars.T())->get_phi();
231 
232  // The temperature shape functions gradients at interior quadrature points.
233  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
234  context.get_element_fe(this->_temp_vars.T())->get_dphi();
235 
236  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
237 
238  unsigned int n_qpoints = context.get_element_qrule().n_points();
239 
240  for (unsigned int qp=0; qp != n_qpoints; qp++)
241  {
242  libMesh::Number u, v;
243  u = context.interior_value(this->_flow_vars.u(), qp);
244  v = context.interior_value(this->_flow_vars.v(), qp);
245 
246  libMesh::Gradient grad_T = context.interior_gradient(this->_temp_vars.T(), qp);
247 
248  libMesh::NumberVectorValue U(u,v);
249  if (this->_flow_vars.dim() == 3)
250  U(2) = context.interior_value(this->_flow_vars.w(), qp); // w
251 
252  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
253  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
254 
255  libMesh::Real mu = this->_mu(T);
256  libMesh::Real k = this->_k(T);
257  libMesh::Real cp = this->_cp(T);
258 
259  libMesh::Number rho_cp = rho*this->_cp(T);
260 
261  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
262 
263  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
264  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
265 
266  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
267  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, this->_is_steady );
268 
269  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
270  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
271 
272  for (unsigned int i=0; i != n_T_dofs; i++)
273  {
274  FT(i) += ( rho_cp*tau_M*RM_s*grad_T*T_phi[i][qp]
275  - rho_cp*tau_E*RE_s*U*T_gradphi[i][qp]
276  + rho_cp*tau_E*RE_s*tau_M*RM_s*T_gradphi[i][qp] )*JxW[qp];
277  }
278 
279  }
280 
281  return;
282  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::Real get_p0_steady(const AssemblyContext &c, 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
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
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 335 of file low_mach_navier_stokes_vms_stab.C.

337  {
338  // The number of local degrees of freedom in each variable.
339  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
340 
341  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
342  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
343  if (this->_flow_vars.dim() == 3)
344  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
345 
346  // Element Jacobian * quadrature weights for interior integration.
347  const std::vector<libMesh::Real> &JxW =
348  context.get_element_fe(this->_flow_vars.u())->get_JxW();
349 
350  // The pressure shape functions at interior quadrature points.
351  const std::vector<std::vector<libMesh::Real> >& u_phi =
352  context.get_element_fe(this->_flow_vars.u())->get_phi();
353 
354  // The velocity shape function gradients at interior quadrature points.
355  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
356  context.get_element_fe(this->_flow_vars.u())->get_dphi();
357 
358  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
359  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
360  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
361 
362  if( this->_flow_vars.dim() == 3 )
363  {
364  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
365  }
366 
367  unsigned int n_qpoints = context.get_element_qrule().n_points();
368  for (unsigned int qp=0; qp != n_qpoints; qp++)
369  {
370  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
371  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
372 
373  libMesh::Real mu = this->_mu(T);
374 
375  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
376  context.fixed_interior_value(this->_flow_vars.v(), qp) );
377 
378  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
379  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
380  libMesh::RealGradient grad_w;
381 
382  if( this->_flow_vars.dim() == 3 )
383  {
384  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
385  grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
386  }
387 
388  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
389 
390  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
391  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
392 
393  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
394  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
395 
396  libMesh::Real RC_t = this->compute_res_continuity_transient( context, qp );
397  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
398  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
399 
400  for (unsigned int i=0; i != n_u_dofs; i++)
401  {
402  Fu(i) -= ( tau_C*RC_t*u_gradphi[i][qp](0)
403  + tau_M*RM_t(0)*rho*U*u_gradphi[i][qp]
404  - rho*tau_M*RM_t*grad_u*u_phi[i][qp]
405  - tau_M*(RM_s(0)+RM_t(0))*rho*tau_M*RM_t*u_gradphi[i][qp]
406  - tau_M*RM_t(0)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
407 
408  Fv(i) -= ( tau_C*RC_t*u_gradphi[i][qp](1)
409  - rho*tau_M*RM_t*grad_v*u_phi[i][qp]
410  + tau_M*RM_t(1)*rho*U*u_gradphi[i][qp]
411  - tau_M*(RM_s(1)+RM_t(1))*rho*tau_M*RM_t*u_gradphi[i][qp]
412  - tau_M*RM_t(1)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
413 
414  if( this->_flow_vars.dim() == 3 )
415  {
416  (*Fw)(i) -= ( tau_C*RC_t*u_gradphi[i][qp](2)
417  - rho*tau_M*RM_t*grad_w*u_phi[i][qp]
418  + tau_M*RM_t(2)*rho*U*u_gradphi[i][qp]
419  - tau_M*(RM_s(2)+RM_t(2))*rho*tau_M*RM_t*u_gradphi[i][qp]
420  - tau_M*RM_t(2)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
421  }
422  }
423 
424  }
425  return;
426  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::Real compute_tau_continuity(libMesh::Real tau_C, libMesh::RealGradient &g) const
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
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
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 128 of file low_mach_navier_stokes_vms_stab.C.

130  {
131  // The number of local degrees of freedom in each variable.
132  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
133 
134  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
135  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
136  if (this->_flow_vars.dim() == 3)
137  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
138 
139  // Element Jacobian * quadrature weights for interior integration.
140  const std::vector<libMesh::Real> &JxW =
141  context.get_element_fe(this->_flow_vars.u())->get_JxW();
142 
143  // The pressure shape functions at interior quadrature points.
144  const std::vector<std::vector<libMesh::Real> >& u_phi =
145  context.get_element_fe(this->_flow_vars.u())->get_phi();
146 
147  // The velocity shape function gradients at interior quadrature points.
148  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
149  context.get_element_fe(this->_flow_vars.u())->get_dphi();
150 
151  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
152  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
153  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
154 
155  if( this->_flow_vars.dim() == 3 )
156  {
157  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
158  }
159 
160  unsigned int n_qpoints = context.get_element_qrule().n_points();
161 
162  for (unsigned int qp=0; qp != n_qpoints; qp++)
163  {
164  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
165  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
166 
167  libMesh::RealGradient U( context.interior_value(this->_flow_vars.u(), qp),
168  context.interior_value(this->_flow_vars.v(), qp) );
169 
170  libMesh::RealGradient grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
171  libMesh::RealGradient grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
172  libMesh::RealGradient grad_w;
173 
174  if( this->_flow_vars.dim() == 3 )
175  {
176  U(2) = context.interior_value(this->_flow_vars.w(), qp);
177  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
178  }
179 
180  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
181 
182  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
183  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
184  libMesh::Real mu = this->_mu(T);
185 
186  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
187  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
188 
189  libMesh::Real RC_s = this->compute_res_continuity_steady( context, qp );
190  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
191 
192  for (unsigned int i=0; i != n_u_dofs; i++)
193  {
194  Fu(i) += ( -tau_C*RC_s*u_gradphi[i][qp](0)
195  - tau_M*RM_s(0)*rho*U*u_gradphi[i][qp]
196  + rho*tau_M*RM_s*grad_u*u_phi[i][qp]
197  + tau_M*RM_s(0)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
198 
199  Fv(i) += ( -tau_C*RC_s*u_gradphi[i][qp](1)
200  - tau_M*RM_s(1)*rho*U*u_gradphi[i][qp]
201  + rho*tau_M*RM_s*grad_v*u_phi[i][qp]
202  + tau_M*RM_s(1)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
203 
204  if( this->_flow_vars.dim() == 3 )
205  {
206  (*Fw)(i) += ( -tau_C*RC_s*u_gradphi[i][qp](2)
207  - tau_M*RM_s(2)*rho*U*u_gradphi[i][qp]
208  + rho*tau_M*RM_s*grad_w*u_phi[i][qp]
209  + tau_M*RM_s(2)*rho*tau_M*RM_s*u_gradphi[i][qp] )*JxW[qp];
210  }
211  }
212 
213  }
214  return;
215  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::Real compute_tau_continuity(libMesh::Real tau_C, libMesh::RealGradient &g) const
libMesh::Real get_p0_steady(const AssemblyContext &c, 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
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
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
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 57 of file low_mach_navier_stokes_vms_stab.C.

59  {
60  this->assemble_continuity_time_deriv( compute_jacobian, context );
61  this->assemble_momentum_time_deriv( compute_jacobian, context );
62  this->assemble_energy_time_deriv( compute_jacobian, context );
63  }
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)
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesVMSStabilization< Mu, SH, TC >::mass_residual ( bool  ,
AssemblyContext  
)
virtual

Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part.

Reimplemented from GRINS::Physics.

Definition at line 67 of file low_mach_navier_stokes_vms_stab.C.

68  {
69  this->assemble_continuity_mass_residual( compute_jacobian, context );
70  this->assemble_momentum_mass_residual( compute_jacobian, context );
71  this->assemble_energy_mass_residual( compute_jacobian, context );
72  }
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)

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