GRINS-0.8.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-2017 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  ( bool compute_jacobian, AssemblyContext & context )
54  {
55  // The number of local degrees of freedom in each variable.
56  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
57 
58  // Element Jacobian * quadrature weights for interior integration.
59  const std::vector<libMesh::Real> &JxW =
60  context.get_element_fe(this->_temp_vars.T())->get_JxW();
61 
62  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
63  context.get_element_fe(this->_temp_vars.T())->get_dphi();
64 
65  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
66 
67  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
68 
69  unsigned int n_qpoints = context.get_element_qrule().n_points();
70 
71  for (unsigned int qp=0; qp != n_qpoints; qp++)
72  {
73  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
74  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
75 
76  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
77  context.interior_value( this->_flow_vars.v(), qp ) );
78  if( this->_flow_vars.dim() == 3 )
79  {
80  U(2) = context.interior_value( this->_flow_vars.w(), qp );
81  }
82 
83  // Compute Conductivity at this qp
84  libMesh::Real _k_qp = this->_k(context, qp);
85 
86  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, G, this->_rho, this->_Cp, _k_qp, U, this->_is_steady );
87 
88  libMesh::Real RE_s = this->_stab_helper.compute_res_energy_steady( context, qp, this->_rho, this->_Cp, _k_qp );
89 
90  for (unsigned int i=0; i != n_T_dofs; i++)
91  {
92  FT(i) += -tau_E*RE_s*this->_rho*this->_Cp*U*T_gradphi[i][qp]*JxW[qp];
93  }
94 
95  if( compute_jacobian )
96  {
97  libmesh_not_implemented();
98  }
99 
100  }
101  }
102 
103  template<class K>
105  AssemblyContext & context )
106  {
107  if( compute_jacobian )
108  libmesh_not_implemented();
109 
110  // The number of local degrees of freedom in each variable.
111  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
112 
113  // Element Jacobian * quadrature weights for interior integration.
114  const std::vector<libMesh::Real> &JxW =
115  context.get_element_fe(this->_temp_vars.T())->get_JxW();
116 
117  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
118  context.get_element_fe(this->_temp_vars.T())->get_dphi();
119 
120  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
121 
122  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
123 
124  unsigned int n_qpoints = context.get_element_qrule().n_points();
125 
126  for (unsigned int qp=0; qp != n_qpoints; qp++)
127  {
128  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
129  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
130 
131  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
132  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
133  if( this->_flow_vars.dim() == 3 )
134  {
135  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
136  }
137 
138  // Compute Conductivity at this qp
139  libMesh::Real _k_qp = this->_k(context, qp);
140 
141  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, G, this->_rho, this->_Cp, _k_qp, U, false );
142 
143  libMesh::Real RE_t = this->_stab_helper.compute_res_energy_transient( context, qp, this->_rho, this->_Cp );
144 
145  for (unsigned int i=0; i != n_T_dofs; i++)
146  {
147  FT(i) -= tau_E*RE_t*this->_rho*this->_Cp*U*T_gradphi[i][qp]*JxW[qp];
148  }
149 
150  }
151  }
152 
153 } // namespace GRINS
154 
155 // Instantiate
156 INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatTransferSPGSMStabilization);
INSTANTIATE_HEAT_TRANSFER_SUBCLASS(HeatTransferSPGSMStabilization)
GRINS namespace.
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 void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.

Generated on Tue Dec 19 2017 12:47:28 for GRINS-0.8.0 by  doxygen 1.8.9.1