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

Physics class for Velocity Drag. More...

#include <velocity_drag_adjoint_stab.h>

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

Public Member Functions

 VelocityDragAdjointStabilization (const std::string &physics_name, const GetPot &input)
 
 ~VelocityDragAdjointStabilization ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. 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...
 
- Public Member Functions inherited from GRINS::VelocityDragBase< Viscosity >
 VelocityDragBase (const std::string &physics_name, const GetPot &input)
 
 ~VelocityDragBase ()
 
bool compute_force (const libMesh::Point &point, const libMesh::Real time, 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 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 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 auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. More...
 
virtual void side_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void side_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void 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 mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void 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

IncompressibleNavierStokesStabilizationHelper _stab_helper
 
- Protected Attributes inherited from GRINS::VelocityDragBase< Viscosity >
libMesh::Real _exponent
 
libMesh::ParsedFunction< libMesh::Number > _coefficient
 
- 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

 VelocityDragAdjointStabilization ()
 

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

Physics class for Velocity Drag.

Definition at line 46 of file velocity_drag_adjoint_stab.h.

Constructor & Destructor Documentation

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

Definition at line 41 of file velocity_drag_adjoint_stab.C.

42  : VelocityDragBase<Mu>(physics_name, input),
43  _stab_helper( physics_name+"StabHelper", input )
44  {}
IncompressibleNavierStokesStabilizationHelper _stab_helper

Definition at line 47 of file velocity_drag_adjoint_stab.C.

48  {
49  return;
50  }
template<class Viscosity >
GRINS::VelocityDragAdjointStabilization< Viscosity >::VelocityDragAdjointStabilization ( )
private

Member Function Documentation

template<class Mu >
void GRINS::VelocityDragAdjointStabilization< 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 237 of file velocity_drag_adjoint_stab.C.

240  {
241 #ifdef GRINS_USE_GRVY_TIMERS
242  this->_timer->BeginTimer("VelocityDragAdjointStabilization::element_constraint");
243 #endif
244 
245  // The number of local degrees of freedom in each variable.
246  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
247  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
248 
249  // Element Jacobian * quadrature weights for interior integration.
250  const std::vector<libMesh::Real> &JxW =
251  context.get_element_fe(this->_flow_vars.u())->get_JxW();
252 
253  const std::vector<libMesh::Point>& u_qpoint =
254  context.get_element_fe(this->_flow_vars.u())->get_xyz();
255 
256  const std::vector<std::vector<libMesh::Real> >& u_phi =
257  context.get_element_fe(this->_flow_vars.u())->get_phi();
258 
259  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
260  context.get_element_fe(this->_press_var.p())->get_dphi();
261 
262  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
263 
264  libMesh::DenseSubMatrix<libMesh::Number> &Kpu =
265  context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.u()); // J_{pu}
266  libMesh::DenseSubMatrix<libMesh::Number> &Kpv =
267  context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.v()); // J_{pv}
268  libMesh::DenseSubMatrix<libMesh::Number> *Kpw = NULL;
269 
270  if(this->_dim == 3)
271  {
272  Kpw = &context.get_elem_jacobian
273  (this->_press_var.p(), this->_flow_vars.w()); // J_{pw}
274  }
275 
276  // Now we will build the element Jacobian and residual.
277  // Constructing the residual requires the solution and its
278  // gradient from the previous timestep. This must be
279  // calculated at each quadrature point by summing the
280  // solution degree-of-freedom values by the appropriate
281  // weight functions.
282  unsigned int n_qpoints = context.get_element_qrule().n_points();
283 
284  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
285 
286  for (unsigned int qp=0; qp != n_qpoints; qp++)
287  {
288  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
289  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
290 
291  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
292  context.interior_value( this->_flow_vars.v(), qp ) );
293  if( this->_dim == 3 )
294  {
295  U(2) = context.interior_value( this->_flow_vars.w(), qp );
296  }
297 
298  // Compute the viscosity at this qp
299  libMesh::Real mu_qp = this->_mu(context, qp);
300 
301  libMesh::Real tau_M;
302  libMesh::Real d_tau_M_d_rho;
303  libMesh::Gradient d_tau_M_dU;
304 
305  if (compute_jacobian)
307  ( context, qp, g, G, this->_rho, U, mu_qp,
308  tau_M, d_tau_M_d_rho, d_tau_M_dU,
309  this->_is_steady );
310  else
311  tau_M = this->_stab_helper.compute_tau_momentum
312  ( context, qp, g, G, this->_rho, U, mu_qp,
313  this->_is_steady );
314 
315  libMesh::NumberVectorValue F;
316  libMesh::NumberTensorValue dFdU;
317  libMesh::NumberTensorValue* dFdU_ptr =
318  compute_jacobian ? &dFdU : NULL;
319  if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
320  continue;
321 
322  // First, an i-loop over the velocity degrees of freedom.
323  // We know that n_u_dofs == n_v_dofs so we can compute contributions
324  // for both at the same time.
325  for (unsigned int i=0; i != n_p_dofs; i++)
326  {
327  Fp(i) += -tau_M*F*p_dphi[i][qp]*JxW[qp];
328 
329  if (compute_jacobian)
330  {
331  const libMesh::Real sol_deriv =
332  context.get_elem_solution_derivative();
333 
334  const libMesh::Real JxWxS = JxW[qp] * sol_deriv;
335 
336  const libMesh::Real fixed_deriv =
337  context.get_fixed_solution_derivative();
338 
339  const libMesh::Real JxWxF = JxW[qp] * fixed_deriv;
340 
341  for (unsigned int j=0; j != n_u_dofs; ++j)
342  {
343  Kpu(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*F*p_dphi[i][qp]*JxWxF;
344  Kpv(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*F*p_dphi[i][qp]*JxWxF;
345  for (unsigned int d=0; d != 3; ++d)
346  {
347  Kpu(i,j) += -tau_M*dFdU(d,0)*u_phi[j][qp]*p_dphi[i][qp](d)*JxWxS;
348  Kpv(i,j) += -tau_M*dFdU(d,1)*u_phi[j][qp]*p_dphi[i][qp](d)*JxWxS;
349  }
350  }
351  if( this->_dim == 3 )
352  for (unsigned int j=0; j != n_u_dofs; ++j)
353  {
354  (*Kpw)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*F*p_dphi[i][qp]*JxWxF;
355  for (unsigned int d=0; d != 3; ++d)
356  {
357  (*Kpw)(i,j) += -tau_M*dFdU(d,2)*u_phi[j][qp]*p_dphi[i][qp](d)*JxWxS;
358  }
359  }
360  }
361  }
362  } // End quadrature loop
363 
364 #ifdef GRINS_USE_GRVY_TIMERS
365  this->_timer->EndTimer("VelocityDragAdjointStabilization::element_constraint");
366 #endif
367 
368  return;
369  }
bool compute_force(const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
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:64
libMesh::Number _rho
Material parameters, read from input.
unsigned int _dim
Physical dimension of problem.
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:277
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
IncompressibleNavierStokesStabilizationHelper _stab_helper
template<class Mu >
void GRINS::VelocityDragAdjointStabilization< 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 67 of file velocity_drag_adjoint_stab.C.

70  {
71 #ifdef GRINS_USE_GRVY_TIMERS
72  this->_timer->BeginTimer("VelocityDragAdjointStabilization::element_time_derivative");
73 #endif
74 
75  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
76 
77  // Element Jacobian * quadrature weights for interior integration
78  const std::vector<libMesh::Real> &JxW = fe->get_JxW();
79 
80  // The shape functions at interior quadrature points.
81  const std::vector<std::vector<libMesh::Real> >& u_phi = fe->get_phi();
82 
83  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
84  context.get_element_fe(this->_flow_vars.u())->get_dphi();
85 
86  const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
87  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
88 
89  const std::vector<libMesh::Point>& u_qpoint = fe->get_xyz();
90 
91  // The number of local degrees of freedom in each variable
92  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
93 
94  // The subvectors and submatrices we need to fill:
95  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // R_{u},{u}
96  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // R_{u},{v}
97  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // R_{v},{u}
98  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // R_{v},{v}
99 
100  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
101  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
102  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
103  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
104  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
105 
106  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
107  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
108  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
109 
110  if( this->_dim == 3 )
111  {
112  Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w()); // R_{u},{w}
113  Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w()); // R_{v},{w}
114 
115  Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u()); // R_{w},{u}
116  Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v()); // R_{w},{v}
117  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
118  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
119  }
120 
121  unsigned int n_qpoints = context.get_element_qrule().n_points();
122 
123  for (unsigned int qp=0; qp != n_qpoints; qp++)
124  {
125  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
126  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
127 
128  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
129  context.interior_value( this->_flow_vars.v(), qp ) );
130  if( this->_dim == 3 )
131  {
132  U(2) = context.interior_value( this->_flow_vars.w(), qp );
133  }
134 
135  // Compute the viscosity at this qp
136  libMesh::Real mu_qp = this->_mu(context, qp);
137 
138  libMesh::Real tau_M;
139  libMesh::Real d_tau_M_d_rho;
140  libMesh::Gradient d_tau_M_dU;
141 
142  if (compute_jacobian)
144  ( context, qp, g, G, this->_rho, U, mu_qp,
145  tau_M, d_tau_M_d_rho, d_tau_M_dU,
146  this->_is_steady );
147  else
148  tau_M = this->_stab_helper.compute_tau_momentum
149  ( context, qp, g, G, this->_rho, U, mu_qp,
150  this->_is_steady );
151 
152  libMesh::NumberVectorValue F;
153  libMesh::NumberTensorValue dFdU;
154  libMesh::NumberTensorValue* dFdU_ptr =
155  compute_jacobian ? &dFdU : NULL;
156  if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
157  continue;
158 
159  for (unsigned int i=0; i != n_u_dofs; i++)
160  {
161  libMesh::Real test_func = this->_rho*U*u_gradphi[i][qp] +
162  mu_qp*( u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2) );
163  Fu(i) += tau_M*F(0)*test_func*JxW[qp];
164 
165  Fv(i) += tau_M*F(1)*test_func*JxW[qp];
166 
167  if (this->_dim == 3)
168  {
169  (*Fw)(i) += tau_M*F(2)*test_func*JxW[qp];
170  }
171 
172  if (compute_jacobian)
173  {
174  const libMesh::Real sol_deriv =
175  context.get_elem_solution_derivative();
176 
177  const libMesh::Real JxWxS = JxW[qp] * sol_deriv;
178 
179  const libMesh::Real fixed_deriv =
180  context.get_fixed_solution_derivative();
181 
182  const libMesh::Real JxWxF = JxW[qp] * fixed_deriv;
183 
184  libMesh::Gradient d_test_func_dU = this->_rho*u_gradphi[i][qp];
185 
186  for (unsigned int j=0; j != n_u_dofs; ++j)
187  {
188  Kuu(i,j) += tau_M*F(0)*d_test_func_dU(0)*u_phi[j][qp]*JxWxS;
189  Kuu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(0)*test_func*JxWxF;
190  Kuu(i,j) += tau_M*dFdU(0,0)*u_phi[j][qp]*test_func*JxWxS;
191  Kuv(i,j) += tau_M*F(0)*d_test_func_dU(1)*u_phi[j][qp]*JxWxS;
192  Kuv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(0)*test_func*JxWxF;
193  Kuv(i,j) += tau_M*dFdU(0,1)*u_phi[j][qp]*test_func*JxWxS;
194  Kvu(i,j) += tau_M*F(1)*d_test_func_dU(0)*u_phi[j][qp]*JxWxS;
195  Kvu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(1)*test_func*JxWxF;
196  Kvu(i,j) += tau_M*dFdU(1,0)*u_phi[j][qp]*test_func*JxWxS;
197  Kvv(i,j) += tau_M*F(1)*d_test_func_dU(1)*u_phi[j][qp]*JxWxS;
198  Kvv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(1)*test_func*JxWxF;
199  Kvv(i,j) += tau_M*dFdU(1,1)*u_phi[j][qp]*test_func*JxWxS;
200  }
201 
202  if (this->_dim == 3)
203  {
204  for (unsigned int j=0; j != n_u_dofs; ++j)
205  {
206  (*Kuw)(i,j) += tau_M*F(0)*d_test_func_dU(2)*u_phi[j][qp]*JxWxS;
207  (*Kuw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(0)*test_func*JxWxF;
208  (*Kuw)(i,j) += tau_M*dFdU(0,2)*u_phi[j][qp]*test_func*JxWxS;
209  (*Kvw)(i,j) += tau_M*F(1)*d_test_func_dU(2)*u_phi[j][qp]*JxWxS;
210  (*Kvw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(1)*test_func*JxWxF;
211  (*Kvw)(i,j) += tau_M*dFdU(1,2)*u_phi[j][qp]*test_func*JxWxS;
212  (*Kwu)(i,j) += tau_M*F(2)*d_test_func_dU(0)*u_phi[j][qp]*JxWxS;
213  (*Kwu)(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(2)*test_func*JxWxF;
214  (*Kwu)(i,j) += tau_M*dFdU(2,0)*u_phi[j][qp]*test_func*JxWxS;
215  (*Kwv)(i,j) += tau_M*F(2)*d_test_func_dU(1)*u_phi[j][qp]*JxWxS;
216  (*Kwv)(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(2)*test_func*JxWxF;
217  (*Kwv)(i,j) += tau_M*dFdU(2,1)*u_phi[j][qp]*test_func*JxWxS;
218  (*Kww)(i,j) += tau_M*F(2)*d_test_func_dU(2)*u_phi[j][qp]*JxWxS;
219  (*Kww)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(2)*test_func*JxWxF;
220  (*Kww)(i,j) += tau_M*dFdU(2,2)*u_phi[j][qp]*test_func*JxWxS;
221  }
222  }
223 
224  } // End compute_jacobian check
225 
226  } // End i dof loop
227  }
228 
229 #ifdef GRINS_USE_GRVY_TIMERS
230  this->_timer->EndTimer("VelocityDragAdjointStabilization::element_time_derivative");
231 #endif
232 
233  return;
234  }
bool compute_force(const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
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:64
libMesh::Number _rho
Material parameters, read from input.
unsigned int _dim
Physical dimension of problem.
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:277
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
IncompressibleNavierStokesStabilizationHelper _stab_helper
template<class Mu >
void GRINS::VelocityDragAdjointStabilization< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

Definition at line 53 of file velocity_drag_adjoint_stab.C.

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

Member Data Documentation

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

Definition at line 69 of file velocity_drag_adjoint_stab.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