www.mooseframework.org
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
INSMassRZ Class Reference

This class computes the mass equation residual and Jacobian contributions for the incompressible Navier-Stokes momentum equation in RZ coordinates. More...

#include <INSMassRZ.h>

Inheritance diagram for INSMassRZ:
[legend]

Public Member Functions

 INSMassRZ (const InputParameters &parameters)
 
virtual ~INSMassRZ ()
 

Protected Member Functions

virtual RealVectorValue strongViscousTermTraction () override
 
virtual RealVectorValue dStrongViscDUCompTraction (unsigned comp) override
 
virtual RealVectorValue strongViscousTermLaplace () override
 
virtual RealVectorValue dStrongViscDUCompLaplace (unsigned comp) override
 
virtual Real computeQpResidual () override
 
virtual Real computeQpOffDiagJacobian (unsigned jvar) override
 
virtual Real computeQpJacobian ()
 
virtual Real computeQpPGResidual ()
 
virtual Real computeQpPGJacobian ()
 
virtual Real computeQpPGOffDiagJacobian (unsigned comp)
 
virtual RealVectorValue convectiveTerm ()
 
virtual RealVectorValue dConvecDUComp (unsigned comp)
 
virtual RealVectorValue weakViscousTermLaplace (unsigned comp)
 
virtual RealVectorValue weakViscousTermTraction (unsigned comp)
 
virtual RealVectorValue dWeakViscDUCompLaplace ()
 
virtual RealVectorValue dWeakViscDUCompTraction ()
 
virtual RealVectorValue strongPressureTerm ()
 
virtual Real weakPressureTerm ()
 
virtual RealVectorValue dStrongPressureDPressure ()
 
virtual Real dWeakPressureDPressure ()
 
virtual RealVectorValue gravityTerm ()
 
virtual RealVectorValue timeDerivativeTerm ()
 
virtual RealVectorValue dTimeDerivativeDUComp (unsigned comp)
 
virtual Real tau ()
 
virtual Real dTauDUComp (unsigned comp)
 
virtual Real tauNodal ()
 Provides tau which yields superconvergence for 1D advection-diffusion. More...
 

Protected Attributes

bool _pspg
 
const Function & _x_ffn
 
const Function & _y_ffn
 
const Function & _z_ffn
 
const VariablePhiSecond & _second_phi
 second derivatives of the shape function More...
 
const VariableValue & _u_vel
 
const VariableValue & _v_vel
 
const VariableValue & _w_vel
 
const VariableValue & _p
 
const VariableGradient & _grad_u_vel
 
const VariableGradient & _grad_v_vel
 
const VariableGradient & _grad_w_vel
 
const VariableGradient & _grad_p
 
const VariableSecond & _second_u_vel
 
const VariableSecond & _second_v_vel
 
const VariableSecond & _second_w_vel
 
const VariableValue & _u_vel_dot
 
const VariableValue & _v_vel_dot
 
const VariableValue & _w_vel_dot
 
const VariableValue & _d_u_vel_dot_du
 
const VariableValue & _d_v_vel_dot_dv
 
const VariableValue & _d_w_vel_dot_dw
 
unsigned _u_vel_var_number
 
unsigned _v_vel_var_number
 
unsigned _w_vel_var_number
 
unsigned _p_var_number
 
RealVectorValue _gravity
 
const MaterialProperty< Real > & _mu
 
const MaterialProperty< Real > & _rho
 
const Real & _alpha
 
bool _laplace
 
bool _convective_term
 
bool _transient_term
 

Detailed Description

This class computes the mass equation residual and Jacobian contributions for the incompressible Navier-Stokes momentum equation in RZ coordinates.

Inherits most of its functionality from INSMass, and calls its computeQpXYZ() functions when necessary.

Definition at line 27 of file INSMassRZ.h.

Constructor & Destructor Documentation

◆ INSMassRZ()

INSMassRZ::INSMassRZ ( const InputParameters &  parameters)

Definition at line 25 of file INSMassRZ.C.

25 : INSMass(parameters) {}

◆ ~INSMassRZ()

virtual INSMassRZ::~INSMassRZ ( )
inlinevirtual

Definition at line 31 of file INSMassRZ.h.

31 {}

Member Function Documentation

◆ computeQpJacobian()

Real INSMass::computeQpJacobian ( )
protectedvirtualinherited

Implements INSBase.

Definition at line 74 of file INSMass.C.

75 {
76  // Derivative wrt to p is zero
77  Real r = 0;
78 
79  // Unless we are doing GLS stabilization
80  if (_pspg)
81  r += computeQpPGJacobian();
82 
83  return r;
84 }

◆ computeQpOffDiagJacobian()

Real INSMassRZ::computeQpOffDiagJacobian ( unsigned  jvar)
overrideprotectedvirtual

Reimplemented from INSMass.

Definition at line 92 of file INSMassRZ.C.

93 {
94  // Base class jacobian contribution
95  Real jac_base = INSMass::computeQpOffDiagJacobian(jvar);
96 
97  // The radial coordinate value.
98  const Real r = _q_point[_qp](0);
99 
100  if (jvar == _u_vel_var_number)
101  jac_base -= _phi[_j][_qp] / r * _test[_i][_qp];
102 
103  return jac_base;
104 }

◆ computeQpPGJacobian()

Real INSMass::computeQpPGJacobian ( )
protectedvirtualinherited

Definition at line 87 of file INSMass.C.

88 {
89  return -1. / _rho[_qp] * tau() * _grad_test[_i][_qp] * dStrongPressureDPressure();
90 }

Referenced by INSMass::computeQpJacobian().

◆ computeQpPGOffDiagJacobian()

Real INSMass::computeQpPGOffDiagJacobian ( unsigned  comp)
protectedvirtualinherited

Definition at line 124 of file INSMass.C.

125 {
126  RealVectorValue convective_term = _convective_term ? convectiveTerm() : RealVectorValue(0, 0, 0);
127  RealVectorValue d_convective_term_d_u_comp =
128  _convective_term ? dConvecDUComp(comp) : RealVectorValue(0, 0, 0);
129  RealVectorValue viscous_term =
131  RealVectorValue d_viscous_term_d_u_comp =
133  RealVectorValue transient_term =
134  _transient_term ? timeDerivativeTerm() : RealVectorValue(0, 0, 0);
135  RealVectorValue d_transient_term_d_u_comp =
136  _transient_term ? dTimeDerivativeDUComp(comp) : RealVectorValue(0, 0, 0);
137 
138  return -1. / _rho[_qp] * tau() * _grad_test[_i][_qp] *
139  (d_convective_term_d_u_comp + d_viscous_term_d_u_comp + d_transient_term_d_u_comp) -
140  1. / _rho[_qp] * dTauDUComp(comp) * _grad_test[_i][_qp] *
141  (convective_term + viscous_term + transient_term + strongPressureTerm() +
142  gravityTerm() - RealVectorValue(_x_ffn.value(_t, _q_point[_qp]),
143  _y_ffn.value(_t, _q_point[_qp]),
144  _z_ffn.value(_t, _q_point[_qp])));
145 }

Referenced by INSMass::computeQpOffDiagJacobian().

◆ computeQpPGResidual()

Real INSMass::computeQpPGResidual ( )
protectedvirtualinherited

Definition at line 57 of file INSMass.C.

58 {
59  RealVectorValue viscous_term =
61  RealVectorValue transient_term =
62  _transient_term ? timeDerivativeTerm() : RealVectorValue(0, 0, 0);
63  RealVectorValue convective_term = _convective_term ? convectiveTerm() : RealVectorValue(0, 0, 0);
64  Real r = -1. / _rho[_qp] * tau() * _grad_test[_i][_qp] *
65  (strongPressureTerm() + gravityTerm() + viscous_term + convective_term + transient_term -
66  RealVectorValue(_x_ffn.value(_t, _q_point[_qp]),
67  _y_ffn.value(_t, _q_point[_qp]),
68  _z_ffn.value(_t, _q_point[_qp])));
69 
70  return r;
71 }

Referenced by INSMass::computeQpResidual().

◆ computeQpResidual()

Real INSMassRZ::computeQpResidual ( )
overrideprotectedvirtual

Reimplemented from INSMass.

Definition at line 77 of file INSMassRZ.C.

78 {
79  // Base class residual contribution
80  Real res_base = INSMass::computeQpResidual();
81 
82  // The radial coordinate value.
83  const Real r = _q_point[_qp](0);
84 
85  // The sign of this term is multiplied by -1 here.
86  res_base -= _u_vel[_qp] / r * _test[_i][_qp];
87 
88  return res_base;
89 }

◆ convectiveTerm()

RealVectorValue INSBase::convectiveTerm ( )
protectedvirtualinherited

Definition at line 96 of file INSBase.C.

97 {
98  RealVectorValue U(_u_vel[_qp], _v_vel[_qp], _w_vel[_qp]);
99  return _rho[_qp] *
100  RealVectorValue(U * _grad_u_vel[_qp], U * _grad_v_vel[_qp], U * _grad_w_vel[_qp]);
101 }

Referenced by INSMomentumBase::computeQpPGJacobian(), INSMass::computeQpPGOffDiagJacobian(), INSMass::computeQpPGResidual(), INSMomentumBase::computeQpPGResidual(), and INSMomentumBase::computeQpResidual().

◆ dConvecDUComp()

RealVectorValue INSBase::dConvecDUComp ( unsigned  comp)
protectedvirtualinherited

Definition at line 104 of file INSBase.C.

105 {
106  RealVectorValue U(_u_vel[_qp], _v_vel[_qp], _w_vel[_qp]);
107  RealVectorValue d_U_d_comp(0, 0, 0);
108  d_U_d_comp(comp) = _phi[_j][_qp];
109 
110  RealVectorValue convective_term = _rho[_qp] * RealVectorValue(d_U_d_comp * _grad_u_vel[_qp],
111  d_U_d_comp * _grad_v_vel[_qp],
112  d_U_d_comp * _grad_w_vel[_qp]);
113  convective_term(comp) += _rho[_qp] * U * _grad_phi[_j][_qp];
114 
115  return convective_term;
116 }

Referenced by INSMomentumBase::computeQpJacobian(), INSMomentumBase::computeQpOffDiagJacobian(), INSMomentumBase::computeQpPGJacobian(), and INSMass::computeQpPGOffDiagJacobian().

◆ dStrongPressureDPressure()

RealVectorValue INSBase::dStrongPressureDPressure ( )
protectedvirtualinherited

Definition at line 226 of file INSBase.C.

227 {
228  return _grad_phi[_j][_qp];
229 }

Referenced by INSMomentumBase::computeQpOffDiagJacobian(), and INSMass::computeQpPGJacobian().

◆ dStrongViscDUCompLaplace()

RealVectorValue INSMassRZ::dStrongViscDUCompLaplace ( unsigned  comp)
overrideprotectedvirtual

Reimplemented from INSBase.

Definition at line 38 of file INSMassRZ.C.

39 {
40  const Real & r = _q_point[_qp](0);
41  RealVectorValue add_jac(0, 0, 0);
42  if (comp == 0)
43  add_jac(0) = _mu[_qp] * (_phi[_j][_qp] / (r * r) - _grad_phi[_j][_qp](0) / r);
44  else if (comp == 1)
45  add_jac(1) = -_mu[_qp] * _grad_phi[_j][_qp](0) / r;
46 
47  return INSBase::dStrongViscDUCompLaplace(comp) + add_jac;
48 }

◆ dStrongViscDUCompTraction()

RealVectorValue INSMassRZ::dStrongViscDUCompTraction ( unsigned  comp)
overrideprotectedvirtual

Reimplemented from INSBase.

Definition at line 61 of file INSMassRZ.C.

62 {
63  const Real & r = _q_point[_qp](0);
64  RealVectorValue add_jac(0, 0, 0);
65  if (comp == 0)
66  {
67  add_jac(0) = 2. * _mu[_qp] * (_phi[_j][_qp] / (r * r) - _grad_phi[_j][_qp](0) / r);
68  add_jac(1) = -_mu[_qp] / r * _grad_phi[_j][_qp](1);
69  }
70  else if (comp == 1)
71  add_jac(1) = -_mu[_qp] * _grad_phi[_j][_qp](0) / r;
72 
73  return INSBase::dStrongViscDUCompTraction(comp) + add_jac;
74 }

◆ dTauDUComp()

Real INSBase::dTauDUComp ( unsigned  comp)
protectedvirtualinherited

Definition at line 298 of file INSBase.C.

299 {
300  Real nu = _mu[_qp] / _rho[_qp];
301  RealVectorValue U(_u_vel[_qp], _v_vel[_qp], _w_vel[_qp]);
302  Real h = _current_elem->hmax();
303  Real transient_part = _transient_term ? 4. / (_dt * _dt) : 0.;
304  return -_alpha / 2. * std::pow(transient_part + (2. * U.norm() / h) * (2. * U.norm() / h) +
305  9. * (4. * nu / (h * h)) * (4. * nu / (h * h)),
306  -1.5) *
307  2. * (2. * U.norm() / h) * 2. / h * U(comp) * _phi[_j][_qp] /
308  (U.norm() + std::numeric_limits<double>::epsilon());
309 }

Referenced by INSMomentumBase::computeQpPGJacobian(), and INSMass::computeQpPGOffDiagJacobian().

◆ dTimeDerivativeDUComp()

RealVectorValue INSBase::dTimeDerivativeDUComp ( unsigned  comp)
protectedvirtualinherited

Definition at line 250 of file INSBase.C.

251 {
252  Real base = _rho[_qp] * _phi[_j][_qp];
253  switch (comp)
254  {
255  case 0:
256  return RealVectorValue(base * _d_u_vel_dot_du[_qp], 0, 0);
257 
258  case 1:
259  return RealVectorValue(0, base * _d_v_vel_dot_dv[_qp], 0);
260 
261  case 2:
262  return RealVectorValue(0, 0, base * _d_w_vel_dot_dw[_qp]);
263 
264  default:
265  mooseError("comp must be 0, 1, or 2");
266  }
267 }

Referenced by INSMomentumBase::computeQpPGJacobian(), and INSMass::computeQpPGOffDiagJacobian().

◆ dWeakPressureDPressure()

Real INSBase::dWeakPressureDPressure ( )
protectedvirtualinherited

Definition at line 232 of file INSBase.C.

233 {
234  return -_phi[_j][_qp];
235 }

Referenced by INSMomentumBase::computeQpOffDiagJacobian().

◆ dWeakViscDUCompLaplace()

RealVectorValue INSBase::dWeakViscDUCompLaplace ( )
protectedvirtualinherited

Definition at line 202 of file INSBase.C.

203 {
204  return _mu[_qp] * _grad_phi[_j][_qp];
205 }

◆ dWeakViscDUCompTraction()

RealVectorValue INSBase::dWeakViscDUCompTraction ( )
protectedvirtualinherited

Definition at line 208 of file INSBase.C.

209 {
210  return _mu[_qp] * _grad_phi[_j][_qp];
211 }

◆ gravityTerm()

RealVectorValue INSBase::gravityTerm ( )
protectedvirtualinherited

◆ strongPressureTerm()

RealVectorValue INSBase::strongPressureTerm ( )
protectedvirtualinherited

◆ strongViscousTermLaplace()

RealVectorValue INSMassRZ::strongViscousTermLaplace ( )
overrideprotectedvirtual

Reimplemented from INSBase.

Definition at line 28 of file INSMassRZ.C.

29 {
30  const Real & r = _q_point[_qp](0);
32  RealVectorValue(_mu[_qp] * (_u_vel[_qp] / (r * r) - _grad_u_vel[_qp](0) / r),
33  -_mu[_qp] * _grad_v_vel[_qp](0) / r,
34  0);
35 }

◆ strongViscousTermTraction()

RealVectorValue INSMassRZ::strongViscousTermTraction ( )
overrideprotectedvirtual

Reimplemented from INSBase.

Definition at line 51 of file INSMassRZ.C.

52 {
53  const Real & r = _q_point[_qp](0);
55  RealVectorValue(2. * _mu[_qp] * (_u_vel[_qp] / (r * r) - _grad_u_vel[_qp](0) / r),
56  -_mu[_qp] / r * (_grad_v_vel[_qp](0) + _grad_u_vel[_qp](1)),
57  0);
58 }

◆ tau()

Real INSBase::tau ( )
protectedvirtualinherited

Definition at line 270 of file INSBase.C.

271 {
272  Real nu = _mu[_qp] / _rho[_qp];
273  RealVectorValue U(_u_vel[_qp], _v_vel[_qp], _w_vel[_qp]);
274  Real h = _current_elem->hmax();
275  Real transient_part = _transient_term ? 4. / (_dt * _dt) : 0.;
276  return _alpha / std::sqrt(transient_part + (2. * U.norm() / h) * (2. * U.norm() / h) +
277  9. * (4. * nu / (h * h)) * (4. * nu / (h * h)));
278 }

Referenced by Advection::computeQpJacobian(), INSMomentumBase::computeQpOffDiagJacobian(), INSMass::computeQpPGJacobian(), INSMomentumBase::computeQpPGJacobian(), INSMass::computeQpPGOffDiagJacobian(), INSMass::computeQpPGResidual(), INSMomentumBase::computeQpPGResidual(), and Advection::computeQpResidual().

◆ tauNodal()

Real INSBase::tauNodal ( )
protectedvirtualinherited

Provides tau which yields superconvergence for 1D advection-diffusion.

Definition at line 281 of file INSBase.C.

282 {
283  Real nu = _mu[_qp] / _rho[_qp];
284  RealVectorValue U(_u_vel[_qp], _v_vel[_qp], _w_vel[_qp]);
285  Real h = _current_elem->hmax();
286  Real xi;
287  if (nu < std::numeric_limits<Real>::epsilon())
288  xi = 1;
289  else
290  {
291  Real alpha = U.norm() * h / (2 * nu);
292  xi = 1. / std::tanh(alpha) - 1. / alpha;
293  }
294  return h / (2 * U.norm()) * xi;
295 }

Referenced by Advection::computeQpJacobian(), and Advection::computeQpResidual().

◆ timeDerivativeTerm()

RealVectorValue INSBase::timeDerivativeTerm ( )
protectedvirtualinherited

Definition at line 244 of file INSBase.C.

245 {
246  return _rho[_qp] * RealVectorValue(_u_vel_dot[_qp], _v_vel_dot[_qp], _w_vel_dot[_qp]);
247 }

Referenced by INSMomentumBase::computeQpPGJacobian(), INSMass::computeQpPGOffDiagJacobian(), INSMass::computeQpPGResidual(), and INSMomentumBase::computeQpPGResidual().

◆ weakPressureTerm()

Real INSBase::weakPressureTerm ( )
protectedvirtualinherited

Definition at line 220 of file INSBase.C.

221 {
222  return -_p[_qp];
223 }

Referenced by INSMomentumBase::computeQpResidual().

◆ weakViscousTermLaplace()

RealVectorValue INSBase::weakViscousTermLaplace ( unsigned  comp)
protectedvirtualinherited

Definition at line 155 of file INSBase.C.

156 {
157  switch (comp)
158  {
159  case 0:
160  return _mu[_qp] * _grad_u_vel[_qp];
161 
162  case 1:
163  return _mu[_qp] * _grad_v_vel[_qp];
164 
165  case 2:
166  return _mu[_qp] * _grad_w_vel[_qp];
167 
168  default:
169  return _zero[_qp];
170  }
171 }

◆ weakViscousTermTraction()

RealVectorValue INSBase::weakViscousTermTraction ( unsigned  comp)
protectedvirtualinherited

Definition at line 174 of file INSBase.C.

175 {
176  switch (comp)
177  {
178  case 0:
179  {
180  RealVectorValue transpose(_grad_u_vel[_qp](0), _grad_v_vel[_qp](0), _grad_w_vel[_qp](0));
181  return _mu[_qp] * _grad_u_vel[_qp] + _mu[_qp] * transpose;
182  }
183 
184  case 1:
185  {
186  RealVectorValue transpose(_grad_u_vel[_qp](1), _grad_v_vel[_qp](1), _grad_w_vel[_qp](1));
187  return _mu[_qp] * _grad_v_vel[_qp] + _mu[_qp] * transpose;
188  }
189 
190  case 2:
191  {
192  RealVectorValue transpose(_grad_u_vel[_qp](2), _grad_v_vel[_qp](2), _grad_w_vel[_qp](2));
193  return _mu[_qp] * _grad_w_vel[_qp] + _mu[_qp] * transpose;
194  }
195 
196  default:
197  return _zero[_qp];
198  }
199 }

Member Data Documentation

◆ _alpha

const Real& INSBase::_alpha
protectedinherited

Definition at line 107 of file INSBase.h.

Referenced by INSBase::dTauDUComp(), and INSBase::tau().

◆ _convective_term

bool INSBase::_convective_term
protectedinherited

◆ _d_u_vel_dot_du

const VariableValue& INSBase::_d_u_vel_dot_du
protectedinherited

Definition at line 91 of file INSBase.h.

Referenced by INSBase::dTimeDerivativeDUComp().

◆ _d_v_vel_dot_dv

const VariableValue& INSBase::_d_v_vel_dot_dv
protectedinherited

Definition at line 92 of file INSBase.h.

Referenced by INSBase::dTimeDerivativeDUComp().

◆ _d_w_vel_dot_dw

const VariableValue& INSBase::_d_w_vel_dot_dw
protectedinherited

Definition at line 93 of file INSBase.h.

Referenced by INSBase::dTimeDerivativeDUComp().

◆ _grad_p

const VariableGradient& INSBase::_grad_p
protectedinherited

Definition at line 78 of file INSBase.h.

Referenced by INSBase::strongPressureTerm().

◆ _grad_u_vel

const VariableGradient& INSBase::_grad_u_vel
protectedinherited

◆ _grad_v_vel

const VariableGradient& INSBase::_grad_v_vel
protectedinherited

◆ _grad_w_vel

const VariableGradient& INSBase::_grad_w_vel
protectedinherited

◆ _gravity

RealVectorValue INSBase::_gravity
protectedinherited

Definition at line 101 of file INSBase.h.

Referenced by INSBase::gravityTerm().

◆ _laplace

bool INSBase::_laplace
protectedinherited

◆ _mu

const MaterialProperty<Real>& INSBase::_mu
protectedinherited

◆ _p

const VariableValue& INSBase::_p
protectedinherited

◆ _p_var_number

unsigned INSBase::_p_var_number
protectedinherited

◆ _pspg

bool INSMass::_pspg
protectedinherited

◆ _rho

const MaterialProperty<Real>& INSBase::_rho
protectedinherited

◆ _second_phi

const VariablePhiSecond& INSBase::_second_phi
protectedinherited

second derivatives of the shape function

Definition at line 66 of file INSBase.h.

Referenced by INSBase::dStrongViscDUCompLaplace(), and INSBase::dStrongViscDUCompTraction().

◆ _second_u_vel

const VariableSecond& INSBase::_second_u_vel
protectedinherited

◆ _second_v_vel

const VariableSecond& INSBase::_second_v_vel
protectedinherited

◆ _second_w_vel

const VariableSecond& INSBase::_second_w_vel
protectedinherited

◆ _transient_term

bool INSBase::_transient_term
protectedinherited

◆ _u_vel

const VariableValue& INSBase::_u_vel
protectedinherited

◆ _u_vel_dot

const VariableValue& INSBase::_u_vel_dot
protectedinherited

Definition at line 86 of file INSBase.h.

Referenced by INSBase::timeDerivativeTerm().

◆ _u_vel_var_number

unsigned INSBase::_u_vel_var_number
protectedinherited

◆ _v_vel

const VariableValue& INSBase::_v_vel
protectedinherited

◆ _v_vel_dot

const VariableValue& INSBase::_v_vel_dot
protectedinherited

Definition at line 87 of file INSBase.h.

Referenced by INSBase::timeDerivativeTerm().

◆ _v_vel_var_number

unsigned INSBase::_v_vel_var_number
protectedinherited

◆ _w_vel

const VariableValue& INSBase::_w_vel
protectedinherited

◆ _w_vel_dot

const VariableValue& INSBase::_w_vel_dot
protectedinherited

Definition at line 88 of file INSBase.h.

Referenced by INSBase::timeDerivativeTerm().

◆ _w_vel_var_number

unsigned INSBase::_w_vel_var_number
protectedinherited

◆ _x_ffn

const Function& INSMass::_x_ffn
protectedinherited

Definition at line 42 of file INSMass.h.

Referenced by INSMass::computeQpPGOffDiagJacobian(), and INSMass::computeQpPGResidual().

◆ _y_ffn

const Function& INSMass::_y_ffn
protectedinherited

Definition at line 43 of file INSMass.h.

Referenced by INSMass::computeQpPGOffDiagJacobian(), and INSMass::computeQpPGResidual().

◆ _z_ffn

const Function& INSMass::_z_ffn
protectedinherited

Definition at line 44 of file INSMass.h.

Referenced by INSMass::computeQpPGOffDiagJacobian(), and INSMass::computeQpPGResidual().


The documentation for this class was generated from the following files:
INSBase::gravityTerm
virtual RealVectorValue gravityTerm()
Definition: INSBase.C:238
INSBase::strongPressureTerm
virtual RealVectorValue strongPressureTerm()
Definition: INSBase.C:214
INSMass::INSMass
INSMass(const InputParameters &parameters)
Definition: INSMass.C:32
INSMass::_z_ffn
const Function & _z_ffn
Definition: INSMass.h:44
INSBase::_d_v_vel_dot_dv
const VariableValue & _d_v_vel_dot_dv
Definition: INSBase.h:92
INSBase::_v_vel
const VariableValue & _v_vel
Definition: INSBase.h:70
INSBase::_convective_term
bool _convective_term
Definition: INSBase.h:109
INSBase::_gravity
RealVectorValue _gravity
Definition: INSBase.h:101
pow
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Definition: ExpressionBuilder.h:673
INSBase::_w_vel_dot
const VariableValue & _w_vel_dot
Definition: INSBase.h:88
INSBase::dStrongPressureDPressure
virtual RealVectorValue dStrongPressureDPressure()
Definition: INSBase.C:226
INSBase::_d_w_vel_dot_dw
const VariableValue & _d_w_vel_dot_dw
Definition: INSBase.h:93
INSBase::dStrongViscDUCompLaplace
virtual RealVectorValue dStrongViscDUCompLaplace(unsigned comp)
Definition: INSBase.C:134
INSBase::dTimeDerivativeDUComp
virtual RealVectorValue dTimeDerivativeDUComp(unsigned comp)
Definition: INSBase.C:250
INSBase::strongViscousTermLaplace
virtual RealVectorValue strongViscousTermLaplace()
Definition: INSBase.C:119
INSBase::_v_vel_dot
const VariableValue & _v_vel_dot
Definition: INSBase.h:87
INSBase::strongViscousTermTraction
virtual RealVectorValue strongViscousTermTraction()
Definition: INSBase.C:126
INSMass::computeQpResidual
virtual Real computeQpResidual()
Definition: INSMass.C:43
INSMass::computeQpPGJacobian
virtual Real computeQpPGJacobian()
Definition: INSMass.C:87
INSBase::tau
virtual Real tau()
Definition: INSBase.C:270
INSBase::_u_vel_var_number
unsigned _u_vel_var_number
Definition: INSBase.h:96
INSBase::_d_u_vel_dot_du
const VariableValue & _d_u_vel_dot_du
Definition: INSBase.h:91
INSBase::_grad_w_vel
const VariableGradient & _grad_w_vel
Definition: INSBase.h:77
INSBase::_p
const VariableValue & _p
Definition: INSBase.h:72
INSBase::convectiveTerm
virtual RealVectorValue convectiveTerm()
Definition: INSBase.C:96
INSMass::computeQpOffDiagJacobian
virtual Real computeQpOffDiagJacobian(unsigned jvar)
Definition: INSMass.C:93
INSBase::_rho
const MaterialProperty< Real > & _rho
Definition: INSBase.h:105
INSBase::dConvecDUComp
virtual RealVectorValue dConvecDUComp(unsigned comp)
Definition: INSBase.C:104
INSBase::_grad_p
const VariableGradient & _grad_p
Definition: INSBase.h:78
INSBase::_grad_u_vel
const VariableGradient & _grad_u_vel
Definition: INSBase.h:75
INSBase::_alpha
const Real & _alpha
Definition: INSBase.h:107
INSBase::_u_vel
const VariableValue & _u_vel
Definition: INSBase.h:69
INSMass::_pspg
bool _pspg
Definition: INSMass.h:41
INSBase::dStrongViscDUCompTraction
virtual RealVectorValue dStrongViscDUCompTraction(unsigned comp)
Definition: INSBase.C:143
INSBase::_grad_v_vel
const VariableGradient & _grad_v_vel
Definition: INSBase.h:76
INSBase::_laplace
bool _laplace
Definition: INSBase.h:108
INSBase::_w_vel
const VariableValue & _w_vel
Definition: INSBase.h:71
INSMass::_x_ffn
const Function & _x_ffn
Definition: INSMass.h:42
INSBase::_transient_term
bool _transient_term
Definition: INSBase.h:110
INSMass::_y_ffn
const Function & _y_ffn
Definition: INSMass.h:43
INSBase::_mu
const MaterialProperty< Real > & _mu
Definition: INSBase.h:104
INSBase::_u_vel_dot
const VariableValue & _u_vel_dot
Definition: INSBase.h:86
INSBase::timeDerivativeTerm
virtual RealVectorValue timeDerivativeTerm()
Definition: INSBase.C:244
INSBase::dTauDUComp
virtual Real dTauDUComp(unsigned comp)
Definition: INSBase.C:298