GRINS-0.8.0
inc_navier_stokes_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 // This class
27 
28 // GRINS
29 #include "grins/assembly_context.h"
31 
32 //libMesh
33 #include "libmesh/quadrature.h"
34 
35 namespace GRINS
36 {
37 
38  template<class Mu>
40  const GetPot& input )
41  : IncompressibleNavierStokesStabilizationBase<Mu>(physics_name,input)
42  {}
43 
44  template<class Mu>
46  ( bool compute_jacobian,
47  AssemblyContext & context )
48  {
49  // The number of local degrees of freedom in each variable.
50  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
51 
52  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
53  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
54 
55 
56  if (this->_flow_vars.dim() == 3)
57  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
58 
59  // Element Jacobian * quadrature weights for interior integration.
60  const std::vector<libMesh::Real> &JxW =
61  context.get_element_fe(this->_flow_vars.u())->get_JxW();
62 
63  // The velocity shape function gradients at interior quadrature points.
64  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
65  context.get_element_fe(this->_flow_vars.u())->get_dphi();
66 
67  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
68  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
69  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
70  if(this->_flow_vars.dim() == 3)
71  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
72 
73  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
74 
75  unsigned int n_qpoints = context.get_element_qrule().n_points();
76 
77  for (unsigned int qp=0; qp != n_qpoints; qp++)
78  {
79  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
80  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
81 
82  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
83  context.interior_value( this->_flow_vars.v(), qp ) );
84  if( this->_flow_vars.dim() == 3 )
85  {
86  U(2) = context.interior_value( this->_flow_vars.w(), qp );
87  }
88 
89  // Compute the viscosity at this qp
90  libMesh::Real _mu_qp = this->_mu(context, qp);
91 
92  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
93  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
94 
95  libMesh::RealGradient RM_s = this->_stab_helper.compute_res_momentum_steady( context, qp, this->_rho, _mu_qp );
96  libMesh::Real RC = this->_stab_helper.compute_res_continuity( context, qp );
97 
98  for (unsigned int i=0; i != n_u_dofs; i++)
99  {
100  Fu(i) += ( - tau_C*RC*u_gradphi[i][qp](0)
101  - tau_M*RM_s(0)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
102 
103  Fv(i) += ( - tau_C*RC*u_gradphi[i][qp](1)
104  - tau_M*RM_s(1)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
105 
106  if( this->_flow_vars.dim() == 3 )
107  {
108  (*Fw)(i) += ( - tau_C*RC*u_gradphi[i][qp](2)
109  - tau_M*RM_s(2)*this->_rho*U*u_gradphi[i][qp] )*JxW[qp];
110  }
111  }
112 
113  if( compute_jacobian )
114  libmesh_not_implemented();
115  }
116  }
117 
118  template<class Mu>
120  ( bool compute_jacobian, AssemblyContext & context )
121  {
122  // The number of local degrees of freedom in each variable.
123  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
124 
125  // Element Jacobian * quadrature weights for interior integration.
126  const std::vector<libMesh::Real> &JxW =
127  context.get_element_fe(this->_flow_vars.u())->get_JxW();
128 
129  // The pressure shape functions at interior quadrature points.
130  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
131  context.get_element_fe(this->_press_var.p())->get_dphi();
132 
133  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
134 
135  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
136 
137  unsigned int n_qpoints = context.get_element_qrule().n_points();
138 
139  for (unsigned int qp=0; qp != n_qpoints; qp++)
140  {
141  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
142  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
143 
144  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
145  context.interior_value( this->_flow_vars.v(), qp ) );
146  if( this->_flow_vars.dim() == 3 )
147  U(2) = context.interior_value( this->_flow_vars.w(), qp );
148 
149  // Compute the viscosity at this qp
150  libMesh::Real _mu_qp = this->_mu(context, qp);
151 
152  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, this->_is_steady );
153 
154  libMesh::RealGradient RM_s = this->_stab_helper.compute_res_momentum_steady( context, qp, this->_rho, _mu_qp );
155 
156  for (unsigned int i=0; i != n_p_dofs; i++)
157  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*JxW[qp];
158 
159  if( compute_jacobian )
160  libmesh_not_implemented();
161  }
162  }
163 
164  template<class Mu>
166  ( bool compute_jacobian, AssemblyContext & context )
167  {
168  // The number of local degrees of freedom in each variable.
169  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
170  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
171 
172  // Element Jacobian * quadrature weights for interior integration.
173  const std::vector<libMesh::Real> &JxW =
174  context.get_element_fe(this->_flow_vars.u())->get_JxW();
175 
176  // The pressure shape functions at interior quadrature points.
177  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
178  context.get_element_fe(this->_press_var.p())->get_dphi();
179 
180  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
181  context.get_element_fe(this->_flow_vars.u())->get_dphi();
182 
183  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{p}
184  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{p}
185  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
186 
187 
188  if(this->_flow_vars.dim() == 3)
189  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
190 
191  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
192 
193  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
194 
195  unsigned int n_qpoints = context.get_element_qrule().n_points();
196 
197  for (unsigned int qp=0; qp != n_qpoints; qp++)
198  {
199  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
200  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
201 
202  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
203  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
204  // Compute the viscosity at this qp
205  libMesh::Real _mu_qp = this->_mu(context, qp);
206 
207  if( this->_flow_vars.dim() == 3 )
208  {
209  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
210  }
211 
212  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, _mu_qp, false );
213 
214  libMesh::RealGradient RM_t = this->_stab_helper.compute_res_momentum_transient( context, qp, this->_rho );
215 
216  for (unsigned int i=0; i != n_p_dofs; i++)
217  {
218  Fp(i) -= tau_M*RM_t*p_dphi[i][qp]*JxW[qp];
219  }
220 
221  for (unsigned int i=0; i != n_u_dofs; i++)
222  {
223  Fu(i) -= tau_M*RM_t(0)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
224 
225  Fv(i) -= tau_M*RM_t(1)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
226 
227  if( this->_flow_vars.dim() == 3 )
228  {
229  (*Fw)(i) -= tau_M*RM_t(2)*this->_rho*U*u_gradphi[i][qp]*JxW[qp];
230  }
231  }
232 
233  if( compute_jacobian )
234  {
235  libmesh_not_implemented();
236  }
237 
238  }
239  }
240 
241 } // end namespace GRINS
242 
243 // Instantiate
244 INSTANTIATE_INC_NS_SUBCLASS(IncompressibleNavierStokesSPGSMStabilization);
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.
GRINS namespace.
INSTANTIATE_INC_NS_SUBCLASS(IncompressibleNavierStokesSPGSMStabilization)
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context)
Constraint 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