GRINS-0.6.0
heat_transfer_spgsm_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-2015 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
27 #include "grins/assembly_context.h"
30 
31 // libMesh
32 #include "libmesh/quadrature.h"
33 
34 namespace GRINS
35 {
36 
37  template<class K>
39  const GetPot& input )
40  : HeatTransferStabilizationBase<K>(physics_name,input)
41  {
42  return;
43  }
44 
45  template<class K>
47  {
48  return;
49  }
50 
51  template<class K>
53  AssemblyContext& context,
54  CachedValues& /*cache*/ )
55  {
56 #ifdef GRINS_USE_GRVY_TIMERS
57  this->_timer->BeginTimer("HeatTransferSPGSMStabilization::element_time_derivative");
58 #endif
59 
60  // The number of local degrees of freedom in each variable.
61  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T_var()).size();
62 
63  // Element Jacobian * quadrature weights for interior integration.
64  const std::vector<libMesh::Real> &JxW =
65  context.get_element_fe(this->_temp_vars.T_var())->get_JxW();
66 
67  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
68  context.get_element_fe(this->_temp_vars.T_var())->get_dphi();
69 
70  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T_var()); // R_{T}
71 
72  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T_var());
73 
74  unsigned int n_qpoints = context.get_element_qrule().n_points();
75 
76  for (unsigned int qp=0; qp != n_qpoints; qp++)
77  {
78  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
79  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
80 
81  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u_var(), qp ),
82  context.interior_value( this->_flow_vars.v_var(), qp ) );
83  if( this->_dim == 3 )
84  {
85  U(2) = context.interior_value( this->_flow_vars.w_var(), qp );
86  }
87 
88  // Compute Conductivity at this qp
89  libMesh::Real _k_qp = this->_k(context, qp);
90 
91  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, G, this->_rho, this->_Cp, _k_qp, U, this->_is_steady );
92 
93  libMesh::Real RE_s = this->_stab_helper.compute_res_energy_steady( context, qp, this->_rho, this->_Cp, _k_qp );
94 
95  for (unsigned int i=0; i != n_T_dofs; i++)
96  {
97  FT(i) += -tau_E*RE_s*this->_rho*this->_Cp*U*T_gradphi[i][qp]*JxW[qp];
98  }
99 
100  if( compute_jacobian )
101  {
102  libmesh_not_implemented();
103  }
104 
105  }
106 
107 #ifdef GRINS_USE_GRVY_TIMERS
108  this->_timer->EndTimer("HeatTransferSPGSMStabilization::element_time_derivative");
109 #endif
110  return;
111  }
112 
113  template<class K>
114  void HeatTransferSPGSMStabilization<K>::mass_residual( bool /*compute_jacobian*/,
115  AssemblyContext& context,
116  CachedValues& /*cache*/ )
117  {
118 #ifdef GRINS_USE_GRVY_TIMERS
119  this->_timer->BeginTimer("HeatTransferSPGSMStabilization::mass_residual");
120 #endif
121 
122  // The number of local degrees of freedom in each variable.
123  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T_var()).size();
124 
125  // Element Jacobian * quadrature weights for interior integration.
126  const std::vector<libMesh::Real> &JxW =
127  context.get_element_fe(this->_temp_vars.T_var())->get_JxW();
128 
129  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
130  context.get_element_fe(this->_temp_vars.T_var())->get_dphi();
131 
132  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T_var()); // R_{T}
133 
134  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T_var());
135 
136  unsigned int n_qpoints = context.get_element_qrule().n_points();
137 
138  for (unsigned int qp=0; qp != n_qpoints; qp++)
139  {
140  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
141  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
142 
143  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u_var(), qp ),
144  context.fixed_interior_value( this->_flow_vars.v_var(), qp ) );
145  if( this->_dim == 3 )
146  {
147  U(2) = context.fixed_interior_value( this->_flow_vars.w_var(), qp );
148  }
149 
150  // Compute Conductivity at this qp
151  libMesh::Real _k_qp = this->_k(context, qp);
152 
153  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, G, this->_rho, this->_Cp, _k_qp, U, false );
154 
155  libMesh::Real RE_t = this->_stab_helper.compute_res_energy_transient( context, qp, this->_rho, this->_Cp );
156 
157  for (unsigned int i=0; i != n_T_dofs; i++)
158  {
159  FT(i) -= tau_E*RE_t*this->_rho*this->_Cp*U*T_gradphi[i][qp]*JxW[qp];
160  }
161 
162  }
163 
164 #ifdef GRINS_USE_GRVY_TIMERS
165  this->_timer->EndTimer("HeatTransferSPGSMStabilization::mass_residual");
166 #endif
167  return;
168  }
169 
170 } // namespace GRINS
171 
172 // Instantiate
173 INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatTransferSPGSMStabilization);
INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatTransferSPGSMStabilization)
GRINS namespace.
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...
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.

Generated on Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1