GRINS-0.8.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-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_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 
53  ( bool compute_jacobian,
54  AssemblyContext & context )
55  {
56  // The number of local degrees of freedom in each variable.
57  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
58  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
59 
60  // Element Jacobian * quadrature weights for interior integration.
61  const std::vector<libMesh::Real> &JxW =
62  context.get_element_fe(_flow_vars.u())->get_JxW();
63 
64  // The velocity shape functions at interior quadrature points.
65  const std::vector<std::vector<libMesh::Real> >& vel_phi =
66  context.get_element_fe(_flow_vars.u())->get_phi();
67 
68  // The temperature shape functions at interior quadrature points.
69  const std::vector<std::vector<libMesh::Real> >& T_phi =
70  context.get_element_fe(_temp_vars.T())->get_phi();
71 
72  // Get residuals
73  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u()); // R_{u}
74  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v()); // R_{v}
75  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
76 
77  // Get Jacobians
78  libMesh::DenseSubMatrix<libMesh::Number> &KuT = context.get_elem_jacobian(_flow_vars.u(), _temp_vars.T()); // R_{u},{T}
79  libMesh::DenseSubMatrix<libMesh::Number> &KvT = context.get_elem_jacobian(_flow_vars.v(), _temp_vars.T()); // R_{v},{T}
80  libMesh::DenseSubMatrix<libMesh::Number>* KwT = NULL;
81 
82 
83 
84  if( this->_flow_vars.dim() == 3 )
85  {
86  Fw = &context.get_elem_residual(_flow_vars.w()); // R_{w}
87  KwT = &context.get_elem_jacobian(_flow_vars.w(), _temp_vars.T()); // R_{w},{T}
88  }
89 
90  // Now we will build the element Jacobian and residual.
91  // Constructing the residual requires the solution and its
92  // gradient from the previous timestep. This must be
93  // calculated at each quadrature point by summing the
94  // solution degree-of-freedom values by the appropriate
95  // weight functions.
96  unsigned int n_qpoints = context.get_element_qrule().n_points();
97 
98  for (unsigned int qp=0; qp != n_qpoints; qp++)
99  {
100  // Compute the solution & its gradient at the old Newton iterate.
101  libMesh::Number T;
102  T = context.interior_value(_temp_vars.T(), qp);
103 
104  // First, an i-loop over the velocity degrees of freedom.
105  // We know that n_u_dofs == n_v_dofs so we can compute contributions
106  // for both at the same time.
107  for (unsigned int i=0; i != n_u_dofs; i++)
108  {
109  Fu(i) += -_rho*_beta_T*(T - _T_ref)*_g(0)*vel_phi[i][qp]*JxW[qp];
110  Fv(i) += -_rho*_beta_T*(T - _T_ref)*_g(1)*vel_phi[i][qp]*JxW[qp];
111 
112  if (this->_flow_vars.dim() == 3)
113  (*Fw)(i) += -_rho*_beta_T*(T - _T_ref)*_g(2)*vel_phi[i][qp]*JxW[qp];
114 
115  if (compute_jacobian)
116  {
117  for (unsigned int j=0; j != n_T_dofs; j++)
118  {
119  KuT(i,j) += context.get_elem_solution_derivative() *
120  -_rho*_beta_T*_g(0)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
121  KvT(i,j) += context.get_elem_solution_derivative() *
122  -_rho*_beta_T*_g(1)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
123 
124  if (this->_flow_vars.dim() == 3)
125  (*KwT)(i,j) += context.get_elem_solution_derivative() *
126  -_rho*_beta_T*_g(2)*vel_phi[i][qp]*T_phi[j][qp]*JxW[qp];
127 
128  } // End j dof loop
129  } // End compute_jacobian check
130 
131  } // End i dof loop
132  } // End quadrature loop
133  }
134 
135 } // namespace GRINS
GRINS namespace.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Source term contribution for BoussinesqBuoyancy.

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