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

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

#include <low_mach_navier_stokes_spgsm_stab.h>

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

Public Member Functions

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

 LowMachNavierStokesSPGSMStabilization ()
 

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

Adds SPGSM-based stabilization to LowMachNavierStokes physics class.

Definition at line 35 of file low_mach_navier_stokes_spgsm_stab.h.

Constructor & Destructor Documentation

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

Definition at line 41 of file low_mach_navier_stokes_spgsm_stab.C.

43  : LowMachNavierStokesStabilizationBase<Mu,SH,TC>(physics_name,input)
44  {
45  }
template<class Mu , class SH , class TC >
GRINS::LowMachNavierStokesSPGSMStabilization< Mu, SH, TC >::~LowMachNavierStokesSPGSMStabilization ( )
virtual

Definition at line 48 of file low_mach_navier_stokes_spgsm_stab.C.

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

Member Function Documentation

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

Definition at line 264 of file low_mach_navier_stokes_spgsm_stab.C.

266  {
267  // The number of local degrees of freedom in each variable.
268  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
269 
270  // Element Jacobian * quadrature weights for interior integration.
271  const std::vector<libMesh::Real> &JxW =
272  context.get_element_fe(this->_flow_vars.u())->get_JxW();
273 
274  // The pressure shape functions at interior quadrature points.
275  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
276  context.get_element_fe(this->_press_var.p())->get_dphi();
277 
278  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
279 
280  unsigned int n_qpoints = context.get_element_qrule().n_points();
281 
282  for (unsigned int qp=0; qp != n_qpoints; qp++)
283  {
284  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
285 
286  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
287  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
288 
289  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
290  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
291 
292  libMesh::Real mu = this->_mu(T);
293 
294  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
295  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
296  if( this->_flow_vars.dim() == 3 )
297  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
298 
299  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
300  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
301 
302  // Now a loop over the pressure degrees of freedom. This
303  // computes the contributions of the continuity equation.
304  for (unsigned int i=0; i != n_p_dofs; i++)
305  {
306  Fp(i) += tau_M*RM_t*p_dphi[i][qp]*JxW[qp];
307  }
308  }
309 
310  return;
311  }
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::LowMachNavierStokesSPGSMStabilization< Mu, SH, TC >::assemble_continuity_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 73 of file low_mach_navier_stokes_spgsm_stab.C.

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

Definition at line 395 of file low_mach_navier_stokes_spgsm_stab.C.

397  {
398  // The number of local degrees of freedom in each variable.
399  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
400 
401  // Element Jacobian * quadrature weights for interior integration.
402  const std::vector<libMesh::Real> &JxW =
403  context.get_element_fe(this->_temp_vars.T())->get_JxW();
404 
405  // The temperature shape functions gradients at interior quadrature points.
406  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
407  context.get_element_fe(this->_temp_vars.T())->get_dphi();
408 
409  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
410 
411  unsigned int n_qpoints = context.get_element_qrule().n_points();
412 
413  for (unsigned int qp=0; qp != n_qpoints; qp++)
414  {
415  libMesh::Number u, v;
416  u = context.fixed_interior_value(this->_flow_vars.u(), qp);
417  v = context.fixed_interior_value(this->_flow_vars.v(), qp);
418 
419  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
420 
421  libMesh::NumberVectorValue U(u,v);
422  if (this->_flow_vars.dim() == 3)
423  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp); // w
424 
425  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
426  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
427 
428  libMesh::Real k = this->_k(T);
429  libMesh::Real cp = this->_cp(T);
430 
431  libMesh::Number rho_cp = rho*this->_cp(T);
432 
433  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
434 
435  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
436  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
437 
438  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
439 
440  libMesh::Real RE_t = this->compute_res_energy_transient( context, qp );
441 
442  for (unsigned int i=0; i != n_T_dofs; i++)
443  {
444  FT(i) -= rho_cp*tau_E*RE_t*U*T_gradphi[i][qp]*JxW[qp];
445  }
446 
447  }
448 
449  return;
450  }
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::LowMachNavierStokesSPGSMStabilization< Mu, SH, TC >::assemble_energy_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 206 of file low_mach_navier_stokes_spgsm_stab.C.

208  {
209  // The number of local degrees of freedom in each variable.
210  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
211 
212  // Element Jacobian * quadrature weights for interior integration.
213  const std::vector<libMesh::Real> &JxW =
214  context.get_element_fe(this->_temp_vars.T())->get_JxW();
215 
216  // The temperature shape functions gradients at interior quadrature points.
217  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
218  context.get_element_fe(this->_temp_vars.T())->get_dphi();
219 
220  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
221 
222  unsigned int n_qpoints = context.get_element_qrule().n_points();
223 
224  for (unsigned int qp=0; qp != n_qpoints; qp++)
225  {
226  libMesh::Number u, v;
227  u = context.interior_value(this->_flow_vars.u(), qp);
228  v = context.interior_value(this->_flow_vars.v(), qp);
229 
230  libMesh::Gradient grad_T = context.interior_gradient(this->_temp_vars.T(), qp);
231 
232  libMesh::NumberVectorValue U(u,v);
233  if (this->_flow_vars.dim() == 3)
234  U(2) = context.interior_value(this->_flow_vars.w(), qp); // w
235 
236  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
237  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
238 
239  libMesh::Real k = this->_k(T);
240  libMesh::Real cp = this->_cp(T);
241 
242  libMesh::Number rho_cp = rho*this->_cp(T);
243 
244  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
245 
246  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
247  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
248 
249  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, this->_is_steady );
250 
251  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
252 
253  for (unsigned int i=0; i != n_T_dofs; i++)
254  {
255  FT(i) -= rho_cp*tau_E*RE_s*U*T_gradphi[i][qp]*JxW[qp];
256  }
257 
258  }
259 
260  return;
261  }
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::LowMachNavierStokesSPGSMStabilization< Mu, SH, TC >::assemble_momentum_mass_residual ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 314 of file low_mach_navier_stokes_spgsm_stab.C.

316  {
317  // The number of local degrees of freedom in each variable.
318  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
319 
320  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
321  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
322  if (this->_flow_vars.dim() == 3)
323  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
324 
325  // Element Jacobian * quadrature weights for interior integration.
326  const std::vector<libMesh::Real> &JxW =
327  context.get_element_fe(this->_flow_vars.u())->get_JxW();
328 
329  // The velocity shape function gradients at interior quadrature points.
330  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
331  context.get_element_fe(this->_flow_vars.u())->get_dphi();
332 
333  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
334  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
335  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
336 
337  if( this->_flow_vars.dim() == 3 )
338  {
339  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
340  }
341 
342  unsigned int n_qpoints = context.get_element_qrule().n_points();
343  for (unsigned int qp=0; qp != n_qpoints; qp++)
344  {
345  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
346  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
347 
348  libMesh::Real mu = this->_mu(T);
349 
350  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
351  context.fixed_interior_value(this->_flow_vars.v(), qp) );
352 
353  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
354  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
355  libMesh::RealGradient grad_w;
356 
357  if( this->_flow_vars.dim() == 3 )
358  {
359  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
360  grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
361  }
362 
363  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
364 
365  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
366  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
367 
368  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
369  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
370 
371  libMesh::Real RC_t = this->compute_res_continuity_transient( context, qp );
372  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
373  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
374 
375  for (unsigned int i=0; i != n_u_dofs; i++)
376  {
377  Fu(i) -= ( tau_C*RC_t*u_gradphi[i][qp](0)
378  + tau_M*RM_t(0)*rho*U*u_gradphi[i][qp] )*JxW[qp];
379 
380  Fv(i) -= ( tau_C*RC_t*u_gradphi[i][qp](1)
381  + tau_M*RM_t(1)*rho*U*u_gradphi[i][qp] )*JxW[qp];
382 
383  if( this->_flow_vars.dim() == 3 )
384  {
385  (*Fw)(i) -= ( tau_C*RC_t*u_gradphi[i][qp](2)
386  + tau_M*RM_t(2)*rho*U*u_gradphi[i][qp] )*JxW[qp];
387  }
388  }
389 
390  }
391  return;
392  }
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::LowMachNavierStokesSPGSMStabilization< Mu, SH, TC >::assemble_momentum_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 126 of file low_mach_navier_stokes_spgsm_stab.C.

128  {
129  // The number of local degrees of freedom in each variable.
130  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
131 
132  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
133  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
134  if (this->_flow_vars.dim() == 3)
135  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
136 
137  // Element Jacobian * quadrature weights for interior integration.
138  const std::vector<libMesh::Real> &JxW =
139  context.get_element_fe(this->_flow_vars.u())->get_JxW();
140 
141  // The velocity shape function gradients at interior quadrature points.
142  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
143  context.get_element_fe(this->_flow_vars.u())->get_dphi();
144 
145  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
146  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
147  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
148 
149  if( this->_flow_vars.dim() == 3 )
150  {
151  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
152  }
153 
154  unsigned int n_qpoints = context.get_element_qrule().n_points();
155 
156  for (unsigned int qp=0; qp != n_qpoints; qp++)
157  {
158  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
159  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
160 
161  libMesh::RealGradient U( context.interior_value(this->_flow_vars.u(), qp),
162  context.interior_value(this->_flow_vars.v(), qp) );
163 
164  libMesh::RealGradient grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
165  libMesh::RealGradient grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
166  libMesh::RealGradient grad_w;
167 
168  if( this->_flow_vars.dim() == 3 )
169  {
170  U(2) = context.interior_value(this->_flow_vars.w(), qp);
171  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
172  }
173 
174  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
175 
176  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
177  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
178  libMesh::Real mu = this->_mu(T);
179 
180  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
181  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
182 
183  libMesh::Real RC_s = this->compute_res_continuity_steady( context, qp );
184  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
185 
186  for (unsigned int i=0; i != n_u_dofs; i++)
187  {
188  Fu(i) += ( - tau_C*RC_s*u_gradphi[i][qp](0)
189  - tau_M*RM_s(0)*rho*U*u_gradphi[i][qp] )*JxW[qp];
190 
191  Fv(i) += ( - tau_C*RC_s*u_gradphi[i][qp](1)
192  - tau_M*RM_s(1)*rho*U*u_gradphi[i][qp] )*JxW[qp];
193 
194  if( this->_flow_vars.dim() == 3 )
195  {
196  (*Fw)(i) += ( - tau_C*RC_s*u_gradphi[i][qp](2)
197  - tau_M*RM_s(2)*rho*U*u_gradphi[i][qp] )*JxW[qp];
198  }
199  }
200 
201  }
202  return;
203  }
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::LowMachNavierStokesSPGSMStabilization< 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 55 of file low_mach_navier_stokes_spgsm_stab.C.

57  {
58  this->assemble_continuity_time_deriv( compute_jacobian, context );
59  this->assemble_momentum_time_deriv( compute_jacobian, context );
60  this->assemble_energy_time_deriv( compute_jacobian, context );
61  }
void assemble_momentum_time_deriv(bool compute_jacobian, AssemblyContext &context)
void assemble_continuity_time_deriv(bool compute_jacobian, AssemblyContext &context)
void assemble_energy_time_deriv(bool compute_jacobian, AssemblyContext &context)
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesSPGSMStabilization< 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 65 of file low_mach_navier_stokes_spgsm_stab.C.

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

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