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

Adds Velocity penalty adjoint stabilization source term. More...

#include <velocity_penalty_adjoint_stab.h>

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

Public Member Functions

 VelocityPenaltyAdjointStabilization (const std::string &physics_name, const GetPot &input)
 
 ~VelocityPenaltyAdjointStabilization ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. 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...
 
- Public Member Functions inherited from GRINS::VelocityPenaltyBase< Viscosity >
 VelocityPenaltyBase (const std::string &physics_name, const GetPot &input)
 
 ~VelocityPenaltyBase ()
 
void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
bool compute_force (const libMesh::Point &point, const AssemblyContext &context, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
 
- 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 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 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 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 mass_residual (bool, AssemblyContext &)
 Mass 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

IncompressibleNavierStokesStabilizationHelper _stab_helper
 
- Protected Attributes inherited from GRINS::VelocityPenaltyBase< Viscosity >
std::string base_physics_name
 
bool _quadratic_scaling
 
libMesh::UniquePtr< libMesh::FEMFunctionBase< libMesh::Number > > normal_vector_function
 
libMesh::UniquePtr< libMesh::FEMFunctionBase< libMesh::Number > > base_velocity_function
 
- 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

 VelocityPenaltyAdjointStabilization ()
 

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::VelocityPenaltyAdjointStabilization< Viscosity >

Adds Velocity penalty adjoint stabilization source term.

This class implements the adjoint stabilization term for the VelocityPenalty Physics. Intended to be used with IncompressibleNavierStokesAdjointStabilization.

Definition at line 42 of file velocity_penalty_adjoint_stab.h.

Constructor & Destructor Documentation

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

Definition at line 42 of file velocity_penalty_adjoint_stab.C.

43  : VelocityPenaltyBase<Mu>(physics_name,input),
44  _stab_helper( physics_name+"StabHelper", input )
45  {
46  }
IncompressibleNavierStokesStabilizationHelper _stab_helper

Definition at line 49 of file velocity_penalty_adjoint_stab.C.

50  {
51  return;
52  }
template<class Viscosity >
GRINS::VelocityPenaltyAdjointStabilization< Viscosity >::VelocityPenaltyAdjointStabilization ( )
private

Member Function Documentation

template<class Mu >
void GRINS::VelocityPenaltyAdjointStabilization< Mu >::element_constraint ( bool  ,
AssemblyContext  
)
virtual

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 235 of file velocity_penalty_adjoint_stab.C.

236  {
237  // The number of local degrees of freedom in each variable.
238  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
239  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
240 
241  // Element Jacobian * quadrature weights for interior integration.
242  const std::vector<libMesh::Real> &JxW =
243  context.get_element_fe(this->_flow_vars.u())->get_JxW();
244 
245  const std::vector<libMesh::Point>& u_qpoint =
246  context.get_element_fe(this->_flow_vars.u())->get_xyz();
247 
248  const std::vector<std::vector<libMesh::Real> >& u_phi =
249  context.get_element_fe(this->_flow_vars.u())->get_phi();
250 
251  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
252  context.get_element_fe(this->_press_var.p())->get_dphi();
253 
254  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
255 
256  libMesh::DenseSubMatrix<libMesh::Number> &Kpu =
257  context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.u()); // J_{pu}
258  libMesh::DenseSubMatrix<libMesh::Number> &Kpv =
259  context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.v()); // J_{pv}
260  libMesh::DenseSubMatrix<libMesh::Number> *Kpw = NULL;
261 
262  if(this->_flow_vars.dim() == 3)
263  {
264  Kpw = &context.get_elem_jacobian
265  (this->_press_var.p(), this->_flow_vars.w()); // J_{pw}
266  }
267 
268  // Now we will build the element Jacobian and residual.
269  // Constructing the residual requires the solution and its
270  // gradient from the previous timestep. This must be
271  // calculated at each quadrature point by summing the
272  // solution degree-of-freedom values by the appropriate
273  // weight functions.
274  unsigned int n_qpoints = context.get_element_qrule().n_points();
275 
276  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
277 
278  for (unsigned int qp=0; qp != n_qpoints; qp++)
279  {
280  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
281  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
282 
283  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
284  context.interior_value( this->_flow_vars.v(), qp ) );
285  if( this->_flow_vars.dim() == 3 )
286  {
287  U(2) = context.interior_value( this->_flow_vars.w(), qp );
288  }
289 
290  // Compute the viscosity at this qp
291  libMesh::Real mu_qp = this->_mu(context, qp);
292 
293  libMesh::Real tau_M;
294  libMesh::Real d_tau_M_d_rho;
295  libMesh::Gradient d_tau_M_dU;
296 
297  if (compute_jacobian)
299  ( context, qp, g, G, this->_rho, U, mu_qp,
300  tau_M, d_tau_M_d_rho, d_tau_M_dU,
301  this->_is_steady );
302  else
303  tau_M = this->_stab_helper.compute_tau_momentum
304  ( context, qp, g, G, this->_rho, U, mu_qp,
305  this->_is_steady );
306 
307  libMesh::NumberVectorValue F;
308  libMesh::NumberTensorValue dFdU;
309  libMesh::NumberTensorValue* dFdU_ptr =
310  compute_jacobian ? &dFdU : NULL;
311  if (!this->compute_force(u_qpoint[qp], context, U, F, dFdU_ptr))
312  continue;
313 
314  // First, an i-loop over the velocity degrees of freedom.
315  // We know that n_u_dofs == n_v_dofs so we can compute contributions
316  // for both at the same time.
317  for (unsigned int i=0; i != n_p_dofs; i++)
318  {
319  Fp(i) += -tau_M*F*p_dphi[i][qp]*JxW[qp];
320 
321  if (compute_jacobian)
322  {
323  for (unsigned int j=0; j != n_u_dofs; ++j)
324  {
325  Kpu(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
326  Kpv(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
327  for (unsigned int d=0; d != 3; ++d)
328  {
329  Kpu(i,j) += -tau_M*dFdU(d,0)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
330  Kpv(i,j) += -tau_M*dFdU(d,1)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
331  }
332  }
333  if( this->_flow_vars.dim() == 3 )
334  for (unsigned int j=0; j != n_u_dofs; ++j)
335  {
336  (*Kpw)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
337  for (unsigned int d=0; d != 3; ++d)
338  {
339  (*Kpw)(i,j) += -tau_M*dFdU(d,2)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
340  }
341  }
342  }
343  }
344  } // End quadrature loop
345  }
void compute_tau_momentum_and_derivs(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real T, libMesh::Real &tau_M, libMesh::Real &d_tau_M_d_rho, libMesh::Gradient &d_tau_M_d_U, bool is_steady) const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::Number _rho
Material parameters, read from input.
bool compute_force(const libMesh::Point &point, const AssemblyContext &context, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
VariableIndex p() const
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
template<class Mu >
void GRINS::VelocityPenaltyAdjointStabilization< 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 velocity_penalty_adjoint_stab.C.

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()).size();
74 
75  // Element Jacobian * quadrature weights for interior integration.
76  const std::vector<libMesh::Real> &JxW =
77  context.get_element_fe(this->_flow_vars.u())->get_JxW();
78 
79  const std::vector<libMesh::Point>& u_qpoint =
80  context.get_element_fe(this->_flow_vars.u())->get_xyz();
81 
82  const std::vector<std::vector<libMesh::Real> >& u_phi =
83  context.get_element_fe(this->_flow_vars.u())->get_phi();
84 
85  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
86  context.get_element_fe(this->_flow_vars.u())->get_dphi();
87 
88  const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
89  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
90 
91  // Get residuals and jacobians
92  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
93  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
94  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
95 
96  libMesh::DenseSubMatrix<libMesh::Number> &Kuu =
97  context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // J_{uu}
98  libMesh::DenseSubMatrix<libMesh::Number> &Kuv =
99  context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // J_{uv}
100  libMesh::DenseSubMatrix<libMesh::Number> &Kvu =
101  context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // J_{vu}
102  libMesh::DenseSubMatrix<libMesh::Number> &Kvv =
103  context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // J_{vv}
104 
105  libMesh::DenseSubMatrix<libMesh::Number> *Kuw = NULL;
106  libMesh::DenseSubMatrix<libMesh::Number> *Kvw = NULL;
107  libMesh::DenseSubMatrix<libMesh::Number> *Kwu = NULL;
108  libMesh::DenseSubMatrix<libMesh::Number> *Kwv = NULL;
109  libMesh::DenseSubMatrix<libMesh::Number> *Kww = NULL;
110 
111  if(this->_flow_vars.dim() == 3)
112  {
113  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
114  Kuw = &context.get_elem_jacobian
115  (this->_flow_vars.u(), this->_flow_vars.w()); // J_{uw}
116  Kvw = &context.get_elem_jacobian
117  (this->_flow_vars.v(), this->_flow_vars.w()); // J_{vw}
118  Kwu = &context.get_elem_jacobian
119  (this->_flow_vars.w(), this->_flow_vars.u()); // J_{wu}
120  Kwv = &context.get_elem_jacobian
121  (this->_flow_vars.w(), this->_flow_vars.v()); // J_{wv}
122  Kww = &context.get_elem_jacobian
123  (this->_flow_vars.w(), this->_flow_vars.w()); // J_{ww}
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  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
135 
136  for (unsigned int qp=0; qp != n_qpoints; qp++)
137  {
138  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
139  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
140 
141  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
142  context.interior_value( this->_flow_vars.v(), qp ) );
143  if( this->_flow_vars.dim() == 3 )
144  {
145  U(2) = context.interior_value( this->_flow_vars.w(), qp );
146  }
147 
148  // Compute the viscosity at this qp
149  libMesh::Real mu_qp = this->_mu(context, qp);
150 
151  libMesh::Real tau_M;
152  libMesh::Real d_tau_M_d_rho;
153  libMesh::Gradient d_tau_M_dU;
154 
155  if (compute_jacobian)
157  ( context, qp, g, G, this->_rho, U, mu_qp,
158  tau_M, d_tau_M_d_rho, d_tau_M_dU,
159  this->_is_steady );
160  else
161  tau_M = this->_stab_helper.compute_tau_momentum
162  ( context, qp, g, G, this->_rho, U, mu_qp,
163  this->_is_steady );
164 
165  libMesh::NumberVectorValue F;
166  libMesh::NumberTensorValue dFdU;
167  libMesh::NumberTensorValue* dFdU_ptr =
168  compute_jacobian ? &dFdU : NULL;
169  if (!this->compute_force(u_qpoint[qp], context, U, F, dFdU_ptr))
170  continue;
171 
172  for (unsigned int i=0; i != n_u_dofs; i++)
173  {
174  libMesh::Real test_func = this->_rho*U*u_gradphi[i][qp] +
175  mu_qp*( u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2) );
176  Fu(i) += tau_M*F(0)*test_func*JxW[qp];
177 
178  Fv(i) += tau_M*F(1)*test_func*JxW[qp];
179 
180  if (this->_flow_vars.dim() == 3)
181  {
182  (*Fw)(i) += tau_M*F(2)*test_func*JxW[qp];
183  }
184 
185  if (compute_jacobian)
186  {
187  libMesh::Gradient d_test_func_dU = this->_rho*u_gradphi[i][qp];
188 
189  for (unsigned int j=0; j != n_u_dofs; ++j)
190  {
191  Kuu(i,j) += tau_M*F(0)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
192  Kuu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
193  Kuu(i,j) += tau_M*dFdU(0,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
194  Kuv(i,j) += tau_M*F(0)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
195  Kuv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
196  Kuv(i,j) += tau_M*dFdU(0,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
197  Kvu(i,j) += tau_M*F(1)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
198  Kvu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
199  Kvu(i,j) += tau_M*dFdU(1,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
200  Kvv(i,j) += tau_M*F(1)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
201  Kvv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
202  Kvv(i,j) += tau_M*dFdU(1,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
203  }
204 
205  if (this->_flow_vars.dim() == 3)
206  {
207  for (unsigned int j=0; j != n_u_dofs; ++j)
208  {
209  (*Kuw)(i,j) += tau_M*F(0)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
210  (*Kuw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
211  (*Kuw)(i,j) += tau_M*dFdU(0,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
212  (*Kvw)(i,j) += tau_M*F(1)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
213  (*Kvw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
214  (*Kvw)(i,j) += tau_M*dFdU(1,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
215  (*Kwu)(i,j) += tau_M*F(2)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
216  (*Kwu)(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
217  (*Kwu)(i,j) += tau_M*dFdU(2,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
218  (*Kwv)(i,j) += tau_M*F(2)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
219  (*Kwv)(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
220  (*Kwv)(i,j) += tau_M*dFdU(2,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
221  (*Kww)(i,j) += tau_M*F(2)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
222  (*Kww)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
223  (*Kww)(i,j) += tau_M*dFdU(2,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
224  }
225  }
226 
227  } // End compute_jacobian check
228 
229  } // End i dof loop
230  } // End quadrature loop
231  }
void compute_tau_momentum_and_derivs(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real T, libMesh::Real &tau_M, libMesh::Real &d_tau_M_d_rho, libMesh::Gradient &d_tau_M_d_U, bool is_steady) const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::Number _rho
Material parameters, read from input.
bool compute_force(const libMesh::Point &point, const AssemblyContext &context, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
unsigned int dim() const
Number of components.
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
template<class Mu >
void GRINS::VelocityPenaltyAdjointStabilization< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

Definition at line 55 of file velocity_penalty_adjoint_stab.C.

56  {
57  context.get_element_fe(this->_press_var.p())->get_dphi();
58 
59  context.get_element_fe(this->_flow_vars.u())->get_xyz();
60  context.get_element_fe(this->_flow_vars.u())->get_phi();
61  context.get_element_fe(this->_flow_vars.u())->get_dphi();
62  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
63 
64  return;
65  }
VariableIndex p() const

Member Data Documentation

template<class Viscosity >
IncompressibleNavierStokesStabilizationHelper GRINS::VelocityPenaltyAdjointStabilization< Viscosity >::_stab_helper
protected

Definition at line 60 of file velocity_penalty_adjoint_stab.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