GRINS-0.8.0
boussinesq_buoyancy_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
26 #include "grins_config.h"
28 
29 // GRINS
30 #include "grins/assembly_context.h"
33 
34 // libMesh
35 #include "libmesh/getpot.h"
36 #include "libmesh/fem_system.h"
37 #include "libmesh/quadrature.h"
38 
39 namespace GRINS
40 {
41  template<class Mu>
42  BoussinesqBuoyancySPGSMStabilization<Mu>::BoussinesqBuoyancySPGSMStabilization( const std::string& physics_name, const GetPot& input )
43  : BoussinesqBuoyancyBase(physics_name,input),
44  _flow_stab_helper(physics_name+"FlowStabHelper", input),
45  _temp_stab_helper(physics_name+"TempStabHelper", input),
46  _Cp(0.0),
47  _k(1.0),
48  _mu(input,MaterialsParsing::material_name(input,PhysicsNaming::boussinesq_buoyancy()))
49  {
51 
52  this->set_parameter
53  (_k, input, "Physics/"+PhysicsNaming::heat_transfer()+"/k", _k);
54  }
55 
56  template<class Mu>
58  {
59  return;
60  }
61 
62  template<class Mu>
64  ( bool compute_jacobian, AssemblyContext & context )
65  {
66  // The number of local degrees of freedom in each variable.
67  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
68 
69  // Element Jacobian * quadrature weights for interior integration.
70  const std::vector<libMesh::Real> &JxW =
71  context.get_element_fe(_flow_vars.u())->get_JxW();
72 
73  /*
74  const std::vector<std::vector<libMesh::Real> >& u_phi =
75  context.get_element_fe(this->_flow_vars.u())->get_phi();
76  */
77 
78  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
79  context.get_element_fe(this->_flow_vars.u())->get_dphi();
80 
81  // Get residuals
82  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u()); // R_{u}
83  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v()); // R_{v}
84  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
85 
86  if(this->_flow_vars.dim() == 3)
87  {
88  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
89  }
90 
91  // Now we will build the element Jacobian and residual.
92  // Constructing the residual requires the solution and its
93  // gradient from the previous timestep. This must be
94  // calculated at each quadrature point by summing the
95  // solution degree-of-freedom values by the appropriate
96  // weight functions.
97  unsigned int n_qpoints = context.get_element_qrule().n_points();
98 
99  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
100 
101  for (unsigned int qp=0; qp != n_qpoints; qp++)
102  {
103  libMesh::RealGradient g = this->_flow_stab_helper.compute_g( fe, context, qp );
104  libMesh::RealTensor G = this->_flow_stab_helper.compute_G( fe, context, qp );
105 
106  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
107  context.interior_value( this->_flow_vars.v(), qp ) );
108  if( this->_flow_vars.dim() == 3 )
109  {
110  U(2) = context.interior_value( this->_flow_vars.w(), qp );
111  }
112 
113  // Compute the viscosity at this qp
114  libMesh::Real mu_qp = this->_mu(context, qp);
115 
116  libMesh::Real tau_M = this->_flow_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, mu_qp, this->_is_steady );
117 
118  //libMesh::Real tau_E = this->_temp_stab_helper.compute_tau_energy( context, G, _rho, _Cp, _k, U, this->_is_steady );
119 
120  //libMesh::Real RE = this->_temp_stab_helper.compute_res_energy_steady( context, qp, _rho, _Cp, _k );
121 
122  // Compute the solution & its gradient at the old Newton iterate.
123  libMesh::Number T;
124  T = context.interior_value(_temp_vars.T(), qp);
125 
126  libMesh::RealGradient residual = _rho*_beta_T*(T-_T_ref)*_g;
127 
128  for (unsigned int i=0; i != n_u_dofs; i++)
129  {
130  Fu(i) += ( -tau_M*residual(0)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
131  // + _rho*_beta_T*tau_E*RE*_g(0)*u_phi[i][qp] )*JxW[qp];
132 
133  Fv(i) += ( -tau_M*residual(1)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
134  // + _rho*_beta_T*tau_E*RE*_g(1)*u_phi[i][qp] )*JxW[qp];
135 
136  if (this->_flow_vars.dim() == 3)
137  {
138  (*Fw)(i) += ( -tau_M*residual(2)*_rho*U*u_gradphi[i][qp] )*JxW[qp];
139  // + _rho*_beta_T*tau_E*RE*_g(2)*u_phi[i][qp] )*JxW[qp];
140  }
141 
142  if (compute_jacobian)
143  {
144  libmesh_not_implemented();
145  } // End compute_jacobian check
146 
147  } // End i dof loop
148  } // End quadrature loop
149  }
150 
151  template<class Mu>
153  ( bool compute_jacobian,
154  AssemblyContext & context )
155  {
156  // The number of local degrees of freedom in each variable.
157  const unsigned int n_p_dofs = context.get_dof_indices(_press_var.p()).size();
158 
159  // Element Jacobian * quadrature weights for interior integration.
160  const std::vector<libMesh::Real> &JxW =
161  context.get_element_fe(_flow_vars.u())->get_JxW();
162 
163  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
164  context.get_element_fe(this->_press_var.p())->get_dphi();
165 
166  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
167 
168  // Now we will build the element Jacobian and residual.
169  // Constructing the residual requires the solution and its
170  // gradient from the previous timestep. This must be
171  // calculated at each quadrature point by summing the
172  // solution degree-of-freedom values by the appropriate
173  // weight functions.
174  unsigned int n_qpoints = context.get_element_qrule().n_points();
175 
176  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
177 
178  for (unsigned int qp=0; qp != n_qpoints; qp++)
179  {
180  libMesh::RealGradient g = this->_flow_stab_helper.compute_g( fe, context, qp );
181  libMesh::RealTensor G = this->_flow_stab_helper.compute_G( fe, context, qp );
182 
183  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
184  context.interior_value( this->_flow_vars.v(), qp ) );
185  if( this->_flow_vars.dim() == 3 )
186  U(2) = context.interior_value( this->_flow_vars.w(), qp );
187 
188  // Compute the viscosity at this qp
189  libMesh::Real mu_qp = this->_mu(context, qp);
190 
191  libMesh::Real tau_M = this->_flow_stab_helper.compute_tau_momentum( context, qp, g, G, this->_rho, U, mu_qp, this->_is_steady );
192 
193  // Compute the solution & its gradient at the old Newton iterate.
194  libMesh::Number T;
195  T = context.interior_value(_temp_vars.T(), qp);
196 
197  libMesh::RealGradient residual = _rho*_beta_T*(T-_T_ref)*_g;
198 
199  // First, an i-loop over the velocity degrees of freedom.
200  // We know that n_u_dofs == n_v_dofs so we can compute contributions
201  // for both at the same time.
202  for (unsigned int i=0; i != n_p_dofs; i++)
203  {
204  Fp(i) += -tau_M*residual*p_dphi[i][qp]*JxW[qp];
205 
206  if (compute_jacobian)
207  {
208  libmesh_not_implemented();
209  } // End compute_jacobian check
210 
211  } // End i dof loop
212 
213  } // End quadrature loop
214  }
215 
216  template<class Mu>
218  AssemblyContext& /*context*/ )
219  {
220  /*
221  // The number of local degrees of freedom in each variable.
222  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
223 
224  // Element Jacobian * quadrature weights for interior integration.
225  const std::vector<libMesh::Real> &JxW =
226  context.get_element_fe(_flow_vars.u())->get_JxW();
227 
228  const std::vector<std::vector<libMesh::Real> >& u_phi =
229  context.get_element_fe(this->_flow_vars.u())->get_phi();
230 
231  // Get residuals
232  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u()); // R_{u}
233  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v()); // R_{v}
234  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
235  if(this->_flow_vars.dim() == 3)
236  {
237  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
238  }
239 
240  // Now we will build the element Jacobian and residual.
241  // Constructing the residual requires the solution and its
242  // gradient from the previous timestep. This must be
243  // calculated at each quadrature point by summing the
244  // solution degree-of-freedom values by the appropriate
245  // weight functions.
246  unsigned int n_qpoints = context.get_element_qrule().n_points();
247 
248  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
249 
250  for (unsigned int qp=0; qp != n_qpoints; qp++)
251  {
252  libMesh::RealGradient g = this->_flow_stab_helper.compute_g( fe, context, qp );
253  libMesh::RealTensor G = this->_flow_stab_helper.compute_G( fe, context, qp );
254 
255  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
256  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
257  if( this->_flow_vars.dim() == 3 )
258  {
259  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
260  }
261 
262  libMesh::Real tau_E = this->_temp_stab_helper.compute_tau_energy( context, G, _rho, _Cp, _k, U, false );
263 
264  libMesh::Real RE = this->_temp_stab_helper.compute_res_energy_transient( context, qp, _rho, _Cp );
265 
266 
267  for (unsigned int i=0; i != n_u_dofs; i++)
268  {
269  Fu(i) += -_rho*_beta_T*tau_E*RE*_g(0)*u_phi[i][qp]*JxW[qp];
270 
271  Fv(i) += -_rho*_beta_T*tau_E*RE*_g(1)*u_phi[i][qp]*JxW[qp];
272 
273  if (this->_flow_vars.dim() == 3)
274  {
275  (*Fw)(i) += -_rho*_beta_T*tau_E*RE*_g(2)*u_phi[i][qp]*JxW[qp];
276  }
277 
278  if (compute_jacobian)
279  {
280  libmesh_not_implemented();
281  } // End compute_jacobian check
282 
283  } // End i dof loop
284  } // End quadrature loop
285  */
286  }
287 
288  template<class Mu>
290  ( const std::string & param_name,
292  const
293  {
294  ParameterUser::register_parameter(param_name, param_pointer);
295  _mu.register_parameter(param_name, param_pointer);
296  }
297 
298 
299 } // namespace GRINS
300 
301 // Instantiate
302 INSTANTIATE_INC_NS_SUBCLASS(BoussinesqBuoyancySPGSMStabilization);
virtual void set_parameter(libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
Each subclass can simultaneously read a parameter value from.
static PhysicsName heat_transfer()
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...
INSTANTIATE_INC_NS_SUBCLASS(BoussinesqBuoyancySPGSMStabilization)
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context)
Constraint part(s) of physics for element interiors.
GRINS namespace.
static void read_specific_heat(const std::string &core_physics_name, const GetPot &input, ParameterUser &params, libMesh::Real &cp)
Helper function to reading scalar specific heat from input.
Helper functions for parsing material properties.
static PhysicsName boussinesq_buoyancy()
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.

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