GRINS-0.6.0
Public Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | List of all members
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...
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
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::ParameterMultiPointer< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 
virtual void read_input_options (const GetPot &input)
 Read options from GetPot input file. By default, nothing is read. 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 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 nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
void init_bcs (libMesh::FEMSystem *system)
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
void attach_neumann_bound_func (GRINS::NBCContainer &neumann_bcs)
 
void attach_dirichlet_bound_func (const GRINS::DBCContainer &dirichlet_bc)
 
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_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)
 
BCHandlingBaseget_bc_handler ()
 
ICHandlingBaseget_ic_handler ()
 
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...
 

Protected Attributes

IncompressibleNavierStokesStabilizationHelper _stab_helper
 
unsigned int _dim
 Physical dimension of problem. More...
 
PrimitiveFlowFEVariables _flow_vars
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
Viscosity _mu
 Viscosity object. More...
 
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::BCHandlingBase_bc_handler
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 
bool _is_axisymmetric
 

Static Protected Attributes

static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 

Private Member Functions

 IncompressibleNavierStokesSPGSMStabilization ()
 

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.

References GRINS::Physics::read_input_options().

41  : IncompressibleNavierStokesStabilizationBase<Mu>(physics_name,input)
42  {
43  this->read_input_options(input);
44 
45  return;
46  }
virtual void read_input_options(const GetPot &input)
Read options from GetPot input file. By default, nothing is read.
Definition: physics.C:70

Definition at line 49 of file inc_navier_stokes_spgsm_stab.C.

50  {
51  return;
52  }

Member Function Documentation

void GRINS::Physics::attach_dirichlet_bound_func ( const GRINS::DBCContainer dirichlet_bc)
inherited

Definition at line 150 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_dirichlet_bound_func().

151  {
152  _bc_handler->attach_dirichlet_bound_func( dirichlet_bc );
153  return;
154  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_dirichlet_bound_func(const GRINS::DBCContainer &dirichlet_bc)
void GRINS::Physics::attach_neumann_bound_func ( GRINS::NBCContainer neumann_bcs)
inherited

Definition at line 144 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_neumann_bound_func().

145  {
146  _bc_handler->attach_neumann_bound_func( neumann_bcs );
147  return;
148  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_neumann_bound_func(GRINS::NBCContainer &neumann_bcs)
void GRINS::Physics::auxiliary_init ( MultiphysicsSystem system)
virtualinherited

Any auxillary initialization a Physics class may need.

This is called after all variables are added, so this method can safely query the MultiphysicsSystem about variable information.

Definition at line 113 of file physics.C.

114  {
115  return;
116  }
void GRINS::Physics::compute_element_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 185 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::element_constraint().

187  {
188  return;
189  }
void GRINS::Physics::compute_element_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited
void GRINS::Physics::compute_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 203 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::mass_residual().

205  {
206  return;
207  }
void GRINS::Physics::compute_nonlocal_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 197 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

199  {
200  return;
201  }
void GRINS::Physics::compute_nonlocal_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 209 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

211  {
212  return;
213  }
void GRINS::Physics::compute_nonlocal_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 179 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

181  {
182  return;
183  }
void GRINS::Physics::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtualinherited
void GRINS::Physics::compute_side_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 191 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

193  {
194  return;
195  }
void GRINS::Physics::compute_side_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Reimplemented in GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 173 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_time_derivative().

175  {
176  return;
177  }
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 144 of file inc_navier_stokes_spgsm_stab.C.

147  {
148 #ifdef GRINS_USE_GRVY_TIMERS
149  this->_timer->BeginTimer("IncompressibleNavierStokesSPGSMStabilization::element_constraint");
150 #endif
151 
152  // The number of local degrees of freedom in each variable.
153  const unsigned int n_p_dofs = context.get_dof_indices(this->_flow_vars.p_var()).size();
154 
155  // Element Jacobian * quadrature weights for interior integration.
156  const std::vector<libMesh::Real> &JxW =
157  context.get_element_fe(this->_flow_vars.u_var())->get_JxW();
158 
159  // The pressure shape functions at interior quadrature points.
160  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
161  context.get_element_fe(this->_flow_vars.p_var())->get_dphi();
162 
163  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_flow_vars.p_var()); // R_{p}
164 
165  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u_var());
166 
167  unsigned int n_qpoints = context.get_element_qrule().n_points();
168 
169  for (unsigned int qp=0; qp != n_qpoints; qp++)
170  {
171  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
172  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
173 
174  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u_var(), qp ),
175  context.interior_value( this->_flow_vars.v_var(), qp ) );
176  if( this->_dim == 3 )
177  {
178  U(2) = context.interior_value( this->_flow_vars.w_var(), qp );
179  }
180 
181  // Compute the viscosity at this qp
182  libMesh::Real _mu_qp = this->_mu(context, qp);
183 
184  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
185 
186  libMesh::RealGradient RM_s = this->_stab_helper.compute_res_momentum_steady( context, qp, this->_rho, _mu_qp );
187 
188  for (unsigned int i=0; i != n_p_dofs; i++)
189  {
190  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*JxW[qp];
191  }
192 
193  if( compute_jacobian )
194  {
195  libmesh_not_implemented();
196  }
197 
198  }
199 
200 #ifdef GRINS_USE_GRVY_TIMERS
201  this->_timer->EndTimer("IncompressibleNavierStokesSPGSMStabilization::element_constraint");
202 #endif
203 
204  return;
205  }
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
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
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 55 of file inc_navier_stokes_spgsm_stab.C.

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

Find if current physics is active on supplied element.

Definition at line 83 of file physics.C.

References GRINS::Physics::_enabled_subdomains.

84  {
85  // Check if enabled_subdomains flag has been set and if we're
86  // looking at a real element (rather than a nonlocal evaluation)
87  if( !elem || _enabled_subdomains.empty() )
88  return true;
89 
90  // Check if current physics is enabled on elem
91  if( _enabled_subdomains.find( elem->subdomain_id() ) == _enabled_subdomains.end() )
92  return false;
93 
94  return true;
95  }
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
BCHandlingBase * GRINS::Physics::get_bc_handler ( )
inlineinherited

Definition at line 282 of file physics.h.

References GRINS::Physics::_bc_handler.

283  {
284  return _bc_handler;
285  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
ICHandlingBase * GRINS::Physics::get_ic_handler ( )
inlineinherited

Definition at line 288 of file physics.h.

References GRINS::Physics::_ic_handler.

289  {
290  return _ic_handler;
291  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
void GRINS::Physics::init_bcs ( libMesh::FEMSystem *  system)
inherited

Definition at line 118 of file physics.C.

References GRINS::Physics::_bc_handler, GRINS::BCHandlingBase::init_bc_data(), GRINS::BCHandlingBase::init_dirichlet_bc_func_objs(), GRINS::BCHandlingBase::init_dirichlet_bcs(), and GRINS::BCHandlingBase::init_periodic_bcs().

119  {
120  // Only need to init BC's if the physics actually created a handler
121  if( _bc_handler )
122  {
123  _bc_handler->init_bc_data( *system );
124  _bc_handler->init_dirichlet_bcs( system );
126  _bc_handler->init_periodic_bcs( system );
127  }
128 
129  return;
130  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
virtual void init_dirichlet_bcs(libMesh::FEMSystem *system) const
virtual void init_periodic_bcs(libMesh::FEMSystem *system) const
virtual void init_dirichlet_bc_func_objs(libMesh::FEMSystem *system) const
virtual void init_bc_data(const libMesh::FEMSystem &system)
Override this method to initialize any system-dependent data.
template<class Mu >
void GRINS::IncompressibleNavierStokesStabilizationBase< Mu >::init_context ( AssemblyContext context)
virtualinherited

Initialize context for added physics variables.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

Definition at line 54 of file inc_navier_stokes_stab_base.C.

References GRINS::IncompressibleNavierStokesBase< Viscosity >::init_context().

55  {
56  // First call base class
58 
59  // We need pressure derivatives
60  context.get_element_fe(this->_flow_vars.p_var())->get_dphi();
61 
62  // We also need second derivatives, so initialize those.
63  context.get_element_fe(this->_flow_vars.u_var())->get_d2phi();
64 
65  return;
66  }
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
void GRINS::Physics::init_ics ( libMesh::FEMSystem *  system,
libMesh::CompositeFunction< libMesh::Number > &  all_ics 
)
inherited

Definition at line 133 of file physics.C.

References GRINS::Physics::_ic_handler, and GRINS::ICHandlingBase::init_ic_data().

135  {
136  if( _ic_handler )
137  {
138  _ic_handler->init_ic_data( *system, all_ics );
139  }
140 
141  return;
142  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
virtual void init_ic_data(const libMesh::FEMSystem &system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
Override this method to initialize any system-dependent data.
template<class Mu >
void GRINS::IncompressibleNavierStokesStabilizationBase< Mu >::init_variables ( libMesh::FEMSystem *  system)
virtualinherited

Initialization of Navier-Stokes variables.

Add velocity and pressure variables to system.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

Definition at line 69 of file inc_navier_stokes_stab_base.C.

References GRINS::IncompressibleNavierStokesBase< Viscosity >::init_variables().

70  {
71  // First call base class
73 
74  _stab_helper.init(*system);
75 
76  return;
77  }
IncompressibleNavierStokesStabilizationHelper _stab_helper
virtual void init_variables(libMesh::FEMSystem *system)
Initialization of Navier-Stokes variables.
bool GRINS::Physics::is_steady ( ) const
inherited

Returns whether or not this physics is being solved with a steady solver.

Definition at line 103 of file physics.C.

References GRINS::Physics::_is_steady.

Referenced by GRINS::Physics::set_is_steady().

104  {
105  return _is_steady;
106  }
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
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 208 of file inc_navier_stokes_spgsm_stab.C.

211  {
212 #ifdef GRINS_USE_GRVY_TIMERS
213  this->_timer->BeginTimer("IncompressibleNavierStokesSPGSMStabilization::mass_residual");
214 #endif
215 
216  // The number of local degrees of freedom in each variable.
217  const unsigned int n_p_dofs = context.get_dof_indices(this->_flow_vars.p_var()).size();
218  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u_var()).size();
219 
220  // Element Jacobian * quadrature weights for interior integration.
221  const std::vector<libMesh::Real> &JxW =
222  context.get_element_fe(this->_flow_vars.u_var())->get_JxW();
223 
224  // The pressure shape functions at interior quadrature points.
225  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
226  context.get_element_fe(this->_flow_vars.p_var())->get_dphi();
227 
228  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
229  context.get_element_fe(this->_flow_vars.u_var())->get_dphi();
230 
231  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u_var()); // R_{p}
232  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v_var()); // R_{p}
233  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
234  if(this->_dim == 3)
235  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
236 
237  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_flow_vars.p_var()); // R_{p}
238 
239  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u_var());
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::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
246  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
247 
248  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u_var(), qp ),
249  context.fixed_interior_value( this->_flow_vars.v_var(), qp ) );
250  // Compute the viscosity at this qp
251  libMesh::Real _mu_qp = this->_mu(context, qp);
252 
253  if( this->_dim == 3 )
254  {
255  U(2) = context.fixed_interior_value( this->_flow_vars.w_var(), qp );
256  }
257 
258  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, false );
259 
260  libMesh::RealGradient RM_t = this->_stab_helper.compute_res_momentum_transient( context, qp, this->_rho );
261 
262  for (unsigned int i=0; i != n_p_dofs; i++)
263  {
264  Fp(i) -= tau_M*RM_t*p_dphi[i][qp]*JxW[qp];
265  }
266 
267  for (unsigned int i=0; i != n_u_dofs; i++)
268  {
269  Fu(i) -= tau_M*RM_t(0)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
270 
271  Fv(i) -= tau_M*RM_t(1)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
272 
273  if( this->_dim == 3 )
274  {
275  (*Fw)(i) -= tau_M*RM_t(2)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
276  }
277  }
278 
279  if( compute_jacobian )
280  {
281  libmesh_not_implemented();
282  }
283 
284  }
285 
286 #ifdef GRINS_USE_GRVY_TIMERS
287  this->_timer->EndTimer("IncompressibleNavierStokesSPGSMStabilization::mass_residual");
288 #endif
289 
290  return;
291  }
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
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
void GRINS::Physics::nonlocal_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for scalar variables.

Reimplemented in GRINS::ScalarODE.

Definition at line 250 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

253  {
254  return;
255  }
void GRINS::Physics::nonlocal_mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Mass matrix part(s) for scalar variables.

Reimplemented in GRINS::ScalarODE, and GRINS::AveragedTurbine< Viscosity >.

Definition at line 264 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

267  {
268  return;
269  }
void GRINS::Physics::nonlocal_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Time dependent part(s) of physics for scalar variables.

Reimplemented in GRINS::AveragedTurbine< Viscosity >, and GRINS::ScalarODE.

Definition at line 229 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

232  {
233  return;
234  }
void GRINS::Physics::read_input_options ( const GetPot &  input)
virtualinherited

Read options from GetPot input file. By default, nothing is read.

Reimplemented in GRINS::ScalarODE, GRINS::AveragedTurbineBase< Viscosity >, GRINS::AxisymmetricBoussinesqBuoyancy, GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::AxisymmetricHeatTransfer< Conductivity >, GRINS::AveragedFanBase< Viscosity >, GRINS::ParsedVelocitySourceBase< Viscosity >, GRINS::VelocityDragBase< Viscosity >, GRINS::VelocityPenaltyBase< Viscosity >, GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::HeatTransfer< Conductivity >, GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 70 of file physics.C.

References GRINS::Physics::_enabled_subdomains, and GRINS::Physics::_physics_name.

Referenced by GRINS::HeatTransferBase< Conductivity >::HeatTransferBase(), GRINS::HeatTransferStabilizationBase< Conductivity >::HeatTransferStabilizationBase(), GRINS::IncompressibleNavierStokesAdjointStabilization< Viscosity >::IncompressibleNavierStokesAdjointStabilization(), GRINS::IncompressibleNavierStokesSPGSMStabilization< Viscosity >::IncompressibleNavierStokesSPGSMStabilization(), GRINS::Physics::Physics(), and GRINS::SpalartAllmarasSPGSMStabilization< Viscosity >::SpalartAllmarasSPGSMStabilization().

71  {
72  int num_ids = input.vector_variable_size( "Physics/"+this->_physics_name+"/enabled_subdomains" );
73 
74  for( int i = 0; i < num_ids; i++ )
75  {
76  libMesh::subdomain_id_type dumvar = input( "Physics/"+this->_physics_name+"/enabled_subdomains", -1, i );
77  _enabled_subdomains.insert( dumvar );
78  }
79 
80  return;
81  }
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:254
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
template<class Mu >
void GRINS::IncompressibleNavierStokesBase< Mu >::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiPointer< libMesh::Number > &  param_pointer 
) const
virtualinherited

Each subclass will register its copy of an independent.

Reimplemented from GRINS::ParameterUser.

Definition at line 114 of file inc_navier_stokes_base.C.

References GRINS::ParameterUser::register_parameter().

117  {
118  ParameterUser::register_parameter(param_name, param_pointer);
119  _mu.register_parameter(param_name, param_pointer);
120  }
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
void GRINS::Physics::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtualinherited

Register name of postprocessed quantity with PostProcessedQuantities.

Each Physics class will need to cache an unsigned int corresponding to each postprocessed quantity. This will be used in computing the values and putting them in the CachedVariables object.

Reimplemented in GRINS::ParsedVelocitySource< Viscosity >, GRINS::VelocityPenalty< Viscosity >, GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::HeatTransfer< Conductivity >, GRINS::ElasticCable< StressStrainLaw >, GRINS::ElasticMembrane< StressStrainLaw >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 161 of file physics.C.

163  {
164  return;
165  }
void GRINS::Physics::set_is_steady ( bool  is_steady)
inherited

Sets whether this physics is to be solved with a steady solver or not.

Since the member variable is static, only needs to be called on a single physics.

Definition at line 97 of file physics.C.

References GRINS::Physics::_is_steady, and GRINS::Physics::is_steady().

98  {
100  return;
101  }
bool is_steady() const
Returns whether or not this physics is being solved with a steady solver.
Definition: physics.C:103
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
void GRINS::ParameterUser::set_parameter ( libMesh::Number &  param_variable,
const GetPot &  input,
const std::string &  param_name,
libMesh::Number  param_default 
)
virtualinherited

Each subclass can simultaneously read a parameter value from.

Definition at line 35 of file parameter_user.C.

References GRINS::ParameterUser::_my_name, and GRINS::ParameterUser::_my_parameters.

Referenced by GRINS::AveragedFanAdjointStabilization< Viscosity >::AveragedFanAdjointStabilization(), GRINS::AveragedTurbineAdjointStabilization< Viscosity >::AveragedTurbineAdjointStabilization(), GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::BoussinesqBuoyancyAdjointStabilization(), GRINS::BoussinesqBuoyancyBase::BoussinesqBuoyancyBase(), GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::BoussinesqBuoyancySPGSMStabilization(), GRINS::ConstantConductivity::ConstantConductivity(), GRINS::ConstantPrandtlConductivity::ConstantPrandtlConductivity(), GRINS::ConstantSourceFunction::ConstantSourceFunction(), GRINS::ConstantSourceTerm::ConstantSourceTerm(), GRINS::ConstantSpecificHeat::ConstantSpecificHeat(), GRINS::ConstantViscosity::ConstantViscosity(), GRINS::ElasticCable< StressStrainLaw >::ElasticCable(), GRINS::ElasticCableConstantGravity::ElasticCableConstantGravity(), GRINS::ElasticMembrane< StressStrainLaw >::ElasticMembrane(), GRINS::ElasticMembraneConstantPressure::ElasticMembraneConstantPressure(), GRINS::HeatConduction< Conductivity >::HeatConduction(), GRINS::HeatTransferBase< Conductivity >::HeatTransferBase(), GRINS::IncompressibleNavierStokesBase< Viscosity >::IncompressibleNavierStokesBase(), GRINS::AverageNusseltNumber::init(), GRINS::MooneyRivlin::MooneyRivlin(), GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >::ReactingLowMachNavierStokesBase(), GRINS::HookesLaw1D::read_input_options(), GRINS::HookesLaw::read_input_options(), GRINS::AxisymmetricBoussinesqBuoyancy::read_input_options(), and GRINS::VelocityDragAdjointStabilization< Viscosity >::VelocityDragAdjointStabilization().

39  {
40  param_variable = input(param_name, param_default);
41 
42  libmesh_assert_msg(!_my_parameters.count(param_name),
43  "ERROR: " << _my_name << " double-registered parameter " <<
44  param_name);
45 
46  _my_parameters[param_name] = &param_variable;
47  }
std::map< std::string, libMesh::Number * > _my_parameters
template<class Mu >
void GRINS::IncompressibleNavierStokesBase< Mu >::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtualinherited

Sets velocity variables to be time-evolving.

Reimplemented from GRINS::Physics.

Reimplemented in GRINS::AveragedTurbineBase< Viscosity >, and GRINS::ParsedVelocitySourceBase< Viscosity >.

Definition at line 75 of file inc_navier_stokes_base.C.

Referenced by GRINS::AveragedTurbineBase< Viscosity >::set_time_evolving_vars().

76  {
77  const unsigned int dim = system->get_mesh().mesh_dimension();
78 
79  // Tell the system to march velocity forward in time, but
80  // leave p as a constraint only
81  system->time_evolving(_flow_vars.u_var());
82  system->time_evolving(_flow_vars.v_var());
83 
84  if (dim == 3)
85  system->time_evolving(_flow_vars.w_var());
86 
87  return;
88  }
void GRINS::Physics::side_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for boundaries of elements on the domain boundary.

Reimplemented in GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 243 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

246  {
247  return;
248  }
void GRINS::Physics::side_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Member Data Documentation

GRINS::BCHandlingBase* GRINS::Physics::_bc_handler
protectedinherited
template<class Viscosity >
unsigned int GRINS::IncompressibleNavierStokesBase< Viscosity >::_dim
protectedinherited

Physical dimension of problem.

Todo:
Do we really need to cache this?

Definition at line 79 of file inc_navier_stokes_base.h.

std::set<libMesh::subdomain_id_type> GRINS::Physics::_enabled_subdomains
protectedinherited

Subdomains on which the current Physics class is enabled.

Definition at line 261 of file physics.h.

Referenced by GRINS::Physics::enabled_on_elem(), and GRINS::Physics::read_input_options().

template<class Viscosity >
PrimitiveFlowFEVariables GRINS::IncompressibleNavierStokesBase< Viscosity >::_flow_vars
protectedinherited

Definition at line 81 of file inc_navier_stokes_base.h.

GRINS::ICHandlingBase* GRINS::Physics::_ic_handler
protectedinherited
bool GRINS::Physics::_is_axisymmetric
protectedinherited
bool GRINS::Physics::_is_steady = false
staticprotectedinherited

Caches whether or not the solver that's being used is steady or not.

This is need, for example, in flow stabilization as the tau terms change depending on whether the solver is steady or unsteady.

Definition at line 266 of file physics.h.

Referenced by GRINS::Physics::is_steady(), and GRINS::Physics::set_is_steady().

template<class Viscosity >
Viscosity GRINS::IncompressibleNavierStokesBase< Viscosity >::_mu
protectedinherited

Viscosity object.

Definition at line 88 of file inc_navier_stokes_base.h.

const PhysicsName GRINS::Physics::_physics_name
protectedinherited

Name of the physics object. Used for reading physics specific inputs.

We use a reference because the physics names are const global objects in GRINS namespace

No, we use a copy, because otherwise as soon as the memory in std::set<std::string> requested_physics gets overwritten we get in trouble.

Definition at line 254 of file physics.h.

Referenced by GRINS::SourceTermBase::parse_var_info(), and GRINS::Physics::read_input_options().

template<class Viscosity >
libMesh::Number GRINS::IncompressibleNavierStokesBase< Viscosity >::_rho
protectedinherited

Material parameters, read from input.

Todo:
Create objects to allow for function specification

Definition at line 85 of file inc_navier_stokes_base.h.

Referenced by GRINS::IncompressibleNavierStokesBase< Viscosity >::IncompressibleNavierStokesBase().

template<class Viscosity >
IncompressibleNavierStokesStabilizationHelper GRINS::IncompressibleNavierStokesStabilizationBase< Viscosity >::_stab_helper
protectedinherited

Definition at line 52 of file inc_navier_stokes_stab_base.h.


The documentation for this class was generated from the following files:

Generated on Mon Jun 22 2015 21:32:23 for GRINS-0.6.0 by  doxygen 1.8.9.1