GRINS-0.8.0
List of all members | Public Member Functions | Private Member Functions
GRINS::BoussinesqBuoyancy Class Reference

Adds Boussinesq bouyancy source term. More...

#include <boussinesq_buoyancy.h>

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

Public Member Functions

 BoussinesqBuoyancy (const std::string &physics_name, const GetPot &input)
 
 ~BoussinesqBuoyancy ()
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Source term contribution for BoussinesqBuoyancy. More...
 
- Public Member Functions inherited from GRINS::BoussinesqBuoyancyBase
 BoussinesqBuoyancyBase (const std::string &physics_name, const GetPot &input)
 
 ~BoussinesqBuoyancyBase ()
 
- Public Member Functions inherited from GRINS::Physics
 Physics (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~Physics ()
 
virtual void init_variables (libMesh::FEMSystem *)
 Initialize variables for this physics. More...
 
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 init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void preassembly (MultiphysicsSystem &)
 Perform any necessary setup before element assembly begins. More...
 
virtual void reinit (MultiphysicsSystem &)
 Any reinitialization that needs to be done. More...
 
virtual void side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void mass_residual (bool, AssemblyContext &)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 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 (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
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...
 

Private Member Functions

 BoussinesqBuoyancy ()
 

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::BoussinesqBuoyancyBase
void read_property (const GetPot &input, const std::string &old_option, const std::string &property, libMesh::Real &value)
 Helper function for parsing/maintaing backward compatibility. More...
 
void no_input_warning (const GetPot &input, const std::string &old_option, const std::string &material, const std::string &property)
 Helper function for parsing/maintaing backward compatibility. 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)
 
void check_var_subdomain_consistency (const FEVariablesBase &var) const
 Check that var is enabled on at least the subdomains this Physics is. More...
 
- Protected Attributes inherited from GRINS::BoussinesqBuoyancyBase
const VelocityVariable_flow_vars
 
const PressureFEVariable_press_var
 
const 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...
 
- 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 Boussinesq bouyancy source term.

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

Definition at line 48 of file boussinesq_buoyancy.h.

Constructor & Destructor Documentation

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

Definition at line 41 of file boussinesq_buoyancy.C.

42  : BoussinesqBuoyancyBase(physics_name,input)
43  {
44  return;
45  }
GRINS::BoussinesqBuoyancy::~BoussinesqBuoyancy ( )

Definition at line 47 of file boussinesq_buoyancy.C.

48  {
49  return;
50  }
GRINS::BoussinesqBuoyancy::BoussinesqBuoyancy ( )
private

Member Function Documentation

void GRINS::BoussinesqBuoyancy::element_time_derivative ( bool  compute_jacobian,
AssemblyContext context 
)
virtual

Source term contribution for BoussinesqBuoyancy.

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

Reimplemented from GRINS::Physics.

Definition at line 53 of file boussinesq_buoyancy.C.

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  }
libMesh::Point _g
Gravitational vector.
VariableIndex T() const
const VelocityVariable & _flow_vars
libMesh::Number _T_ref
reference temperature
const PrimitiveTempFEVariables & _temp_vars
libMesh::Number _beta_T
coefficient of thermal expansion
unsigned int dim() const
Number of components.

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

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