GRINS-0.8.0
axisym_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
28 
29 // GRINS
30 #include "grins/assembly_context.h"
31 #include "grins/grins_enums.h"
34 
35 // libMesh
36 #include "libmesh/utility.h"
37 #include "libmesh/string_to_enum.h"
38 #include "libmesh/getpot.h"
39 #include "libmesh/fem_system.h"
40 #include "libmesh/quadrature.h"
41 
42 namespace GRINS
43 {
44 
46  const GetPot& input )
47  : Physics(physics_name, input),
48  _flow_vars(GRINSPrivate::VariableWarehouse::get_variable_subclass<VelocityVariable>(VariablesParsing::velocity_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
49  _press_var(GRINSPrivate::VariableWarehouse::get_variable_subclass<PressureFEVariable>(VariablesParsing::press_variable_name(input,physics_name,VariablesParsing::PHYSICS))),
50  _temp_vars(GRINSPrivate::VariableWarehouse::get_variable_subclass<PrimitiveTempFEVariables>(VariablesParsing::temp_variable_name(input,physics_name,VariablesParsing::PHYSICS)))
51  {
52  this->read_input_options(input);
53  }
54 
56  {
57  this->set_parameter
58  (_rho, input,
59  "Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/rho_ref", 1.0);
60 
61  this->set_parameter
62  (_T_ref, input, "Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/T_ref", 1.0);
63 
64  this->set_parameter
65  (_beta_T, input, "Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/beta_T", 1.0);
66 
67  _g(0) = input("Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/g", 0.0, 0 );
68  _g(1) = input("Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/g", 0.0, 1 );
69  }
70 
72  {
73  context.get_element_fe(_flow_vars.u())->get_JxW();
74  context.get_element_fe(_flow_vars.u())->get_phi();
75  context.get_element_fe(_flow_vars.u())->get_xyz();
76  context.get_element_fe(_temp_vars.T())->get_phi();
77  }
78 
80  ( bool compute_jacobian,
81  AssemblyContext & context )
82  {
83  // The number of local degrees of freedom in each variable.
84  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
85  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
86 
87  // Element Jacobian * quadrature weights for interior integration.
88  const std::vector<libMesh::Real> &JxW =
89  context.get_element_fe(_flow_vars.u())->get_JxW();
90 
91  // The velocity shape functions at interior quadrature points.
92  const std::vector<std::vector<libMesh::Real> >& vel_phi =
93  context.get_element_fe(_flow_vars.u())->get_phi();
94 
95  // The temperature shape functions at interior quadrature points.
96  const std::vector<std::vector<libMesh::Real> >& T_phi =
97  context.get_element_fe(_temp_vars.T())->get_phi();
98 
99  // Physical location of the quadrature points
100  const std::vector<libMesh::Point>& u_qpoint =
101  context.get_element_fe(_flow_vars.u())->get_xyz();
102 
103  // Get residuals
104  libMesh::DenseSubVector<libMesh::Number> &Fr = context.get_elem_residual(_flow_vars.u()); // R_{r}
105  libMesh::DenseSubVector<libMesh::Number> &Fz = context.get_elem_residual(_flow_vars.v()); // R_{z}
106 
107  // Get Jacobians
108  libMesh::DenseSubMatrix<libMesh::Number> &KrT = context.get_elem_jacobian(_flow_vars.u(), _temp_vars.T()); // R_{r},{T}
109  libMesh::DenseSubMatrix<libMesh::Number> &KzT = context.get_elem_jacobian(_flow_vars.v(), _temp_vars.T()); // R_{z},{T}
110 
111  // Now we will build the element Jacobian and residual.
112  // Constructing the residual requires the solution and its
113  // gradient from the previous timestep. This must be
114  // calculated at each quadrature point by summing the
115  // solution degree-of-freedom values by the appropriate
116  // weight functions.
117  unsigned int n_qpoints = context.get_element_qrule().n_points();
118 
119  for (unsigned int qp=0; qp != n_qpoints; qp++)
120  {
121  const libMesh::Number r = u_qpoint[qp](0);
122 
123  // Compute the solution & its gradient at the old Newton iterate.
124  libMesh::Number T;
125  T = context.interior_value(_temp_vars.T(), qp);
126 
127  // First, an i-loop over the velocity degrees of freedom.
128  // We know that n_u_dofs == n_v_dofs so we can compute contributions
129  // for both at the same time.
130  for (unsigned int i=0; i != n_u_dofs; i++)
131  {
132  Fr(i) += -_rho*_beta_T*(T - _T_ref)*_g(0)*vel_phi[i][qp]*r*JxW[qp];
133  Fz(i) += -_rho*_beta_T*(T - _T_ref)*_g(1)*vel_phi[i][qp]*r*JxW[qp];
134 
135  if (compute_jacobian && context.get_elem_solution_derivative())
136  {
137  for (unsigned int j=0; j != n_T_dofs; j++)
138  {
139  const libMesh::Number val =
140  -_rho*_beta_T*vel_phi[i][qp]*T_phi[j][qp]*r*JxW[qp]
141  * context.get_elem_solution_derivative();
142  KrT(i,j) += val*_g(0);
143  KzT(i,j) += val*_g(1);
144  } // End j dof loop
145  } // End compute_jacobian check
146 
147  } // End i dof loop
148  } // End quadrature loop
149  }
150 
151 } // namespace GRINS
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.
VariableIndex T() const
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
static PhysicsName axisymmetric_boussinesq_buoyancy()
void read_input_options(const GetPot &input)
Read options from GetPot input file.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
GRINS namespace.
const PrimitiveTempFEVariables & _temp_vars
libMesh::Point _g
Gravitational vector.
libMesh::Number _T_ref
reference temperature
libMesh::Number _beta_T
coefficient of thermal expansion
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Source term contribution for AxisymmetricBoussinesqBuoyancy.

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