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

#include <reacting_low_mach_navier_stokes_spgsm_stab.h>

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

Public Member Functions

 ReactingLowMachNavierStokesSPGSMStabilization (const GRINS::PhysicsName &physics_name, const GetPot &input, libMesh::UniquePtr< Mixture > &gas_mix)
 
virtual ~ReactingLowMachNavierStokesSPGSMStabilization ()
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent 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...
 
- Public Member Functions inherited from GRINS::ReactingLowMachNavierStokesStabilizationBase< Mixture, Evaluator >
 ReactingLowMachNavierStokesStabilizationBase (const GRINS::PhysicsName &physics_name, const GetPot &input, libMesh::UniquePtr< Mixture > &gas_mix)
 
virtual ~ReactingLowMachNavierStokesStabilizationBase ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
void compute_res_steady (AssemblyContext &context, unsigned int qp, libMesh::Real &RP_s, libMesh::RealGradient &RM_s, libMesh::Real &RE_s, std::vector< libMesh::Real > &Rs_s)
 
void compute_res_transient (AssemblyContext &context, unsigned int qp, libMesh::Real &RP_t, libMesh::RealGradient &RM_t, libMesh::Real &RE_t, std::vector< libMesh::Real > &Rs_t)
 
- 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 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 preassembly (MultiphysicsSystem &)
 Perform any necessary setup before element assembly begins. More...
 
virtual void reinit (MultiphysicsSystem &)
 Any reinitialization that needs to be done. More...
 
virtual void side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

Private Member Functions

 ReactingLowMachNavierStokesSPGSMStabilization ()
 

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...
 
- Protected Attributes inherited from GRINS::ReactingLowMachNavierStokesStabilizationBase< Mixture, Evaluator >
ReactingLowMachNavierStokesStabilizationHelper _stab_helper
 
- 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...
 
- 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::ReactingLowMachNavierStokesSPGSMStabilization< Mixture, Evaluator >

Definition at line 34 of file reacting_low_mach_navier_stokes_spgsm_stab.h.

Constructor & Destructor Documentation

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

Definition at line 38 of file reacting_low_mach_navier_stokes_spgsm_stab.C.

39  : ReactingLowMachNavierStokesStabilizationBase<Mixture,Evaluator>(physics_name,input,gas_mix)
40  {}
template<typename Mixture , typename Evaluator >
virtual GRINS::ReactingLowMachNavierStokesSPGSMStabilization< Mixture, Evaluator >::~ReactingLowMachNavierStokesSPGSMStabilization ( )
inlinevirtual

Definition at line 40 of file reacting_low_mach_navier_stokes_spgsm_stab.h.

40 {};
template<typename Mixture , typename Evaluator >
GRINS::ReactingLowMachNavierStokesSPGSMStabilization< Mixture, Evaluator >::ReactingLowMachNavierStokesSPGSMStabilization ( )
private

Member Function Documentation

template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokesSPGSMStabilization< Mixture, Evaluator >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 44 of file reacting_low_mach_navier_stokes_spgsm_stab.C.

46  {
47  // The number of local degrees of freedom in each variable.
48  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
49  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
50  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
51  const VariableIndex s0_var = this->_species_vars.species(0);
52  const unsigned int n_s_dofs = context.get_dof_indices(s0_var).size();
53 
54  // Element Jacobian * quadrature weights for interior integration.
55  const std::vector<libMesh::Real> &JxW =
56  context.get_element_fe(this->_flow_vars.u())->get_JxW();
57 
58  // The pressure shape functions at interior quadrature points.
59  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
60  context.get_element_fe(this->_press_var.p())->get_dphi();
61 
62  const std::vector<std::vector<libMesh::Real> >& u_phi =
63  context.get_element_fe(this->_flow_vars.u())->get_phi();
64 
65  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
66  context.get_element_fe(this->_flow_vars.u())->get_dphi();
67 
68  // The temperature shape functions gradients at interior quadrature points.
69  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
70  context.get_element_fe(this->_temp_vars.T())->get_dphi();
71 
72  const std::vector<std::vector<libMesh::Gradient> >& s_gradphi = context.get_element_fe(s0_var)->get_dphi();
73 
74  // We're assuming the quadrature rule is the same for all variables
75  const std::vector<libMesh::Point>& u_qpoint =
76  context.get_element_fe(this->_flow_vars.u())->get_xyz();
77 
78  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
79 
80  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
81 
82  libMesh::DenseSubVector<libMesh::Number>* Fv = NULL;
83  if( this->_flow_vars.dim() > 1 )
84  Fv = &context.get_elem_residual(this->_flow_vars.v()); // R_{v}
85 
86  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
87  if( this->_flow_vars.dim() == 3 )
88  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
89 
90  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
91 
92  unsigned int n_qpoints = context.get_element_qrule().n_points();
93 
94 
95  libMesh::FEBase* u_fe = context.get_element_fe(this->_flow_vars.u());
96  for (unsigned int qp=0; qp != n_qpoints; qp++)
97  {
98  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
99 
100  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ) );
101  if( this->_flow_vars.dim() > 1 )
102  U(1) = context.interior_value( this->_flow_vars.v(), qp );
103  if( this->_flow_vars.dim() == 3 )
104  U(2) = context.interior_value( this->_flow_vars.w(), qp );
105 
106  std::vector<libMesh::Real> ws(this->n_species());
107  for(unsigned int s=0; s < this->_n_species; s++ )
108  {
109  ws[s] = context.fixed_interior_value(this->_species_vars.species(s), qp);
110  }
111 
112  Evaluator gas_evaluator( *(this->_gas_mixture) );
113  const libMesh::Real R_mix = gas_evaluator.R_mix(ws);
114  const libMesh::Real p0 = this->get_p0_steady(context,qp);
115  libMesh::Real rho = this->rho(T, p0, R_mix);
116 
117  const libMesh::Real cp = gas_evaluator.cp(T,p0,ws);
118 
119  std::vector<libMesh::Real> D( this->n_species() );
120  libMesh::Real mu, k;
121 
122  gas_evaluator.mu_and_k_and_D( T, rho, cp, ws, mu, k, D );
123 
124  libMesh::RealGradient g = this->_stab_helper.compute_g( u_fe, context, qp );
125  libMesh::RealTensor G = this->_stab_helper.compute_G( u_fe, context, qp );
126 
127  // Taus
128  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, this->_is_steady );
129 
130  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
131 
132  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, this->_is_steady );
133 
134 
135 
136  // Strong form residuals
137  libMesh::RealGradient RM_s = 0.0;
138  libMesh::Real RC_s = 0.0;
139  libMesh::Real RE_s = 0.0;
140  std::vector<libMesh::Real> Rs_s;
141 
142  this->compute_res_steady( context, qp, RC_s, RM_s, RE_s, Rs_s );
143 
144  const libMesh::Number r = u_qpoint[qp](0);
145 
146  libMesh::Real jac = JxW[qp];
147 
148  if( this->_is_axisymmetric )
149  jac *= r;
150 
151  // Pressure PSPG term
152  for (unsigned int i=0; i != n_p_dofs; i++)
153  Fp(i) += tau_M*RM_s*p_dphi[i][qp]*jac;
154 
155  // Momentum SUPG + div-div terms
156  for (unsigned int i=0; i != n_u_dofs; i++)
157  {
158  Fu(i) += ( - tau_C*RC_s*u_gradphi[i][qp](0)
159  - tau_M*RM_s(0)*rho*U*u_gradphi[i][qp] )*jac;
160 
161  if( this->_is_axisymmetric )
162  Fu(i) += (-tau_C*RC_s/r)*u_phi[i][qp]*jac;
163 
164  if( this->_flow_vars.dim() > 1 )
165  (*Fv)(i) += ( - tau_C*RC_s*u_gradphi[i][qp](1)
166  - tau_M*RM_s(1)*rho*U*u_gradphi[i][qp] )*jac;
167 
168  if( this->_flow_vars.dim() == 3 )
169  (*Fw)(i) += ( - tau_C*RC_s*u_gradphi[i][qp](2)
170  - tau_M*RM_s(2)*rho*U*u_gradphi[i][qp] )*jac;
171  }
172 
173  // Energy SUPG terms
174  for (unsigned int i=0; i != n_T_dofs; i++)
175  FT(i) -= rho*cp*tau_E*RE_s*U*T_gradphi[i][qp]*jac;
176 
177  // Species SUPG terms
178  for(unsigned int s=0; s < this->n_species(); s++)
179  {
180  libMesh::DenseSubVector<libMesh::Number> &Fs =
181  context.get_elem_residual(this->_species_vars.species(s));
182 
183  libMesh::Real tau_s = this->_stab_helper.compute_tau_species( context, qp, g, G, rho, U, D[s], this->_is_steady );
184 
185  for (unsigned int i=0; i != n_s_dofs; i++)
186  Fs(i) -= rho*tau_s*Rs_s[s]*U*s_gradphi[i][qp]*jac;
187  }
188 
189  if(compute_jacobian)
190  libmesh_not_implemented();
191  }
192  }
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::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::Real compute_tau_continuity(libMesh::Real tau_C, libMesh::RealGradient &g) const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
VariableIndex p() const
void compute_res_steady(AssemblyContext &context, unsigned int qp, libMesh::Real &RP_s, libMesh::RealGradient &RM_s, libMesh::Real &RE_s, std::vector< libMesh::Real > &Rs_s)
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
libMesh::Real compute_tau_species(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real D_s, bool is_steady) const
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int dim() const
Number of components.
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
static bool _is_axisymmetric
Caches whether we are solving an axisymmetric problem or not.
Definition: physics.h:269
libMesh::Real compute_tau_energy(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real k, libMesh::Real cp, bool is_steady) const
template<typename Mixture , typename Evaluator >
void GRINS::ReactingLowMachNavierStokesSPGSMStabilization< 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 196 of file reacting_low_mach_navier_stokes_spgsm_stab.C.

197  {
198  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
199  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
200  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
201  const VariableIndex s0_var = this->_species_vars.species(0);
202  const unsigned int n_s_dofs = context.get_dof_indices(s0_var).size();
203 
204  const std::vector<libMesh::Real> &JxW =
205  context.get_element_fe(this->_flow_vars.u())->get_JxW();
206 
207  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
208  context.get_element_fe(this->_press_var.p())->get_dphi();
209 
210  const std::vector<std::vector<libMesh::Real> >& u_phi =
211  context.get_element_fe(this->_flow_vars.u())->get_phi();
212 
213  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
214  context.get_element_fe(this->_flow_vars.u())->get_dphi();
215 
216  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
217  context.get_element_fe(this->_temp_vars.T())->get_dphi();
218 
219  const std::vector<std::vector<libMesh::Gradient> >& s_gradphi = context.get_element_fe(s0_var)->get_dphi();
220 
221  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
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());
233 
234  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
235 
236  unsigned int n_qpoints = context.get_element_qrule().n_points();
237 
238  const std::vector<libMesh::Point>& u_qpoint =
239  context.get_element_fe(this->_flow_vars.u())->get_xyz();
240 
241  for (unsigned int qp=0; qp != n_qpoints; qp++)
242  {
243  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
244  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
245 
246  libMesh::Real T = context.interior_value( this->_temp_vars.T(), qp );
247 
248  libMesh::Number u = context.fixed_interior_value(this->_flow_vars.u(), qp);
249 
250  libMesh::NumberVectorValue U(u);
251  if (this->_flow_vars.dim() > 1)
252  U(1) = context.fixed_interior_value(this->_flow_vars.v(), qp);
253  if (this->_flow_vars.dim() == 3)
254  U(2) = context.fixed_interior_value(this->_flow_vars.w(), qp);
255 
256  std::vector<libMesh::Real> ws(this->n_species());
257  for(unsigned int s=0; s < this->_n_species; s++ )
258  ws[s] = context.fixed_interior_value(this->_species_vars.species(s), qp);
259 
260  Evaluator gas_evaluator( *(this->_gas_mixture) );
261  const libMesh::Real R_mix = gas_evaluator.R_mix(ws);
262  const libMesh::Real p0 = this->get_p0_steady(context,qp);
263  libMesh::Real rho = this->rho(T, p0, R_mix);
264 
265  const libMesh::Real cp = gas_evaluator.cp(T,p0,ws);
266 
267  std::vector<libMesh::Real> D( this->n_species() );
268  libMesh::Real mu, k;
269 
270  gas_evaluator.mu_and_k_and_D( T, rho, cp, ws, mu, k, D );
271 
272  libMesh::Real tau_M = this->_stab_helper.compute_tau_momentum( context, qp, g, G, rho, U, mu, false );
273  libMesh::Real tau_C = this->_stab_helper.compute_tau_continuity( tau_M, g );
274  libMesh::Real tau_E = this->_stab_helper.compute_tau_energy( context, qp, g, G, rho, U, k, cp, false );
275 
276  // Strong residuals
277  libMesh::Real RC_t;
278  libMesh::RealGradient RM_t;
279  libMesh::Real RE_t;
280  std::vector<libMesh::Real> Rs_t(this->n_species());
281 
282  this->compute_res_transient( context, qp, RC_t, RM_t, RE_t, Rs_t );
283 
284  libMesh::Real jac = JxW[qp];
285  const libMesh::Number r = u_qpoint[qp](0);
286 
287  if( this->_is_axisymmetric )
288  jac *= r;
289 
290  for (unsigned int i=0; i != n_p_dofs; i++)
291  Fp(i) -= tau_M*RM_t*p_dphi[i][qp]*jac;
292 
293  for(unsigned int s=0; s < this->n_species(); s++)
294  {
295  libMesh::DenseSubVector<libMesh::Number> &Fs =
296  context.get_elem_residual(this->_species_vars.species(s));
297 
298  libMesh::Real tau_s = this->_stab_helper.compute_tau_species( context, qp, g, G, rho, U, D[s], false );
299 
300  for (unsigned int i=0; i != n_s_dofs; i++)
301  Fs(i) -= rho*tau_s*Rs_t[s]*U*s_gradphi[i][qp]*jac;
302  }
303 
304  for (unsigned int i=0; i != n_u_dofs; i++)
305  {
306  Fu(i) -= ( tau_C*RC_t*u_gradphi[i][qp](0)
307  + tau_M*RM_t(0)*rho*U*u_gradphi[i][qp] )*jac;
308 
309  if( this->_is_axisymmetric )
310  Fu(i) += (-tau_C*RC_t/r)*u_phi[i][qp]*jac;
311 
312  if( this->_flow_vars.dim() > 1 )
313  (*Fv)(i) -= ( tau_C*RC_t*u_gradphi[i][qp](1)
314  + tau_M*RM_t(1)*rho*U*u_gradphi[i][qp] )*jac;
315 
316  if( this->_flow_vars.dim() == 3 )
317  (*Fw)(i) -= ( tau_C*RC_t*u_gradphi[i][qp](2)
318  + tau_M*RM_t(2)*rho*U*u_gradphi[i][qp] )*jac;
319  }
320 
321  for (unsigned int i=0; i != n_T_dofs; i++)
322  FT(i) -= rho*cp*tau_E*RE_t*U*T_gradphi[i][qp]*jac;
323 
324  if(compute_jacobian)
325  libmesh_not_implemented();
326  }
327 
328  return;
329  }
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::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
libMesh::Real compute_tau_continuity(libMesh::Real tau_C, libMesh::RealGradient &g) const
libMesh::Real T(const libMesh::Point &p, const AssemblyContext &c) const
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
libMesh::Real rho(libMesh::Real T, libMesh::Real p0, libMesh::Real R_mix) const
VariableIndex p() const
void compute_res_transient(AssemblyContext &context, unsigned int qp, libMesh::Real &RP_t, libMesh::RealGradient &RM_t, libMesh::Real &RE_t, std::vector< libMesh::Real > &Rs_t)
libMesh::Real compute_tau_species(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real D_s, bool is_steady) const
libMesh::Real get_p0_steady(const AssemblyContext &c, unsigned int qp) const
unsigned int dim() const
Number of components.
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
static bool _is_axisymmetric
Caches whether we are solving an axisymmetric problem or not.
Definition: physics.h:269
libMesh::Real compute_tau_energy(AssemblyContext &c, unsigned int qp, libMesh::RealGradient &g, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Gradient U, libMesh::Real k, libMesh::Real cp, bool is_steady) const

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