GRINS-0.6.0
Public Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | List of all members
GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity > Class Template Reference

Adds Boussinesq bouyancy adjoint stabilization source term. More...

#include <boussinesq_buoyancy_adjoint_stab.h>

Inheritance diagram for GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >:
Collaboration graph
[legend]

Public Member Functions

 BoussinesqBuoyancyAdjointStabilization (const std::string &physics_name, const GetPot &input)
 
 ~BoussinesqBuoyancyAdjointStabilization ()
 
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)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for element interiors. More...
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialization of BoussinesqBuoyancy variables. More...
 
virtual void read_input_options (const GetPot &input)
 Read options from GetPot input file. By default, nothing is read. 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 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 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 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_bcs (libMesh::FEMSystem *system)
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
void attach_neumann_bound_func (GRINS::NBCContainer &neumann_bcs)
 
void attach_dirichlet_bound_func (const GRINS::DBCContainer &dirichlet_bc)
 
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_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)
 
BCHandlingBaseget_bc_handler ()
 
ICHandlingBaseget_ic_handler ()
 
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...
 

Protected Attributes

libMesh::Number _rho
 
Viscosity _mu
 Viscosity object. More...
 
IncompressibleNavierStokesStabilizationHelper _stab_helper
 
PrimitiveFlowFEVariables _flow_vars
 
PrimitiveTempFEVariables _temp_vars
 
libMesh::Number _rho_ref
 $ \rho_0 = $ reference 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...
 
unsigned int _dim
 Physical dimension of problem. More...
 
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::BCHandlingBase_bc_handler
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 
bool _is_axisymmetric
 

Static Protected Attributes

static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 

Private Member Functions

 BoussinesqBuoyancyAdjointStabilization ()
 

Detailed Description

template<class Viscosity>
class GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >

Adds Boussinesq bouyancy adjoint stabilization source term.

This class implements the adjiont stabilization term for the BoussinesqBuoyancy Physics. Intended to be used with IncompressibleNavierStokesAdjointStabilization and HeatTransferStabilization.

Definition at line 42 of file boussinesq_buoyancy_adjoint_stab.h.

Constructor & Destructor Documentation

template<class Mu >
GRINS::BoussinesqBuoyancyAdjointStabilization< Mu >::BoussinesqBuoyancyAdjointStabilization ( const std::string &  physics_name,
const GetPot &  input 
)

Definition at line 43 of file boussinesq_buoyancy_adjoint_stab.C.

References GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::_rho, GRINS::incompressible_navier_stokes, and GRINS::ParameterUser::set_parameter().

44  : BoussinesqBuoyancyBase(physics_name,input),
45  /* \todo Do we want to have these come from a BoussinesqBuoyancyAdjointStabilization section instead? */
46  _rho(1.0),
47  _mu(input),
48  _stab_helper( physics_name+"StabHelper", input )
49  {
50  this->set_parameter
51  (_rho, input, "Physics/"+incompressible_navier_stokes+"/rho", _rho);
52  }
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.
const PhysicsName incompressible_navier_stokes
IncompressibleNavierStokesStabilizationHelper _stab_helper

Definition at line 55 of file boussinesq_buoyancy_adjoint_stab.C.

56  {
57  return;
58  }
template<class Viscosity >
GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::BoussinesqBuoyancyAdjointStabilization ( )
private

Member Function Documentation

void GRINS::Physics::attach_dirichlet_bound_func ( const GRINS::DBCContainer dirichlet_bc)
inherited

Definition at line 150 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_dirichlet_bound_func().

151  {
152  _bc_handler->attach_dirichlet_bound_func( dirichlet_bc );
153  return;
154  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_dirichlet_bound_func(const GRINS::DBCContainer &dirichlet_bc)
void GRINS::Physics::attach_neumann_bound_func ( GRINS::NBCContainer neumann_bcs)
inherited

Definition at line 144 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_neumann_bound_func().

145  {
146  _bc_handler->attach_neumann_bound_func( neumann_bcs );
147  return;
148  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_neumann_bound_func(GRINS::NBCContainer &neumann_bcs)
void GRINS::Physics::auxiliary_init ( MultiphysicsSystem system)
virtualinherited

Any auxillary initialization a Physics class may need.

This is called after all variables are added, so this method can safely query the MultiphysicsSystem about variable information.

Definition at line 113 of file physics.C.

114  {
115  return;
116  }
void GRINS::Physics::compute_element_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 185 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::element_constraint().

187  {
188  return;
189  }
void GRINS::Physics::compute_element_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited
void GRINS::Physics::compute_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 203 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::mass_residual().

205  {
206  return;
207  }
void GRINS::Physics::compute_nonlocal_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 197 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

199  {
200  return;
201  }
void GRINS::Physics::compute_nonlocal_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 209 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

211  {
212  return;
213  }
void GRINS::Physics::compute_nonlocal_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 179 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

181  {
182  return;
183  }
void GRINS::Physics::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtualinherited
void GRINS::Physics::compute_side_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 191 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

193  {
194  return;
195  }
void GRINS::Physics::compute_side_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Reimplemented in GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 173 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_time_derivative().

175  {
176  return;
177  }
template<class Mu >
void GRINS::BoussinesqBuoyancyAdjointStabilization< Mu >::element_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 262 of file boussinesq_buoyancy_adjoint_stab.C.

265  {
266 #ifdef GRINS_USE_GRVY_TIMERS
267  this->_timer->BeginTimer("BoussinesqBuoyancyAdjointStabilization::element_constraint");
268 #endif
269 
270  // The number of local degrees of freedom in each variable.
271  const unsigned int n_p_dofs = context.get_dof_indices(_flow_vars.p_var()).size();
272  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u_var()).size();
273  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T_var()).size();
274 
275  // Element Jacobian * quadrature weights for interior integration.
276  const std::vector<libMesh::Real> &JxW =
277  context.get_element_fe(_flow_vars.u_var())->get_JxW();
278 
279  const std::vector<std::vector<libMesh::Real> >& T_phi =
280  context.get_element_fe(this->_temp_vars.T_var())->get_phi();
281 
282  const std::vector<std::vector<libMesh::Real> >& u_phi =
283  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
284 
285  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
286  context.get_element_fe(this->_flow_vars.p_var())->get_dphi();
287 
288  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_flow_vars.p_var()); // R_{p}
289 
290  libMesh::DenseSubMatrix<libMesh::Number> &KpT =
291  context.get_elem_jacobian(_flow_vars.p_var(), _temp_vars.T_var()); // J_{pT}
292  libMesh::DenseSubMatrix<libMesh::Number> &Kpu =
293  context.get_elem_jacobian(_flow_vars.p_var(), _flow_vars.u_var()); // J_{pu}
294  libMesh::DenseSubMatrix<libMesh::Number> &Kpv =
295  context.get_elem_jacobian(_flow_vars.p_var(), _flow_vars.v_var()); // J_{pv}
296  libMesh::DenseSubMatrix<libMesh::Number> *Kpw = NULL;
297 
298  if(this->_dim == 3)
299  {
300  Kpw = &context.get_elem_jacobian
301  (_flow_vars.p_var(), _flow_vars.w_var()); // J_{pw}
302  }
303 
304  // Now we will build the element Jacobian and residual.
305  // Constructing the residual requires the solution and its
306  // gradient from the previous timestep. This must be
307  // calculated at each quadrature point by summing the
308  // solution degree-of-freedom values by the appropriate
309  // weight functions.
310  unsigned int n_qpoints = context.get_element_qrule().n_points();
311 
312  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u_var());
313 
314  for (unsigned int qp=0; qp != n_qpoints; qp++)
315  {
316  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
317  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
318 
319  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u_var(), qp ),
320  context.interior_value( this->_flow_vars.v_var(), qp ) );
321  if( this->_dim == 3 )
322  {
323  U(2) = context.interior_value( this->_flow_vars.w_var(), qp );
324  }
325 
326  // Compute the viscosity at this qp
327  libMesh::Real mu_qp = this->_mu(context, qp);
328 
329  libMesh::Real tau_M;
330  libMesh::Real d_tau_M_d_rho;
331  libMesh::Gradient d_tau_M_dU;
332 
333  if (compute_jacobian)
335  ( context, qp, g, G, this->_rho, U, mu_qp,
336  tau_M, d_tau_M_d_rho, d_tau_M_dU,
337  this->_is_steady );
338  else
339  tau_M = this->_stab_helper.compute_tau_momentum
340  ( context, qp, g, G, this->_rho, U, mu_qp,
341  this->_is_steady );
342 
343  // Compute the solution & its gradient at the old Newton iterate.
344  libMesh::Number T;
345  T = context.interior_value(_temp_vars.T_var(), qp);
346 
347  libMesh::RealGradient d_residual_dT = _rho_ref*_beta_T*_g;
348  // d_residual_dU = 0
349  libMesh::RealGradient residual = (T-_T_ref)*d_residual_dT;
350 
351  // First, an i-loop over the velocity degrees of freedom.
352  // We know that n_u_dofs == n_v_dofs so we can compute contributions
353  // for both at the same time.
354  for (unsigned int i=0; i != n_p_dofs; i++)
355  {
356  Fp(i) += tau_M*residual*p_dphi[i][qp]*JxW[qp];
357 
358  if (compute_jacobian)
359  {
360  for (unsigned int j=0; j != n_T_dofs; ++j)
361  {
362  KpT(i,j) += tau_M*d_residual_dT*T_phi[j][qp]*p_dphi[i][qp]*JxW[qp] * context.get_elem_solution_derivative();
363  }
364 
365  for (unsigned int j=0; j != n_u_dofs; ++j)
366  {
367  Kpu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*residual*p_dphi[i][qp]*JxW[qp] * context.get_elem_solution_derivative();
368  Kpv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*residual*p_dphi[i][qp]*JxW[qp] * context.get_elem_solution_derivative();
369  }
370  if( this->_dim == 3 )
371  for (unsigned int j=0; j != n_u_dofs; ++j)
372  {
373  (*Kpw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*residual*p_dphi[i][qp]*JxW[qp] * context.get_elem_solution_derivative();
374  }
375  }
376  }
377  } // End quadrature loop
378 
379 #ifdef GRINS_USE_GRVY_TIMERS
380  this->_timer->EndTimer("BoussinesqBuoyancyAdjointStabilization::element_constraint");
381 #endif
382 
383  return;
384  }
libMesh::Point _g
Gravitational vector.
void compute_tau_momentum_and_derivs(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real T, libMesh::Real &tau_M, libMesh::Real &d_tau_M_d_rho, libMesh::Gradient &d_tau_M_d_U, bool is_steady) const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
unsigned int _dim
Physical dimension of problem.
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Number _T_ref
reference temperature
PrimitiveFlowFEVariables _flow_vars
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
libMesh::Number _beta_T
coefficient of thermal expansion
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
libMesh::Number _rho_ref
reference density
PrimitiveTempFEVariables _temp_vars
template<class Mu >
void GRINS::BoussinesqBuoyancyAdjointStabilization< Mu >::element_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Time dependent part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 72 of file boussinesq_buoyancy_adjoint_stab.C.

75  {
76 #ifdef GRINS_USE_GRVY_TIMERS
77  this->_timer->BeginTimer("BoussinesqBuoyancyAdjointStabilization::element_time_derivative");
78 #endif
79 
80  // The number of local degrees of freedom in each variable.
81  const unsigned int n_u_dofs = context.get_dof_indices(_flow_vars.u_var()).size();
82  const unsigned int n_T_dofs = context.get_dof_indices(_temp_vars.T_var()).size();
83 
84  // Element Jacobian * quadrature weights for interior integration.
85  const std::vector<libMesh::Real> &JxW =
86  context.get_element_fe(_flow_vars.u_var())->get_JxW();
87 
88  const std::vector<std::vector<libMesh::Real> >& T_phi =
89  context.get_element_fe(this->_temp_vars.T_var())->get_phi();
90 
91  const std::vector<std::vector<libMesh::Real> >& u_phi =
92  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
93 
94  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
95  context.get_element_fe(this->_flow_vars.u_var())->get_dphi();
96 
97  const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
98  context.get_element_fe(this->_flow_vars.u_var())->get_d2phi();
99 
100  // Get residuals and jacobians
101  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_flow_vars.u_var()); // R_{u}
102  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_flow_vars.v_var()); // R_{v}
103  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
104 
105  libMesh::DenseSubMatrix<libMesh::Number> &KuT =
106  context.get_elem_jacobian(_flow_vars.u_var(), _temp_vars.T_var()); // J_{uT}
107  libMesh::DenseSubMatrix<libMesh::Number> &KvT =
108  context.get_elem_jacobian(_flow_vars.v_var(), _temp_vars.T_var()); // J_{vT}
109  libMesh::DenseSubMatrix<libMesh::Number> &Kuu =
110  context.get_elem_jacobian(_flow_vars.u_var(), _flow_vars.u_var()); // J_{uu}
111  libMesh::DenseSubMatrix<libMesh::Number> &Kuv =
112  context.get_elem_jacobian(_flow_vars.u_var(), _flow_vars.v_var()); // J_{uv}
113  libMesh::DenseSubMatrix<libMesh::Number> &Kvu =
114  context.get_elem_jacobian(_flow_vars.v_var(), _flow_vars.u_var()); // J_{vu}
115  libMesh::DenseSubMatrix<libMesh::Number> &Kvv =
116  context.get_elem_jacobian(_flow_vars.v_var(), _flow_vars.v_var()); // J_{vv}
117 
118  libMesh::DenseSubMatrix<libMesh::Number> *KwT = NULL;
119  libMesh::DenseSubMatrix<libMesh::Number> *Kuw = NULL;
120  libMesh::DenseSubMatrix<libMesh::Number> *Kvw = NULL;
121  libMesh::DenseSubMatrix<libMesh::Number> *Kwu = NULL;
122  libMesh::DenseSubMatrix<libMesh::Number> *Kwv = NULL;
123  libMesh::DenseSubMatrix<libMesh::Number> *Kww = NULL;
124 
125  if(this->_dim == 3)
126  {
127  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
128  KwT = &context.get_elem_jacobian
129  (_flow_vars.w_var(), _temp_vars.T_var()); // J_{wT}
130  Kuw = &context.get_elem_jacobian
131  (_flow_vars.u_var(), _flow_vars.w_var()); // J_{uw}
132  Kvw = &context.get_elem_jacobian
133  (_flow_vars.v_var(), _flow_vars.w_var()); // J_{vw}
134  Kwu = &context.get_elem_jacobian
135  (_flow_vars.w_var(), _flow_vars.u_var()); // J_{wu}
136  Kwv = &context.get_elem_jacobian
137  (_flow_vars.w_var(), _flow_vars.v_var()); // J_{wv}
138  Kww = &context.get_elem_jacobian
139  (_flow_vars.w_var(), _flow_vars.w_var()); // J_{ww}
140  }
141 
142  // Now we will build the element Jacobian and residual.
143  // Constructing the residual requires the solution and its
144  // gradient from the previous timestep. This must be
145  // calculated at each quadrature point by summing the
146  // solution degree-of-freedom values by the appropriate
147  // weight functions.
148  unsigned int n_qpoints = context.get_element_qrule().n_points();
149 
150  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u_var());
151 
152  for (unsigned int qp=0; qp != n_qpoints; qp++)
153  {
154  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
155  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
156 
157  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u_var(), qp ),
158  context.interior_value( this->_flow_vars.v_var(), qp ) );
159  if( this->_dim == 3 )
160  {
161  U(2) = context.interior_value( this->_flow_vars.w_var(), qp );
162  }
163 
164  // Compute the viscosity at this qp
165  libMesh::Real mu_qp = this->_mu(context, qp);
166 
167  libMesh::Real tau_M;
168  libMesh::Real d_tau_M_d_rho;
169  libMesh::Gradient d_tau_M_dU;
170 
171  if (compute_jacobian)
173  ( context, qp, g, G, this->_rho, U, mu_qp,
174  tau_M, d_tau_M_d_rho, d_tau_M_dU,
175  this->_is_steady );
176  else
177  tau_M = this->_stab_helper.compute_tau_momentum
178  ( context, qp, g, G, this->_rho, U, mu_qp,
179  this->_is_steady );
180 
181  // Compute the solution & its gradient at the old Newton iterate.
182  libMesh::Number T;
183  T = context.interior_value(_temp_vars.T_var(), qp);
184 
185  libMesh::RealGradient d_residual_dT = _rho_ref*_beta_T*_g;
186  // d_residual_dU = 0
187  libMesh::RealGradient residual = (T-_T_ref)*d_residual_dT;
188 
189  for (unsigned int i=0; i != n_u_dofs; i++)
190  {
191  libMesh::Real test_func = this->_rho*U*u_gradphi[i][qp] +
192  mu_qp*( u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2) );
193  Fu(i) += -tau_M*residual(0)*test_func*JxW[qp];
194 
195  Fv(i) += -tau_M*residual(1)*test_func*JxW[qp];
196 
197  if (_dim == 3)
198  {
199  (*Fw)(i) += -tau_M*residual(2)*test_func*JxW[qp];
200  }
201 
202  if (compute_jacobian)
203  {
204  libMesh::Gradient d_test_func_dU = this->_rho*u_gradphi[i][qp];
205  // d_test_func_dT = 0
206 
207  for (unsigned int j=0; j != n_u_dofs; ++j)
208  {
209  Kuu(i,j) += -tau_M*residual(0)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
210  Kuu(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*residual(0)*test_func*JxW[qp] * context.get_elem_solution_derivative();
211  Kuv(i,j) += -tau_M*residual(0)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
212  Kuv(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*residual(0)*test_func*JxW[qp] * context.get_elem_solution_derivative();
213  Kvu(i,j) += -tau_M*residual(1)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
214  Kvu(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*residual(1)*test_func*JxW[qp] * context.get_elem_solution_derivative();
215  Kvv(i,j) += -tau_M*residual(1)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
216  Kvv(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*residual(1)*test_func*JxW[qp] * context.get_elem_solution_derivative();
217  }
218 
219  for (unsigned int j=0; j != n_T_dofs; ++j)
220  {
221  // KuT(i,j) += -tau_M*residual(0)*dtest_func_dT[j]*JxW[qp] * context.get_elem_solution_derivative();
222  KuT(i,j) += -tau_M*d_residual_dT(0)*T_phi[j][qp]*test_func*JxW[qp] * context.get_elem_solution_derivative();
223  // KvT(i,j) += -tau_M*residual(1)*dtest_func_dT[j]*JxW[qp] * context.get_elem_solution_derivative();
224  KvT(i,j) += -tau_M*d_residual_dT(1)*T_phi[j][qp]*test_func*JxW[qp] * context.get_elem_solution_derivative();
225  }
226  if (_dim == 3)
227  {
228  for (unsigned int j=0; j != n_T_dofs; ++j)
229  {
230  // KwT(i,j) += -tau_M*residual(2)*dtest_func_dT[j]*JxW[qp] * context.get_elem_solution_derivative();
231  (*KwT)(i,j) += -tau_M*d_residual_dT(2)*T_phi[j][qp]*test_func*JxW[qp] * context.get_elem_solution_derivative();
232  }
233 
234  for (unsigned int j=0; j != n_u_dofs; ++j)
235  {
236  (*Kuw)(i,j) += -tau_M*residual(0)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
237  (*Kuw)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*residual(0)*test_func*JxW[qp] * context.get_elem_solution_derivative();
238  (*Kvw)(i,j) += -tau_M*residual(1)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
239  (*Kvw)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*residual(1)*test_func*JxW[qp] * context.get_elem_solution_derivative();
240  (*Kwu)(i,j) += -tau_M*residual(2)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
241  (*Kwu)(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*residual(2)*test_func*JxW[qp] * context.get_elem_solution_derivative();
242  (*Kwv)(i,j) += -tau_M*residual(2)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
243  (*Kwv)(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*residual(2)*test_func*JxW[qp] * context.get_elem_solution_derivative();
244  (*Kww)(i,j) += -tau_M*residual(2)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp] * context.get_elem_solution_derivative();
245  (*Kww)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*residual(2)*test_func*JxW[qp] * context.get_elem_solution_derivative();
246  }
247  }
248 
249  } // End compute_jacobian check
250 
251  } // End i dof loop
252  } // End quadrature loop
253 
254 #ifdef GRINS_USE_GRVY_TIMERS
255  this->_timer->EndTimer("BoussinesqBuoyancyAdjointStabilization::element_time_derivative");
256 #endif
257 
258  return;
259  }
libMesh::Point _g
Gravitational vector.
void compute_tau_momentum_and_derivs(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real T, libMesh::Real &tau_M, libMesh::Real &d_tau_M_d_rho, libMesh::Gradient &d_tau_M_d_U, bool is_steady) const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:64
unsigned int _dim
Physical dimension of problem.
IncompressibleNavierStokesStabilizationHelper _stab_helper
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Number _T_ref
reference temperature
PrimitiveFlowFEVariables _flow_vars
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
libMesh::Number _beta_T
coefficient of thermal expansion
libMesh::Real compute_tau_momentum(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real mu, bool is_steady) const
libMesh::Number _rho_ref
reference density
PrimitiveTempFEVariables _temp_vars
bool GRINS::Physics::enabled_on_elem ( const libMesh::Elem *  elem)
virtualinherited

Find if current physics is active on supplied element.

Definition at line 83 of file physics.C.

References GRINS::Physics::_enabled_subdomains.

84  {
85  // Check if enabled_subdomains flag has been set and if we're
86  // looking at a real element (rather than a nonlocal evaluation)
87  if( !elem || _enabled_subdomains.empty() )
88  return true;
89 
90  // Check if current physics is enabled on elem
91  if( _enabled_subdomains.find( elem->subdomain_id() ) == _enabled_subdomains.end() )
92  return false;
93 
94  return true;
95  }
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
BCHandlingBase * GRINS::Physics::get_bc_handler ( )
inlineinherited

Definition at line 282 of file physics.h.

References GRINS::Physics::_bc_handler.

283  {
284  return _bc_handler;
285  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
ICHandlingBase * GRINS::Physics::get_ic_handler ( )
inlineinherited

Definition at line 288 of file physics.h.

References GRINS::Physics::_ic_handler.

289  {
290  return _ic_handler;
291  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
void GRINS::Physics::init_bcs ( libMesh::FEMSystem *  system)
inherited

Definition at line 118 of file physics.C.

References GRINS::Physics::_bc_handler, GRINS::BCHandlingBase::init_bc_data(), GRINS::BCHandlingBase::init_dirichlet_bc_func_objs(), GRINS::BCHandlingBase::init_dirichlet_bcs(), and GRINS::BCHandlingBase::init_periodic_bcs().

119  {
120  // Only need to init BC's if the physics actually created a handler
121  if( _bc_handler )
122  {
123  _bc_handler->init_bc_data( *system );
124  _bc_handler->init_dirichlet_bcs( system );
126  _bc_handler->init_periodic_bcs( system );
127  }
128 
129  return;
130  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
virtual void init_dirichlet_bcs(libMesh::FEMSystem *system) const
virtual void init_periodic_bcs(libMesh::FEMSystem *system) const
virtual void init_dirichlet_bc_func_objs(libMesh::FEMSystem *system) const
virtual void init_bc_data(const libMesh::FEMSystem &system)
Override this method to initialize any system-dependent data.
template<class Mu >
void GRINS::BoussinesqBuoyancyAdjointStabilization< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::Physics.

Definition at line 61 of file boussinesq_buoyancy_adjoint_stab.C.

62  {
63  context.get_element_fe(this->_flow_vars.p_var())->get_dphi();
64 
65  context.get_element_fe(this->_flow_vars.u_var())->get_dphi();
66  context.get_element_fe(this->_flow_vars.u_var())->get_d2phi();
67 
68  return;
69  }
PrimitiveFlowFEVariables _flow_vars
void GRINS::Physics::init_ics ( libMesh::FEMSystem *  system,
libMesh::CompositeFunction< libMesh::Number > &  all_ics 
)
inherited

Definition at line 133 of file physics.C.

References GRINS::Physics::_ic_handler, and GRINS::ICHandlingBase::init_ic_data().

135  {
136  if( _ic_handler )
137  {
138  _ic_handler->init_ic_data( *system, all_ics );
139  }
140 
141  return;
142  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
virtual void init_ic_data(const libMesh::FEMSystem &system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
Override this method to initialize any system-dependent data.
void GRINS::BoussinesqBuoyancyBase::init_variables ( libMesh::FEMSystem *  system)
virtualinherited

Initialization of BoussinesqBuoyancy variables.

Implements GRINS::Physics.

Definition at line 74 of file boussinesq_buoyancy_base.C.

References GRINS::BoussinesqBuoyancyBase::_dim, GRINS::BoussinesqBuoyancyBase::_flow_vars, GRINS::BoussinesqBuoyancyBase::_temp_vars, GRINS::PrimitiveFlowFEVariables::init(), and GRINS::PrimitiveTempFEVariables::init().

75  {
76  // Get libMesh to assign an index for each variable
77  this->_dim = system->get_mesh().mesh_dimension();
78 
79  _temp_vars.init(system);
80  _flow_vars.init(system);
81 
82  return;
83  }
virtual void init(libMesh::FEMSystem *system)
unsigned int _dim
Physical dimension of problem.
PrimitiveFlowFEVariables _flow_vars
PrimitiveTempFEVariables _temp_vars
virtual void init(libMesh::FEMSystem *system)
bool GRINS::Physics::is_steady ( ) const
inherited

Returns whether or not this physics is being solved with a steady solver.

Definition at line 103 of file physics.C.

References GRINS::Physics::_is_steady.

Referenced by GRINS::Physics::set_is_steady().

104  {
105  return _is_steady;
106  }
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
void GRINS::Physics::mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited
void GRINS::Physics::nonlocal_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for scalar variables.

Reimplemented in GRINS::ScalarODE.

Definition at line 250 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

253  {
254  return;
255  }
void GRINS::Physics::nonlocal_mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Mass matrix part(s) for scalar variables.

Reimplemented in GRINS::ScalarODE, and GRINS::AveragedTurbine< Viscosity >.

Definition at line 264 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

267  {
268  return;
269  }
void GRINS::Physics::nonlocal_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Time dependent part(s) of physics for scalar variables.

Reimplemented in GRINS::AveragedTurbine< Viscosity >, and GRINS::ScalarODE.

Definition at line 229 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

232  {
233  return;
234  }
void GRINS::Physics::read_input_options ( const GetPot &  input)
virtualinherited

Read options from GetPot input file. By default, nothing is read.

Reimplemented in GRINS::ScalarODE, GRINS::AveragedTurbineBase< Viscosity >, GRINS::AxisymmetricBoussinesqBuoyancy, GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::AxisymmetricHeatTransfer< Conductivity >, GRINS::AveragedFanBase< Viscosity >, GRINS::ParsedVelocitySourceBase< Viscosity >, GRINS::VelocityDragBase< Viscosity >, GRINS::VelocityPenaltyBase< Viscosity >, GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::HeatTransfer< Conductivity >, GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 70 of file physics.C.

References GRINS::Physics::_enabled_subdomains, and GRINS::Physics::_physics_name.

Referenced by GRINS::HeatTransferBase< Conductivity >::HeatTransferBase(), GRINS::HeatTransferStabilizationBase< Conductivity >::HeatTransferStabilizationBase(), GRINS::IncompressibleNavierStokesAdjointStabilization< Viscosity >::IncompressibleNavierStokesAdjointStabilization(), GRINS::IncompressibleNavierStokesSPGSMStabilization< Viscosity >::IncompressibleNavierStokesSPGSMStabilization(), GRINS::Physics::Physics(), and GRINS::SpalartAllmarasSPGSMStabilization< Viscosity >::SpalartAllmarasSPGSMStabilization().

71  {
72  int num_ids = input.vector_variable_size( "Physics/"+this->_physics_name+"/enabled_subdomains" );
73 
74  for( int i = 0; i < num_ids; i++ )
75  {
76  libMesh::subdomain_id_type dumvar = input( "Physics/"+this->_physics_name+"/enabled_subdomains", -1, i );
77  _enabled_subdomains.insert( dumvar );
78  }
79 
80  return;
81  }
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:254
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
template<class Mu >
void GRINS::BoussinesqBuoyancyAdjointStabilization< Mu >::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiPointer< libMesh::Number > &  param_pointer 
) const
virtual

Each subclass will register its copy of an independent.

Reimplemented from GRINS::ParameterUser.

Definition at line 388 of file boussinesq_buoyancy_adjoint_stab.C.

References GRINS::ParameterUser::register_parameter().

391  {
392  ParameterUser::register_parameter(param_name, param_pointer);
393  _mu.register_parameter(param_name, param_pointer);
394  }
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
void GRINS::Physics::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtualinherited

Register name of postprocessed quantity with PostProcessedQuantities.

Each Physics class will need to cache an unsigned int corresponding to each postprocessed quantity. This will be used in computing the values and putting them in the CachedVariables object.

Reimplemented in GRINS::ParsedVelocitySource< Viscosity >, GRINS::VelocityPenalty< Viscosity >, GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::HeatTransfer< Conductivity >, GRINS::ElasticCable< StressStrainLaw >, GRINS::ElasticMembrane< StressStrainLaw >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 161 of file physics.C.

163  {
164  return;
165  }
void GRINS::Physics::set_is_steady ( bool  is_steady)
inherited

Sets whether this physics is to be solved with a steady solver or not.

Since the member variable is static, only needs to be called on a single physics.

Definition at line 97 of file physics.C.

References GRINS::Physics::_is_steady, and GRINS::Physics::is_steady().

98  {
100  return;
101  }
bool is_steady() const
Returns whether or not this physics is being solved with a steady solver.
Definition: physics.C:103
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
void GRINS::ParameterUser::set_parameter ( libMesh::Number &  param_variable,
const GetPot &  input,
const std::string &  param_name,
libMesh::Number  param_default 
)
virtualinherited

Each subclass can simultaneously read a parameter value from.

Definition at line 35 of file parameter_user.C.

References GRINS::ParameterUser::_my_name, and GRINS::ParameterUser::_my_parameters.

Referenced by GRINS::AveragedFanAdjointStabilization< Viscosity >::AveragedFanAdjointStabilization(), GRINS::AveragedTurbineAdjointStabilization< Viscosity >::AveragedTurbineAdjointStabilization(), GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::BoussinesqBuoyancyAdjointStabilization(), GRINS::BoussinesqBuoyancyBase::BoussinesqBuoyancyBase(), GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::BoussinesqBuoyancySPGSMStabilization(), GRINS::ConstantConductivity::ConstantConductivity(), GRINS::ConstantPrandtlConductivity::ConstantPrandtlConductivity(), GRINS::ConstantSourceFunction::ConstantSourceFunction(), GRINS::ConstantSourceTerm::ConstantSourceTerm(), GRINS::ConstantSpecificHeat::ConstantSpecificHeat(), GRINS::ConstantViscosity::ConstantViscosity(), GRINS::ElasticCable< StressStrainLaw >::ElasticCable(), GRINS::ElasticCableConstantGravity::ElasticCableConstantGravity(), GRINS::ElasticMembrane< StressStrainLaw >::ElasticMembrane(), GRINS::ElasticMembraneConstantPressure::ElasticMembraneConstantPressure(), GRINS::HeatConduction< Conductivity >::HeatConduction(), GRINS::HeatTransferBase< Conductivity >::HeatTransferBase(), GRINS::IncompressibleNavierStokesBase< Viscosity >::IncompressibleNavierStokesBase(), GRINS::AverageNusseltNumber::init(), GRINS::MooneyRivlin::MooneyRivlin(), GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >::ReactingLowMachNavierStokesBase(), GRINS::HookesLaw1D::read_input_options(), GRINS::HookesLaw::read_input_options(), GRINS::AxisymmetricBoussinesqBuoyancy::read_input_options(), and GRINS::VelocityDragAdjointStabilization< Viscosity >::VelocityDragAdjointStabilization().

39  {
40  param_variable = input(param_name, param_default);
41 
42  libmesh_assert_msg(!_my_parameters.count(param_name),
43  "ERROR: " << _my_name << " double-registered parameter " <<
44  param_name);
45 
46  _my_parameters[param_name] = &param_variable;
47  }
std::map< std::string, libMesh::Number * > _my_parameters
void GRINS::Physics::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtualinherited

Set which variables are time evolving.

Set those variables which evolve in time (as opposed to variables that behave like constraints). This is done separately from init_variables() because the MultiphysicsSystem must initialize its base class before time_evolving variables can be set. Default implementation is no time evolving variables.

Reimplemented in GRINS::SpalartAllmaras< Viscosity >, GRINS::ScalarODE, GRINS::AxisymmetricHeatTransfer< Conductivity >, GRINS::IncompressibleNavierStokesBase< Viscosity >, GRINS::AveragedTurbineBase< Viscosity >, GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::HeatTransferBase< Conductivity >, GRINS::ParsedVelocitySourceBase< Viscosity >, GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >, GRINS::ElasticCableBase, GRINS::ElasticMembraneBase, and GRINS::HeatConduction< Conductivity >.

Definition at line 108 of file physics.C.

109  {
110  return;
111  }
void GRINS::Physics::side_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for boundaries of elements on the domain boundary.

Reimplemented in GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 243 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

246  {
247  return;
248  }
void GRINS::Physics::side_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Member Data Documentation

GRINS::BCHandlingBase* GRINS::Physics::_bc_handler
protectedinherited
libMesh::Number GRINS::BoussinesqBuoyancyBase::_beta_T
protectedinherited
unsigned int GRINS::BoussinesqBuoyancyBase::_dim
protectedinherited
std::set<libMesh::subdomain_id_type> GRINS::Physics::_enabled_subdomains
protectedinherited

Subdomains on which the current Physics class is enabled.

Definition at line 261 of file physics.h.

Referenced by GRINS::Physics::enabled_on_elem(), and GRINS::Physics::read_input_options().

PrimitiveFlowFEVariables GRINS::BoussinesqBuoyancyBase::_flow_vars
protectedinherited
libMesh::Point GRINS::BoussinesqBuoyancyBase::_g
protectedinherited
GRINS::ICHandlingBase* GRINS::Physics::_ic_handler
protectedinherited
bool GRINS::Physics::_is_axisymmetric
protectedinherited
bool GRINS::Physics::_is_steady = false
staticprotectedinherited

Caches whether or not the solver that's being used is steady or not.

This is need, for example, in flow stabilization as the tau terms change depending on whether the solver is steady or unsteady.

Definition at line 266 of file physics.h.

Referenced by GRINS::Physics::is_steady(), and GRINS::Physics::set_is_steady().

template<class Viscosity >
Viscosity GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::_mu
protected

Viscosity object.

Definition at line 72 of file boussinesq_buoyancy_adjoint_stab.h.

const PhysicsName GRINS::Physics::_physics_name
protectedinherited

Name of the physics object. Used for reading physics specific inputs.

We use a reference because the physics names are const global objects in GRINS namespace

No, we use a copy, because otherwise as soon as the memory in std::set<std::string> requested_physics gets overwritten we get in trouble.

Definition at line 254 of file physics.h.

Referenced by GRINS::SourceTermBase::parse_var_info(), and GRINS::Physics::read_input_options().

template<class Viscosity >
libMesh::Number GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::_rho
protected
libMesh::Number GRINS::BoussinesqBuoyancyBase::_rho_ref
protectedinherited
template<class Viscosity >
IncompressibleNavierStokesStabilizationHelper GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::_stab_helper
protected

Definition at line 74 of file boussinesq_buoyancy_adjoint_stab.h.

libMesh::Number GRINS::BoussinesqBuoyancyBase::_T_ref
protectedinherited
PrimitiveTempFEVariables GRINS::BoussinesqBuoyancyBase::_temp_vars
protectedinherited

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

Generated on Mon Jun 22 2015 21:32:22 for GRINS-0.6.0 by  doxygen 1.8.9.1