GRINS-0.8.0
List of all members | Public Member Functions | Protected Attributes | Private Member Functions
GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator > Class Template Reference

#include <reacting_low_mach_navier_stokes.h>

Inheritance diagram for GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >:
Collaboration graph
[legend]

Public Member Functions

 ReactingLowMachNavierStokes (const PhysicsName &physics_name, const GetPot &input, libMesh::UniquePtr< Mixture > &gas_mix)
 
virtual ~ReactingLowMachNavierStokes ()
 
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 postprocessing variables for ReactingLowMachNavierStokes. More...
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context)
 Constraint part(s) of physics for element interiors. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void compute_element_time_derivative_cache (AssemblyContext &context)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
- Public Member Functions inherited from GRINS::ReactingLowMachNavierStokesBase< Mixture >
 ReactingLowMachNavierStokesBase (const PhysicsName &physics_name, const GetPot &input, libMesh::UniquePtr< Mixture > &gas_mix)
 
virtual ~ReactingLowMachNavierStokesBase ()
 
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...
 
const Mixture & gas_mixture () const
 
- Public Member Functions inherited from GRINS::ReactingLowMachNavierStokesAbstract
 ReactingLowMachNavierStokesAbstract (const PhysicsName &physics_name, const GetPot &input)
 
virtual ~ReactingLowMachNavierStokesAbstract ()
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
unsigned int n_species () const
 
libMesh::Real T (const libMesh::Point &p, const AssemblyContext &c) const
 
void mass_fractions (const libMesh::Point &p, const AssemblyContext &c, std::vector< libMesh::Real > &mass_fracs) const
 
libMesh::Real rho (libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
 
libMesh::Real get_p0_steady (const AssemblyContext &c, unsigned int qp) const
 
libMesh::Real get_p0_steady_side (const AssemblyContext &c, unsigned int qp) const
 
libMesh::Real get_p0_steady (const AssemblyContext &c, const libMesh::Point &p) const
 
libMesh::Real get_p0_transient (const AssemblyContext &c, unsigned int qp) const
 
- 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 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 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 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_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 &)
 
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...
 

Protected Attributes

bool _pin_pressure
 Enable pressure pinning. More...
 
PressurePinning _p_pinning
 
unsigned int _rho_index
 Index from registering this quantity. More...
 
std::vector< unsigned int > _species_viscosity
 Index from registering this quantity. Each species will have it's own index. More...
 
unsigned int _mu_index
 Index from registering this quantity. More...
 
unsigned int _k_index
 Index from registering this quantity. More...
 
unsigned int _cp_index
 Index from registering this quantity. More...
 
std::vector< unsigned int > _mole_fractions_index
 Index from registering this quantity. Each species will have it's own index. More...
 
std::vector< unsigned int > _h_s_index
 Index from registering this quantity. Each species will have it's own index. More...
 
std::vector< unsigned int > _omega_dot_index
 Index from registering this quantity. Each species will have it's own index. More...
 
std::vector< unsigned int > _Ds_index
 Index from registering this quantity. Each species will have it's own index. More...
 
- Protected Attributes inherited from GRINS::ReactingLowMachNavierStokesBase< Mixture >
libMesh::UniquePtr< Mixture > _gas_mixture
 
- Protected Attributes inherited from GRINS::ReactingLowMachNavierStokesAbstract
libMesh::Number _p0
 
VelocityVariable_flow_vars
 
PressureFEVariable_press_var
 
PrimitiveTempFEVariables_temp_vars
 
SpeciesMassFractionsVariable_species_vars
 
ThermoPressureVariable_p0_var
 
unsigned int _n_species
 Number of species. More...
 
libMesh::Point _g
 Gravity vector. More...
 
bool _enable_thermo_press_calc
 Flag to enable thermodynamic pressure calculation. More...
 
bool _fixed_density
 
libMesh::Real _fixed_rho_value
 
- 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

 ReactingLowMachNavierStokes ()
 

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)
 
void check_var_subdomain_consistency (const FEVariablesBase &var) const
 Check that var is enabled on at least the subdomains this Physics is. 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

template<typename Mixture, typename Evaluator>
class GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >

Definition at line 35 of file reacting_low_mach_navier_stokes.h.

Constructor & Destructor Documentation

template<typename Mixture , typename Evaluator >
GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::ReactingLowMachNavierStokes ( const PhysicsName physics_name,
const GetPot &  input,
libMesh::UniquePtr< Mixture > &  gas_mix 
)

Definition at line 42 of file reacting_low_mach_navier_stokes.C.

References GRINS::Physics::_ic_handler, GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_pin_pressure, and GRINS::PhysicsNaming::reacting_low_mach_navier_stokes().

45  : ReactingLowMachNavierStokesBase<Mixture>(physics_name,input,gas_mix),
46  _p_pinning(input,physics_name),
47  _rho_index(0),
48  _mu_index(0),
49  _k_index(0),
50  _cp_index(0)
51  {
52  this->_pin_pressure = input("Physics/"+PhysicsNaming::reacting_low_mach_navier_stokes()+"/pin_pressure", false );
53 
54  this->_ic_handler = new GenericICHandler( physics_name, input );
55  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
static PhysicsName reacting_low_mach_navier_stokes()
unsigned int _k_index
Index from registering this quantity.
unsigned int _mu_index
Index from registering this quantity.
unsigned int _rho_index
Index from registering this quantity.
unsigned int _cp_index
Index from registering this quantity.
template<typename Mixture , typename Evaluator >
virtual GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::~ReactingLowMachNavierStokes ( )
inlinevirtual

Definition at line 42 of file reacting_low_mach_navier_stokes.h.

42 {};
template<typename Mixture , typename Evaluator >
GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::ReactingLowMachNavierStokes ( )
private

Member Function Documentation

template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::auxiliary_init ( MultiphysicsSystem system)
virtual

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.

Reimplemented from GRINS::Physics.

Definition at line 58 of file reacting_low_mach_navier_stokes.C.

59  {
60  if( _pin_pressure )
61  _p_pinning.check_pin_location(system.get_mesh());
62  }
void check_pin_location(const libMesh::MeshBase &mesh)
Check the mesh to ensure pin location is found.
template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::compute_element_time_derivative_cache ( AssemblyContext context)
virtual
Todo:
Need to figure out something smarter for controling species that go slightly negative.

Reimplemented from GRINS::Physics.

Definition at line 546 of file reacting_low_mach_navier_stokes.C.

References GRINS::Cache::DIFFUSION_COEFFS, GRINS::AssemblyContext::get_cached_values(), GRINS::Cache::MASS_FRACTIONS, GRINS::Cache::MASS_FRACTIONS_GRAD, GRINS::Cache::MIXTURE_DENSITY, GRINS::Cache::MIXTURE_GAS_CONSTANT, GRINS::Cache::MIXTURE_SPECIFIC_HEAT_P, GRINS::Cache::MIXTURE_THERMAL_CONDUCTIVITY, GRINS::Cache::MIXTURE_VISCOSITY, GRINS::Cache::MOLAR_MASS, GRINS::Cache::OMEGA_DOT, GRINS::Cache::PRESSURE, GRINS::CachedValues::set_gradient_values(), GRINS::CachedValues::set_values(), GRINS::CachedValues::set_vector_gradient_values(), GRINS::CachedValues::set_vector_values(), GRINS::Cache::SPECIES_ENTHALPY, GRINS::Cache::TEMPERATURE, GRINS::Cache::TEMPERATURE_GRAD, GRINS::Cache::THERMO_PRESSURE, GRINS::Cache::X_VELOCITY, GRINS::Cache::X_VELOCITY_GRAD, GRINS::Cache::Y_VELOCITY, GRINS::Cache::Y_VELOCITY_GRAD, GRINS::Cache::Z_VELOCITY, and GRINS::Cache::Z_VELOCITY_GRAD.

547  {
548  CachedValues & cache = context.get_cached_values();
549 
550  Evaluator gas_evaluator( *(this->_gas_mixture) );
551 
552  const unsigned int n_qpoints = context.get_element_qrule().n_points();
553 
554  std::vector<libMesh::Real> u, v, w, T, p, p0;
555  u.resize(n_qpoints);
556  v.resize(n_qpoints);
557  if( this->_flow_vars.dim() > 2 )
558  w.resize(n_qpoints);
559 
560  T.resize(n_qpoints);
561  p.resize(n_qpoints);
562  p0.resize(n_qpoints);
563 
564  std::vector<libMesh::Gradient> grad_u, grad_v, grad_w, grad_T;
565  grad_u.resize(n_qpoints);
566  grad_v.resize(n_qpoints);
567  if( this->_flow_vars.dim() > 2 )
568  grad_w.resize(n_qpoints);
569 
570  grad_T.resize(n_qpoints);
571 
572  std::vector<std::vector<libMesh::Real> > mass_fractions;
573  std::vector<std::vector<libMesh::Gradient> > grad_mass_fractions;
574  mass_fractions.resize(n_qpoints);
575  grad_mass_fractions.resize(n_qpoints);
576 
577  std::vector<libMesh::Real> M;
578  M.resize(n_qpoints);
579 
580  std::vector<libMesh::Real> R;
581  R.resize(n_qpoints);
582 
583  std::vector<libMesh::Real> rho;
584  rho.resize(n_qpoints);
585 
586  std::vector<libMesh::Real> cp;
587  cp.resize(n_qpoints);
588 
589  std::vector<libMesh::Real> mu;
590  mu.resize(n_qpoints);
591 
592  std::vector<libMesh::Real> k;
593  k.resize(n_qpoints);
594 
595  std::vector<std::vector<libMesh::Real> > h_s;
596  h_s.resize(n_qpoints);
597 
598  std::vector<std::vector<libMesh::Real> > D_s;
599  D_s.resize(n_qpoints);
600 
601  std::vector<std::vector<libMesh::Real> > omega_dot_s;
602  omega_dot_s.resize(n_qpoints);
603 
604  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
605  {
606  u[qp] = context.interior_value(this->_flow_vars.u(), qp);
607  v[qp] = context.interior_value(this->_flow_vars.v(), qp);
608 
609  grad_u[qp] = context.interior_gradient(this->_flow_vars.u(), qp);
610  grad_v[qp] = context.interior_gradient(this->_flow_vars.v(), qp);
611  if( this->_flow_vars.dim() > 2 )
612  {
613  w[qp] = context.interior_value(this->_flow_vars.w(), qp);
614  grad_w[qp] = context.interior_gradient(this->_flow_vars.w(), qp);
615  }
616 
617  T[qp] = context.interior_value(this->_temp_vars.T(), qp);
618  grad_T[qp] = context.interior_gradient(this->_temp_vars.T(), qp);
619 
620  p[qp] = context.interior_value(this->_press_var.p(), qp);
621  p0[qp] = this->get_p0_steady(context, qp);
622 
623  mass_fractions[qp].resize(this->_n_species);
624  grad_mass_fractions[qp].resize(this->_n_species);
625  h_s[qp].resize(this->_n_species);
626 
627  for( unsigned int s = 0; s < this->_n_species; s++ )
628  {
631  mass_fractions[qp][s] = std::max( context.interior_value(this->_species_vars.species(s),qp), 0.0 );
632  grad_mass_fractions[qp][s] = context.interior_gradient(this->_species_vars.species(s),qp);
633  h_s[qp][s] = gas_evaluator.h_s( T[qp], s );
634  }
635 
636  M[qp] = gas_evaluator.M_mix( mass_fractions[qp] );
637 
638  R[qp] = gas_evaluator.R_mix( mass_fractions[qp] );
639 
640  rho[qp] = this->rho( T[qp], p0[qp], R[qp] );
641 
642  cp[qp] = gas_evaluator.cp(T[qp], p0[qp], mass_fractions[qp]);
643 
644  D_s[qp].resize(this->_n_species);
645 
646  gas_evaluator.mu_and_k_and_D( T[qp], rho[qp], cp[qp], mass_fractions[qp],
647  mu[qp], k[qp], D_s[qp] );
648 
649  omega_dot_s[qp].resize(this->_n_species);
650  gas_evaluator.omega_dot( T[qp], rho[qp], mass_fractions[qp], omega_dot_s[qp] );
651  }
652 
653  cache.set_values(Cache::X_VELOCITY, u);
654  cache.set_values(Cache::Y_VELOCITY, v);
655  cache.set_gradient_values(Cache::X_VELOCITY_GRAD, grad_u);
656  cache.set_gradient_values(Cache::Y_VELOCITY_GRAD, grad_v);
657 
658  if(this->_flow_vars.dim() > 2)
659  {
660  cache.set_values(Cache::Z_VELOCITY, w);
661  cache.set_gradient_values(Cache::Z_VELOCITY_GRAD, grad_w);
662  }
663 
664  cache.set_values(Cache::TEMPERATURE, T);
665  cache.set_gradient_values(Cache::TEMPERATURE_GRAD, grad_T);
666  cache.set_values(Cache::PRESSURE, p);
667  cache.set_values(Cache::THERMO_PRESSURE, p0);
668  cache.set_vector_values(Cache::MASS_FRACTIONS, mass_fractions);
669  cache.set_vector_gradient_values(Cache::MASS_FRACTIONS_GRAD, grad_mass_fractions);
670  cache.set_values(Cache::MOLAR_MASS, M);
671  cache.set_values(Cache::MIXTURE_GAS_CONSTANT, R);
672  cache.set_values(Cache::MIXTURE_DENSITY, rho);
673  cache.set_values(Cache::MIXTURE_SPECIFIC_HEAT_P, cp);
674  cache.set_values(Cache::MIXTURE_VISCOSITY, mu);
675  cache.set_values(Cache::MIXTURE_THERMAL_CONDUCTIVITY, k);
676  cache.set_vector_values(Cache::DIFFUSION_COEFFS, D_s);
677  cache.set_vector_values(Cache::SPECIES_ENTHALPY, h_s);
678  cache.set_vector_values(Cache::OMEGA_DOT, omega_dot_s);
679  }
VariableIndex species(unsigned int species) const
VariableIndex T() const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
void mass_fractions(const libMesh::Point &p, const AssemblyContext &c, std::vector< libMesh::Real > &mass_fracs) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
VariableIndex p() const
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int dim() const
Number of components.
template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtual

Reimplemented from GRINS::Physics.

Definition at line 682 of file reacting_low_mach_navier_stokes.C.

686  {
687  Evaluator gas_evaluator( *(this->_gas_mixture) );
688 
689  if( quantity_index == this->_rho_index )
690  {
691  std::vector<libMesh::Real> Y( this->_n_species );
692  libMesh::Real T = this->T(point,context);
693  libMesh::Real p0 = this->get_p0_steady(context,point);
694  this->mass_fractions( point, context, Y );
695 
696  value = this->rho( T, p0, gas_evaluator.R_mix(Y) );
697  }
698  else if( quantity_index == this->_mu_index )
699  {
700  std::vector<libMesh::Real> Y( this->_n_species );
701  libMesh::Real T = this->T(point,context);
702  this->mass_fractions( point, context, Y );
703  libMesh::Real p0 = this->get_p0_steady(context,point);
704 
705  value = gas_evaluator.mu( T, p0, Y );
706  }
707  else if( quantity_index == this->_k_index )
708  {
709  std::vector<libMesh::Real> Y( this->_n_species );
710 
711  libMesh::Real T = this->T(point,context);
712  this->mass_fractions( point, context, Y );
713  libMesh::Real p0 = this->get_p0_steady(context,point);
714 
715  libMesh::Real cp = gas_evaluator.cp( T, p0, Y );
716 
717  libMesh::Real rho = this->rho( T, p0, gas_evaluator.R_mix(Y) );
718 
719  libMesh::Real mu, k;
720  std::vector<libMesh::Real> D( this->_n_species );
721 
722  gas_evaluator.mu_and_k_and_D( T, rho, cp, Y, mu, k, D );
723 
724  value = k;
725  return;
726  }
727  else if( quantity_index == this->_cp_index )
728  {
729  std::vector<libMesh::Real> Y( this->_n_species );
730  libMesh::Real T = this->T(point,context);
731  this->mass_fractions( point, context, Y );
732  libMesh::Real p0 = this->get_p0_steady(context,point);
733 
734  value = gas_evaluator.cp( T, p0, Y );
735  }
736  // Now check for species dependent stuff
737  else
738  {
739  if( !this->_h_s_index.empty() )
740  {
741  libmesh_assert_equal_to( _h_s_index.size(), this->n_species() );
742 
743  for( unsigned int s = 0; s < this->n_species(); s++ )
744  {
745  if( quantity_index == this->_h_s_index[s] )
746  {
747  libMesh::Real T = this->T(point,context);
748 
749  value = gas_evaluator.h_s( T, s );
750  return;
751  }
752  }
753  }
754 
755  if( !this->_mole_fractions_index.empty() )
756  {
757  libmesh_assert_equal_to( _mole_fractions_index.size(), this->n_species() );
758 
759  for( unsigned int s = 0; s < this->n_species(); s++ )
760  {
761  if( quantity_index == this->_mole_fractions_index[s] )
762  {
763  std::vector<libMesh::Real> Y( this->_n_species );
764  this->mass_fractions( point, context, Y );
765 
766  libMesh::Real M = gas_evaluator.M_mix(Y);
767 
768  value = gas_evaluator.X( s, M, Y[s] );
769  return;
770  }
771  }
772  }
773 
774  if( !this->_omega_dot_index.empty() )
775  {
776  libmesh_assert_equal_to( _omega_dot_index.size(), this->n_species() );
777 
778  for( unsigned int s = 0; s < this->n_species(); s++ )
779  {
780  if( quantity_index == this->_omega_dot_index[s] )
781  {
782  std::vector<libMesh::Real> Y( this->n_species() );
783  this->mass_fractions( point, context, Y );
784 
785  libMesh::Real T = this->T(point,context);
786 
787  libMesh::Real p0 = this->get_p0_steady(context,point);
788 
789  libMesh::Real rho = this->rho( T, p0, gas_evaluator.R_mix(Y) );
790 
791  std::vector<libMesh::Real> omega_dot( this->n_species() );
792  gas_evaluator.omega_dot( T, rho, Y, omega_dot );
793 
794  value = omega_dot[s];
795  return;
796  }
797  }
798  }
799 
800  if( !this->_Ds_index.empty() )
801  {
802  libmesh_assert_equal_to( _Ds_index.size(), this->n_species() );
803 
804  for( unsigned int s = 0; s < this->n_species(); s++ )
805  {
806  if( quantity_index == this->_Ds_index[s] )
807  {
808  std::vector<libMesh::Real> Y( this->_n_species );
809 
810  libMesh::Real T = this->T(point,context);
811  this->mass_fractions( point, context, Y );
812  libMesh::Real p0 = this->get_p0_steady(context,point);
813 
814  libMesh::Real cp = gas_evaluator.cp( T, p0, Y );
815 
816  libMesh::Real rho = this->rho( T, p0, gas_evaluator.R_mix(Y) );
817 
818  libMesh::Real mu, k;
819  std::vector<libMesh::Real> D( this->_n_species );
820 
821  gas_evaluator.mu_and_k_and_D( T, rho, cp, Y, mu, k, D );
822 
823  value = D[s];
824  return;
825  }
826  }
827  }
828 
829  } // if/else quantity_index
830 
831  return;
832  }
std::vector< unsigned int > _mole_fractions_index
Index from registering this quantity. Each species will have it's own index.
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
void mass_fractions(const libMesh::Point &p, const AssemblyContext &c, std::vector< libMesh::Real > &mass_fracs) const
unsigned int _k_index
Index from registering this quantity.
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
std::vector< unsigned int > _Ds_index
Index from registering this quantity. Each species will have it's own index.
unsigned int _mu_index
Index from registering this quantity.
std::vector< unsigned int > _omega_dot_index
Index from registering this quantity. Each species will have it's own index.
unsigned int _rho_index
Index from registering this quantity.
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
std::vector< unsigned int > _h_s_index
Index from registering this quantity. Each species will have it's own index.
unsigned int _cp_index
Index from registering this quantity.
template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::element_constraint ( bool  ,
AssemblyContext  
)
virtual

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 404 of file reacting_low_mach_navier_stokes.C.

406  {
407  // Pin p = p_value at p_point
408  if( this->_pin_pressure )
409  {
410  _p_pinning.pin_value( context, compute_jacobian, this->_press_var.p() );
411  }
412 
413  return;
414  }
VariableIndex p() const
void pin_value(libMesh::DiffContext &context, const bool request_jacobian, const GRINS::VariableIndex var, const double penalty=1.0)
The idea here is to pin a variable to a particular value if there is a null space - e...
template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

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

Todo:
Need to add SCEBD term.
Todo:
Would it be better to put this in its own DoF loop and do the if check once?

Reimplemented from GRINS::Physics.

Definition at line 164 of file reacting_low_mach_navier_stokes.C.

References GRINS::Cache::DIFFUSION_COEFFS, GRINS::CachedValues::get_cached_gradient_values(), GRINS::AssemblyContext::get_cached_values(), GRINS::CachedValues::get_cached_values(), GRINS::CachedValues::get_cached_vector_gradient_values(), GRINS::CachedValues::get_cached_vector_values(), GRINS::Physics::is_axisymmetric(), GRINS::Cache::MASS_FRACTIONS_GRAD, GRINS::Cache::MIXTURE_DENSITY, GRINS::Cache::MIXTURE_SPECIFIC_HEAT_P, GRINS::Cache::MIXTURE_THERMAL_CONDUCTIVITY, GRINS::Cache::MIXTURE_VISCOSITY, GRINS::Cache::MOLAR_MASS, GRINS::Cache::OMEGA_DOT, GRINS::Cache::PRESSURE, GRINS::Cache::SPECIES_ENTHALPY, GRINS::Cache::TEMPERATURE, GRINS::Cache::TEMPERATURE_GRAD, GRINS::Cache::X_VELOCITY, GRINS::Cache::X_VELOCITY_GRAD, GRINS::Cache::Y_VELOCITY, GRINS::Cache::Y_VELOCITY_GRAD, GRINS::Cache::Z_VELOCITY, and GRINS::Cache::Z_VELOCITY_GRAD.

165  {
166  if( compute_jacobian )
167  libmesh_not_implemented();
168 
169  const CachedValues & cache = context.get_cached_values();
170 
171  // Convenience
172  const VariableIndex s0_var = this->_species_vars.species(0);
173 
174  // The number of local degrees of freedom in each variable.
175  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
176  const unsigned int n_s_dofs = context.get_dof_indices(s0_var).size();
177  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
178  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
179 
180  // Check number of dofs is same for _flow_vars.u(), v_var and w_var.
181  if (this->_flow_vars.dim() > 1)
182  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.v()).size());
183 
184  if (this->_flow_vars.dim() == 3)
185  libmesh_assert (n_u_dofs == context.get_dof_indices(this->_flow_vars.w()).size());
186 
187  // Element Jacobian * quadrature weights for interior integration.
188  const std::vector<libMesh::Real>& JxW =
189  context.get_element_fe(this->_flow_vars.u())->get_JxW();
190 
191  // The pressure shape functions at interior quadrature points.
192  const std::vector<std::vector<libMesh::Real> >& p_phi =
193  context.get_element_fe(this->_press_var.p())->get_phi();
194 
195  // The species shape functions at interior quadrature points.
196  const std::vector<std::vector<libMesh::Real> >& s_phi = context.get_element_fe(s0_var)->get_phi();
197 
198  // The species shape function gradients at interior quadrature points.
199  const std::vector<std::vector<libMesh::Gradient> >& s_grad_phi = context.get_element_fe(s0_var)->get_dphi();
200 
201  // The pressure shape functions at interior quadrature points.
202  const std::vector<std::vector<libMesh::Real> >& u_phi =
203  context.get_element_fe(this->_flow_vars.u())->get_phi();
204 
205  // The velocity shape function gradients at interior quadrature points.
206  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
207  context.get_element_fe(this->_flow_vars.u())->get_dphi();
208 
209  // The temperature shape functions at interior quadrature points.
210  const std::vector<std::vector<libMesh::Real> >& T_phi =
211  context.get_element_fe(this->_temp_vars.T())->get_phi();
212 
213  // The temperature shape functions gradients at interior quadrature points.
214  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
215  context.get_element_fe(this->_temp_vars.T())->get_dphi();
216 
217  const std::vector<libMesh::Point>& u_qpoint =
218  context.get_element_fe(this->_flow_vars.u())->get_xyz();
219 
220  libMesh::DenseSubVector<libMesh::Number>& Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
221 
222  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
223 
224  libMesh::DenseSubVector<libMesh::Number>* Fv = NULL;
225  if( this->_flow_vars.dim() > 1 )
226  Fv = &context.get_elem_residual(this->_flow_vars.v()); // R_{v}
227 
228  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
229  if( this->_flow_vars.dim() == 3 )
230  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
231 
232  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
233 
234  unsigned int n_qpoints = context.get_element_qrule().n_points();
235  for (unsigned int qp=0; qp != n_qpoints; qp++)
236  {
237  libMesh::Number rho = cache.get_cached_values(Cache::MIXTURE_DENSITY)[qp];
238 
239  libMesh::Number u, T;
240  u = cache.get_cached_values(Cache::X_VELOCITY)[qp];
241 
242  T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
243 
244  const libMesh::Gradient& grad_T =
245  cache.get_cached_gradient_values(Cache::TEMPERATURE_GRAD)[qp];
246 
247  libMesh::NumberVectorValue U(u);
248  if (this->_flow_vars.dim() > 1)
249  U(1) = cache.get_cached_values(Cache::Y_VELOCITY)[qp];
250  if (this->_flow_vars.dim() == 3)
251  U(2) = cache.get_cached_values(Cache::Z_VELOCITY)[qp]; // w
252 
253  libMesh::Gradient grad_u = cache.get_cached_gradient_values(Cache::X_VELOCITY_GRAD)[qp];
254  libMesh::Gradient grad_v;
255 
256  if (this->_flow_vars.dim() > 1)
257  grad_v = cache.get_cached_gradient_values(Cache::Y_VELOCITY_GRAD)[qp];
258 
259  libMesh::Gradient grad_w;
260  if (this->_flow_vars.dim() == 3)
261  grad_w = cache.get_cached_gradient_values(Cache::Z_VELOCITY_GRAD)[qp];
262 
263  libMesh::Number divU = grad_u(0);
264  if (this->_flow_vars.dim() > 1)
265  divU += grad_v(1);
266  if (this->_flow_vars.dim() == 3)
267  divU += grad_w(2);
268 
269  libMesh::NumberVectorValue grad_uT( grad_u(0) );
270  libMesh::NumberVectorValue grad_vT;
271  if( this->_flow_vars.dim() > 1 )
272  {
273  grad_uT(1) = grad_v(0);
274  grad_vT = libMesh::NumberVectorValue( grad_u(1), grad_v(1) );
275  }
276  libMesh::NumberVectorValue grad_wT;
277  if( this->_flow_vars.dim() == 3 )
278  {
279  grad_uT(2) = grad_w(0);
280  grad_vT(2) = grad_w(1);
281  grad_wT = libMesh::NumberVectorValue( grad_u(2), grad_v(2), grad_w(2) );
282  }
283 
284  libMesh::Number mu = cache.get_cached_values(Cache::MIXTURE_VISCOSITY)[qp];
285  libMesh::Number p = cache.get_cached_values(Cache::PRESSURE)[qp];
286 
287  const std::vector<libMesh::Real>& D =
288  cache.get_cached_vector_values(Cache::DIFFUSION_COEFFS)[qp];
289 
290  libMesh::Number cp =
291  cache.get_cached_values(Cache::MIXTURE_SPECIFIC_HEAT_P)[qp];
292 
293  libMesh::Number k =
294  cache.get_cached_values(Cache::MIXTURE_THERMAL_CONDUCTIVITY)[qp];
295 
296  const std::vector<libMesh::Real>& omega_dot =
297  cache.get_cached_vector_values(Cache::OMEGA_DOT)[qp];
298 
299  const std::vector<libMesh::Real>& h =
300  cache.get_cached_vector_values(Cache::SPECIES_ENTHALPY)[qp];
301 
302  const libMesh::Number r = u_qpoint[qp](0);
303 
304  libMesh::Real jac = JxW[qp];
305 
307  {
308  divU += U(0)/r;
309  jac *= r;
310  }
311 
312  libMesh::Real M = cache.get_cached_values(Cache::MOLAR_MASS)[qp];
313 
314  std::vector<libMesh::Gradient> grad_ws = cache.get_cached_vector_gradient_values(Cache::MASS_FRACTIONS_GRAD)[qp];
315  libmesh_assert_equal_to( grad_ws.size(), this->_n_species );
316 
317  // Continuity Residual
318  libMesh::Gradient mass_term(0.0,0.0,0.0);
319  for(unsigned int s=0; s < this->_n_species; s++ )
320  {
321  mass_term += grad_ws[s]/this->_gas_mixture->M(s);
322  }
323  mass_term *= M;
324 
325  for (unsigned int i=0; i != n_p_dofs; i++)
326  {
327  Fp(i) += (-U*(mass_term + grad_T/T) + divU)*jac*p_phi[i][qp];
328  }
329 
330  // Species residuals
331  for(unsigned int s=0; s < this->_n_species; s++ )
332  {
333  libMesh::DenseSubVector<libMesh::Number> &Fs =
334  context.get_elem_residual(this->_species_vars.species(s)); // R_{s}
335 
336  const libMesh::Real term1 = -rho*(U*grad_ws[s]) + omega_dot[s];
337  const libMesh::Gradient term2 = -rho*D[s]*grad_ws[s];
338 
339  for (unsigned int i=0; i != n_s_dofs; i++)
340  {
342  Fs(i) += ( term1*s_phi[i][qp] + term2*s_grad_phi[i][qp] )*jac;
343  }
344  }
345 
346  // Momentum residuals
347  for (unsigned int i=0; i != n_u_dofs; i++)
348  {
349  Fu(i) += ( -rho*U*grad_u*u_phi[i][qp]
350  + p*u_gradphi[i][qp](0)
351  - mu*(u_gradphi[i][qp]*grad_u + u_gradphi[i][qp]*grad_uT
352  - 2.0/3.0*divU*u_gradphi[i][qp](0) )
353  + rho*this->_g(0)*u_phi[i][qp]
354  )*jac;
355 
359  {
360  Fu(i) += u_phi[i][qp]*( p/r - 2*mu*U(0)/(r*r) - 2.0/3.0*mu*divU/r )*jac;
361  }
362 
363  if (this->_flow_vars.dim() > 1)
364  {
365  (*Fv)(i) += ( -rho*U*grad_v*u_phi[i][qp]
366  + p*u_gradphi[i][qp](1)
367  - mu*(u_gradphi[i][qp]*grad_v + u_gradphi[i][qp]*grad_vT
368  - 2.0/3.0*divU*u_gradphi[i][qp](1) )
369  + rho*this->_g(1)*u_phi[i][qp]
370  )*jac;
371  }
372 
373  if (this->_flow_vars.dim() == 3)
374  {
375  (*Fw)(i) += ( -rho*U*grad_w*u_phi[i][qp]
376  + p*u_gradphi[i][qp](2)
377  - mu*(u_gradphi[i][qp]*grad_w + u_gradphi[i][qp]*grad_wT
378  - 2.0/3.0*divU*u_gradphi[i][qp](2) )
379  + rho*this->_g(2)*u_phi[i][qp]
380  )*jac;
381  }
382  }
383 
384  // Energy residual
385  libMesh::Real chem_term = 0.0;
386  for(unsigned int s=0; s < this->_n_species; s++ )
387  {
388  chem_term += h[s]*omega_dot[s];
389  }
390 
391  for (unsigned int i=0; i != n_T_dofs; i++)
392  {
393  FT(i) += ( ( -rho*cp*U*grad_T - chem_term )*T_phi[i][qp]
394  - k*grad_T*T_gradphi[i][qp] )*jac;
395  }
396 
397  } // quadrature loop
398 
399  return;
400  }
static bool is_axisymmetric()
Definition: physics.h:132
unsigned int VariableIndex
More descriptive name of the type used for variable indices.
Definition: var_typedefs.h:42
VariableIndex species(unsigned int species) const
VariableIndex T() const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
VariableIndex p() const
unsigned int dim() const
Number of components.
template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::ReactingLowMachNavierStokesAbstract.

Definition at line 145 of file reacting_low_mach_navier_stokes.C.

References GRINS::ReactingLowMachNavierStokesAbstract::init_context().

146  {
147  // First call base class
149 
150  // We also need the side shape functions, etc.
151  context.get_side_fe(this->_flow_vars.u())->get_JxW();
152  context.get_side_fe(this->_flow_vars.u())->get_phi();
153  context.get_side_fe(this->_flow_vars.u())->get_dphi();
154  context.get_side_fe(this->_flow_vars.u())->get_xyz();
155 
156  context.get_side_fe(this->_temp_vars.T())->get_JxW();
157  context.get_side_fe(this->_temp_vars.T())->get_phi();
158  context.get_side_fe(this->_temp_vars.T())->get_dphi();
159  context.get_side_fe(this->_temp_vars.T())->get_xyz();
160  }
VariableIndex T() const
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::mass_residual ( bool  ,
AssemblyContext  
)
virtual

Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part.

Reimplemented from GRINS::Physics.

Definition at line 418 of file reacting_low_mach_navier_stokes.C.

419  {
420  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
421  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
422  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
423  const VariableIndex s0_var = this->_species_vars.species(0);
424  const unsigned int n_s_dofs = context.get_dof_indices(s0_var).size();
425 
426  const std::vector<libMesh::Real> &JxW =
427  context.get_element_fe(this->_flow_vars.u())->get_JxW();
428 
429  const std::vector<std::vector<libMesh::Real> >& p_phi =
430  context.get_element_fe(this->_press_var.p())->get_phi();
431 
432  const std::vector<std::vector<libMesh::Real> >& u_phi =
433  context.get_element_fe(this->_flow_vars.u())->get_phi();
434 
435  const std::vector<std::vector<libMesh::Real> >& T_phi =
436  context.get_element_fe(this->_temp_vars.T())->get_phi();
437 
438  const std::vector<std::vector<libMesh::Real> >& s_phi =
439  context.get_element_fe(s0_var)->get_phi();
440 
441 
442  // The subvectors and submatrices we need to fill:
443  libMesh::DenseSubVector<libMesh::Real> &F_p = context.get_elem_residual(this->_press_var.p());
444 
445  // The subvectors and submatrices we need to fill:
446  libMesh::DenseSubVector<libMesh::Real> &F_u = context.get_elem_residual(this->_flow_vars.u());
447 
448  libMesh::DenseSubVector<libMesh::Real>* F_v = NULL;
449  if( this->_flow_vars.dim() > 1 )
450  F_v = &context.get_elem_residual(this->_flow_vars.v());
451 
452  libMesh::DenseSubVector<libMesh::Real>* F_w = NULL;
453  if( this->_flow_vars.dim() == 3 )
454  F_w = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
455 
456  libMesh::DenseSubVector<libMesh::Real> &F_T = context.get_elem_residual(this->_temp_vars.T());
457 
458  unsigned int n_qpoints = context.get_element_qrule().n_points();
459 
460  const std::vector<libMesh::Point>& u_qpoint =
461  context.get_element_fe(this->_flow_vars.u())->get_xyz();
462 
463  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
464  {
465  libMesh::Real u_dot, v_dot = 0.0, w_dot = 0.0;
466  context.interior_rate(this->_flow_vars.u(), qp, u_dot);
467 
468  if( this->_flow_vars.dim() > 1 )
469  context.interior_rate(this->_flow_vars.v(), qp, v_dot);
470 
471  if( this->_flow_vars.dim() == 3 )
472  context.interior_rate(this->_flow_vars.w(), qp, w_dot);
473 
474  libMesh::Real T_dot;
475  context.interior_rate(this->_temp_vars.T(), qp, T_dot);
476 
477  libMesh::Real T = context.interior_value(this->_temp_vars.T(), qp);
478 
479  std::vector<libMesh::Real> ws(this->n_species());
480  for(unsigned int s=0; s < this->_n_species; s++ )
481  ws[s] = context.interior_value(this->_species_vars.species(s), qp);
482 
483  Evaluator gas_evaluator( *(this->_gas_mixture) );
484  const libMesh::Real R_mix = gas_evaluator.R_mix(ws);
485  const libMesh::Real p0 = this->get_p0_steady(context,qp);
486  const libMesh::Real rho = this->rho(T, p0, R_mix);
487  const libMesh::Real cp = gas_evaluator.cp(T,p0,ws);
488  const libMesh::Real M = gas_evaluator.M_mix( ws );
489 
490  libMesh::Real jac = JxW[qp];
491  const libMesh::Number r = u_qpoint[qp](0);
492 
493  if( this->_is_axisymmetric )
494  jac *= r;
495 
496  // M_dot = -M^2 \sum_s w_dot[s]/Ms
497  libMesh::Real M_dot = 0.0;
498 
499  // Species residual
500  for(unsigned int s=0; s < this->n_species(); s++)
501  {
502  libMesh::DenseSubVector<libMesh::Number> &F_s =
503  context.get_elem_residual(this->_species_vars.species(s));
504 
505  libMesh::Real ws_dot;
506  context.interior_rate(this->_species_vars.species(s), qp, ws_dot);
507 
508  for (unsigned int i = 0; i != n_s_dofs; ++i)
509  F_s(i) -= rho*ws_dot*s_phi[i][qp]*jac;
510 
511  // Start accumulating M_dot
512  M_dot += ws_dot/this->_gas_mixture->M(s);
513  }
514 
515  // Continuity residual
516  // M_dot = -M^2 \sum_s w_dot[s]/Ms
517  libMesh::Real M_dot_over_M = M_dot*(-M);
518 
519  for (unsigned int i = 0; i != n_p_dofs; ++i)
520  F_p(i) -= (T_dot/T-M_dot_over_M)*p_phi[i][qp]*jac;
521 
522  // Momentum residual
523  for (unsigned int i = 0; i != n_u_dofs; ++i)
524  {
525  F_u(i) -= rho*u_dot*u_phi[i][qp]*jac;
526 
527  if( this->_flow_vars.dim() > 1 )
528  (*F_v)(i) -= rho*v_dot*u_phi[i][qp]*jac;
529 
530  if( this->_flow_vars.dim() == 3 )
531  (*F_w)(i) -= rho*w_dot*u_phi[i][qp]*jac;
532  }
533 
534  // Energy residual
535  for (unsigned int i = 0; i != n_T_dofs; ++i)
536  F_T(i) -= rho*cp*T_dot*T_phi[i][qp]*jac;
537 
538  if( compute_jacobian )
539  libmesh_not_implemented();
540 
541  }
542  }
unsigned int VariableIndex
More descriptive name of the type used for variable indices.
Definition: var_typedefs.h:42
VariableIndex species(unsigned int species) const
VariableIndex T() const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
VariableIndex p() const
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int dim() const
Number of components.
static bool _is_axisymmetric
Caches whether we are solving an axisymmetric problem or not.
Definition: physics.h:269
template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtual

Register postprocessing variables for ReactingLowMachNavierStokes.

Reimplemented from GRINS::Physics.

Definition at line 65 of file reacting_low_mach_navier_stokes.C.

References GRINS::PhysicsNaming::reacting_low_mach_navier_stokes(), and GRINS::PostProcessedQuantities< NumericType >::register_quantity().

67  {
68  std::string section = "Physics/"+PhysicsNaming::reacting_low_mach_navier_stokes()+"/output_vars";
69 
70  if( input.have_variable(section) )
71  {
72  unsigned int n_vars = input.vector_variable_size(section);
73 
74  for( unsigned int v = 0; v < n_vars; v++ )
75  {
76  std::string name = input(section,"DIE!",v);
77 
78  if( name == std::string("rho") )
79  {
80  this->_rho_index = postprocessing.register_quantity( name );
81  }
82  else if( name == std::string("mu") )
83  {
84  this->_mu_index = postprocessing.register_quantity( name );
85  }
86  else if( name == std::string("k") )
87  {
88  this->_k_index = postprocessing.register_quantity( name );
89  }
90  else if( name == std::string("cp") )
91  {
92  this->_cp_index = postprocessing.register_quantity( name );
93  }
94  else if( name == std::string("mole_fractions") )
95  {
96  this->_mole_fractions_index.resize(this->n_species());
97 
98  for( unsigned int s = 0; s < this->n_species(); s++ )
99  {
100  this->_mole_fractions_index[s] = postprocessing.register_quantity( "X_"+this->_gas_mixture->species_name(s) );
101  }
102  }
103  else if( name == std::string("h_s") )
104  {
105  this->_h_s_index.resize(this->n_species());
106 
107  for( unsigned int s = 0; s < this->n_species(); s++ )
108  {
109  this->_h_s_index[s] = postprocessing.register_quantity( "h_"+this->_gas_mixture->species_name(s) );
110  }
111  }
112  else if( name == std::string("omega_dot") )
113  {
114  this->_omega_dot_index.resize(this->n_species());
115 
116  for( unsigned int s = 0; s < this->n_species(); s++ )
117  this->_omega_dot_index[s] = postprocessing.register_quantity( "omega_dot_"+this->_gas_mixture->species_name(s) );
118  }
119  else if( name == std::string("D_s") )
120  {
121  this->_Ds_index.resize(this->n_species());
122 
123  for( unsigned int s = 0; s < this->n_species(); s++ )
124  this->_Ds_index[s] = postprocessing.register_quantity( "D_"+this->_gas_mixture->species_name(s) );
125  }
126  else
127  {
128  std::cerr << "Error: Invalue output_vars value for "+PhysicsNaming::reacting_low_mach_navier_stokes() << std::endl
129  << " Found " << name << std::endl
130  << " Acceptable values are: rho" << std::endl
131  << " mu" << std::endl
132  << " k" << std::endl
133  << " cp" << std::endl
134  << " mole_fractions" << std::endl
135  << " omega_dot" << std::endl;
136  libmesh_error();
137  }
138  }
139  }
140 
141  return;
142  }
static PhysicsName reacting_low_mach_navier_stokes()
std::vector< unsigned int > _mole_fractions_index
Index from registering this quantity. Each species will have it's own index.
unsigned int _k_index
Index from registering this quantity.
std::vector< unsigned int > _Ds_index
Index from registering this quantity. Each species will have it's own index.
unsigned int _mu_index
Index from registering this quantity.
std::vector< unsigned int > _omega_dot_index
Index from registering this quantity. Each species will have it's own index.
unsigned int _rho_index
Index from registering this quantity.
std::vector< unsigned int > _h_s_index
Index from registering this quantity. Each species will have it's own index.
unsigned int _cp_index
Index from registering this quantity.

Member Data Documentation

template<typename Mixture , typename Evaluator >
unsigned int GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_cp_index
protected

Index from registering this quantity.

Definition at line 93 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
std::vector<unsigned int> GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_Ds_index
protected

Index from registering this quantity. Each species will have it's own index.

Definition at line 105 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
std::vector<unsigned int> GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_h_s_index
protected

Index from registering this quantity. Each species will have it's own index.

Definition at line 99 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
unsigned int GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_k_index
protected

Index from registering this quantity.

Definition at line 90 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
std::vector<unsigned int> GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_mole_fractions_index
protected

Index from registering this quantity. Each species will have it's own index.

Definition at line 96 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
unsigned int GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_mu_index
protected

Index from registering this quantity.

Definition at line 87 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
std::vector<unsigned int> GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_omega_dot_index
protected

Index from registering this quantity. Each species will have it's own index.

Definition at line 102 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
PressurePinning GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_p_pinning
protected

Definition at line 78 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
bool GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_pin_pressure
protected
template<typename Mixture , typename Evaluator >
unsigned int GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_rho_index
protected

Index from registering this quantity.

Definition at line 81 of file reacting_low_mach_navier_stokes.h.

template<typename Mixture , typename Evaluator >
std::vector<unsigned int> GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >::_species_viscosity
protected

Index from registering this quantity. Each species will have it's own index.

Definition at line 84 of file reacting_low_mach_navier_stokes.h.


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

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