GRINS-0.7.0
low_mach_navier_stokes_braack_stab.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2016 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 // This class
28 
29 // GRINS
30 #include "grins/assembly_context.h"
34 
35 // libMesh
36 #include "libmesh/quadrature.h"
37 
38 namespace GRINS
39 {
40 
41  template<class Mu, class SH, class TC>
43  const GetPot& input )
44  : LowMachNavierStokesStabilizationBase<Mu,SH,TC>(physics_name,input)
45  {
46  return;
47  }
48 
49  template<class Mu, class SH, class TC>
51  {
52  return;
53  }
54 
55  template<class Mu, class SH, class TC>
57  AssemblyContext& context,
58  CachedValues& /*cache*/ )
59  {
60 #ifdef GRINS_USE_GRVY_TIMERS
61  this->_timer->BeginTimer("LowMachNavierStokesBraackStabilization::element_time_derivative");
62 #endif
63 
64  this->assemble_continuity_time_deriv( compute_jacobian, context );
65  this->assemble_momentum_time_deriv( compute_jacobian, context );
66  this->assemble_energy_time_deriv( compute_jacobian, context );
67 
68 #ifdef GRINS_USE_GRVY_TIMERS
69  this->_timer->EndTimer("LowMachNavierStokesBraackStabilization::element_time_derivative");
70 #endif
71  return;
72  }
73 
74  template<class Mu, class SH, class TC>
76  AssemblyContext& context,
77  CachedValues& /*cache*/ )
78  {
79 #ifdef GRINS_USE_GRVY_TIMERS
80  this->_timer->BeginTimer("LowMachNavierStokesBraackStabilization::mass_residual");
81 #endif
82 
83  this->assemble_continuity_mass_residual( compute_jacobian, context );
84  this->assemble_momentum_mass_residual( compute_jacobian, context );
85  this->assemble_energy_mass_residual( compute_jacobian, context );
86 
87 #ifdef GRINS_USE_GRVY_TIMERS
88  this->_timer->EndTimer("LowMachNavierStokesBraackStabilization::mass_residual");
89 #endif
90  return;
91  }
92 
93  template<class Mu, class SH, class TC>
95  AssemblyContext& context )
96  {
97  // The number of local degrees of freedom in each variable.
98  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
99 
100  // Element Jacobian * quadrature weights for interior integration.
101  const std::vector<libMesh::Real> &JxW =
102  context.get_element_fe(this->_flow_vars.u())->get_JxW();
103 
104  // The pressure shape functions at interior quadrature points.
105  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
106  context.get_element_fe(this->_press_var.p())->get_dphi();
107 
108  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
109 
110  unsigned int n_qpoints = context.get_element_qrule().n_points();
111 
112  for (unsigned int qp=0; qp != n_qpoints; qp++)
113  {
114  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
115 
116  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
117  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
118 
119  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
120  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
121 
122  libMesh::Real mu = this->_mu(T);
123  libMesh::Real k = this->_k(T);
124  libMesh::Real cp = this->_cp(T);
125 
126  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
127  context.interior_value( this->_flow_vars.v(), qp ) );
128  if( this->_dim == 3 )
129  U(2) = context.interior_value( this->_flow_vars.w(), qp ); // w
130 
131  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
132  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, this->_is_steady );
133 
134  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
135  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
136 
137  // Now a loop over the pressure degrees of freedom. This
138  // computes the contributions of the continuity equation.
139  for (unsigned int i=0; i != n_p_dofs; i++)
140  {
141  Fp(i) += ( tau_M*RM_s*p_dphi[i][qp]
142  + tau_E*RE_s*(U*p_dphi[i][qp])/T )*JxW[qp];
143  }
144 
145  }
146 
147  return;
148  }
149 
150  template<class Mu, class SH, class TC>
152  AssemblyContext& context )
153  {
154  // The number of local degrees of freedom in each variable.
155  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
156 
157  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
158  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
159  if (this->_dim == 3)
160  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
161 
162  // Element Jacobian * quadrature weights for interior integration.
163  const std::vector<libMesh::Real> &JxW =
164  context.get_element_fe(this->_flow_vars.u())->get_JxW();
165 
166  // The velocity shape function gradients at interior quadrature points.
167  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
168  context.get_element_fe(this->_flow_vars.u())->get_dphi();
169 
170  const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
171  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
172 
173  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
174  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
175  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
176 
177  if( this->_dim == 3 )
178  {
179  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
180  }
181 
182  unsigned int n_qpoints = context.get_element_qrule().n_points();
183 
184  for (unsigned int qp=0; qp != n_qpoints; qp++)
185  {
186  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
187  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
188 
189  libMesh::Real mu = this->_mu(T);
190 
191  libMesh::RealGradient U( context.interior_value(this->_flow_vars.u(), qp),
192  context.interior_value(this->_flow_vars.v(), qp) );
193 
194  libMesh::RealGradient grad_u = context.interior_gradient(this->_flow_vars.u(), qp);
195  libMesh::RealGradient grad_v = context.interior_gradient(this->_flow_vars.v(), qp);
196  libMesh::RealGradient grad_w;
197 
198  if( this->_dim == 3 )
199  {
200  U(2) = context.interior_value(this->_flow_vars.w(), qp);
201  grad_w = context.interior_gradient(this->_flow_vars.w(), qp);
202  }
203 
204  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
205 
206  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
207  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
208 
209  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
210  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
211 
212  libMesh::Real RC_s = this->compute_res_continuity_steady( context, qp );
213  libMesh::RealGradient RM_s = this->compute_res_momentum_steady( context, qp );
214 
215  for (unsigned int i=0; i != n_u_dofs; i++)
216  {
217  Fu(i) += ( tau_C*RC_s*u_gradphi[i][qp](0)
218  + tau_M*RM_s(0)*rho*U*u_gradphi[i][qp]
219  + mu*tau_M*RM_s(0)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1)
220  + u_hessphi[i][qp](0,0) + u_hessphi[i][qp](0,1)
221  - 2.0/3.0*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,0)) )
222  )*JxW[qp];
223 
224  Fv(i) += ( tau_C*RC_s*u_gradphi[i][qp](1)
225  + tau_M*RM_s(1)*rho*U*u_gradphi[i][qp]
226  + mu*tau_M*RM_s(1)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1)
227  + u_hessphi[i][qp](1,0) + u_hessphi[i][qp](1,1)
228  - 2.0/3.0*(u_hessphi[i][qp](0,1) + u_hessphi[i][qp](1,1)) )
229  )*JxW[qp];
230 
231  if( this->_dim == 3 )
232  {
233  Fu(i) += mu*tau_M*RM_s(0)*(u_hessphi[i][qp](2,2) + u_hessphi[i][qp](0,2)
234  - 2.0/3.0*u_hessphi[i][qp](2,0))*JxW[qp];
235 
236  Fv(i) += mu*tau_M*RM_s(1)*(u_hessphi[i][qp](2,2) + u_hessphi[i][qp](1,2)
237  - 2.0/3.0*u_hessphi[i][qp](2,1))*JxW[qp];
238 
239  (*Fw)(i) += ( tau_C*RC_s*u_gradphi[i][qp](2)
240  + tau_M*RM_s(2)*rho*U*u_gradphi[i][qp]
241  + mu*tau_M*RM_s(2)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2)
242  + u_hessphi[i][qp](2,0) + u_hessphi[i][qp](2,1) + u_hessphi[i][qp](2,2)
243  - 2.0/3.0*(u_hessphi[i][qp](0,2) + u_hessphi[i][qp](1,2)
244  + u_hessphi[i][qp](2,2))
245  )
246  )*JxW[qp];
247  }
248  }
249 
250  }
251  return;
252  }
253 
254  template<class Mu, class SH, class TC>
256  AssemblyContext& context )
257  {
258  // The number of local degrees of freedom in each variable.
259  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
260 
261  // Element Jacobian * quadrature weights for interior integration.
262  const std::vector<libMesh::Real> &JxW =
263  context.get_element_fe(this->_temp_vars.T())->get_JxW();
264 
265  // The temperature shape functions gradients at interior quadrature points.
266  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
267  context.get_element_fe(this->_temp_vars.T())->get_dphi();
268 
269  const std::vector<std::vector<libMesh::RealTensor> >& T_hessphi =
270  context.get_element_fe(this->_temp_vars.T())->get_d2phi();
271 
272  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
273 
274  unsigned int n_qpoints = context.get_element_qrule().n_points();
275 
276  for (unsigned int qp=0; qp != n_qpoints; qp++)
277  {
278  libMesh::Number u, v;
279  u = context.interior_value(this->_flow_vars.u(), qp);
280  v = context.interior_value(this->_flow_vars.v(), qp);
281 
282  libMesh::Gradient grad_T = context.interior_gradient(this->_temp_vars.T(), qp);
283 
284  libMesh::NumberVectorValue U(u,v);
285  if (this->_dim == 3)
286  U(2) = context.interior_value(this->_flow_vars.w(), qp);
287 
288  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
289  libMesh::Real rho = this->rho( T, this->get_p0_steady( context, qp ) );
290 
291  libMesh::Real k = this->_k(T);
292  libMesh::Real cp = this->_cp(T);
293 
294  libMesh::Number rho_cp = rho*this->_cp(T);
295 
296  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
297 
298  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
299  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
300 
301  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, this->_is_steady );
302 
303  libMesh::Real RE_s = this->compute_res_energy_steady( context, qp );
304 
305  for (unsigned int i=0; i != n_T_dofs; i++)
306  {
307  FT(i) += ( rho_cp*tau_E*RE_s*U*T_gradphi[i][qp]
308  + tau_E*RE_s*k*(T_hessphi[i][qp](0,0) + T_hessphi[i][qp](1,1) + T_hessphi[i][qp](2,2) )
309  )*JxW[qp];
310  }
311 
312  }
313 
314  return;
315  }
316 
317  template<class Mu, class SH, class TC>
319  AssemblyContext& context )
320  {
321  // The number of local degrees of freedom in each variable.
322  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
323 
324  // Element Jacobian * quadrature weights for interior integration.
325  const std::vector<libMesh::Real> &JxW =
326  context.get_element_fe(this->_flow_vars.u())->get_JxW();
327 
328  // The pressure shape functions at interior quadrature points.
329  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
330  context.get_element_fe(this->_press_var.p())->get_dphi();
331 
332  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
333 
334  unsigned int n_qpoints = context.get_element_qrule().n_points();
335 
336  for (unsigned int qp=0; qp != n_qpoints; qp++)
337  {
338  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
339 
340  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
341  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
342 
343  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
344  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
345 
346  libMesh::Real mu = this->_mu(T);
347  libMesh::Real k = this->_k(T);
348  libMesh::Real cp = this->_cp(T);
349 
350  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
351  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
352  if( this->_dim == 3 )
353  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
354 
355  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
356  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
357 
358  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
359  libMesh::Real RE_t = this->compute_res_energy_transient( context, qp );
360 
361  // Now a loop over the pressure degrees of freedom. This
362  // computes the contributions of the continuity equation.
363  for (unsigned int i=0; i != n_p_dofs; i++)
364  {
365  Fp(i) += ( tau_M*RM_t*p_dphi[i][qp]
366  + tau_E*RE_t*(U*p_dphi[i][qp])/T
367  )*JxW[qp];
368  }
369  }
370 
371  return;
372  }
373 
374  template<class Mu, class SH, class TC>
376  AssemblyContext& context )
377  {
378  // The number of local degrees of freedom in each variable.
379  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
380 
381  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
382  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
383  if (this->_dim == 3)
384  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
385 
386  // Element Jacobian * quadrature weights for interior integration.
387  const std::vector<libMesh::Real> &JxW =
388  context.get_element_fe(this->_flow_vars.u())->get_JxW();
389 
390  // The velocity shape function gradients at interior quadrature points.
391  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
392  context.get_element_fe(this->_flow_vars.u())->get_dphi();
393 
394  const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
395  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
396 
397  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
398  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
399  libMesh::DenseSubVector<libMesh::Real>* Fw = NULL;
400 
401  if( this->_dim == 3 )
402  {
403  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
404  }
405 
406  unsigned int n_qpoints = context.get_element_qrule().n_points();
407  for (unsigned int qp=0; qp != n_qpoints; qp++)
408  {
409  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
410  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
411 
412  libMesh::Real mu = this->_mu(T);
413 
414  libMesh::RealGradient U( context.fixed_interior_value(this->_flow_vars.u(), qp),
415  context.fixed_interior_value(this->_flow_vars.v(), qp) );
416 
417  libMesh::RealGradient grad_u = context.fixed_interior_gradient(this->_flow_vars.u(), qp);
418  libMesh::RealGradient grad_v = context.fixed_interior_gradient(this->_flow_vars.v(), qp);
419  libMesh::RealGradient grad_w;
420 
421  if( this->_dim == 3 )
422  {
423  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
424  grad_w = context.fixed_interior_gradient(this->_flow_vars.w(), qp);
425  }
426 
427  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
428 
429  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
430  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
431 
432  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
433  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
434 
435  libMesh::Real RC_t = this->compute_res_continuity_transient( context, qp );
436  libMesh::RealGradient RM_t = this->compute_res_momentum_transient( context, qp );
437 
438  for (unsigned int i=0; i != n_u_dofs; i++)
439  {
440  Fu(i) += ( tau_C*RC_t*u_gradphi[i][qp](0)
441  + tau_M*RM_t(0)*rho*U*u_gradphi[i][qp]
442  + mu*tau_M*RM_t(0)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1)
443  + u_hessphi[i][qp](0,0) + u_hessphi[i][qp](0,1)
444  - 2.0/3.0*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,0)) )
445  )*JxW[qp];
446 
447  Fv(i) += ( tau_C*RC_t*u_gradphi[i][qp](1)
448  + tau_M*RM_t(1)*rho*U*u_gradphi[i][qp]
449  + mu*tau_M*RM_t(1)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1)
450  + u_hessphi[i][qp](1,0) + u_hessphi[i][qp](1,1)
451  - 2.0/3.0*(u_hessphi[i][qp](0,1) + u_hessphi[i][qp](1,1)) )
452  )*JxW[qp];
453 
454  if( this->_dim == 3 )
455  {
456  Fu(i) -= mu*tau_M*RM_t(0)*(u_hessphi[i][qp](2,2) + u_hessphi[i][qp](0,2)
457  - 2.0/3.0*u_hessphi[i][qp](2,0))*JxW[qp];
458 
459  Fv(i) -= mu*tau_M*RM_t(1)*(u_hessphi[i][qp](2,2) + u_hessphi[i][qp](1,2)
460  - 2.0/3.0*u_hessphi[i][qp](2,1))*JxW[qp];
461 
462  (*Fw)(i) -= ( tau_C*RC_t*u_gradphi[i][qp](2)
463  + tau_M*RM_t(2)*rho*U*u_gradphi[i][qp]
464  + mu*tau_M*RM_t(2)*(u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2)
465  + u_hessphi[i][qp](2,0) + u_hessphi[i][qp](2,1) + u_hessphi[i][qp](2,2)
466  - 2.0/3.0*(u_hessphi[i][qp](0,2) + u_hessphi[i][qp](1,2)
467  + u_hessphi[i][qp](2,2))
468  )
469  )*JxW[qp];
470  }
471  }
472 
473  }
474  return;
475  }
476 
477  template<class Mu, class SH, class TC>
479  AssemblyContext& context )
480  {
481  // The number of local degrees of freedom in each variable.
482  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
483 
484  // Element Jacobian * quadrature weights for interior integration.
485  const std::vector<libMesh::Real> &JxW =
486  context.get_element_fe(this->_temp_vars.T())->get_JxW();
487 
488  // The temperature shape functions gradients at interior quadrature points.
489  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
490  context.get_element_fe(this->_temp_vars.T())->get_dphi();
491 
492  const std::vector<std::vector<libMesh::RealTensor> >& T_hessphi =
493  context.get_element_fe(this->_temp_vars.T())->get_d2phi();
494 
495  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
496 
497  unsigned int n_qpoints = context.get_element_qrule().n_points();
498 
499  for (unsigned int qp=0; qp != n_qpoints; qp++)
500  {
501  libMesh::Number u, v;
502  u = context.fixed_interior_value(this->_flow_vars.u(), qp);
503  v = context.fixed_interior_value(this->_flow_vars.v(), qp);
504 
505  libMesh::Gradient grad_T = context.fixed_interior_gradient(this->_temp_vars.T(), qp);
506 
507  libMesh::NumberVectorValue U(u,v);
508  if (this->_dim == 3)
509  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp); // w
510 
511  libMesh::Real T = context.fixed_interior_value( this->_temp_vars.T(), qp );
512  libMesh::Real rho = this->rho( T, this->get_p0_transient( context, qp ) );
513 
514  libMesh::Real k = this->_k(T);
515  libMesh::Real cp = this->_cp(T);
516 
517  libMesh::Number rho_cp = rho*cp;
518 
519  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
520 
521  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
522  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
523 
524  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
525 
526  libMesh::Real RE_t = this->compute_res_energy_transient( context, qp );
527 
528  for (unsigned int i=0; i != n_T_dofs; i++)
529  {
530  FT(i) += ( rho_cp*tau_E*RE_t*U*T_gradphi[i][qp]
531  + tau_E*RE_t*k*(T_hessphi[i][qp](0,0) + T_hessphi[i][qp](1,1) + T_hessphi[i][qp](2,2))
532  )*JxW[qp];
533  }
534 
535  }
536 
537  return;
538  }
539 
540 
541 } // namespace GRINS
542 
543 // Instantiate
void assemble_energy_mass_residual(bool compute_jacobian, AssemblyContext &context)
Adds VMS-based stabilization to LowMachNavierStokes physics class.
void assemble_energy_time_deriv(bool compute_jacobian, AssemblyContext &context)
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
GRINS namespace.
void assemble_momentum_mass_residual(bool compute_jacobian, AssemblyContext &context)
void assemble_continuity_time_deriv(bool compute_jacobian, AssemblyContext &context)
Adds VMS-based stabilization to LowMachNavierStokes physics class.
void assemble_continuity_mass_residual(bool compute_jacobian, AssemblyContext &context)
void assemble_momentum_time_deriv(bool compute_jacobian, AssemblyContext &context)
virtual void mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...

Generated on Thu Jun 2 2016 21:52:28 for GRINS-0.7.0 by  doxygen 1.8.10