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

Adds Boussinesq bouyancy adjoint stabilization source term. More...

#include <boussinesq_buoyancy_adjoint_stab.h>

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

Public Member Functions

 BoussinesqBuoyancyAdjointStabilization (const std::string &physics_name, const GetPot &input)
 
 ~BoussinesqBuoyancyAdjointStabilization ()
 
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...
 
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...
 
- Public Member Functions inherited from GRINS::BoussinesqBuoyancyBase
 BoussinesqBuoyancyBase (const std::string &physics_name, const GetPot &input)
 
 ~BoussinesqBuoyancyBase ()
 
- 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 set_time_evolving_vars (libMesh::FEMSystem *system)
 Set which variables are time evolving. 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

Viscosity _mu
 Viscosity object. More...
 
IncompressibleNavierStokesStabilizationHelper _stab_helper
 
- Protected Attributes inherited from GRINS::BoussinesqBuoyancyBase
const VelocityVariable_flow_vars
 
const PressureFEVariable_press_var
 
const PrimitiveTempFEVariables_temp_vars
 
libMesh::Number _rho
 $ \rho = $ density More...
 
libMesh::Number _T_ref
 $ T_0 = $ reference temperature More...
 
libMesh::Number _beta_T
 $ \beta_T = $ coefficient of thermal expansion More...
 
libMesh::Point _g
 Gravitational vector. 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

 BoussinesqBuoyancyAdjointStabilization ()
 

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::BoussinesqBuoyancyBase
void read_property (const GetPot &input, const std::string &old_option, const std::string &property, libMesh::Real &value)
 Helper function for parsing/maintaing backward compatibility. More...
 
void no_input_warning (const GetPot &input, const std::string &old_option, const std::string &material, const std::string &property)
 Helper function for parsing/maintaing backward compatibility. 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::BoussinesqBuoyancyAdjointStabilization< Viscosity >

Adds Boussinesq bouyancy adjoint stabilization source term.

This class implements the adjiont stabilization term for the BoussinesqBuoyancy Physics. Intended to be used with IncompressibleNavierStokesAdjointStabilization and HeatTransferStabilization.

Definition at line 42 of file boussinesq_buoyancy_adjoint_stab.h.

Constructor & Destructor Documentation

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

Definition at line 44 of file boussinesq_buoyancy_adjoint_stab.C.

45  : BoussinesqBuoyancyBase(physics_name,input),
47  _stab_helper( physics_name+"StabHelper", input )
48  {}
IncompressibleNavierStokesStabilizationHelper _stab_helper
static PhysicsName boussinesq_buoyancy()
static std::string material_name(const GetPot &input, const std::string &physics)
Get the name of the material in the Physics/physics section.

Definition at line 51 of file boussinesq_buoyancy_adjoint_stab.C.

52  {
53  return;
54  }
template<class Viscosity >
GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::BoussinesqBuoyancyAdjointStabilization ( )
private

Member Function Documentation

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

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 249 of file boussinesq_buoyancy_adjoint_stab.C.

250  {
251  // The number of local degrees of freedom in each variable.
252  const unsigned int n_p_dofs = context.get_dof_indices(_press_var.p()).size();
253  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
254  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
255 
256  // Element Jacobian * quadrature weights for interior integration.
257  const std::vector<libMesh::Real> &JxW =
258  context.get_element_fe(_flow_vars.u())->get_JxW();
259 
260  const std::vector<std::vector<libMesh::Real> >& T_phi =
261  context.get_element_fe(this->_temp_vars.T())->get_phi();
262 
263  const std::vector<std::vector<libMesh::Real> >& u_phi =
264  context.get_element_fe(this->_flow_vars.u())->get_phi();
265 
266  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
267  context.get_element_fe(this->_press_var.p())->get_dphi();
268 
269  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
270 
271  libMesh::DenseSubMatrix<libMesh::Number> &KpT =
272  context.get_elem_jacobian(_press_var.p(), _temp_vars.T()); // J_{pT}
273  libMesh::DenseSubMatrix<libMesh::Number> &Kpu =
274  context.get_elem_jacobian(_press_var.p(), _flow_vars.u()); // J_{pu}
275  libMesh::DenseSubMatrix<libMesh::Number> &Kpv =
276  context.get_elem_jacobian(_press_var.p(), _flow_vars.v()); // J_{pv}
277  libMesh::DenseSubMatrix<libMesh::Number> *Kpw = NULL;
278 
279  if(this->_flow_vars.dim() == 3)
280  {
281  Kpw = &context.get_elem_jacobian
282  (_press_var.p(), _flow_vars.w()); // J_{pw}
283  }
284 
285  // Now we will build the element Jacobian and residual.
286  // Constructing the residual requires the solution and its
287  // gradient from the previous timestep. This must be
288  // calculated at each quadrature point by summing the
289  // solution degree-of-freedom values by the appropriate
290  // weight functions.
291  unsigned int n_qpoints = context.get_element_qrule().n_points();
292 
293  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
294 
295  for (unsigned int qp=0; qp != n_qpoints; qp++)
296  {
297  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
298  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
299 
300  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
301  context.interior_value( this->_flow_vars.v(), qp ) );
302  if( this->_flow_vars.dim() == 3 )
303  {
304  U(2) = context.interior_value( this->_flow_vars.w(), qp );
305  }
306 
307  // Compute the viscosity at this qp
308  libMesh::Real mu_qp = this->_mu(context, qp);
309 
310  libMesh::Real tau_M;
311  libMesh::Real d_tau_M_d_rho;
312  libMesh::Gradient d_tau_M_dU;
313 
314  if (compute_jacobian)
316  ( context, qp, g, G, this->_rho, U, mu_qp,
317  tau_M, d_tau_M_d_rho, d_tau_M_dU,
318  this->_is_steady );
319  else
320  tau_M = this->_stab_helper.compute_tau_momentum
321  ( context, qp, g, G, this->_rho, U, mu_qp,
322  this->_is_steady );
323 
324  // Compute the solution & its gradient at the old Newton iterate.
325  libMesh::Number T;
326  T = context.interior_value(_temp_vars.T(), qp);
327 
328  libMesh::RealGradient d_residual_dT = _rho*_beta_T*_g;
329  // d_residual_dU = 0
330  libMesh::RealGradient residual = (T-_T_ref)*d_residual_dT;
331 
332  // First, an i-loop over the velocity degrees of freedom.
333  // We know that n_u_dofs == n_v_dofs so we can compute contributions
334  // for both at the same time.
335  for (unsigned int i=0; i != n_p_dofs; i++)
336  {
337  Fp(i) += tau_M*residual*p_dphi[i][qp]*JxW[qp];
338 
339  if (compute_jacobian)
340  {
341  for (unsigned int j=0; j != n_T_dofs; ++j)
342  {
343  KpT(i,j) += tau_M*d_residual_dT*T_phi[j][qp]*p_dphi[i][qp]*JxW[qp] * context.get_elem_solution_derivative();
344  }
345 
346  for (unsigned int j=0; j != n_u_dofs; ++j)
347  {
348  Kpu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*residual*p_dphi[i][qp]*JxW[qp] * context.get_elem_solution_derivative();
349  Kpv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*residual*p_dphi[i][qp]*JxW[qp] * context.get_elem_solution_derivative();
350  }
351  if( this->_flow_vars.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]*residual*p_dphi[i][qp]*JxW[qp] * context.get_elem_solution_derivative();
355  }
356  }
357  }
358  } // End quadrature loop
359  }
libMesh::Point _g
Gravitational vector.
const PressureFEVariable & _press_var
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
VariableIndex T() const
const VelocityVariable & _flow_vars
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
VariableIndex p() const
libMesh::Number _T_ref
reference temperature
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
const PrimitiveTempFEVariables & _temp_vars
libMesh::Number _beta_T
coefficient of thermal expansion
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::BoussinesqBuoyancyAdjointStabilization< 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 boussinesq_buoyancy_adjoint_stab.C.

70  {
71  // The number of local degrees of freedom in each variable.
72  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
73  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
74 
75  // Element Jacobian * quadrature weights for interior integration.
76  const std::vector<libMesh::Real> &JxW =
77  context.get_element_fe(_flow_vars.u())->get_JxW();
78 
79  const std::vector<std::vector<libMesh::Real> >& T_phi =
80  context.get_element_fe(this->_temp_vars.T())->get_phi();
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(_flow_vars.u()); // R_{u}
93  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v()); // R_{v}
94  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
95 
96  libMesh::DenseSubMatrix<libMesh::Number> &KuT =
97  context.get_elem_jacobian(_flow_vars.u(), _temp_vars.T()); // J_{uT}
98  libMesh::DenseSubMatrix<libMesh::Number> &KvT =
99  context.get_elem_jacobian(_flow_vars.v(), _temp_vars.T()); // J_{vT}
100  libMesh::DenseSubMatrix<libMesh::Number> &Kuu =
101  context.get_elem_jacobian(_flow_vars.u(), _flow_vars.u()); // J_{uu}
102  libMesh::DenseSubMatrix<libMesh::Number> &Kuv =
103  context.get_elem_jacobian(_flow_vars.u(), _flow_vars.v()); // J_{uv}
104  libMesh::DenseSubMatrix<libMesh::Number> &Kvu =
105  context.get_elem_jacobian(_flow_vars.v(), _flow_vars.u()); // J_{vu}
106  libMesh::DenseSubMatrix<libMesh::Number> &Kvv =
107  context.get_elem_jacobian(_flow_vars.v(), _flow_vars.v()); // J_{vv}
108 
109  libMesh::DenseSubMatrix<libMesh::Number> *KwT = NULL;
110  libMesh::DenseSubMatrix<libMesh::Number> *Kuw = NULL;
111  libMesh::DenseSubMatrix<libMesh::Number> *Kvw = NULL;
112  libMesh::DenseSubMatrix<libMesh::Number> *Kwu = NULL;
113  libMesh::DenseSubMatrix<libMesh::Number> *Kwv = NULL;
114  libMesh::DenseSubMatrix<libMesh::Number> *Kww = NULL;
115 
116 
117  if(this->_flow_vars.dim() == 3)
118  {
119  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
120  KwT = &context.get_elem_jacobian
121  (_flow_vars.w(), _temp_vars.T()); // J_{wT}
122  Kuw = &context.get_elem_jacobian
123  (_flow_vars.u(), _flow_vars.w()); // J_{uw}
124  Kvw = &context.get_elem_jacobian
125  (_flow_vars.v(), _flow_vars.w()); // J_{vw}
126  Kwu = &context.get_elem_jacobian
127  (_flow_vars.w(), _flow_vars.u()); // J_{wu}
128  Kwv = &context.get_elem_jacobian
129  (_flow_vars.w(), _flow_vars.v()); // J_{wv}
130  Kww = &context.get_elem_jacobian
131  (_flow_vars.w(), _flow_vars.w()); // J_{ww}
132  }
133 
134  // Now we will build the element Jacobian and residual.
135  // Constructing the residual requires the solution and its
136  // gradient from the previous timestep. This must be
137  // calculated at each quadrature point by summing the
138  // solution degree-of-freedom values by the appropriate
139  // weight functions.
140  unsigned int n_qpoints = context.get_element_qrule().n_points();
141 
142  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
143 
144  for (unsigned int qp=0; qp != n_qpoints; qp++)
145  {
146  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
147  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
148 
149  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
150  context.interior_value( this->_flow_vars.v(), qp ) );
151  if( this->_flow_vars.dim() == 3 )
152  {
153  U(2) = context.interior_value( this->_flow_vars.w(), qp );
154  }
155 
156  // Compute the viscosity at this qp
157  libMesh::Real mu_qp = this->_mu(context, qp);
158 
159  libMesh::Real tau_M;
160  libMesh::Real d_tau_M_d_rho;
161  libMesh::Gradient d_tau_M_dU;
162 
163  if (compute_jacobian)
165  ( context, qp, g, G, this->_rho, U, mu_qp,
166  tau_M, d_tau_M_d_rho, d_tau_M_dU,
167  this->_is_steady );
168  else
169  tau_M = this->_stab_helper.compute_tau_momentum
170  ( context, qp, g, G, this->_rho, U, mu_qp,
171  this->_is_steady );
172 
173  // Compute the solution & its gradient at the old Newton iterate.
174  libMesh::Number T;
175  T = context.interior_value(_temp_vars.T(), qp);
176 
177  libMesh::RealGradient d_residual_dT = _rho*_beta_T*_g;
178  // d_residual_dU = 0
179  libMesh::RealGradient residual = (T-_T_ref)*d_residual_dT;
180 
181  for (unsigned int i=0; i != n_u_dofs; i++)
182  {
183  libMesh::Real test_func = this->_rho*U*u_gradphi[i][qp] +
184  mu_qp*( u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2) );
185  Fu(i) += -tau_M*residual(0)*test_func*JxW[qp];
186 
187  Fv(i) += -tau_M*residual(1)*test_func*JxW[qp];
188 
189  if (this->_flow_vars.dim() == 3)
190  {
191  (*Fw)(i) += -tau_M*residual(2)*test_func*JxW[qp];
192  }
193 
194  if (compute_jacobian)
195  {
196  libMesh::Gradient d_test_func_dU = this->_rho*u_gradphi[i][qp];
197  // d_test_func_dT = 0
198 
199  for (unsigned int j=0; j != n_u_dofs; ++j)
200  {
201  Kuu(i,j) += -tau_M*residual(0)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
202  Kuu(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*residual(0)*test_func*JxW[qp] * context.get_elem_solution_derivative();
203  Kuv(i,j) += -tau_M*residual(0)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
204  Kuv(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*residual(0)*test_func*JxW[qp] * context.get_elem_solution_derivative();
205  Kvu(i,j) += -tau_M*residual(1)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
206  Kvu(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*residual(1)*test_func*JxW[qp] * context.get_elem_solution_derivative();
207  Kvv(i,j) += -tau_M*residual(1)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
208  Kvv(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*residual(1)*test_func*JxW[qp] * context.get_elem_solution_derivative();
209  }
210 
211  for (unsigned int j=0; j != n_T_dofs; ++j)
212  {
213  // KuT(i,j) += -tau_M*residual(0)*dtest_func_dT[j]*JxW[qp] * context.get_elem_solution_derivative();
214  KuT(i,j) += -tau_M*d_residual_dT(0)*T_phi[j][qp]*test_func*JxW[qp] * context.get_elem_solution_derivative();
215  // KvT(i,j) += -tau_M*residual(1)*dtest_func_dT[j]*JxW[qp] * context.get_elem_solution_derivative();
216  KvT(i,j) += -tau_M*d_residual_dT(1)*T_phi[j][qp]*test_func*JxW[qp] * context.get_elem_solution_derivative();
217  }
218  if (this->_flow_vars.dim() == 3)
219  {
220  for (unsigned int j=0; j != n_T_dofs; ++j)
221  {
222  // KwT(i,j) += -tau_M*residual(2)*dtest_func_dT[j]*JxW[qp] * context.get_elem_solution_derivative();
223  (*KwT)(i,j) += -tau_M*d_residual_dT(2)*T_phi[j][qp]*test_func*JxW[qp] * context.get_elem_solution_derivative();
224  }
225 
226  for (unsigned int j=0; j != n_u_dofs; ++j)
227  {
228  (*Kuw)(i,j) += -tau_M*residual(0)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
229  (*Kuw)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*residual(0)*test_func*JxW[qp] * context.get_elem_solution_derivative();
230  (*Kvw)(i,j) += -tau_M*residual(1)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
231  (*Kvw)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*residual(1)*test_func*JxW[qp] * context.get_elem_solution_derivative();
232  (*Kwu)(i,j) += -tau_M*residual(2)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
233  (*Kwu)(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*residual(2)*test_func*JxW[qp] * context.get_elem_solution_derivative();
234  (*Kwv)(i,j) += -tau_M*residual(2)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
235  (*Kwv)(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*residual(2)*test_func*JxW[qp] * context.get_elem_solution_derivative();
236  (*Kww)(i,j) += -tau_M*residual(2)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
237  (*Kww)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*residual(2)*test_func*JxW[qp] * context.get_elem_solution_derivative();
238  }
239  }
240 
241  } // End compute_jacobian check
242 
243  } // End i dof loop
244  } // End quadrature loop
245  }
libMesh::Point _g
Gravitational vector.
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
VariableIndex T() const
const VelocityVariable & _flow_vars
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Number _T_ref
reference temperature
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
const PrimitiveTempFEVariables & _temp_vars
libMesh::Number _beta_T
coefficient of thermal expansion
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::BoussinesqBuoyancyAdjointStabilization< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::Physics.

Definition at line 57 of file boussinesq_buoyancy_adjoint_stab.C.

58  {
59  context.get_element_fe(this->_press_var.p())->get_dphi();
60 
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  }
const PressureFEVariable & _press_var
const VelocityVariable & _flow_vars
VariableIndex p() const
template<class Mu >
void GRINS::BoussinesqBuoyancyAdjointStabilization< Mu >::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiAccessor< libMesh::Number > &  param_pointer 
) const
virtual

Each subclass will register its copy of an independent.

Reimplemented from GRINS::ParameterUser.

Definition at line 363 of file boussinesq_buoyancy_adjoint_stab.C.

References GRINS::ParameterUser::register_parameter().

366  {
367  ParameterUser::register_parameter(param_name, param_pointer);
368  _mu.register_parameter(param_name, param_pointer);
369  }
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.

Member Data Documentation

template<class Viscosity >
Viscosity GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::_mu
protected

Viscosity object.

Definition at line 68 of file boussinesq_buoyancy_adjoint_stab.h.

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

Definition at line 70 of file boussinesq_buoyancy_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