GRINS-0.7.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, CachedValues &cache)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
- 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 init_variables (libMesh::FEMSystem *system)
 Initialize variables for this physics. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
libMesh::Real T (const libMesh::Point &p, const AssemblyContext &c) const
 
libMesh::Real rho (libMesh::Real T, libMesh::Real p0) const
 
libMesh::Real d_rho_dT (libMesh::Real T, libMesh::Real p0) const
 
libMesh::Real get_p0_steady (const AssemblyContext &c, unsigned int qp) const
 
libMesh::Real get_p0_steady_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 bool enabled_on_elem (const libMesh::Elem *elem)
 Find if current physics is active on supplied element. More...
 
void set_is_steady (bool is_steady)
 Sets whether this physics is to be solved with a steady solver or not. More...
 
bool is_steady () const
 Returns whether or not this physics is being solved with a steady solver. More...
 
virtual void auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. More...
 
virtual void side_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 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 (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_element_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_damping_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
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)
 

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
 
unsigned int _dim
 Physical dimension of problem. More...
 
VelocityFEVariables _flow_vars
 
PressureFEVariable _press_var
 
PrimitiveTempFEVariables _temp_vars
 
libMesh::UniquePtr< ThermoPressureFEVariable_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 283 of file low_mach_navier_stokes_spgsm_stab.C.

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

Definition at line 92 of file low_mach_navier_stokes_spgsm_stab.C.

94  {
95  // The number of local degrees of freedom in each variable.
96  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
97 
98  // Element Jacobian * quadrature weights for interior integration.
99  const std::vector<libMesh::Real> &JxW =
100  context.get_element_fe(this->_flow_vars.u())->get_JxW();
101 
102  // The pressure shape functions at interior quadrature points.
103  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
104  context.get_element_fe(this->_press_var.p())->get_dphi();
105 
106  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
107 
108  unsigned int n_qpoints = context.get_element_qrule().n_points();
109 
110  for (unsigned int qp=0; qp != n_qpoints; qp++)
111  {
112  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
113 
114  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
115  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
116 
117  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
118 
119  libMesh::Real mu = this->_mu(T);
120 
121  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
122 
123  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
124  context.interior_value( this->_flow_vars.v(), qp ) );
125  if( this->_dim == 3 )
126  U(2) = context.interior_value( this->_flow_vars.w(), qp );
127 
128  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
129 
130  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
131 
132  // Now a loop over the pressure degrees of freedom. This
133  // computes the contributions of the continuity equation.
134  for (unsigned int i=0; i != n_p_dofs; i++)
135  {
136  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*JxW[qp];
137  }
138 
139  }
140 
141  return;
142  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
VariableIndex p() const
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:277
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 414 of file low_mach_navier_stokes_spgsm_stab.C.

416  {
417  // The number of local degrees of freedom in each variable.
418  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
419 
420  // Element Jacobian * quadrature weights for interior integration.
421  const std::vector<libMesh::Real> &JxW =
422  context.get_element_fe(this->_temp_vars.T())->get_JxW();
423 
424  // The temperature shape functions gradients at interior quadrature points.
425  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
426  context.get_element_fe(this->_temp_vars.T())->get_dphi();
427 
428  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
429 
430  unsigned int n_qpoints = context.get_element_qrule().n_points();
431 
432  for (unsigned int qp=0; qp != n_qpoints; qp++)
433  {
434  libMesh::Number u, v;
435  u = context.fixed_interior_value(this->_flow_vars.u(), qp);
436  v = context.fixed_interior_value(this->_flow_vars.v(), qp);
437 
438  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
439 
440  libMesh::NumberVectorValue U(u,v);
441  if (this->_dim == 3)
442  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp); // w
443 
444  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
445  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
446 
447  libMesh::Real k = this->_k(T);
448  libMesh::Real cp = this->_cp(T);
449 
450  libMesh::Number rho_cp = rho*this->_cp(T);
451 
452  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
453 
454  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
455  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
456 
457  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
458 
459  libMesh::Real RE_t = this->compute_res_energy_transient( context, qp );
460 
461  for (unsigned int i=0; i != n_T_dofs; i++)
462  {
463  FT(i) -= rho_cp*tau_E*RE_t*U*T_gradphi[i][qp]*JxW[qp];
464  }
465 
466  }
467 
468  return;
469  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
unsigned int _dim
Physical dimension of problem.
SpecificHeat _cp
Specific heat object.
ThermalConductivity _k
Thermal conductivity object.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real compute_res_energy_transient(AssemblyContext &context, unsigned int qp) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
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 225 of file low_mach_navier_stokes_spgsm_stab.C.

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

335  {
336  // The number of local degrees of freedom in each variable.
337  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
338 
339  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
340  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
341  if (this->_dim == 3)
342  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
343 
344  // Element Jacobian * quadrature weights for interior integration.
345  const std::vector<libMesh::Real> &JxW =
346  context.get_element_fe(this->_flow_vars.u())->get_JxW();
347 
348  // The velocity shape function gradients at interior quadrature points.
349  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
350  context.get_element_fe(this->_flow_vars.u())->get_dphi();
351 
352  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
353  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
354  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
355 
356  if( this->_dim == 3 )
357  {
358  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
359  }
360 
361  unsigned int n_qpoints = context.get_element_qrule().n_points();
362  for (unsigned int qp=0; qp != n_qpoints; qp++)
363  {
364  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
365  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
366 
367  libMesh::Real mu = this->_mu(T);
368 
369  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
370  context.fixed_interior_value(this->_flow_vars.v(), qp) );
371 
372  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
373  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
374  libMesh::RealGradient grad_w;
375 
376  if( this->_dim == 3 )
377  {
378  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
379  grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
380  }
381 
382  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
383 
384  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
385  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
386 
387  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
388  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
389 
390  libMesh::Real RC_t = this->compute_res_continuity_transient( context, qp );
391  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
392  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
393 
394  for (unsigned int i=0; i != n_u_dofs; i++)
395  {
396  Fu(i) -= ( tau_C*RC_t*u_gradphi[i][qp](0)
397  + tau_M*RM_t(0)*rho*U*u_gradphi[i][qp] )*JxW[qp];
398 
399  Fv(i) -= ( tau_C*RC_t*u_gradphi[i][qp](1)
400  + tau_M*RM_t(1)*rho*U*u_gradphi[i][qp] )*JxW[qp];
401 
402  if( this->_dim == 3 )
403  {
404  (*Fw)(i) -= ( tau_C*RC_t*u_gradphi[i][qp](2)
405  + tau_M*RM_t(2)*rho*U*u_gradphi[i][qp] )*JxW[qp];
406  }
407  }
408 
409  }
410  return;
411  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
libMesh::Real compute_tau_continuity(libMesh::Real tau_C, libMesh::RealGradient &g) const
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient compute_res_momentum_transient(AssemblyContext &context, unsigned int qp) const
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real get_p0_transient(AssemblyContext &c, unsigned int qp) const
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
libMesh::RealGradient compute_res_momentum_steady(AssemblyContext &context, unsigned int qp) const
LowMachNavierStokesStabilizationHelper _stab_helper
libMesh::Real compute_res_continuity_transient(AssemblyContext &context, unsigned int qp) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0) const
template<class Mu , class SH , class TC >
void GRINS::LowMachNavierStokesSPGSMStabilization< Mu, SH, TC >::assemble_momentum_time_deriv ( bool  compute_jacobian,
AssemblyContext context 
)
protected

Definition at line 145 of file low_mach_navier_stokes_spgsm_stab.C.

147  {
148  // The number of local degrees of freedom in each variable.
149  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
150 
151  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
152  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
153  if (this->_dim == 3)
154  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
155 
156  // Element Jacobian * quadrature weights for interior integration.
157  const std::vector<libMesh::Real> &JxW =
158  context.get_element_fe(this->_flow_vars.u())->get_JxW();
159 
160  // The velocity shape function gradients at interior quadrature points.
161  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
162  context.get_element_fe(this->_flow_vars.u())->get_dphi();
163 
164  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
165  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
166  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
167 
168  if( this->_dim == 3 )
169  {
170  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
171  }
172 
173  unsigned int n_qpoints = context.get_element_qrule().n_points();
174 
175  for (unsigned int qp=0; qp != n_qpoints; qp++)
176  {
177  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
178  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
179 
180  libMesh::RealGradient U( context.interior_value(this->_flow_vars.u(), qp),
181  context.interior_value(this->_flow_vars.v(), qp) );
182 
183  libMesh::RealGradient grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
184  libMesh::RealGradient grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
185  libMesh::RealGradient grad_w;
186 
187  if( this->_dim == 3 )
188  {
189  U(2) = context.interior_value(this->_flow_vars.w(), qp);
190  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
191  }
192 
193  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
194 
195  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
196  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
197  libMesh::Real mu = this->_mu(T);
198 
199  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
200  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
201 
202  libMesh::Real RC_s = this->compute_res_continuity_steady( context, qp );
203  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
204 
205  for (unsigned int i=0; i != n_u_dofs; i++)
206  {
207  Fu(i) += ( - tau_C*RC_s*u_gradphi[i][qp](0)
208  - tau_M*RM_s(0)*rho*U*u_gradphi[i][qp] )*JxW[qp];
209 
210  Fv(i) += ( - tau_C*RC_s*u_gradphi[i][qp](1)
211  - tau_M*RM_s(1)*rho*U*u_gradphi[i][qp] )*JxW[qp];
212 
213  if( this->_dim == 3 )
214  {
215  (*Fw)(i) += ( - tau_C*RC_s*u_gradphi[i][qp](2)
216  - tau_M*RM_s(2)*rho*U*u_gradphi[i][qp] )*JxW[qp];
217  }
218  }
219 
220  }
221  return;
222  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
libMesh::Real compute_tau_continuity(libMesh::Real tau_C, libMesh::RealGradient &g) const
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int _dim
Physical dimension of problem.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:277
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  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 54 of file low_mach_navier_stokes_spgsm_stab.C.

57  {
58 #ifdef GRINS_USE_GRVY_TIMERS
59  this->_timer->BeginTimer("LowMachNavierStokesSPGSMStabilization::element_time_derivative");
60 #endif
61 
62  this->assemble_continuity_time_deriv( compute_jacobian, context );
63  this->assemble_momentum_time_deriv( compute_jacobian, context );
64  this->assemble_energy_time_deriv( compute_jacobian, context );
65 
66 #ifdef GRINS_USE_GRVY_TIMERS
67  this->_timer->EndTimer("LowMachNavierStokesSPGSMStabilization::element_time_derivative");
68 #endif
69  return;
70  }
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  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 73 of file low_mach_navier_stokes_spgsm_stab.C.

76  {
77 #ifdef GRINS_USE_GRVY_TIMERS
78  this->_timer->BeginTimer("LowMachNavierStokesSPGSMStabilization::mass_residual");
79 #endif
80 
81  this->assemble_continuity_mass_residual( compute_jacobian, context );
82  this->assemble_momentum_mass_residual( compute_jacobian, context );
83  this->assemble_energy_mass_residual( compute_jacobian, context );
84 
85 #ifdef GRINS_USE_GRVY_TIMERS
86  this->_timer->EndTimer("LowMachNavierStokesSPGSMStabilization::mass_residual");
87 #endif
88  return;
89  }
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 Thu Jun 2 2016 21:52:31 for GRINS-0.7.0 by  doxygen 1.8.10