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

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

#include <averaged_fan_adjoint_stab.h>

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

Public Member Functions

 AveragedFanAdjointStabilization (const std::string &physics_name, const GetPot &input)
 
 ~AveragedFanAdjointStabilization ()
 
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::AveragedFanBase< Viscosity >
 AveragedFanBase (const std::string &physics_name, const GetPot &input)
 
 ~AveragedFanBase ()
 
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 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 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::AveragedFanBase< 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 > chord_function
 
libMesh::ParsedFunction< libMesh::Number > area_swept_function
 
libMesh::ParsedFunction< libMesh::Number > aoa_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

 AveragedFanAdjointStabilization ()
 

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

Physics class for spatially-averaged fan.

Definition at line 45 of file averaged_fan_adjoint_stab.h.

Constructor & Destructor Documentation

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

Definition at line 41 of file averaged_fan_adjoint_stab.C.

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

Definition at line 47 of file averaged_fan_adjoint_stab.C.

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

Member Function Documentation

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

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 227 of file averaged_fan_adjoint_stab.C.

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

Time dependent part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 67 of file averaged_fan_adjoint_stab.C.

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

Initialize context for added physics variables.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

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

Definition at line 66 of file averaged_fan_adjoint_stab.h.


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

Generated on Tue Dec 19 2017 12:47:30 for GRINS-0.8.0 by  doxygen 1.8.9.1