GRINS-0.7.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, CachedValues &cache)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for element interiors. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
- Public Member Functions inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
 IncompressibleNavierStokesBase (const std::string &my_physics_name, const std::string &core_physics_name, const GetPot &input)
 
 ~IncompressibleNavierStokesBase ()
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialization of Navier-Stokes variables. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
virtual void 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 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 side_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void side_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_element_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_damping_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

Protected Attributes

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

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)
 
- 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:269
bool _pin_pressure
Enable pressure pinning.
Definition: stokes.h:76
static PhysicsName stokes()
PressurePinning _p_pinning
Definition: stokes.h:73
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:76
void check_pin_location(const libMesh::MeshBase &mesh)
Check the mesh to ensure pin location is found.
PressurePinning _p_pinning
Definition: stokes.h:73
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 220 of file stokes.C.

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

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

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

Member Data Documentation

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

Definition at line 73 of file stokes.h.

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

Enable pressure pinning.

Definition at line 76 of file stokes.h.


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

Generated on Thu Jun 2 2016 21:52:32 for GRINS-0.7.0 by  doxygen 1.8.10