32 #include "libmesh/quadrature.h"
53 (
bool compute_jacobian,
57 const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
58 const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
61 const std::vector<libMesh::Real> &JxW =
62 context.get_element_fe(this->_temp_vars.T())->get_JxW();
64 const std::vector<std::vector<libMesh::Real> >& u_phi =
65 context.get_element_fe(this->_flow_vars.u())->get_phi();
67 const std::vector<std::vector<libMesh::Real> >& T_phi =
68 context.get_element_fe(this->_temp_vars.T())->get_phi();
70 const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
71 context.get_element_fe(this->_temp_vars.T())->get_dphi();
73 const std::vector<std::vector<libMesh::RealTensor> >& T_hessphi =
74 context.get_element_fe(this->_temp_vars.T())->get_d2phi();
91 libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T());
92 libMesh::DenseSubMatrix<libMesh::Number> &KTT =
93 context.get_elem_jacobian(this->_temp_vars.T(), this->_temp_vars.T());
94 libMesh::DenseSubMatrix<libMesh::Number> &KTu =
95 context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.u());
96 libMesh::DenseSubMatrix<libMesh::Number> &KTv =
97 context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.v());
98 libMesh::DenseSubMatrix<libMesh::Number> *KTw = NULL;
100 if(this->_flow_vars.dim() == 3)
102 KTw = &context.get_elem_jacobian
103 (this->_temp_vars.T(), this->_flow_vars.w());
106 unsigned int n_qpoints = context.get_element_qrule().n_points();
108 for (
unsigned int qp=0; qp != n_qpoints; qp++)
110 libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
112 libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
113 libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
115 libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
116 context.interior_value( this->_flow_vars.v(), qp ) );
117 if( this->_flow_vars.dim() == 3 )
119 U(2) = context.interior_value( this->_flow_vars.w(), qp );
124 libMesh::Real tau_E, RE_s;
125 libMesh::Real d_tau_E_dT_rho, d_RE_s_dT;
126 libMesh::Gradient d_tau_E_dU, d_RE_s_dgradT, d_RE_s_dU;
127 libMesh::Tensor d_RE_s_dhessT;
130 libMesh::Real _k_qp = this->_k(context, qp);
132 if (compute_jacobian)
134 this->_stab_helper.compute_tau_energy_and_derivs
135 ( context, G, this->_rho, this->_Cp, _k_qp, U,
136 tau_E, d_tau_E_dT_rho, d_tau_E_dU, this->_is_steady );
137 this->_stab_helper.compute_res_energy_steady_and_derivs
138 ( context, qp, this->_rho, this->_Cp, _k_qp,
139 RE_s, d_RE_s_dT, d_RE_s_dgradT, d_RE_s_dhessT,
144 tau_E = this->_stab_helper.compute_tau_energy
145 ( context, G, this->_rho, this->_Cp, _k_qp, U, this->_is_steady );
147 RE_s = this->_stab_helper.compute_res_energy_steady
148 ( context, qp, this->_rho, this->_Cp, _k_qp );
163 for (
unsigned int i=0; i != n_T_dofs; i++)
165 FT(i) += -tau_E*RE_s*( this->_rho*this->_Cp*U*T_gradphi[i][qp]
166 + _k_qp*(T_hessphi[i][qp](0,0) + T_hessphi[i][qp](1,1) + T_hessphi[i][qp](2,2))
168 if (compute_jacobian)
170 for (
unsigned int j=0; j != n_T_dofs; ++j)
173 (d_RE_s_dT*T_phi[j][qp] +
174 d_RE_s_dgradT*T_gradphi[j][qp] +
175 d_RE_s_dhessT.contract(T_hessphi[j][qp])
177 ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
178 + _k_qp*(T_hessphi[i][qp](0,0) +
179 T_hessphi[i][qp](1,1) +
180 T_hessphi[i][qp](2,2))
182 * context.get_fixed_solution_derivative();
184 for (
unsigned int j=0; j != n_u_dofs; ++j)
186 KTu(i,j) += -tau_E*RE_s*
187 ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](0) )*JxW[qp]
188 * context.get_fixed_solution_derivative();
190 -(tau_E*d_RE_s_dU(0)+d_tau_E_dU(0)*RE_s)*u_phi[j][qp]*
191 ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
192 + _k_qp*(T_hessphi[i][qp](0,0) +
193 T_hessphi[i][qp](1,1) +
194 T_hessphi[i][qp](2,2))
196 * context.get_fixed_solution_derivative();
197 KTv(i,j) += -tau_E*RE_s*
198 ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](1) )*JxW[qp]
199 * context.get_fixed_solution_derivative();
201 -(tau_E*d_RE_s_dU(1)+d_tau_E_dU(1)*RE_s)*u_phi[j][qp]*
202 ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
203 + _k_qp*(T_hessphi[i][qp](0,0) +
204 T_hessphi[i][qp](1,1) +
205 T_hessphi[i][qp](2,2))
207 * context.get_fixed_solution_derivative();
209 if(this->_flow_vars.dim() == 3)
211 for (
unsigned int j=0; j != n_u_dofs; ++j)
213 (*KTw)(i,j) += -tau_E*RE_s*
214 ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](2) )*JxW[qp]
215 * context.get_fixed_solution_derivative();
217 -(tau_E*d_RE_s_dU(2)+d_tau_E_dU(2)*RE_s)*u_phi[j][qp]*
218 ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
219 + _k_qp*(T_hessphi[i][qp](0,0) +
220 T_hessphi[i][qp](1,1) +
221 T_hessphi[i][qp](2,2))
223 * context.get_fixed_solution_derivative();
236 const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
237 const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
240 const std::vector<libMesh::Real> &JxW =
241 context.get_element_fe(this->_temp_vars.T())->get_JxW();
243 const std::vector<std::vector<libMesh::Real> >& u_phi =
244 context.get_element_fe(this->_flow_vars.u())->get_phi();
246 const std::vector<std::vector<libMesh::Real> >& T_phi =
247 context.get_element_fe(this->_temp_vars.T())->get_phi();
249 const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
250 context.get_element_fe(this->_temp_vars.T())->get_dphi();
252 const std::vector<std::vector<libMesh::RealTensor> >& T_hessphi =
253 context.get_element_fe(this->_temp_vars.T())->get_d2phi();
270 libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T());
271 libMesh::DenseSubMatrix<libMesh::Number> &KTT =
272 context.get_elem_jacobian(this->_temp_vars.T(), this->_temp_vars.T());
273 libMesh::DenseSubMatrix<libMesh::Number> &KTu =
274 context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.u());
275 libMesh::DenseSubMatrix<libMesh::Number> &KTv =
276 context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.v());
277 libMesh::DenseSubMatrix<libMesh::Number> *KTw = NULL;
279 if(this->_flow_vars.dim() == 3)
281 KTw = &context.get_elem_jacobian
282 (this->_temp_vars.T(), this->_flow_vars.w());
285 unsigned int n_qpoints = context.get_element_qrule().n_points();
287 for (
unsigned int qp=0; qp != n_qpoints; qp++)
289 libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
291 libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
292 libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
294 libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
295 context.fixed_interior_value( this->_flow_vars.v(), qp ) );
296 if( this->_flow_vars.dim() == 3 )
297 U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
301 libMesh::Real tau_E, RE_t;
302 libMesh::Real d_tau_E_d_rho, d_RE_t_dT;
303 libMesh::Gradient d_tau_E_dU;
306 libMesh::Real _k_qp = this->_k(context, qp);
308 if (compute_jacobian)
310 this->_stab_helper.compute_tau_energy_and_derivs
311 ( context, G, this->_rho, this->_Cp, _k_qp, U,
312 tau_E, d_tau_E_d_rho, d_tau_E_dU,
false );
313 this->_stab_helper.compute_res_energy_transient_and_derivs
314 ( context, qp, this->_rho, this->_Cp,
319 tau_E = this->_stab_helper.compute_tau_energy
320 ( context, G, this->_rho, this->_Cp, _k_qp, U,
false );
322 RE_t = this->_stab_helper.compute_res_energy_transient
323 ( context, qp, this->_rho, this->_Cp );
339 for (
unsigned int i=0; i != n_T_dofs; i++)
341 FT(i) -= tau_E*RE_t*( this->_rho*this->_Cp*U*T_gradphi[i][qp]
342 + _k_qp*(T_hessphi[i][qp](0,0) + T_hessphi[i][qp](1,1) + T_hessphi[i][qp](2,2))
344 if (compute_jacobian)
346 const libMesh::Real fixed_deriv =
347 context.get_fixed_solution_derivative();
349 for (
unsigned int j=0; j != n_T_dofs; ++j)
352 (tau_E*d_RE_t_dT)*T_phi[j][qp]*
353 ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
354 + _k_qp*(T_hessphi[i][qp](0,0) +
355 T_hessphi[i][qp](1,1) +
356 T_hessphi[i][qp](2,2))
359 for (
unsigned int j=0; j != n_u_dofs; ++j)
362 d_tau_E_dU(0)*u_phi[j][qp]*RE_t*
363 ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
364 + _k_qp*(T_hessphi[i][qp](0,0) +
365 T_hessphi[i][qp](1,1) +
366 T_hessphi[i][qp](2,2))
367 )*fixed_deriv*JxW[qp];
370 ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](0)*fixed_deriv
373 d_tau_E_dU(1)*u_phi[j][qp]*RE_t*
374 ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
375 + _k_qp*(T_hessphi[i][qp](0,0) +
376 T_hessphi[i][qp](1,1) +
377 T_hessphi[i][qp](2,2))
378 )*fixed_deriv*JxW[qp];
381 ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](1)*fixed_deriv
384 if(this->_flow_vars.dim() == 3)
386 for (
unsigned int j=0; j != n_u_dofs; ++j)
389 d_tau_E_dU(2)*u_phi[j][qp]*RE_t*
390 ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
391 + _k_qp*(T_hessphi[i][qp](0,0) +
392 T_hessphi[i][qp](1,1) +
393 T_hessphi[i][qp](2,2))
394 )*fixed_deriv*JxW[qp];
397 ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](2)*fixed_deriv
INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatTransferAdjointStabilization)
HeatTransferAdjointStabilization()
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...
virtual ~HeatTransferAdjointStabilization()
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.