GRINS-0.7.0
List of all members | Public Member Functions | Protected Attributes | Private Member Functions
GRINS::AxisymmetricBoussinesqBuoyancy Class Reference

Adds Axisymmetric Boussinesq bouyancy source term. More...

#include <axisym_boussinesq_buoyancy.h>

Inheritance diagram for GRINS::AxisymmetricBoussinesqBuoyancy:
Inheritance graph
[legend]
Collaboration diagram for GRINS::AxisymmetricBoussinesqBuoyancy:
Collaboration graph
[legend]

Public Member Functions

 AxisymmetricBoussinesqBuoyancy (const std::string &physics_name, const GetPot &input)
 
 ~AxisymmetricBoussinesqBuoyancy ()
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialization of AxisymmetricBoussinesqBuoyancy variables. More...
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Source term contribution for AxisymmetricBoussinesqBuoyancy. More...
 
- Public Member Functions inherited from GRINS::Physics
 Physics (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~Physics ()
 
virtual bool enabled_on_elem (const libMesh::Elem *elem)
 Find if current physics is active on supplied element. More...
 
void set_is_steady (bool is_steady)
 Sets whether this physics is to be solved with a steady solver or not. More...
 
bool is_steady () const
 Returns whether or not this physics is being solved with a steady solver. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Set which variables are time evolving. More...
 
virtual void auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. More...
 
virtual void side_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_element_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_damping_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
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. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
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. More...
 

Protected Attributes

unsigned int _dim
 Physical dimension of problem. More...
 
VelocityFEVariables _flow_vars
 
PressureFEVariable _press_var
 
PrimitiveTempFEVariables _temp_vars
 
libMesh::Number _rho
 $ \rho = $ density More...
 
libMesh::Number _T_ref
 $ T_0 = $ reference temperature More...
 
libMesh::Number _beta_T
 $ \beta_T = $ coefficient of thermal expansion More...
 
libMesh::Point _g
 Gravitational vector. More...
 
- Protected Attributes inherited from GRINS::Physics
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 

Private Member Functions

 AxisymmetricBoussinesqBuoyancy ()
 
void read_input_options (const GetPot &input)
 Read options from GetPot input file. More...
 
void register_variables ()
 

Additional Inherited Members

- Static Public Member Functions inherited from GRINS::Physics
static void set_is_axisymmetric (bool is_axisymmetric)
 Set whether we should treat the problem as axisymmetric. More...
 
static bool is_axisymmetric ()
 
- Static Public Attributes inherited from GRINS::ParameterUser
static std::string zero_vector_function = std::string("{0}")
 A parseable function string with LIBMESH_DIM components, all 0. More...
 
- Protected Member Functions inherited from GRINS::Physics
libMesh::UniquePtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe (const libMesh::Elem *elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
 
void parse_enabled_subdomains (const GetPot &input, const std::string &physics_name)
 
- Static Protected Attributes inherited from GRINS::Physics
static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 
static bool _is_axisymmetric = false
 Caches whether we are solving an axisymmetric problem or not. More...
 

Detailed Description

Adds Axisymmetric Boussinesq bouyancy source term.

This class implements the Axisymmetric Boussinesq approximation for thermal buoyancy. Namely: $ \mathbf{F} = -\rho_0 \beta_T \left( T - T_0 \right) \mathbf{g} $ where $ \rho = $ density, $ T_0 = $ reference temperature, $ \beta_T = $ coefficient of thermal expansion, and $ \mathbf{g} = $ the gravitional vector. This source term is added to the governing flow equations through the element_time_derivative routine. This class requires an axisymmetric flow physics enabled and the AxisymmetricHeatTransfer physics class enabled.

Definition at line 58 of file axisym_boussinesq_buoyancy.h.

Constructor & Destructor Documentation

GRINS::AxisymmetricBoussinesqBuoyancy::AxisymmetricBoussinesqBuoyancy ( const std::string &  physics_name,
const GetPot &  input 
)

Definition at line 45 of file axisym_boussinesq_buoyancy.C.

References read_input_options(), and register_variables().

47  : Physics(physics_name, input),
49  _press_var(input,PhysicsNaming::incompressible_navier_stokes(), true /*is_constraint_var*/),
51  {
52  this->read_input_options(input);
53  this->register_variables();
54  }
void read_input_options(const GetPot &input)
Read options from GetPot input file.
static PhysicsName incompressible_navier_stokes()
static PhysicsName axisymmetric_heat_transfer()
GRINS::AxisymmetricBoussinesqBuoyancy::~AxisymmetricBoussinesqBuoyancy ( )
inline

Definition at line 64 of file axisym_boussinesq_buoyancy.h.

64 {};
GRINS::AxisymmetricBoussinesqBuoyancy::AxisymmetricBoussinesqBuoyancy ( )
private

Member Function Documentation

void GRINS::AxisymmetricBoussinesqBuoyancy::element_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Source term contribution for AxisymmetricBoussinesqBuoyancy.

This is the main part of the class. This will add the source term to the AxisymmetricIncompNavierStokes class.

Reimplemented from GRINS::Physics.

Definition at line 101 of file axisym_boussinesq_buoyancy.C.

References _beta_T, _flow_vars, _g, _rho, _T_ref, _temp_vars, GRINS::PrimitiveTempFEVariables::T(), GRINS::VelocityFEVariables::u(), and GRINS::VelocityFEVariables::v().

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  }
libMesh::Point _g
Gravitational vector.
libMesh::Number _T_ref
reference temperature
libMesh::Number _beta_T
coefficient of thermal expansion
void GRINS::AxisymmetricBoussinesqBuoyancy::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::Physics.

Definition at line 93 of file axisym_boussinesq_buoyancy.C.

References _flow_vars, _temp_vars, GRINS::PrimitiveTempFEVariables::T(), and GRINS::VelocityFEVariables::u().

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  }
void GRINS::AxisymmetricBoussinesqBuoyancy::init_variables ( libMesh::FEMSystem *  system)
virtual

Initialization of AxisymmetricBoussinesqBuoyancy variables.

Implements GRINS::Physics.

Definition at line 84 of file axisym_boussinesq_buoyancy.C.

References _dim, _flow_vars, _press_var, _temp_vars, GRINS::VelocityFEVariables::init(), and GRINS::SingleFETypeVariable::init().

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  }
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
unsigned int _dim
Physical dimension of problem.
void GRINS::AxisymmetricBoussinesqBuoyancy::read_input_options ( const GetPot &  input)
private

Read options from GetPot input file.

Definition at line 56 of file axisym_boussinesq_buoyancy.C.

References _beta_T, _g, _rho, _T_ref, GRINS::PhysicsNaming::axisymmetric_boussinesq_buoyancy(), and GRINS::ParameterUser::set_parameter().

Referenced by AxisymmetricBoussinesqBuoyancy().

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  }
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 axisymmetric_boussinesq_buoyancy()
libMesh::Point _g
Gravitational vector.
libMesh::Number _T_ref
reference temperature
libMesh::Number _beta_T
coefficient of thermal expansion
void GRINS::AxisymmetricBoussinesqBuoyancy::register_variables ( )
private

Definition at line 74 of file axisym_boussinesq_buoyancy.C.

References _flow_vars, _press_var, _temp_vars, GRINS::GRINSPrivate::VariableWarehouse::check_and_register_variable(), GRINS::VariablesParsing::pressure_section(), GRINS::VariablesParsing::temperature_section(), and GRINS::VariablesParsing::velocity_section().

Referenced by AxisymmetricBoussinesqBuoyancy().

Member Data Documentation

libMesh::Number GRINS::AxisymmetricBoussinesqBuoyancy::_beta_T
protected

$ \beta_T = $ coefficient of thermal expansion

Definition at line 95 of file axisym_boussinesq_buoyancy.h.

Referenced by element_time_derivative(), and read_input_options().

unsigned int GRINS::AxisymmetricBoussinesqBuoyancy::_dim
protected

Physical dimension of problem.

Definition at line 82 of file axisym_boussinesq_buoyancy.h.

Referenced by init_variables().

VelocityFEVariables GRINS::AxisymmetricBoussinesqBuoyancy::_flow_vars
protected
libMesh::Point GRINS::AxisymmetricBoussinesqBuoyancy::_g
protected

Gravitational vector.

Definition at line 99 of file axisym_boussinesq_buoyancy.h.

Referenced by element_time_derivative(), and read_input_options().

PressureFEVariable GRINS::AxisymmetricBoussinesqBuoyancy::_press_var
protected

Definition at line 85 of file axisym_boussinesq_buoyancy.h.

Referenced by init_variables(), and register_variables().

libMesh::Number GRINS::AxisymmetricBoussinesqBuoyancy::_rho
protected

$ \rho = $ density

Definition at line 89 of file axisym_boussinesq_buoyancy.h.

Referenced by element_time_derivative(), and read_input_options().

libMesh::Number GRINS::AxisymmetricBoussinesqBuoyancy::_T_ref
protected

$ T_0 = $ reference temperature

Definition at line 92 of file axisym_boussinesq_buoyancy.h.

Referenced by element_time_derivative(), and read_input_options().

PrimitiveTempFEVariables GRINS::AxisymmetricBoussinesqBuoyancy::_temp_vars
protected

The documentation for this class was generated from the following files:

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