GRINS-0.6.0
Public Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | List of all members
GRINS::Stokes< Viscosity > Class Template Reference

Physics class for Stokes. More...

#include <stokes.h>

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

Public Member Functions

 Stokes (const std::string &physics_name, const GetPot &input)
 
 ~Stokes ()
 
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_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 init_context (AssemblyContext &context)
 Initialize context for added physics variables. 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

PressurePinning _p_pinning
 
bool _pin_pressure
 Enable pressure pinning. More...
 
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

 Stokes ()
 

Detailed Description

template<class Viscosity>
class GRINS::Stokes< Viscosity >

Physics class for Stokes.

This physics class implements the classical Stokes equations.

Definition at line 42 of file stokes.h.

Constructor & Destructor Documentation

template<class Mu >
GRINS::Stokes< Mu >::Stokes ( const std::string &  physics_name,
const GetPot &  input 
)

Definition at line 43 of file stokes.C.

References GRINS::Physics::_bc_handler, and GRINS::Physics::_ic_handler.

44  : IncompressibleNavierStokesBase<Mu>(physics_name,
45  stokes, /* "core" Physics name */
46  input),
47  _p_pinning(input,physics_name),
48  _pin_pressure( input("Physics/"+stokes+"/pin_pressure", false ) )
49  {
50  // This is deleted in the base class
51  this->_bc_handler = new IncompressibleNavierStokesBCHandling( physics_name, input );
52  this->_ic_handler = new GenericICHandler( physics_name, input );
53 
54  return;
55  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
bool _pin_pressure
Enable pressure pinning.
Definition: stokes.h:73
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
const PhysicsName stokes
PressurePinning _p_pinning
Definition: stokes.h:70
template<class Mu >
GRINS::Stokes< Mu >::~Stokes ( )

Definition at line 58 of file stokes.C.

59  {
60  return;
61  }
template<class Viscosity >
GRINS::Stokes< Viscosity >::Stokes ( )
private

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::Stokes< 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 216 of file stokes.C.

219  {
220 #ifdef GRINS_USE_GRVY_TIMERS
221  this->_timer->BeginTimer("Stokes::element_constraint");
222 #endif
223 
224  // The number of local degrees of freedom in each variable.
225  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u_var()).size();
226  const unsigned int n_p_dofs = context.get_dof_indices(this->_flow_vars.p_var()).size();
227 
228  // We get some references to cell-specific data that
229  // will be used to assemble the linear system.
230 
231  // Element Jacobian * quadrature weights for interior integration.
232  const std::vector<libMesh::Real> &JxW =
233  context.get_element_fe(this->_flow_vars.u_var())->get_JxW();
234 
235  // The velocity shape function gradients (in global coords.)
236  // at interior quadrature points.
237  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
238  context.get_element_fe(this->_flow_vars.u_var())->get_dphi();
239 
240  // The pressure shape functions at interior quadrature points.
241  const std::vector<std::vector<libMesh::Real> >& p_phi =
242  context.get_element_fe(this->_flow_vars.p_var())->get_phi();
243 
244  // The subvectors and submatrices we need to fill:
245  //
246  // Kpu, Kpv, Kpw, Fp
247 
248  libMesh::DenseSubMatrix<libMesh::Number> &Kpu = context.get_elem_jacobian(this->_flow_vars.p_var(), this->_flow_vars.u_var()); // R_{p},{u}
249  libMesh::DenseSubMatrix<libMesh::Number> &Kpv = context.get_elem_jacobian(this->_flow_vars.p_var(), this->_flow_vars.v_var()); // R_{p},{v}
250  libMesh::DenseSubMatrix<libMesh::Number>* Kpw = NULL;
251 
252  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_flow_vars.p_var()); // R_{p}
253 
254  if( this->_dim == 3 )
255  {
256  Kpw = &context.get_elem_jacobian(this->_flow_vars.p_var(), this->_flow_vars.w_var()); // R_{p},{w}
257  }
258 
259  // Add the constraint given by the continuity equation.
260  unsigned int n_qpoints = context.get_element_qrule().n_points();
261  for (unsigned int qp=0; qp != n_qpoints; qp++)
262  {
263  // Compute the velocity gradient at the old Newton iterate.
264  libMesh::Gradient grad_u, grad_v, grad_w;
265  grad_u = context.interior_gradient(this->_flow_vars.u_var(), qp);
266  grad_v = context.interior_gradient(this->_flow_vars.v_var(), qp);
267  if (this->_dim == 3)
268  grad_w = context.interior_gradient(this->_flow_vars.w_var(), qp);
269 
270  // Now a loop over the pressure degrees of freedom. This
271  // computes the contributions of the continuity equation.
272  for (unsigned int i=0; i != n_p_dofs; i++)
273  {
274  Fp(i) += JxW[qp] * p_phi[i][qp] *
275  (grad_u(0) + grad_v(1));
276  if (this->_dim == 3)
277  Fp(i) += JxW[qp] * p_phi[i][qp] *
278  (grad_w(2));
279 
280  if (compute_jacobian)
281  {
282  for (unsigned int j=0; j != n_u_dofs; j++)
283  {
284  Kpu(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](0);
285  Kpv(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](1);
286  if (this->_dim == 3)
287  (*Kpw)(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](2);
288  } // end of the inner dof (j) loop
289 
290  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
291 
292  } // end of the outer dof (i) loop
293  } // end of the quadrature point (qp) loop
294 
295 
296  // Pin p = p_value at p_point
297  if( _pin_pressure )
298  {
299  _p_pinning.pin_value( context, compute_jacobian, this->_flow_vars.p_var() );
300  }
301 
302 
303 #ifdef GRINS_USE_GRVY_TIMERS
304  this->_timer->EndTimer("Stokes::element_constraint");
305 #endif
306 
307  return;
308  }
bool _pin_pressure
Enable pressure pinning.
Definition: stokes.h:73
unsigned int _dim
Physical dimension of problem.
void pin_value(libMesh::DiffContext &context, const bool request_jacobian, const GRINS::VariableIndex var, const double penalty=1.0)
The idea here is to pin a variable to a particular value if there is a null space - e...
PressurePinning _p_pinning
Definition: stokes.h:70
template<class Mu >
void GRINS::Stokes< 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 64 of file stokes.C.

67  {
68 #ifdef GRINS_USE_GRVY_TIMERS
69  this->_timer->BeginTimer("Stokes::element_time_derivative");
70 #endif
71 
72  // The number of local degrees of freedom in each variable.
73  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u_var()).size();
74  const unsigned int n_p_dofs = context.get_dof_indices(this->_flow_vars.p_var()).size();
75 
76  // Check number of dofs is same for this->_flow_vars.u_var(), v_var and w_var.
77  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v_var()).size());
78  if (this->_dim == 3)
79  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w_var()).size());
80 
81  // We get some references to cell-specific data that
82  // will be used to assemble the linear system.
83 
84  // Element Jacobian * quadrature weights for interior integration.
85  const std::vector<libMesh::Real> &JxW =
86  context.get_element_fe(this->_flow_vars.u_var())->get_JxW();
87 
88  // The velocity shape function gradients (in global coords.)
89  // at interior quadrature points.
90  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
91  context.get_element_fe(this->_flow_vars.u_var())->get_dphi();
92 
93  // The pressure shape functions at interior quadrature points.
94  const std::vector<std::vector<libMesh::Real> >& p_phi =
95  context.get_element_fe(this->_flow_vars.p_var())->get_phi();
96 
97  // The subvectors and submatrices we need to fill:
98  //
99  // K_{\alpha \beta} = R_{\alpha},{\beta} = \partial{ R_{\alpha} } / \partial{ {\beta} } (where R denotes residual)
100  // e.g., for \alpha = v and \beta = u we get: K{vu} = R_{v},{u}
101  // Note that Kpu, Kpv, Kpw and Fp comes as constraint.
102 
103  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.u_var()); // R_{u},{u}
104  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.v_var()); // R_{v},{v}
105  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
106 
107  libMesh::DenseSubMatrix<libMesh::Number> &Kup = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.p_var()); // R_{u},{p}
108  libMesh::DenseSubMatrix<libMesh::Number> &Kvp = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.p_var()); // R_{v},{p}
109  libMesh::DenseSubMatrix<libMesh::Number>* Kwp = NULL;
110 
111  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u_var()); // R_{u}
112  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v_var()); // R_{v}
113  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
114 
115  if( this->_dim == 3 )
116  {
117  Kww = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.w_var()); // R_{w},{w}
118  Kwp = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.p_var()); // R_{w},{p}
119  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
120  }
121 
122  // Now we will build the element Jacobian and residual.
123  // Constructing the residual requires the solution and its
124  // gradient from the previous timestep. This must be
125  // calculated at each quadrature point by summing the
126  // solution degree-of-freedom values by the appropriate
127  // weight functions.
128  unsigned int n_qpoints = context.get_element_qrule().n_points();
129 
130  for (unsigned int qp=0; qp != n_qpoints; qp++)
131  {
132  // Compute the solution & its gradient at the old Newton iterate.
133  libMesh::Number p, u, v, w;
134  p = context.interior_value(this->_flow_vars.p_var(), qp);
135  u = context.interior_value(this->_flow_vars.u_var(), qp);
136  v = context.interior_value(this->_flow_vars.v_var(), qp);
137  if (this->_dim == 3)
138  w = context.interior_value(this->_flow_vars.w_var(), qp);
139 
140  libMesh::Gradient grad_u, grad_v, grad_w;
141  grad_u = context.interior_gradient(this->_flow_vars.u_var(), qp);
142  grad_v = context.interior_gradient(this->_flow_vars.v_var(), qp);
143  if (this->_dim == 3)
144  grad_w = context.interior_gradient(this->_flow_vars.w_var(), qp);
145 
146  libMesh::NumberVectorValue Uvec (u,v);
147  if (this->_dim == 3)
148  Uvec(2) = w;
149 
150  // Compute the viscosity at this qp
151  libMesh::Real _mu_qp = this->_mu(context, qp);
152 
153  // First, an i-loop over the velocity degrees of freedom.
154  // We know that n_u_dofs == n_v_dofs so we can compute contributions
155  // for both at the same time.
156  for (unsigned int i=0; i != n_u_dofs; i++)
157  {
158  Fu(i) += JxW[qp] *
159  ( p*u_gradphi[i][qp](0) // pressure term
160  -_mu_qp*(u_gradphi[i][qp]*grad_u) ); // diffusion term
161 
162  Fv(i) += JxW[qp] *
163  ( p*u_gradphi[i][qp](1) // pressure term
164  -_mu_qp*(u_gradphi[i][qp]*grad_v) ); // diffusion term
165  if (this->_dim == 3)
166  {
167  (*Fw)(i) += JxW[qp] *
168  ( p*u_gradphi[i][qp](2) // pressure term
169  -_mu_qp*(u_gradphi[i][qp]*grad_w) ); // diffusion term
170  }
171 
172  if (compute_jacobian)
173  {
174  for (unsigned int j=0; j != n_u_dofs; j++)
175  {
176  // TODO: precompute some terms like:
177  // (Uvec*vel_gblgradphivec[j][qp]),
178  // vel_phi[i][qp]*vel_phi[j][qp],
179  // (vel_gblgradphivec[i][qp]*vel_gblgradphivec[j][qp])
180 
181  Kuu(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
182  (-_mu_qp*(u_gradphi[i][qp]*u_gradphi[j][qp])); // diffusion term
183 
184  Kvv(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
185  (-_mu_qp*(u_gradphi[i][qp]*u_gradphi[j][qp])); // diffusion term
186 
187  if (this->_dim == 3)
188  {
189  (*Kww)(i,j) += JxW[qp] * context.get_elem_solution_derivative() *
190  (-_mu_qp*(u_gradphi[i][qp]*u_gradphi[j][qp])); // diffusion term
191  }
192  } // end of the inner dof (j) loop
193 
194  // Matrix contributions for the up, vp and wp couplings
195  for (unsigned int j=0; j != n_p_dofs; j++)
196  {
197  Kup(i,j) += context.get_elem_solution_derivative() * JxW[qp]*u_gradphi[i][qp](0)*p_phi[j][qp];
198  Kvp(i,j) += context.get_elem_solution_derivative() * JxW[qp]*u_gradphi[i][qp](1)*p_phi[j][qp];
199  if (this->_dim == 3)
200  (*Kwp)(i,j) += context.get_elem_solution_derivative() * JxW[qp]*u_gradphi[i][qp](2)*p_phi[j][qp];
201  } // end of the inner dof (j) loop
202 
203  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
204 
205  } // end of the outer dof (i) loop
206  } // end of the quadrature point (qp) loop
207 
208 #ifdef GRINS_USE_GRVY_TIMERS
209  this->_timer->EndTimer("Stokes::element_time_derivative");
210 #endif
211 
212  return;
213  }
unsigned int _dim
Physical dimension of problem.
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::IncompressibleNavierStokesBase< Mu >::init_context ( AssemblyContext context)
virtualinherited

Initialize context for added physics variables.

Reimplemented from GRINS::Physics.

Reimplemented in GRINS::AveragedTurbine< Viscosity >, GRINS::VelocityDrag< Viscosity >, GRINS::AveragedFan< Viscosity >, GRINS::ParsedVelocitySource< Viscosity >, GRINS::VelocityPenalty< Viscosity >, GRINS::AveragedTurbineAdjointStabilization< Viscosity >, GRINS::VelocityDragAdjointStabilization< Viscosity >, GRINS::AveragedFanAdjointStabilization< Viscosity >, GRINS::ParsedVelocitySourceAdjointStabilization< Viscosity >, GRINS::VelocityPenaltyAdjointStabilization< Viscosity >, and GRINS::IncompressibleNavierStokesStabilizationBase< Viscosity >.

Definition at line 91 of file inc_navier_stokes_base.C.

Referenced by GRINS::IncompressibleNavierStokesStabilizationBase< Viscosity >::init_context().

92  {
93  // We should prerequest all the data
94  // we will need to build the linear system
95  // or evaluate a quantity of interest.
96  context.get_element_fe(_flow_vars.u_var())->get_JxW();
97  context.get_element_fe(_flow_vars.u_var())->get_phi();
98  context.get_element_fe(_flow_vars.u_var())->get_dphi();
99  context.get_element_fe(_flow_vars.u_var())->get_xyz();
100 
101  context.get_element_fe(_flow_vars.p_var())->get_phi();
102  context.get_element_fe(_flow_vars.p_var())->get_xyz();
103 
104  context.get_side_fe(_flow_vars.u_var())->get_JxW();
105  context.get_side_fe(_flow_vars.u_var())->get_phi();
106  context.get_side_fe(_flow_vars.u_var())->get_dphi();
107  context.get_side_fe(_flow_vars.u_var())->get_xyz();
108 
109  return;
110  }
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::IncompressibleNavierStokesBase< Mu >::init_variables ( libMesh::FEMSystem *  system)
virtualinherited

Initialization of Navier-Stokes variables.

Add velocity and pressure variables to system.

Implements GRINS::Physics.

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

Definition at line 63 of file inc_navier_stokes_base.C.

Referenced by GRINS::IncompressibleNavierStokesStabilizationBase< Viscosity >::init_variables(), and GRINS::AveragedTurbineBase< Viscosity >::init_variables().

64  {
65  this->_dim = system->get_mesh().mesh_dimension();
66 
67  this->_flow_vars.init(system);
68 
69  this->_mu.init(system);
70 
71  return;
72  }
virtual void init(libMesh::FEMSystem *system)
unsigned int _dim
Physical dimension of problem.
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::Stokes< 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 311 of file stokes.C.

314  {
315  // Element Jacobian * quadrature weights for interior integration
316  // We assume the same for each flow variable
317  const std::vector<libMesh::Real> &JxW =
318  context.get_element_fe(this->_flow_vars.u_var())->get_JxW();
319 
320  // The shape functions at interior quadrature points.
321  // We assume the same for each flow variable
322  const std::vector<std::vector<libMesh::Real> >& u_phi =
323  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
324 
325  // The number of local degrees of freedom in each variable
326  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u_var()).size();
327 
328  // The subvectors and submatrices we need to fill:
329  libMesh::DenseSubVector<libMesh::Real> &F_u = context.get_elem_residual(this->_flow_vars.u_var());
330  libMesh::DenseSubVector<libMesh::Real> &F_v = context.get_elem_residual(this->_flow_vars.v_var());
331  libMesh::DenseSubVector<libMesh::Real>* F_w = NULL;
332 
333  libMesh::DenseSubMatrix<libMesh::Real> &M_uu = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.u_var());
334  libMesh::DenseSubMatrix<libMesh::Real> &M_vv = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.v_var());
335  libMesh::DenseSubMatrix<libMesh::Real>* M_ww = NULL;
336 
337  if( this->_dim == 3 )
338  {
339  F_w = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
340  M_ww = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.w_var());
341  }
342 
343  unsigned int n_qpoints = context.get_element_qrule().n_points();
344 
345  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
346  {
347  // For the mass residual, we need to be a little careful.
348  // The time integrator is handling the time-discretization
349  // for us so we need to supply M(u_fixed)*u' for the residual.
350  // u_fixed will be given by the fixed_interior_value function
351  // while u' will be given by the interior_rate function.
352  libMesh::Real u_dot, v_dot, w_dot = 0.0;
353  context.interior_rate(this->_flow_vars.u_var(), qp, u_dot);
354  context.interior_rate(this->_flow_vars.v_var(), qp, v_dot);
355 
356  if( this->_dim == 3 )
357  context.interior_rate(this->_flow_vars.w_var(), qp, w_dot);
358 
359  for (unsigned int i = 0; i != n_u_dofs; ++i)
360  {
361  F_u(i) -= JxW[qp]*this->_rho*u_dot*u_phi[i][qp];
362  F_v(i) -= JxW[qp]*this->_rho*v_dot*u_phi[i][qp];
363 
364  if( this->_dim == 3 )
365  (*F_w)(i) -= JxW[qp]*this->_rho*w_dot*u_phi[i][qp];
366 
367  if( compute_jacobian )
368  {
369  for (unsigned int j=0; j != n_u_dofs; j++)
370  {
371  // Assuming rho is constant w.r.t. u, v, w
372  // and T (if Boussinesq added).
373  libMesh::Real value = context.get_elem_solution_derivative() * JxW[qp]*this->_rho*u_phi[i][qp]*u_phi[j][qp];
374 
375  M_uu(i,j) -= value;
376  M_vv(i,j) -= value;
377 
378  if( this->_dim == 3)
379  {
380  (*M_ww)(i,j) -= value;
381  }
382 
383  } // End dof loop
384  } // End Jacobian check
385  } // End dof loop
386  } // End quadrature loop
387 
388  return;
389  }
libMesh::Number _rho
Material parameters, read from input.
unsigned int _dim
Physical dimension of problem.
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.

template<class Viscosity >
PressurePinning GRINS::Stokes< Viscosity >::_p_pinning
protected

Definition at line 70 of file stokes.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 >
bool GRINS::Stokes< Viscosity >::_pin_pressure
protected

Enable pressure pinning.

Definition at line 73 of file stokes.h.

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().


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

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