26 #include "grins_config.h"
35 #include "libmesh/getpot.h"
36 #include "libmesh/fem_system.h"
37 #include "libmesh/quadrature.h"
44 _flow_stab_helper(physics_name+
"FlowStabHelper", input),
45 _temp_stab_helper(physics_name+
"TempStabHelper", input),
67 const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
70 const std::vector<libMesh::Real> &JxW =
71 context.get_element_fe(_flow_vars.u())->get_JxW();
78 const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
79 context.get_element_fe(this->_flow_vars.u())->get_dphi();
82 libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u());
83 libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v());
84 libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
86 if(this->_flow_vars.dim() == 3)
88 Fw = &context.get_elem_residual(this->_flow_vars.w());
97 unsigned int n_qpoints = context.get_element_qrule().n_points();
99 libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
101 for (
unsigned int qp=0; qp != n_qpoints; qp++)
103 libMesh::RealGradient g = this->_flow_stab_helper.compute_g( fe, context, qp );
104 libMesh::RealTensor G = this->_flow_stab_helper.compute_G( fe, context, qp );
106 libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
107 context.interior_value( this->_flow_vars.v(), qp ) );
108 if( this->_flow_vars.dim() == 3 )
110 U(2) = context.interior_value( this->_flow_vars.w(), qp );
114 libMesh::Real mu_qp = this->_mu(context, qp);
116 libMesh::Real tau_M = this->_flow_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, mu_qp, this->_is_steady );
124 T = context.interior_value(_temp_vars.T(), qp);
126 libMesh::RealGradient residual = _rho*_beta_T*(T-_T_ref)*_g;
128 for (
unsigned int i=0; i != n_u_dofs; i++)
130 Fu(i) += ( -tau_M*residual(0)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
133 Fv(i) += ( -tau_M*residual(1)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
136 if (this->_flow_vars.dim() == 3)
138 (*Fw)(i) += ( -tau_M*residual(2)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
142 if (compute_jacobian)
144 libmesh_not_implemented();
153 (
bool compute_jacobian,
157 const unsigned int n_p_dofs = context.get_dof_indices(_press_var.p()).size();
160 const std::vector<libMesh::Real> &JxW =
161 context.get_element_fe(_flow_vars.u())->get_JxW();
163 const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
164 context.get_element_fe(this->_press_var.p())->get_dphi();
166 libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p());
174 unsigned int n_qpoints = context.get_element_qrule().n_points();
176 libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
178 for (
unsigned int qp=0; qp != n_qpoints; qp++)
180 libMesh::RealGradient g = this->_flow_stab_helper.compute_g( fe, context, qp );
181 libMesh::RealTensor G = this->_flow_stab_helper.compute_G( fe, context, qp );
183 libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
184 context.interior_value( this->_flow_vars.v(), qp ) );
185 if( this->_flow_vars.dim() == 3 )
186 U(2) = context.interior_value( this->_flow_vars.w(), qp );
189 libMesh::Real mu_qp = this->_mu(context, qp);
191 libMesh::Real tau_M = this->_flow_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, mu_qp, this->_is_steady );
195 T = context.interior_value(_temp_vars.T(), qp);
197 libMesh::RealGradient residual = _rho*_beta_T*(T-_T_ref)*_g;
202 for (
unsigned int i=0; i != n_p_dofs; i++)
204 Fp(i) += -tau_M*residual*p_dphi[i][qp]*JxW[qp];
206 if (compute_jacobian)
208 libmesh_not_implemented();
290 (
const std::string & param_name,
295 _mu.register_parameter(param_name, param_pointer);
virtual void set_parameter(libMesh::Number ¶m_variable, const GetPot &input, const std::string ¶m_name, libMesh::Number param_default)
Each subclass can simultaneously read a parameter value from.
static PhysicsName heat_transfer()
~BoussinesqBuoyancySPGSMStabilization()
virtual void mass_residual(bool compute_jacobian, AssemblyContext &context)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
INSTANTIATE_INC_NS_SUBCLASS(BoussinesqBuoyancySPGSMStabilization)
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context)
Constraint part(s) of physics for element interiors.
static void read_specific_heat(const std::string &core_physics_name, const GetPot &input, ParameterUser ¶ms, libMesh::Real &cp)
Helper function to reading scalar specific heat from input.
Helper functions for parsing material properties.
static PhysicsName boussinesq_buoyancy()
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
BoussinesqBuoyancySPGSMStabilization()
virtual void register_parameter(const std::string ¶m_name, libMesh::ParameterMultiAccessor< libMesh::Number > ¶m_pointer) const
Each subclass will register its copy of an independent.
virtual void register_parameter(const std::string ¶m_name, libMesh::ParameterMultiAccessor< libMesh::Number > ¶m_pointer) const
Each subclass will register its copy of an independent.