GRINS-0.8.0
List of all members | Public Member Functions | Private Member Functions
GRINS::IncompressibleNavierStokesSPGSMStabilization< Viscosity > Class Template Reference

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

#include <inc_navier_stokes_spgsm_stab.h>

Inheritance diagram for GRINS::IncompressibleNavierStokesSPGSMStabilization< Viscosity >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::IncompressibleNavierStokesSPGSMStabilization< Viscosity >:
Collaboration graph
[legend]

Public Member Functions

 IncompressibleNavierStokesSPGSMStabilization (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~IncompressibleNavierStokesSPGSMStabilization ()
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context)
 Constraint 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::IncompressibleNavierStokesStabilizationBase< Viscosity >
 IncompressibleNavierStokesStabilizationBase (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~IncompressibleNavierStokesStabilizationBase ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
- Public Member Functions inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
 IncompressibleNavierStokesBase (const std::string &my_physics_name, const std::string &core_physics_name, const GetPot &input)
 
 ~IncompressibleNavierStokesBase ()
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
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...
 
libMesh::Real get_viscosity_value (AssemblyContext &context, unsigned int qp) const
 
- 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 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...
 

Private Member Functions

 IncompressibleNavierStokesSPGSMStabilization ()
 

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 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...
 
- Protected Attributes inherited from GRINS::IncompressibleNavierStokesStabilizationBase< Viscosity >
IncompressibleNavierStokesStabilizationHelper _stab_helper
 
- Protected Attributes inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
VelocityVariable_flow_vars
 
PressureFEVariable_press_var
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
Viscosity _mu
 Viscosity object. 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 GRINS::IncompressibleNavierStokesSPGSMStabilization< Viscosity >

Adds VMS-based stabilization to LowMachNavierStokes physics class.

Definition at line 35 of file inc_navier_stokes_spgsm_stab.h.

Constructor & Destructor Documentation

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

Definition at line 39 of file inc_navier_stokes_spgsm_stab.C.

41  : IncompressibleNavierStokesStabilizationBase<Mu>(physics_name,input)
42  {}
template<class Viscosity >
virtual GRINS::IncompressibleNavierStokesSPGSMStabilization< Viscosity >::~IncompressibleNavierStokesSPGSMStabilization ( )
inlinevirtual

Definition at line 41 of file inc_navier_stokes_spgsm_stab.h.

41 {};

Member Function Documentation

template<class Mu >
void GRINS::IncompressibleNavierStokesSPGSMStabilization< Mu >::element_constraint ( bool  ,
AssemblyContext  
)
virtual

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 120 of file inc_navier_stokes_spgsm_stab.C.

121  {
122  // The number of local degrees of freedom in each variable.
123  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
124 
125  // Element Jacobian * quadrature weights for interior integration.
126  const std::vector<libMesh::Real> &JxW =
127  context.get_element_fe(this->_flow_vars.u())->get_JxW();
128 
129  // The pressure shape functions at interior quadrature points.
130  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
131  context.get_element_fe(this->_press_var.p())->get_dphi();
132 
133  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
134 
135  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
136 
137  unsigned int n_qpoints = context.get_element_qrule().n_points();
138 
139  for (unsigned int qp=0; qp != n_qpoints; qp++)
140  {
141  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
142  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
143 
144  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
145  context.interior_value( this->_flow_vars.v(), qp ) );
146  if( this->_flow_vars.dim() == 3 )
147  U(2) = context.interior_value( this->_flow_vars.w(), qp );
148 
149  // Compute the viscosity at this qp
150  libMesh::Real _mu_qp = this->_mu(context, qp);
151 
152  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
153 
154  libMesh::RealGradient RM_s = this->_stab_helper.compute_res_momentum_steady( context, qp, this->_rho, _mu_qp );
155 
156  for (unsigned int i=0; i != n_p_dofs; i++)
157  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*JxW[qp];
158 
159  if( compute_jacobian )
160  libmesh_not_implemented();
161  }
162  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::RealGradient compute_res_momentum_steady(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real mu) const
libMesh::Number _rho
Material parameters, read from input.
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
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
template<class Mu >
void GRINS::IncompressibleNavierStokesSPGSMStabilization< Mu >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 46 of file inc_navier_stokes_spgsm_stab.C.

48  {
49  // The number of local degrees of freedom in each variable.
50  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
51 
52  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
53  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
54 
55 
56  if (this->_flow_vars.dim() == 3)
57  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
58 
59  // Element Jacobian * quadrature weights for interior integration.
60  const std::vector<libMesh::Real> &JxW =
61  context.get_element_fe(this->_flow_vars.u())->get_JxW();
62 
63  // The velocity shape function gradients at interior quadrature points.
64  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
65  context.get_element_fe(this->_flow_vars.u())->get_dphi();
66 
67  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
68  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
69  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
70  if(this->_flow_vars.dim() == 3)
71  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
72 
73  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
74 
75  unsigned int n_qpoints = context.get_element_qrule().n_points();
76 
77  for (unsigned int qp=0; qp != n_qpoints; qp++)
78  {
79  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
80  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
81 
82  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
83  context.interior_value( this->_flow_vars.v(), qp ) );
84  if( this->_flow_vars.dim() == 3 )
85  {
86  U(2) = context.interior_value( this->_flow_vars.w(), qp );
87  }
88 
89  // Compute the viscosity at this qp
90  libMesh::Real _mu_qp = this->_mu(context, qp);
91 
92  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
93  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
94 
95  libMesh::RealGradient RM_s = this->_stab_helper.compute_res_momentum_steady( context, qp, this->_rho, _mu_qp );
96  libMesh::Real RC = this->_stab_helper.compute_res_continuity( context, qp );
97 
98  for (unsigned int i=0; i != n_u_dofs; i++)
99  {
100  Fu(i) += ( - tau_C*RC*u_gradphi[i][qp](0)
101  - tau_M*RM_s(0)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
102 
103  Fv(i) += ( - tau_C*RC*u_gradphi[i][qp](1)
104  - tau_M*RM_s(1)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
105 
106  if( this->_flow_vars.dim() == 3 )
107  {
108  (*Fw)(i) += ( - tau_C*RC*u_gradphi[i][qp](2)
109  - tau_M*RM_s(2)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
110  }
111  }
112 
113  if( compute_jacobian )
114  libmesh_not_implemented();
115  }
116  }
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_steady(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real mu) const
libMesh::Number _rho
Material parameters, read from input.
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
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::Real compute_res_continuity(AssemblyContext &context, unsigned int qp) const
template<class Mu >
void GRINS::IncompressibleNavierStokesSPGSMStabilization< Mu >::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 166 of file inc_navier_stokes_spgsm_stab.C.

167  {
168  // The number of local degrees of freedom in each variable.
169  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
170  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
171 
172  // Element Jacobian * quadrature weights for interior integration.
173  const std::vector<libMesh::Real> &JxW =
174  context.get_element_fe(this->_flow_vars.u())->get_JxW();
175 
176  // The pressure shape functions at interior quadrature points.
177  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
178  context.get_element_fe(this->_press_var.p())->get_dphi();
179 
180  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
181  context.get_element_fe(this->_flow_vars.u())->get_dphi();
182 
183  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{p}
184  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{p}
185  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
186 
187 
188  if(this->_flow_vars.dim() == 3)
189  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
190 
191  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
192 
193  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
194 
195  unsigned int n_qpoints = context.get_element_qrule().n_points();
196 
197  for (unsigned int qp=0; qp != n_qpoints; qp++)
198  {
199  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
200  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
201 
202  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
203  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
204  // Compute the viscosity at this qp
205  libMesh::Real _mu_qp = this->_mu(context, qp);
206 
207  if( this->_flow_vars.dim() == 3 )
208  {
209  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
210  }
211 
212  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, false );
213 
214  libMesh::RealGradient RM_t = this->_stab_helper.compute_res_momentum_transient( context, qp, this->_rho );
215 
216  for (unsigned int i=0; i != n_p_dofs; i++)
217  {
218  Fp(i) -= tau_M*RM_t*p_dphi[i][qp]*JxW[qp];
219  }
220 
221  for (unsigned int i=0; i != n_u_dofs; i++)
222  {
223  Fu(i) -= tau_M*RM_t(0)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
224 
225  Fv(i) -= tau_M*RM_t(1)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
226 
227  if( this->_flow_vars.dim() == 3 )
228  {
229  (*Fw)(i) -= tau_M*RM_t(2)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
230  }
231  }
232 
233  if( compute_jacobian )
234  {
235  libmesh_not_implemented();
236  }
237 
238  }
239  }
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::Real rho) const
libMesh::Number _rho
Material parameters, read from input.
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
VariableIndex p() 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

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