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

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

#include <low_mach_navier_stokes_braack_stab.h>

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

Public Member Functions

 LowMachNavierStokesBraackStabilization (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~LowMachNavierStokesBraackStabilization ()
 
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

 LowMachNavierStokesBraackStabilization ()
 

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::LowMachNavierStokesBraackStabilization< Viscosity, SpecificHeat, ThermalConductivity >

Adds VMS-based stabilization to LowMachNavierStokes physics class.

Definition at line 36 of file low_mach_navier_stokes_braack_stab.h.

Constructor & Destructor Documentation

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

Definition at line 42 of file low_mach_navier_stokes_braack_stab.C.

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

Definition at line 50 of file low_mach_navier_stokes_braack_stab.C.

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

Member Function Documentation

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

Definition at line 299 of file low_mach_navier_stokes_braack_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->_press_var.p()).size();
304 
305  // Element Jacobian * quadrature weights for interior integration.
306  const std::vector<libMesh::Real> &JxW =
307  context.get_element_fe(this->_flow_vars.u())->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->_press_var.p())->get_dphi();
312 
313  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // 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->_flow_vars.u());
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->_temp_vars.T(), qp );
325  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
326 
327  libMesh::Real mu = this->_mu(T);
328  libMesh::Real k = this->_k(T);
329  libMesh::Real cp = this->_cp(T);
330 
331  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
332  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
333  if( this->_flow_vars.dim() == 3 )
334  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
335 
336  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
337  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
338 
339  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
340  libMesh::Real RE_t = this->compute_res_energy_transient( context, qp );
341 
342  // Now a loop over the pressure degrees of freedom. This
343  // computes the contributions of the continuity equation.
344  for (unsigned int i=0; i != n_p_dofs; i++)
345  {
346  Fp(i) += ( tau_M*RM_t*p_dphi[i][qp]
347  + tau_E*RE_t*(U*p_dphi[i][qp])/T
348  )*JxW[qp];
349  }
350  }
351 
352  return;
353  }
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
VariableIndex p() 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
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::LowMachNavierStokesBraackStabilization< Mu, SH, TC >::assemble_continuity_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 75 of file low_mach_navier_stokes_braack_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  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
102 
103  libMesh::Real mu = this->_mu(T);
104  libMesh::Real k = this->_k(T);
105  libMesh::Real cp = this->_cp(T);
106 
107  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
108  context.interior_value( this->_flow_vars.v(), qp ) );
109  if( this->_flow_vars.dim() == 3 )
110  U(2) = context.interior_value( this->_flow_vars.w(), qp ); // w
111 
112  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
113  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, this->_is_steady );
114 
115  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
116  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
117 
118  // Now a loop over the pressure degrees of freedom. This
119  // computes the contributions of the continuity equation.
120  for (unsigned int i=0; i != n_p_dofs; i++)
121  {
122  Fp(i) += ( tau_M*RM_s*p_dphi[i][qp]
123  + tau_E*RE_s*(U*p_dphi[i][qp])/T )*JxW[qp];
124  }
125 
126  }
127 
128  return;
129  }
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
VariableIndex p() 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::LowMachNavierStokesBraackStabilization< Mu, SH, TC >::assemble_energy_mass_residual ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 459 of file low_mach_navier_stokes_braack_stab.C.

461  {
462  // The number of local degrees of freedom in each variable.
463  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
464 
465  // Element Jacobian * quadrature weights for interior integration.
466  const std::vector<libMesh::Real> &JxW =
467  context.get_element_fe(this->_temp_vars.T())->get_JxW();
468 
469  // The temperature shape functions gradients at interior quadrature points.
470  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
471  context.get_element_fe(this->_temp_vars.T())->get_dphi();
472 
473  const std::vector<std::vector<libMesh::RealTensor> >& T_hessphi =
474  context.get_element_fe(this->_temp_vars.T())->get_d2phi();
475 
476  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
477 
478  unsigned int n_qpoints = context.get_element_qrule().n_points();
479 
480  for (unsigned int qp=0; qp != n_qpoints; qp++)
481  {
482  libMesh::Number u, v;
483  u = context.fixed_interior_value(this->_flow_vars.u(), qp);
484  v = context.fixed_interior_value(this->_flow_vars.v(), qp);
485 
486  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
487 
488  libMesh::NumberVectorValue U(u,v);
489  if (this->_flow_vars.dim() == 3)
490  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp); // w
491 
492  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
493  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
494 
495  libMesh::Real k = this->_k(T);
496  libMesh::Real cp = this->_cp(T);
497 
498  libMesh::Number rho_cp = rho*cp;
499 
500  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
501 
502  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
503  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
504 
505  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
506 
507  libMesh::Real RE_t = this->compute_res_energy_transient( context, qp );
508 
509  for (unsigned int i=0; i != n_T_dofs; i++)
510  {
511  FT(i) += ( rho_cp*tau_E*RE_t*U*T_gradphi[i][qp]
512  + tau_E*RE_t*k*(T_hessphi[i][qp](0,0) + T_hessphi[i][qp](1,1) + T_hessphi[i][qp](2,2))
513  )*JxW[qp];
514  }
515 
516  }
517 
518  return;
519  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
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_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.
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::LowMachNavierStokesBraackStabilization< Mu, SH, TC >::assemble_energy_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 236 of file low_mach_navier_stokes_braack_stab.C.

238  {
239  // The number of local degrees of freedom in each variable.
240  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
241 
242  // Element Jacobian * quadrature weights for interior integration.
243  const std::vector<libMesh::Real> &JxW =
244  context.get_element_fe(this->_temp_vars.T())->get_JxW();
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->_temp_vars.T())->get_dphi();
249 
250  const std::vector<std::vector<libMesh::RealTensor> >& T_hessphi =
251  context.get_element_fe(this->_temp_vars.T())->get_d2phi();
252 
253  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
254 
255  unsigned int n_qpoints = context.get_element_qrule().n_points();
256 
257  for (unsigned int qp=0; qp != n_qpoints; qp++)
258  {
259  libMesh::Number u, v;
260  u = context.interior_value(this->_flow_vars.u(), qp);
261  v = context.interior_value(this->_flow_vars.v(), qp);
262 
263  libMesh::Gradient grad_T = context.interior_gradient(this->_temp_vars.T(), qp);
264 
265  libMesh::NumberVectorValue U(u,v);
266  if (this->_flow_vars.dim() == 3)
267  U(2) = context.interior_value(this->_flow_vars.w(), qp);
268 
269  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
270  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
271 
272  libMesh::Real k = this->_k(T);
273  libMesh::Real cp = this->_cp(T);
274 
275  libMesh::Number rho_cp = rho*this->_cp(T);
276 
277  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
278 
279  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
280  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
281 
282  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, this->_is_steady );
283 
284  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
285 
286  for (unsigned int i=0; i != n_T_dofs; i++)
287  {
288  FT(i) += ( rho_cp*tau_E*RE_s*U*T_gradphi[i][qp]
289  + tau_E*RE_s*k*(T_hessphi[i][qp](0,0) + T_hessphi[i][qp](1,1) + T_hessphi[i][qp](2,2) )
290  )*JxW[qp];
291  }
292 
293  }
294 
295  return;
296  }
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.
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::LowMachNavierStokesBraackStabilization< Mu, SH, TC >::assemble_momentum_mass_residual ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 356 of file low_mach_navier_stokes_braack_stab.C.

358  {
359  // The number of local degrees of freedom in each variable.
360  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
361 
362  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
363  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
364  if (this->_flow_vars.dim() == 3)
365  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
366 
367  // Element Jacobian * quadrature weights for interior integration.
368  const std::vector<libMesh::Real> &JxW =
369  context.get_element_fe(this->_flow_vars.u())->get_JxW();
370 
371  // The velocity shape function gradients at interior quadrature points.
372  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
373  context.get_element_fe(this->_flow_vars.u())->get_dphi();
374 
375  const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
376  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
377 
378  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
379  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
380  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
381 
382  if( this->_flow_vars.dim() == 3 )
383  {
384  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
385  }
386 
387  unsigned int n_qpoints = context.get_element_qrule().n_points();
388  for (unsigned int qp=0; qp != n_qpoints; qp++)
389  {
390  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
391  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
392 
393  libMesh::Real mu = this->_mu(T);
394 
395  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
396  context.fixed_interior_value(this->_flow_vars.v(), qp) );
397 
398  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
399  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
400  libMesh::RealGradient grad_w;
401 
402  if( this->_flow_vars.dim() == 3 )
403  {
404  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
405  grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
406  }
407 
408  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
409 
410  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
411  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
412 
413  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
414  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
415 
416  libMesh::Real RC_t = this->compute_res_continuity_transient( context, qp );
417  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
418 
419  for (unsigned int i=0; i != n_u_dofs; i++)
420  {
421  Fu(i) += ( tau_C*RC_t*u_gradphi[i][qp](0)
422  + tau_M*RM_t(0)*rho*U*u_gradphi[i][qp]
423  + mu*tau_M*RM_t(0)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1)
424  + u_hessphi[i][qp](0,0) + u_hessphi[i][qp](0,1)
425  - 2.0/3.0*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,0)) )
426  )*JxW[qp];
427 
428  Fv(i) += ( tau_C*RC_t*u_gradphi[i][qp](1)
429  + tau_M*RM_t(1)*rho*U*u_gradphi[i][qp]
430  + mu*tau_M*RM_t(1)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1)
431  + u_hessphi[i][qp](1,0) + u_hessphi[i][qp](1,1)
432  - 2.0/3.0*(u_hessphi[i][qp](0,1) + u_hessphi[i][qp](1,1)) )
433  )*JxW[qp];
434 
435  if( this->_flow_vars.dim() == 3 )
436  {
437  Fu(i) -= mu*tau_M*RM_t(0)*(u_hessphi[i][qp](2,2) + u_hessphi[i][qp](0,2)
438  - 2.0/3.0*u_hessphi[i][qp](2,0))*JxW[qp];
439 
440  Fv(i) -= mu*tau_M*RM_t(1)*(u_hessphi[i][qp](2,2) + u_hessphi[i][qp](1,2)
441  - 2.0/3.0*u_hessphi[i][qp](2,1))*JxW[qp];
442 
443  (*Fw)(i) -= ( tau_C*RC_t*u_gradphi[i][qp](2)
444  + tau_M*RM_t(2)*rho*U*u_gradphi[i][qp]
445  + mu*tau_M*RM_t(2)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2)
446  + u_hessphi[i][qp](2,0) + u_hessphi[i][qp](2,1) + u_hessphi[i][qp](2,2)
447  - 2.0/3.0*(u_hessphi[i][qp](0,2) + u_hessphi[i][qp](1,2)
448  + u_hessphi[i][qp](2,2))
449  )
450  )*JxW[qp];
451  }
452  }
453 
454  }
455  return;
456  }
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
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::LowMachNavierStokesBraackStabilization< Mu, SH, TC >::assemble_momentum_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 132 of file low_mach_navier_stokes_braack_stab.C.

134  {
135  // The number of local degrees of freedom in each variable.
136  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
137 
138  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
139  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
140  if (this->_flow_vars.dim() == 3)
141  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
142 
143  // Element Jacobian * quadrature weights for interior integration.
144  const std::vector<libMesh::Real> &JxW =
145  context.get_element_fe(this->_flow_vars.u())->get_JxW();
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  const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
152  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
153 
154  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
155  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
156  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
157 
158  if( this->_flow_vars.dim() == 3 )
159  {
160  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
161  }
162 
163  unsigned int n_qpoints = context.get_element_qrule().n_points();
164 
165  for (unsigned int qp=0; qp != n_qpoints; qp++)
166  {
167  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
168  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
169 
170  libMesh::Real mu = this->_mu(T);
171 
172  libMesh::RealGradient U( context.interior_value(this->_flow_vars.u(), qp),
173  context.interior_value(this->_flow_vars.v(), qp) );
174 
175  libMesh::RealGradient grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
176  libMesh::RealGradient grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
177  libMesh::RealGradient grad_w;
178 
179  if( this->_flow_vars.dim() == 3 )
180  {
181  U(2) = context.interior_value(this->_flow_vars.w(), qp);
182  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
183  }
184 
185  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
186 
187  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
188  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
189 
190  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
191  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
192 
193  libMesh::Real RC_s = this->compute_res_continuity_steady( context, qp );
194  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
195 
196  for (unsigned int i=0; i != n_u_dofs; i++)
197  {
198  Fu(i) += ( tau_C*RC_s*u_gradphi[i][qp](0)
199  + tau_M*RM_s(0)*rho*U*u_gradphi[i][qp]
200  + mu*tau_M*RM_s(0)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1)
201  + u_hessphi[i][qp](0,0) + u_hessphi[i][qp](0,1)
202  - 2.0/3.0*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,0)) )
203  )*JxW[qp];
204 
205  Fv(i) += ( tau_C*RC_s*u_gradphi[i][qp](1)
206  + tau_M*RM_s(1)*rho*U*u_gradphi[i][qp]
207  + mu*tau_M*RM_s(1)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1)
208  + u_hessphi[i][qp](1,0) + u_hessphi[i][qp](1,1)
209  - 2.0/3.0*(u_hessphi[i][qp](0,1) + u_hessphi[i][qp](1,1)) )
210  )*JxW[qp];
211 
212  if( this->_flow_vars.dim() == 3 )
213  {
214  Fu(i) += mu*tau_M*RM_s(0)*(u_hessphi[i][qp](2,2) + u_hessphi[i][qp](0,2)
215  - 2.0/3.0*u_hessphi[i][qp](2,0))*JxW[qp];
216 
217  Fv(i) += mu*tau_M*RM_s(1)*(u_hessphi[i][qp](2,2) + u_hessphi[i][qp](1,2)
218  - 2.0/3.0*u_hessphi[i][qp](2,1))*JxW[qp];
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  + mu*tau_M*RM_s(2)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2)
223  + u_hessphi[i][qp](2,0) + u_hessphi[i][qp](2,1) + u_hessphi[i][qp](2,2)
224  - 2.0/3.0*(u_hessphi[i][qp](0,2) + u_hessphi[i][qp](1,2)
225  + u_hessphi[i][qp](2,2))
226  )
227  )*JxW[qp];
228  }
229  }
230 
231  }
232  return;
233  }
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::LowMachNavierStokesBraackStabilization< 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_braack_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_continuity_time_deriv(bool compute_jacobian, AssemblyContext &context)
void assemble_momentum_time_deriv(bool compute_jacobian, AssemblyContext &context)
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesBraackStabilization< 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_braack_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_energy_mass_residual(bool compute_jacobian, AssemblyContext &context)
void assemble_momentum_mass_residual(bool compute_jacobian, AssemblyContext &context)
void assemble_continuity_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