GRINS-0.7.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-2016 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(input, PhysicsNaming::incompressible_navier_stokes()),
49  _press_var(input,PhysicsNaming::incompressible_navier_stokes(), true /*is_constraint_var*/),
50  _temp_vars(input, PhysicsNaming::axisymmetric_heat_transfer())
51  {
52  this->read_input_options(input);
53  this->register_variables();
54  }
55 
57  {
58  this->set_parameter
59  (_rho, input,
60  "Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/rho_ref", 1.0);
61 
62  this->set_parameter
63  (_T_ref, input, "Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/T_ref", 1.0);
64 
65  this->set_parameter
66  (_beta_T, input, "Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/beta_T", 1.0);
67 
68  _g(0) = input("Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/g", 0.0, 0 );
69  _g(1) = input("Physics/"+PhysicsNaming::axisymmetric_boussinesq_buoyancy()+"/g", 0.0, 1 );
70 
71  return;
72  }
73 
75  {
77  this->_press_var);
79  this->_flow_vars);
81  this->_temp_vars);
82  }
83 
84  void AxisymmetricBoussinesqBuoyancy::init_variables( libMesh::FEMSystem* system )
85  {
86  this->_dim = system->get_mesh().mesh_dimension();
87 
88  this->_temp_vars.init(system);
89  this->_flow_vars.init(system);
90  this->_press_var.init(system);
91  }
92 
94  {
95  context.get_element_fe(_flow_vars.u())->get_JxW();
96  context.get_element_fe(_flow_vars.u())->get_phi();
97  context.get_element_fe(_flow_vars.u())->get_xyz();
98  context.get_element_fe(_temp_vars.T())->get_phi();
99  }
100 
102  AssemblyContext& context,
103  CachedValues& /*cache*/ )
104  {
105 #ifdef GRINS_USE_GRVY_TIMERS
106  this->_timer->BeginTimer("AxisymmetricBoussinesqBuoyancy::element_time_derivative");
107 #endif
108 
109  // The number of local degrees of freedom in each variable.
110  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u()).size();
111  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T()).size();
112 
113  // Element Jacobian * quadrature weights for interior integration.
114  const std::vector<libMesh::Real> &JxW =
115  context.get_element_fe(_flow_vars.u())->get_JxW();
116 
117  // The velocity shape functions at interior quadrature points.
118  const std::vector<std::vector<libMesh::Real> >& vel_phi =
119  context.get_element_fe(_flow_vars.u())->get_phi();
120 
121  // The temperature shape functions at interior quadrature points.
122  const std::vector<std::vector<libMesh::Real> >& T_phi =
123  context.get_element_fe(_temp_vars.T())->get_phi();
124 
125  // Physical location of the quadrature points
126  const std::vector<libMesh::Point>& u_qpoint =
127  context.get_element_fe(_flow_vars.u())->get_xyz();
128 
129  // Get residuals
130  libMesh::DenseSubVector<libMesh::Number> &Fr = context.get_elem_residual(_flow_vars.u()); // R_{r}
131  libMesh::DenseSubVector<libMesh::Number> &Fz = context.get_elem_residual(_flow_vars.v()); // R_{z}
132 
133  // Get Jacobians
134  libMesh::DenseSubMatrix<libMesh::Number> &KrT = context.get_elem_jacobian(_flow_vars.u(), _temp_vars.T()); // R_{r},{T}
135  libMesh::DenseSubMatrix<libMesh::Number> &KzT = context.get_elem_jacobian(_flow_vars.v(), _temp_vars.T()); // R_{z},{T}
136 
137  // Now we will build the element Jacobian and residual.
138  // Constructing the residual requires the solution and its
139  // gradient from the previous timestep. This must be
140  // calculated at each quadrature point by summing the
141  // solution degree-of-freedom values by the appropriate
142  // weight functions.
143  unsigned int n_qpoints = context.get_element_qrule().n_points();
144 
145  for (unsigned int qp=0; qp != n_qpoints; qp++)
146  {
147  const libMesh::Number r = u_qpoint[qp](0);
148 
149  // Compute the solution & its gradient at the old Newton iterate.
150  libMesh::Number T;
151  T = context.interior_value(_temp_vars.T(), qp);
152 
153  // First, an i-loop over the velocity degrees of freedom.
154  // We know that n_u_dofs == n_v_dofs so we can compute contributions
155  // for both at the same time.
156  for (unsigned int i=0; i != n_u_dofs; i++)
157  {
158  Fr(i) += -_rho*_beta_T*(T - _T_ref)*_g(0)*vel_phi[i][qp]*r*JxW[qp];
159  Fz(i) += -_rho*_beta_T*(T - _T_ref)*_g(1)*vel_phi[i][qp]*r*JxW[qp];
160 
161  if (compute_jacobian && context.get_elem_solution_derivative())
162  {
163  for (unsigned int j=0; j != n_T_dofs; j++)
164  {
165  const libMesh::Number val =
166  -_rho*_beta_T*vel_phi[i][qp]*T_phi[j][qp]*r*JxW[qp]
167  * context.get_elem_solution_derivative();
168  KrT(i,j) += val*_g(0);
169  KzT(i,j) += val*_g(1);
170  } // End j dof loop
171  } // End compute_jacobian check
172 
173  } // End i dof loop
174  } // End quadrature loop
175 
176 #ifdef GRINS_USE_GRVY_TIMERS
177  this->_timer->EndTimer("AxisymmetricBoussinesqBuoyancy::element_time_derivative");
178 #endif
179 
180  return;
181  }
182 
183 } // 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.
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:107
static std::string temperature_section()
static PhysicsName axisymmetric_boussinesq_buoyancy()
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
static void check_and_register_variable(const std::string &var_name, const FEVariablesBase &variable)
First check if var_name is registered and then register.
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.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Source term contribution for AxisymmetricBoussinesqBuoyancy.
static std::string velocity_section()
static std::string pressure_section()
virtual void init_variables(libMesh::FEMSystem *system)
Initialization of AxisymmetricBoussinesqBuoyancy variables.
libMesh::Point _g
Gravitational vector.
libMesh::Number _T_ref
reference temperature
libMesh::Number _beta_T
coefficient of thermal expansion
unsigned int _dim
Physical dimension of problem.

Generated on Thu Jun 2 2016 21:52:27 for GRINS-0.7.0 by  doxygen 1.8.10