27 #include "grins_config.h"
35 #include "libmesh/getpot.h"
36 #include "libmesh/quadrature.h"
44 _stab_helper( physics_name+
"StabHelper", input )
57 context.get_element_fe(this->_press_var.p())->get_dphi();
59 context.get_element_fe(this->_flow_vars.u())->get_xyz();
60 context.get_element_fe(this->_flow_vars.u())->get_phi();
61 context.get_element_fe(this->_flow_vars.u())->get_dphi();
62 context.get_element_fe(this->_flow_vars.u())->get_d2phi();
69 (
bool compute_jacobian,
73 const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
76 const std::vector<libMesh::Real> &JxW =
77 context.get_element_fe(this->_flow_vars.u())->get_JxW();
79 const std::vector<libMesh::Point>& u_qpoint =
80 context.get_element_fe(this->_flow_vars.u())->get_xyz();
82 const std::vector<std::vector<libMesh::Real> >& u_phi =
83 context.get_element_fe(this->_flow_vars.u())->get_phi();
85 const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
86 context.get_element_fe(this->_flow_vars.u())->get_dphi();
88 const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
89 context.get_element_fe(this->_flow_vars.u())->get_d2phi();
92 libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u());
93 libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v());
94 libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
96 libMesh::DenseSubMatrix<libMesh::Number> &Kuu =
97 context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u());
98 libMesh::DenseSubMatrix<libMesh::Number> &Kuv =
99 context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v());
100 libMesh::DenseSubMatrix<libMesh::Number> &Kvu =
101 context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u());
102 libMesh::DenseSubMatrix<libMesh::Number> &Kvv =
103 context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v());
105 libMesh::DenseSubMatrix<libMesh::Number> *Kuw = NULL;
106 libMesh::DenseSubMatrix<libMesh::Number> *Kvw = NULL;
107 libMesh::DenseSubMatrix<libMesh::Number> *Kwu = NULL;
108 libMesh::DenseSubMatrix<libMesh::Number> *Kwv = NULL;
109 libMesh::DenseSubMatrix<libMesh::Number> *Kww = NULL;
111 if(this->_flow_vars.dim() == 3)
113 Fw = &context.get_elem_residual(this->_flow_vars.w());
114 Kuw = &context.get_elem_jacobian
115 (this->_flow_vars.u(), this->_flow_vars.w());
116 Kvw = &context.get_elem_jacobian
117 (this->_flow_vars.v(), this->_flow_vars.w());
118 Kwu = &context.get_elem_jacobian
119 (this->_flow_vars.w(), this->_flow_vars.u());
120 Kwv = &context.get_elem_jacobian
121 (this->_flow_vars.w(), this->_flow_vars.v());
122 Kww = &context.get_elem_jacobian
123 (this->_flow_vars.w(), this->_flow_vars.w());
132 unsigned int n_qpoints = context.get_element_qrule().n_points();
134 libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
136 for (
unsigned int qp=0; qp != n_qpoints; qp++)
138 libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
139 libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
141 libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
142 context.interior_value( this->_flow_vars.v(), qp ) );
143 if( this->_flow_vars.dim() == 3 )
145 U(2) = context.interior_value( this->_flow_vars.w(), qp );
149 libMesh::Real mu_qp = this->_mu(context, qp);
152 libMesh::Real d_tau_M_d_rho;
153 libMesh::Gradient d_tau_M_dU;
155 if (compute_jacobian)
156 this->_stab_helper.compute_tau_momentum_and_derivs
157 ( context, qp, g, G, this->_rho, U, mu_qp,
158 tau_M, d_tau_M_d_rho, d_tau_M_dU,
161 tau_M = this->_stab_helper.compute_tau_momentum
162 ( context, qp, g, G, this->_rho, U, mu_qp,
165 libMesh::NumberVectorValue F;
166 libMesh::NumberTensorValue dFdU;
167 libMesh::NumberTensorValue* dFdU_ptr =
168 compute_jacobian ? &dFdU : NULL;
169 if (!this->compute_force(u_qpoint[qp], context, U, F, dFdU_ptr))
172 for (
unsigned int i=0; i != n_u_dofs; i++)
174 libMesh::Real test_func = this->_rho*U*u_gradphi[i][qp] +
175 mu_qp*( u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2) );
176 Fu(i) += tau_M*F(0)*test_func*JxW[qp];
178 Fv(i) += tau_M*F(1)*test_func*JxW[qp];
180 if (this->_flow_vars.dim() == 3)
182 (*Fw)(i) += tau_M*F(2)*test_func*JxW[qp];
185 if (compute_jacobian)
187 libMesh::Gradient d_test_func_dU = this->_rho*u_gradphi[i][qp];
189 for (
unsigned int j=0; j != n_u_dofs; ++j)
191 Kuu(i,j) += tau_M*F(0)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
192 Kuu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
193 Kuu(i,j) += tau_M*dFdU(0,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
194 Kuv(i,j) += tau_M*F(0)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
195 Kuv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
196 Kuv(i,j) += tau_M*dFdU(0,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
197 Kvu(i,j) += tau_M*F(1)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
198 Kvu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
199 Kvu(i,j) += tau_M*dFdU(1,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
200 Kvv(i,j) += tau_M*F(1)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
201 Kvv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
202 Kvv(i,j) += tau_M*dFdU(1,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
205 if (this->_flow_vars.dim() == 3)
207 for (
unsigned int j=0; j != n_u_dofs; ++j)
209 (*Kuw)(i,j) += tau_M*F(0)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
210 (*Kuw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
211 (*Kuw)(i,j) += tau_M*dFdU(0,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
212 (*Kvw)(i,j) += tau_M*F(1)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
213 (*Kvw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
214 (*Kvw)(i,j) += tau_M*dFdU(1,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
215 (*Kwu)(i,j) += tau_M*F(2)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
216 (*Kwu)(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
217 (*Kwu)(i,j) += tau_M*dFdU(2,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
218 (*Kwv)(i,j) += tau_M*F(2)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
219 (*Kwv)(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
220 (*Kwv)(i,j) += tau_M*dFdU(2,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
221 (*Kww)(i,j) += tau_M*F(2)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
222 (*Kww)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
223 (*Kww)(i,j) += tau_M*dFdU(2,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
238 const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
239 const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
242 const std::vector<libMesh::Real> &JxW =
243 context.get_element_fe(this->_flow_vars.u())->get_JxW();
245 const std::vector<libMesh::Point>& u_qpoint =
246 context.get_element_fe(this->_flow_vars.u())->get_xyz();
248 const std::vector<std::vector<libMesh::Real> >& u_phi =
249 context.get_element_fe(this->_flow_vars.u())->get_phi();
251 const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
252 context.get_element_fe(this->_press_var.p())->get_dphi();
254 libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p());
256 libMesh::DenseSubMatrix<libMesh::Number> &Kpu =
257 context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.u());
258 libMesh::DenseSubMatrix<libMesh::Number> &Kpv =
259 context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.v());
260 libMesh::DenseSubMatrix<libMesh::Number> *Kpw = NULL;
262 if(this->_flow_vars.dim() == 3)
264 Kpw = &context.get_elem_jacobian
265 (this->_press_var.p(), this->_flow_vars.w());
274 unsigned int n_qpoints = context.get_element_qrule().n_points();
276 libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
278 for (
unsigned int qp=0; qp != n_qpoints; qp++)
280 libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
281 libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
283 libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
284 context.interior_value( this->_flow_vars.v(), qp ) );
285 if( this->_flow_vars.dim() == 3 )
287 U(2) = context.interior_value( this->_flow_vars.w(), qp );
291 libMesh::Real mu_qp = this->_mu(context, qp);
294 libMesh::Real d_tau_M_d_rho;
295 libMesh::Gradient d_tau_M_dU;
297 if (compute_jacobian)
298 this->_stab_helper.compute_tau_momentum_and_derivs
299 ( context, qp, g, G, this->_rho, U, mu_qp,
300 tau_M, d_tau_M_d_rho, d_tau_M_dU,
303 tau_M = this->_stab_helper.compute_tau_momentum
304 ( context, qp, g, G, this->_rho, U, mu_qp,
307 libMesh::NumberVectorValue F;
308 libMesh::NumberTensorValue dFdU;
309 libMesh::NumberTensorValue* dFdU_ptr =
310 compute_jacobian ? &dFdU : NULL;
311 if (!this->compute_force(u_qpoint[qp], context, U, F, dFdU_ptr))
317 for (
unsigned int i=0; i != n_p_dofs; i++)
319 Fp(i) += -tau_M*F*p_dphi[i][qp]*JxW[qp];
321 if (compute_jacobian)
323 for (
unsigned int j=0; j != n_u_dofs; ++j)
325 Kpu(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
326 Kpv(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
327 for (
unsigned int d=0; d != 3; ++d)
329 Kpu(i,j) += -tau_M*dFdU(d,0)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
330 Kpv(i,j) += -tau_M*dFdU(d,1)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
333 if( this->_flow_vars.dim() == 3 )
334 for (
unsigned int j=0; j != n_u_dofs; ++j)
336 (*Kpw)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
337 for (
unsigned int d=0; d != 3; ++d)
339 (*Kpw)(i,j) += -tau_M*dFdU(d,2)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
INSTANTIATE_INC_NS_SUBCLASS(VelocityPenaltyAdjointStabilization)
VelocityPenaltyAdjointStabilization()
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context)
Constraint part(s) of physics for element interiors.
~VelocityPenaltyAdjointStabilization()