27 #include "grins_config.h"
34 #include "libmesh/getpot.h"
35 #include "libmesh/fem_system.h"
36 #include "libmesh/quadrature.h"
53 (
bool compute_jacobian,
57 const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
58 const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
61 const std::vector<libMesh::Real> &JxW =
62 context.get_element_fe(_flow_vars.u())->get_JxW();
65 const std::vector<std::vector<libMesh::Real> >& vel_phi =
66 context.get_element_fe(_flow_vars.u())->get_phi();
69 const std::vector<std::vector<libMesh::Real> >& T_phi =
70 context.get_element_fe(_temp_vars.T())->get_phi();
73 libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u());
74 libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v());
75 libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
78 libMesh::DenseSubMatrix<libMesh::Number> &KuT = context.get_elem_jacobian(_flow_vars.u(), _temp_vars.T());
79 libMesh::DenseSubMatrix<libMesh::Number> &KvT = context.get_elem_jacobian(_flow_vars.v(), _temp_vars.T());
80 libMesh::DenseSubMatrix<libMesh::Number>* KwT = NULL;
84 if( this->_flow_vars.dim() == 3 )
86 Fw = &context.get_elem_residual(_flow_vars.w());
87 KwT = &context.get_elem_jacobian(_flow_vars.w(), _temp_vars.T());
96 unsigned int n_qpoints = context.get_element_qrule().n_points();
98 for (
unsigned int qp=0; qp != n_qpoints; qp++)
102 T = context.interior_value(_temp_vars.T(), qp);
107 for (
unsigned int i=0; i != n_u_dofs; i++)
109 Fu(i) += -_rho*_beta_T*(T - _T_ref)*_g(0)*vel_phi[i][qp]*JxW[qp];
110 Fv(i) += -_rho*_beta_T*(T - _T_ref)*_g(1)*vel_phi[i][qp]*JxW[qp];
112 if (this->_flow_vars.dim() == 3)
113 (*Fw)(i) += -_rho*_beta_T*(T - _T_ref)*_g(2)*vel_phi[i][qp]*JxW[qp];
115 if (compute_jacobian)
117 for (
unsigned int j=0; j != n_T_dofs; j++)
119 KuT(i,j) += context.get_elem_solution_derivative() *
120 -_rho*_beta_T*_g(0)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
121 KvT(i,j) += context.get_elem_solution_derivative() *
122 -_rho*_beta_T*_g(1)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
124 if (this->_flow_vars.dim() == 3)
125 (*KwT)(i,j) += context.get_elem_solution_derivative() *
126 -_rho*_beta_T*_g(2)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Source term contribution for BoussinesqBuoyancy.