GRINS-0.8.0
List of all members | Public Member Functions | Protected Attributes | Private Member Functions
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 auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context)
 Constraint part(s) of physics for element interiors. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
- Public Member Functions inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
 IncompressibleNavierStokesBase (const std::string &my_physics_name, const std::string &core_physics_name, const GetPot &input)
 
 ~IncompressibleNavierStokesBase ()
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 
libMesh::Real get_viscosity_value (AssemblyContext &context, unsigned int qp) const
 
- Public Member Functions inherited from GRINS::Physics
 Physics (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~Physics ()
 
virtual void init_variables (libMesh::FEMSystem *)
 Initialize variables for this physics. More...
 
virtual bool enabled_on_elem (const libMesh::Elem *elem)
 Find if current physics is active on supplied element. More...
 
void set_is_steady (bool is_steady)
 Sets whether this physics is to be solved with a steady solver or not. More...
 
bool is_steady () const
 Returns whether or not this physics is being solved with a steady solver. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. More...
 
virtual void preassembly (MultiphysicsSystem &)
 Perform any necessary setup before element assembly begins. More...
 
virtual void reinit (MultiphysicsSystem &)
 Any reinitialization that needs to be done. More...
 
virtual void side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

Protected Attributes

PressurePinning _p_pinning
 
bool _pin_pressure
 Enable pressure pinning. More...
 
- Protected Attributes inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
VelocityVariable_flow_vars
 
PressureFEVariable_press_var
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
Viscosity _mu
 Viscosity object. More...
 
- Protected Attributes inherited from GRINS::Physics
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 

Private Member Functions

 Stokes ()
 

Additional Inherited Members

- Static Public Member Functions inherited from GRINS::Physics
static void set_is_axisymmetric (bool is_axisymmetric)
 Set whether we should treat the problem as axisymmetric. More...
 
static bool is_axisymmetric ()
 
- Static Public Attributes inherited from GRINS::ParameterUser
static std::string zero_vector_function = std::string("{0}")
 A parseable function string with LIBMESH_DIM components, all 0. More...
 
- Protected Member Functions inherited from GRINS::Physics
libMesh::UniquePtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe (const libMesh::Elem *elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
 
void parse_enabled_subdomains (const GetPot &input, const std::string &physics_name)
 
void check_var_subdomain_consistency (const FEVariablesBase &var) const
 Check that var is enabled on at least the subdomains this Physics is. More...
 
- 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::Stokes< Viscosity >

Physics class for Stokes.

This physics class implements the classical Stokes equations.

Definition at line 43 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 44 of file stokes.C.

References GRINS::Physics::_ic_handler.

45  : IncompressibleNavierStokesBase<Mu>(physics_name,
46  PhysicsNaming::stokes(), /* "core" Physics name */
47  input),
48  _p_pinning(input,physics_name),
49  _pin_pressure( input("Physics/"+PhysicsNaming::stokes()+"/pin_pressure", false ) )
50  {
51  this->_ic_handler = new GenericICHandler( physics_name, input );
52  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
bool _pin_pressure
Enable pressure pinning.
Definition: stokes.h:73
static PhysicsName stokes()
PressurePinning _p_pinning
Definition: stokes.h:70
template<class Mu >
GRINS::Stokes< Mu >::~Stokes ( )

Definition at line 55 of file stokes.C.

56  {
57  return;
58  }
template<class Viscosity >
GRINS::Stokes< Viscosity >::Stokes ( )
private

Member Function Documentation

template<class Mu >
void GRINS::Stokes< Mu >::auxiliary_init ( MultiphysicsSystem system)
virtual

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.

Reimplemented from GRINS::Physics.

Definition at line 61 of file stokes.C.

62  {
63  if( _pin_pressure )
64  _p_pinning.check_pin_location(system.get_mesh());
65  }
bool _pin_pressure
Enable pressure pinning.
Definition: stokes.h:73
void check_pin_location(const libMesh::MeshBase &mesh)
Check the mesh to ensure pin location is found.
PressurePinning _p_pinning
Definition: stokes.h:70
template<class Mu >
void GRINS::Stokes< Mu >::element_constraint ( bool  ,
AssemblyContext  
)
virtual

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 211 of file stokes.C.

212  {
213  // The number of local degrees of freedom in each variable.
214  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
215  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
216 
217  // We get some references to cell-specific data that
218  // will be used to assemble the linear system.
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())->get_JxW();
223 
224  // The velocity shape function gradients (in global coords.)
225  // at interior quadrature points.
226  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
227  context.get_element_fe(this->_flow_vars.u())->get_dphi();
228 
229  // The pressure shape functions at interior quadrature points.
230  const std::vector<std::vector<libMesh::Real> >& p_phi =
231  context.get_element_fe(this->_press_var.p())->get_phi();
232 
233  // The subvectors and submatrices we need to fill:
234  //
235  // Kpu, Kpv, Kpw, Fp
236 
237  libMesh::DenseSubMatrix<libMesh::Number> &Kpu = context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.u()); // R_{p},{u}
238  libMesh::DenseSubMatrix<libMesh::Number> &Kpv = context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.v()); // R_{p},{v}
239  libMesh::DenseSubMatrix<libMesh::Number>* Kpw = NULL;
240 
241  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
242 
243 
244  if( this->_flow_vars.dim() == 3 )
245  {
246  Kpw = &context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.w()); // R_{p},{w}
247  }
248 
249  // Add the constraint given by the continuity equation.
250  unsigned int n_qpoints = context.get_element_qrule().n_points();
251  for (unsigned int qp=0; qp != n_qpoints; qp++)
252  {
253  // Compute the velocity gradient at the old Newton iterate.
254  libMesh::Gradient grad_u, grad_v, grad_w;
255  grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
256  grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
257  if (this->_flow_vars.dim() == 3)
258  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
259 
260  // Now a loop over the pressure degrees of freedom. This
261  // computes the contributions of the continuity equation.
262  for (unsigned int i=0; i != n_p_dofs; i++)
263  {
264  Fp(i) += JxW[qp] * p_phi[i][qp] *
265  (grad_u(0) + grad_v(1));
266  if (this->_flow_vars.dim() == 3)
267  Fp(i) += JxW[qp] * p_phi[i][qp] *
268  (grad_w(2));
269 
270  if (compute_jacobian)
271  {
272  for (unsigned int j=0; j != n_u_dofs; j++)
273  {
274  Kpu(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](0);
275  Kpv(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](1);
276  if (this->_flow_vars.dim() == 3)
277  (*Kpw)(i,j) += context.get_elem_solution_derivative() * JxW[qp]*p_phi[i][qp]*u_gradphi[j][qp](2);
278  } // end of the inner dof (j) loop
279 
280  } // end - if (compute_jacobian && context.get_elem_solution_derivative())
281 
282  } // end of the outer dof (i) loop
283  } // end of the quadrature point (qp) loop
284 
285 
286  // Pin p = p_value at p_point
287  if( _pin_pressure )
288  {
289  _p_pinning.pin_value( context, compute_jacobian, this->_press_var.p() );
290  }
291  }
bool _pin_pressure
Enable pressure pinning.
Definition: stokes.h:73
VariableIndex p() const
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...
unsigned int dim() const
Number of components.
PressurePinning _p_pinning
Definition: stokes.h:70
template<class Mu >
void GRINS::Stokes< Mu >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 69 of file stokes.C.

70  {
71  // The number of local degrees of freedom in each variable.
72  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
73  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
74 
75  // Check number of dofs is same for this->_flow_vars.u(), v_var and w_var.
76  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
77 
78  if (this->_flow_vars.dim() == 3)
79  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).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())->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())->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->_press_var.p())->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(), this->_flow_vars.u()); // R_{u},{u}
104  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // 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(), this->_press_var.p()); // R_{u},{p}
108  libMesh::DenseSubMatrix<libMesh::Number> &Kvp = context.get_elem_jacobian(this->_flow_vars.v(), this->_press_var.p()); // 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()); // R_{u}
112  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
113  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
114 
115  if( this->_flow_vars.dim() == 3 )
116  {
117  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
118  Kwp = &context.get_elem_jacobian(this->_flow_vars.w(), this->_press_var.p()); // R_{w},{p}
119  Fw = &context.get_elem_residual(this->_flow_vars.w()); // 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->_press_var.p(), qp);
135  u = context.interior_value(this->_flow_vars.u(), qp);
136  v = context.interior_value(this->_flow_vars.v(), qp);
137  if (this->_flow_vars.dim() == 3)
138  w = context.interior_value(this->_flow_vars.w(), qp);
139 
140  libMesh::Gradient grad_u, grad_v, grad_w;
141  grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
142  grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
143  if (this->_flow_vars.dim() == 3)
144  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
145 
146  libMesh::NumberVectorValue Uvec (u,v);
147  if (this->_flow_vars.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->_flow_vars.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->_flow_vars.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->_flow_vars.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  }
VariableIndex p() const
unsigned int dim() const
Number of components.
template<class Mu >
void GRINS::Stokes< Mu >::mass_residual ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 295 of file stokes.C.

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

Member Data Documentation

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

Definition at line 70 of file stokes.h.

template<class Viscosity >
bool GRINS::Stokes< Viscosity >::_pin_pressure
protected

Enable pressure pinning.

Definition at line 73 of file stokes.h.


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

Generated on Tue Dec 19 2017 12:47:32 for GRINS-0.8.0 by  doxygen 1.8.9.1