GRINS-0.6.0
boussinesq_buoyancy.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_config.h"
29 
30 // GRINS
31 #include "grins/assembly_context.h"
32 
33 // libMesh
34 #include "libmesh/getpot.h"
35 #include "libmesh/fem_system.h"
36 #include "libmesh/quadrature.h"
37 
38 namespace GRINS
39 {
40 
41  BoussinesqBuoyancy::BoussinesqBuoyancy( const std::string& physics_name, const GetPot& input )
42  : BoussinesqBuoyancyBase(physics_name,input)
43  {
44  return;
45  }
46 
48  {
49  return;
50  }
51 
52  void BoussinesqBuoyancy::element_time_derivative( bool compute_jacobian,
53  AssemblyContext& context,
54  CachedValues& /*cache*/ )
55  {
56 #ifdef GRINS_USE_GRVY_TIMERS
57  this->_timer->BeginTimer("BoussinesqBuoyancy::element_time_derivative");
58 #endif
59 
60  // The number of local degrees of freedom in each variable.
61  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u_var()).size();
62  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T_var()).size();
63 
64  // Element Jacobian * quadrature weights for interior integration.
65  const std::vector<libMesh::Real> &JxW =
66  context.get_element_fe(_flow_vars.u_var())->get_JxW();
67 
68  // The velocity shape functions at interior quadrature points.
69  const std::vector<std::vector<libMesh::Real> >& vel_phi =
70  context.get_element_fe(_flow_vars.u_var())->get_phi();
71 
72  // The temperature shape functions at interior quadrature points.
73  const std::vector<std::vector<libMesh::Real> >& T_phi =
74  context.get_element_fe(_temp_vars.T_var())->get_phi();
75 
76  // Get residuals
77  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u_var()); // R_{u}
78  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v_var()); // R_{v}
79  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
80 
81  // Get Jacobians
82  libMesh::DenseSubMatrix<libMesh::Number> &KuT = context.get_elem_jacobian(_flow_vars.u_var(), _temp_vars.T_var()); // R_{u},{T}
83  libMesh::DenseSubMatrix<libMesh::Number> &KvT = context.get_elem_jacobian(_flow_vars.v_var(), _temp_vars.T_var()); // R_{v},{T}
84  libMesh::DenseSubMatrix<libMesh::Number>* KwT = NULL;
85 
86  if( _dim == 3 )
87  {
88  Fw = &context.get_elem_residual(_flow_vars.w_var()); // R_{w}
89  KwT = &context.get_elem_jacobian(_flow_vars.w_var(), _temp_vars.T_var()); // R_{w},{T}
90  }
91 
92  // Now we will build the element Jacobian and residual.
93  // Constructing the residual requires the solution and its
94  // gradient from the previous timestep. This must be
95  // calculated at each quadrature point by summing the
96  // solution degree-of-freedom values by the appropriate
97  // weight functions.
98  unsigned int n_qpoints = context.get_element_qrule().n_points();
99 
100  for (unsigned int qp=0; qp != n_qpoints; qp++)
101  {
102  // Compute the solution & its gradient at the old Newton iterate.
103  libMesh::Number T;
104  T = context.interior_value(_temp_vars.T_var(), qp);
105 
106  // First, an i-loop over the velocity degrees of freedom.
107  // We know that n_u_dofs == n_v_dofs so we can compute contributions
108  // for both at the same time.
109  for (unsigned int i=0; i != n_u_dofs; i++)
110  {
111  Fu(i) += -_rho_ref*_beta_T*(T - _T_ref)*_g(0)*vel_phi[i][qp]*JxW[qp];
112  Fv(i) += -_rho_ref*_beta_T*(T - _T_ref)*_g(1)*vel_phi[i][qp]*JxW[qp];
113 
114  if (_dim == 3)
115  (*Fw)(i) += -_rho_ref*_beta_T*(T - _T_ref)*_g(2)*vel_phi[i][qp]*JxW[qp];
116 
117  if (compute_jacobian)
118  {
119  for (unsigned int j=0; j != n_T_dofs; j++)
120  {
121  KuT(i,j) += context.get_elem_solution_derivative() *
122  -_rho_ref*_beta_T*_g(0)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
123  KvT(i,j) += context.get_elem_solution_derivative() *
124  -_rho_ref*_beta_T*_g(1)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
125 
126  if (_dim == 3)
127  (*KwT)(i,j) += context.get_elem_solution_derivative() *
128  -_rho_ref*_beta_T*_g(2)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
129 
130  } // End j dof loop
131  } // End compute_jacobian check
132 
133  } // End i dof loop
134  } // End quadrature loop
135 
136 #ifdef GRINS_USE_GRVY_TIMERS
137  this->_timer->EndTimer("BoussinesqBuoyancy::element_time_derivative");
138 #endif
139 
140  return;
141  }
142 
143 } // namespace GRINS
libMesh::Point _g
Gravitational vector.
unsigned int _dim
Physical dimension of problem.
GRINS namespace.
libMesh::Number _T_ref
reference temperature
PrimitiveFlowFEVariables _flow_vars
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Source term contribution for BoussinesqBuoyancy.
libMesh::Number _beta_T
coefficient of thermal expansion
libMesh::Number _rho_ref
reference density
PrimitiveTempFEVariables _temp_vars

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