34 #include "libmesh/quadrature.h"
35 #include "libmesh/boundary_info.h"
48 context.get_element_fe(this->_flow_vars.u())->get_xyz();
49 context.get_element_fe(this->_flow_vars.u())->get_phi();
60 #ifdef GRINS_USE_GRVY_TIMERS
61 this->_timer->BeginTimer(
"VelocityDrag::element_time_derivative");
65 const std::vector<libMesh::Real> &JxW =
66 context.get_element_fe(this->_flow_vars.u())->get_JxW();
69 const std::vector<std::vector<libMesh::Real> >& u_phi =
70 context.get_element_fe(this->_flow_vars.u())->get_phi();
72 const std::vector<libMesh::Point>& u_qpoint =
73 context.get_element_fe(this->_flow_vars.u())->get_xyz();
76 const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
79 libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u());
80 libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v());
81 libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u());
82 libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v());
84 libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
85 libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
86 libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
87 libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
88 libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
90 libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u());
91 libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v());
92 libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
96 Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w());
97 Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w());
99 Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u());
100 Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v());
101 Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w());
102 Fw = &context.get_elem_residual(this->_flow_vars.w());
105 unsigned int n_qpoints = context.get_element_qrule().n_points();
107 for (
unsigned int qp=0; qp != n_qpoints; qp++)
110 libMesh::Number u, v;
111 u = context.interior_value(this->_flow_vars.u(), qp);
112 v = context.interior_value(this->_flow_vars.v(), qp);
114 libMesh::NumberVectorValue U(u,v);
116 U(2) = context.interior_value(this->_flow_vars.w(), qp);
118 libMesh::NumberVectorValue F;
119 libMesh::NumberTensorValue dFdU;
120 libMesh::NumberTensorValue* dFdU_ptr =
121 compute_jacobian ? &dFdU : NULL;
122 if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
125 libMesh::Real jac = JxW[qp];
127 for (
unsigned int i=0; i != n_u_dofs; i++)
129 const libMesh::Number jac_i = jac * u_phi[i][qp];
134 if( this->_dim == 3 )
135 (*Fw)(i) += F(2)*jac_i;
137 if( compute_jacobian )
139 for (
unsigned int j=0; j != n_u_dofs; j++)
141 const libMesh::Number jac_ij =
142 jac_i * context.get_elem_solution_derivative() *
144 Kuu(i,j) += jac_ij * dFdU(0,0);
145 Kuv(i,j) += jac_ij * dFdU(0,1);
146 Kvu(i,j) += jac_ij * dFdU(1,0);
147 Kvv(i,j) += jac_ij * dFdU(1,1);
149 if( this->_dim == 3 )
151 (*Kuw)(i,j) += jac_ij * dFdU(0,2);
152 (*Kvw)(i,j) += jac_ij * dFdU(1,2);
153 (*Kwu)(i,j) += jac_ij * dFdU(2,0);
154 (*Kwv)(i,j) += jac_ij * dFdU(2,1);
155 (*Kww)(i,j) += jac_ij * dFdU(2,2);
163 #ifdef GRINS_USE_GRVY_TIMERS
164 this->_timer->EndTimer(
"VelocityDrag::element_time_derivative");
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
INSTANTIATE_INC_NS_SUBCLASS(VelocityDrag)