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

Physics class for spatially-averaged turbine. More...

#include <averaged_turbine_adjoint_stab.h>

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

Public Member Functions

 AveragedTurbineAdjointStabilization (const std::string &physics_name, const GetPot &input)
 
 ~AveragedTurbineAdjointStabilization ()
 
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::AveragedTurbineBase< Viscosity >
 AveragedTurbineBase (const std::string &physics_name, const GetPot &input)
 
 ~AveragedTurbineBase ()
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialization of variables. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets turbine_speed and velocity variables to be time-evolving. More...
 
bool compute_force (const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::Number s, libMesh::NumberVectorValue &U_B_1, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL, libMesh::NumberVectorValue *dFds=NULL)
 
VariableIndex fan_speed_var () const
 
- 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 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::AveragedTurbineBase< Viscosity >
libMesh::ParsedFunction< libMesh::Number > base_velocity_function
 
libMesh::ParsedFunction< libMesh::Number > local_vertical_function
 
libMesh::ParsedFunction< libMesh::Number > lift_function
 
libMesh::ParsedFunction< libMesh::Number > drag_function
 
libMesh::ParsedFunction< libMesh::Number > torque_function
 
libMesh::Number moment_of_inertia
 
libMesh::Number initial_speed
 
libMesh::ParsedFunction< libMesh::Number > chord_function
 
libMesh::ParsedFunction< libMesh::Number > area_swept_function
 
libMesh::ParsedFunction< libMesh::Number > aoa_function
 
VariableIndex _fan_speed_var
 
std::string _fan_speed_var_name
 
- 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

 AveragedTurbineAdjointStabilization ()
 

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

Physics class for spatially-averaged turbine.

Definition at line 46 of file averaged_turbine_adjoint_stab.h.

Constructor & Destructor Documentation

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

Definition at line 41 of file averaged_turbine_adjoint_stab.C.

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

Definition at line 47 of file averaged_turbine_adjoint_stab.C.

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

Member Function Documentation

template<class Mu >
void GRINS::AveragedTurbineAdjointStabilization< 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 252 of file averaged_turbine_adjoint_stab.C.

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

70  {
71 #ifdef GRINS_USE_GRVY_TIMERS
72  this->_timer->BeginTimer("AveragedTurbineAdjointStabilization::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> &Kus =
101  context.get_elem_jacobian(this->_flow_vars.u(),
102  this->fan_speed_var()); // R_{u},{s}
103  libMesh::DenseSubMatrix<libMesh::Number> &Kvs =
104  context.get_elem_jacobian(this->_flow_vars.v(),
105  this->fan_speed_var()); // R_{v},{s}
106 
107  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
108  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
109  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
110  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
111  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
112 
113  libMesh::DenseSubMatrix<libMesh::Number>* Kws = 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  Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w()); // R_{u},{w}
122  Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w()); // R_{v},{w}
123 
124  Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u()); // R_{w},{u}
125  Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v()); // R_{w},{v}
126  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
127 
128  Kws = &context.get_elem_jacobian(this->_flow_vars.w(), this->fan_speed_var()); // R_{w},{s}
129 
130  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
131  }
132 
133  unsigned int n_qpoints = context.get_element_qrule().n_points();
134 
135  for (unsigned int qp=0; qp != n_qpoints; qp++)
136  {
137  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
138  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
139 
140  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
141  context.interior_value( this->_flow_vars.v(), qp ) );
142 
143  if( this->_dim == 3 )
144  {
145  U(2) = context.interior_value( this->_flow_vars.w(), qp );
146  }
147 
148  libMesh::Number s = context.interior_value(this->fan_speed_var(), qp);
149 
150  // Compute the viscosity at this qp
151  libMesh::Real mu_qp = this->_mu(context, qp);
152 
153  libMesh::Real tau_M;
154  libMesh::Real d_tau_M_d_rho;
155  libMesh::Gradient d_tau_M_dU;
156 
157  if (compute_jacobian)
159  ( context, qp, g, G, this->_rho, U, mu_qp,
160  tau_M, d_tau_M_d_rho, d_tau_M_dU,
161  this->_is_steady );
162  else
163  tau_M = this->_stab_helper.compute_tau_momentum
164  ( context, qp, g, G, this->_rho, U, mu_qp,
165  this->_is_steady );
166 
167  libMesh::NumberVectorValue U_B_1;
168  libMesh::NumberVectorValue F;
169  libMesh::NumberTensorValue dFdU;
170  libMesh::NumberTensorValue* dFdU_ptr =
171  compute_jacobian ? &dFdU : NULL;
172  libMesh::NumberVectorValue dFds;
173  libMesh::NumberVectorValue* dFds_ptr =
174  compute_jacobian ? &dFds : NULL;
175  if (!this->compute_force(u_qpoint[qp], context.time, U, s,
176  U_B_1, F, dFdU_ptr, dFds_ptr))
177  continue;
178 
179  for (unsigned int i=0; i != n_u_dofs; i++)
180  {
181  libMesh::Real test_func = this->_rho*U*u_gradphi[i][qp] +
182  mu_qp*( u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2) );
183  Fu(i) += tau_M*F(0)*test_func*JxW[qp];
184 
185  Fv(i) += tau_M*F(1)*test_func*JxW[qp];
186 
187  if (this->_dim == 3)
188  {
189  (*Fw)(i) += tau_M*F(2)*test_func*JxW[qp];
190  }
191 
192  if (compute_jacobian)
193  {
194  libMesh::Gradient d_test_func_dU = this->_rho*u_gradphi[i][qp];
195 
196  for (unsigned int j=0; j != n_u_dofs; ++j)
197  {
198  Kuu(i,j) += tau_M*F(0)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp];
199  Kuu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(0)*test_func*JxW[qp];
200  Kuu(i,j) += tau_M*dFdU(0,0)*u_phi[j][qp]*test_func*JxW[qp];
201  Kuv(i,j) += tau_M*F(0)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp];
202  Kuv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(0)*test_func*JxW[qp];
203  Kuv(i,j) += tau_M*dFdU(0,1)*u_phi[j][qp]*test_func*JxW[qp];
204  Kvu(i,j) += tau_M*F(1)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp];
205  Kvu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(1)*test_func*JxW[qp];
206  Kvu(i,j) += tau_M*dFdU(1,0)*u_phi[j][qp]*test_func*JxW[qp];
207  Kvv(i,j) += tau_M*F(1)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp];
208  Kvv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(1)*test_func*JxW[qp];
209  Kvv(i,j) += tau_M*dFdU(1,1)*u_phi[j][qp]*test_func*JxW[qp];
210 
211  Kus(i,0) += tau_M*dFds(0)*test_func*JxW[qp];
212  Kvs(i,0) += tau_M*dFds(1)*test_func*JxW[qp];
213  }
214 
215  if (this->_dim == 3)
216  {
217  for (unsigned int j=0; j != n_u_dofs; ++j)
218  {
219  (*Kuw)(i,j) += tau_M*F(0)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp];
220  (*Kuw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(0)*test_func*JxW[qp];
221  (*Kuw)(i,j) += tau_M*dFdU(0,2)*u_phi[j][qp]*test_func*JxW[qp];
222  (*Kvw)(i,j) += tau_M*F(1)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp];
223  (*Kvw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(1)*test_func*JxW[qp];
224  (*Kvw)(i,j) += tau_M*dFdU(1,2)*u_phi[j][qp]*test_func*JxW[qp];
225  (*Kwu)(i,j) += tau_M*F(2)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp];
226  (*Kwu)(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(2)*test_func*JxW[qp];
227  (*Kwu)(i,j) += tau_M*dFdU(2,0)*u_phi[j][qp]*test_func*JxW[qp];
228  (*Kwv)(i,j) += tau_M*F(2)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp];
229  (*Kwv)(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(2)*test_func*JxW[qp];
230  (*Kwv)(i,j) += tau_M*dFdU(2,1)*u_phi[j][qp]*test_func*JxW[qp];
231  (*Kww)(i,j) += tau_M*F(2)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp];
232  (*Kww)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(2)*test_func*JxW[qp];
233  (*Kww)(i,j) += tau_M*dFdU(2,2)*u_phi[j][qp]*test_func*JxW[qp];
234 
235  (*Kws)(i,0) += tau_M*dFds(2)*test_func*JxW[qp];
236  }
237  }
238 
239  } // End compute_jacobian check
240 
241  } // End i dof loop
242  }
243 
244 #ifdef GRINS_USE_GRVY_TIMERS
245  this->_timer->EndTimer("AveragedTurbineAdjointStabilization::element_time_derivative");
246 #endif
247 
248  return;
249  }
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
IncompressibleNavierStokesStabilizationHelper _stab_helper
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
bool compute_force(const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::Number s, libMesh::NumberVectorValue &U_B_1, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL, libMesh::NumberVectorValue *dFds=NULL)
VariableIndex fan_speed_var() const
template<class Mu >
void GRINS::AveragedTurbineAdjointStabilization< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

Definition at line 53 of file averaged_turbine_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::AveragedTurbineAdjointStabilization< Viscosity >::_stab_helper
protected

Definition at line 69 of file averaged_turbine_adjoint_stab.h.


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

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