GRINS-0.7.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, CachedValues &cache)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint 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::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 init_variables (libMesh::FEMSystem *system)
 Initialization of Navier-Stokes variables. More...
 
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 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 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...
 

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)
 
- Protected Attributes inherited from GRINS::IncompressibleNavierStokesStabilizationBase< Viscosity >
IncompressibleNavierStokesStabilizationHelper _stab_helper
 
- Protected Attributes inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
unsigned int _dim
 Physical dimension of problem. More...
 
VelocityFEVariables _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  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 134 of file inc_navier_stokes_spgsm_stab.C.

137  {
138 #ifdef GRINS_USE_GRVY_TIMERS
139  this->_timer->BeginTimer("IncompressibleNavierStokesSPGSMStabilization::element_constraint");
140 #endif
141 
142  // The number of local degrees of freedom in each variable.
143  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
144 
145  // Element Jacobian * quadrature weights for interior integration.
146  const std::vector<libMesh::Real> &JxW =
147  context.get_element_fe(this->_flow_vars.u())->get_JxW();
148 
149  // The pressure shape functions at interior quadrature points.
150  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
151  context.get_element_fe(this->_press_var.p())->get_dphi();
152 
153  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
154 
155  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
156 
157  unsigned int n_qpoints = context.get_element_qrule().n_points();
158 
159  for (unsigned int qp=0; qp != n_qpoints; qp++)
160  {
161  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
162  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
163 
164  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
165  context.interior_value( this->_flow_vars.v(), qp ) );
166  if( this->_dim == 3 )
167  {
168  U(2) = context.interior_value( this->_flow_vars.w(), qp );
169  }
170 
171  // Compute the viscosity at this qp
172  libMesh::Real _mu_qp = this->_mu(context, qp);
173 
174  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
175 
176  libMesh::RealGradient RM_s = this->_stab_helper.compute_res_momentum_steady( context, qp, this->_rho, _mu_qp );
177 
178  for (unsigned int i=0; i != n_p_dofs; i++)
179  {
180  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*JxW[qp];
181  }
182 
183  if( compute_jacobian )
184  {
185  libmesh_not_implemented();
186  }
187 
188  }
189 
190 #ifdef GRINS_USE_GRVY_TIMERS
191  this->_timer->EndTimer("IncompressibleNavierStokesSPGSMStabilization::element_constraint");
192 #endif
193 
194  return;
195  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
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.
unsigned int _dim
Physical dimension of problem.
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: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
template<class Mu >
void GRINS::IncompressibleNavierStokesSPGSMStabilization< Mu >::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 45 of file inc_navier_stokes_spgsm_stab.C.

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

201  {
202 #ifdef GRINS_USE_GRVY_TIMERS
203  this->_timer->BeginTimer("IncompressibleNavierStokesSPGSMStabilization::mass_residual");
204 #endif
205 
206  // The number of local degrees of freedom in each variable.
207  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
208  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
209 
210  // Element Jacobian * quadrature weights for interior integration.
211  const std::vector<libMesh::Real> &JxW =
212  context.get_element_fe(this->_flow_vars.u())->get_JxW();
213 
214  // The pressure shape functions at interior quadrature points.
215  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
216  context.get_element_fe(this->_press_var.p())->get_dphi();
217 
218  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
219  context.get_element_fe(this->_flow_vars.u())->get_dphi();
220 
221  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{p}
222  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{p}
223  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
224  if(this->_dim == 3)
225  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
226 
227  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
228 
229  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
230 
231  unsigned int n_qpoints = context.get_element_qrule().n_points();
232 
233  for (unsigned int qp=0; qp != n_qpoints; qp++)
234  {
235  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
236  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
237 
238  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
239  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
240  // Compute the viscosity at this qp
241  libMesh::Real _mu_qp = this->_mu(context, qp);
242 
243  if( this->_dim == 3 )
244  {
245  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
246  }
247 
248  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, false );
249 
250  libMesh::RealGradient RM_t = this->_stab_helper.compute_res_momentum_transient( context, qp, this->_rho );
251 
252  for (unsigned int i=0; i != n_p_dofs; i++)
253  {
254  Fp(i) -= tau_M*RM_t*p_dphi[i][qp]*JxW[qp];
255  }
256 
257  for (unsigned int i=0; i != n_u_dofs; i++)
258  {
259  Fu(i) -= tau_M*RM_t(0)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
260 
261  Fv(i) -= tau_M*RM_t(1)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
262 
263  if( this->_dim == 3 )
264  {
265  (*Fw)(i) -= tau_M*RM_t(2)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
266  }
267  }
268 
269  if( compute_jacobian )
270  {
271  libmesh_not_implemented();
272  }
273 
274  }
275 
276 #ifdef GRINS_USE_GRVY_TIMERS
277  this->_timer->EndTimer("IncompressibleNavierStokesSPGSMStabilization::mass_residual");
278 #endif
279 
280  return;
281  }
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
libMesh::RealGradient compute_res_momentum_transient(AssemblyContext &context, unsigned int qp, const libMesh::Real rho) const
libMesh::Number _rho
Material parameters, read from input.
unsigned int _dim
Physical dimension of problem.
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
VariableIndex p() 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

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