OpenFCST: The open-source Fuel Cell Simulation Toolbox
|
This class implements the multi-component mass transport equations proposed by Kerkhof-Geboers for fluid transport. More...
#include <compressible_multi_component_KG_equations_coupled.h>
Public Member Functions | |
Constructors, destructor, and initialization | |
CompressibleMultiComponentKGEquationsCoupled (FuelCell::SystemManagement &system_management, boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > data=boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >()) | |
Constructor. More... | |
virtual | ~CompressibleMultiComponentKGEquationsCoupled () |
Destructor. More... | |
virtual void | declare_parameters (ParameterHandler ¶m) const |
Declare parameters. More... | |
virtual void | initialize (ParameterHandler ¶m) |
Initialize parameters. More... | |
Local CG FEM based assemblers | |
virtual void | assemble_cell_matrix (FuelCell::ApplicationCore::MatrixVector &cell_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local cell matrix. More... | |
virtual void | assemble_cell_residual (FuelCell::ApplicationCore::FEVector &cell_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local cell residual. More... | |
virtual void | assemble_bdry_matrix (FuelCell::ApplicationCore::MatrixVector &bdry_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local boundary matrix. More... | |
virtual void | assemble_bdry_residual (FuelCell::ApplicationCore::FEVector &bdry_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local boundary residual. More... | |
Accessors and info | |
const std::vector< double > & | get_inlet_outlet_velocity_max () const |
This function returns inlet_outlet_velocity_max . More... | |
const double | get_inlet_outlet_velocity_mixture_max () const |
This function returns inlet_outlet_velocity_mixture_max . More... | |
const bool | get_use_parabolic_profile () const |
This function tells the application if the user wishes to use a parabolic profile for a boundary condition. More... | |
const double | get_inlet_outlet_velocity_channel_base () const |
This function tells the application the height [cm] of the base of the channel with respect to origin. More... | |
const double | get_inlet_outlet_velocity_channel_roof () const |
This function tells the application the height [cm] of the roof of the channel with respect to origin. More... | |
const int | get_inlet_outlet_boundary_ID () const |
This function tells the application the boundary id to apply velocity profile equation to. More... | |
const std::string | get_press_vel_comp_apply_to () const |
Outputs which component to apply variable boundary equation to (i.e. More... | |
const std::string | get_inlet_outlet_velocity_mixture_equation () const |
If get_use_parabolic_profile() is false then user can specify an equation to use for velocity profile at boundary. More... | |
const std::map< unsigned int, std::string > | get_gas_species_map () const |
get_gas_species_map() returns a map relating the different KG species equations to the desired gas species being simulated; e.g 1:oxygen, 2:water vapour, 3:nitrogen means that the first 3 KG equations (Mass Cons., Mom. More... | |
const unsigned int | get_num_of_species () const |
If returns the number of species to use in the parameter file. More... | |
virtual void | print_equation_info () const |
This function prints out the equations info. More... | |
Public Member Functions inherited from FuelCellShop::Equation::EquationBase< dim > | |
const couplings_map & | get_internal_cell_couplings () const |
This function returns internal_cell_couplings of a derived equation class. More... | |
const couplings_map & | get_internal_flux_couplings () const |
This function returns internal_flux_couplings (DG FEM only) of a derived equation class. More... | |
const component_materialID_value_map & | get_component_materialID_value () const |
This function returns component_materialID_value of a derived equation class. More... | |
const component_boundaryID_value_map & | get_component_boundaryID_value () const |
This function returns component_boundaryID_value of a derived equation class. More... | |
const std::vector< BoundaryType > & | get_boundary_types () const |
This function returns boundary_types of a derived equation class. More... | |
const std::vector< std::vector < BoundaryType > > & | get_multi_boundary_types () const |
This function returns multi_boundary_types of a derived equation class. More... | |
const std::vector< OutputType > & | get_output_types () const |
This function returns output_types of a derived equation class. More... | |
const std::vector< std::vector < OutputType > > & | get_multi_output_types () const |
This function returns multi_output_types of a derived equation class. More... | |
const std::string & | get_equation_name () const |
This function returns equation_name of a derived equation class. More... | |
const std::vector< unsigned int > & | get_matrix_block_indices () const |
This function returns matrix_block_indices of a derived equation class. More... | |
const std::vector< unsigned int > & | get_residual_indices () const |
This function returns residual_indices of a derived equation class. More... | |
Protected Member Functions | |
Local CG FEM based assemblers - make_ functions | |
template<typename INFO > | |
void | make_assemblers_generic_constant_data (const INFO &InFo, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
This function computes Computational constants and Local CG FEM based assemblers - constant data (generic). More... | |
virtual void | make_assemblers_cell_constant_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info) |
This function computes Local CG FEM based assemblers - constant data (cell) and allocates the memory for. More... | |
virtual void | make_assemblers_bdry_constant_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info) |
This function computes Local CG FEM based assemblers - constant data (boundary) and allocates the memory for. More... | |
virtual void | make_assemblers_cell_variable_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
This function computes. More... | |
virtual void | make_assemblers_bdry_variable_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
This function computes. More... | |
Other - make_ functions | |
virtual void | make_internal_cell_couplings () |
This function fills out internal_cell_couplings . More... | |
virtual void | make_matrix_block_indices () |
This function fills out matrix_block_indices . More... | |
virtual void | make_residual_indices () |
This function fills out residual_indices . More... | |
Auxiliary classes | |
template<typename LAYER , typename INFO > | |
void | update_porosity (const LAYER *ptr, const INFO info, std::vector< double > &porosity) |
Inline function to obtain update the porosity for a given layer class. More... | |
template<typename LAYER , typename INFO > | |
void | update_permeability (const LAYER *ptr, const INFO info, std::vector< SymmetricTensor< 2, dim > > &permeability_INV, std::vector< SymmetricTensor< 2, dim > > &Forchheimer_permeability) |
Inline member function used to obtain the inverse permeability and the Forchheimer permeability for a given layer. More... | |
template<typename LAYER , typename INFO > | |
void | update_tortuosity (const LAYER *ptr, const INFO info, std::vector< SymmetricTensor< 2, dim > > &tortuosity) |
Inline member function used to obtain the tortuosity for a given layer. More... | |
template<typename LAYER > | |
void | update_gas_properties (const LAYER *ptr) |
Inline member function used to update. More... | |
template<typename LAYER > | |
void | update_inv_diffusion_coefficients (const LAYER *ptr, std::vector< Table< 2, SymmetricTensor< 2, dim > > > &inv_pDeff) const |
This private member function is used in make_assemblers_cell_variable_data in order to compute the inverse of the product of pressure and molecular diffusion coefficient. More... | |
template<typename LAYER > | |
void | update_partial_viscosities (const LAYER *ptr, const unsigned int &quadraturePoints, const std::vector< double > &porosity, const std::vector< std::vector< double > > &density, std::vector< std::vector< std::vector< double > > > ¶mMatrix, std::vector< FullMatrix< double > > &PInv, std::vector< std::vector< double > > &partialViscosity, std::vector< std::vector< double > > &bulkViscosity) const |
This private member function is used to calculate the partial viscosity of the mixture and bulk viscosity, from dynamicViscosity, density, molar_mass. More... | |
template<typename LAYER > | |
void | update_delta_partial_viscosities (const LAYER *ptr, const std::vector< std::vector< double > > &partialViscosity, const std::vector< std::vector< std::vector< double > > > ¶mMatrix, std::vector< FullMatrix< double > > &PInv, const std::vector< double > &porosity, const std::vector< std::vector< std::vector< double > > > &deltaDensity, const std::vector< std::vector< double > > &density, std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity) const |
This private member function is used to calculate the variation in partial viscosity of the mixture. More... | |
template<typename LAYER > | |
void | update_delta_bulk_viscosities (const LAYER *ptr, const std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity, std::vector< std::vector< std::vector< double > > > &deltaBulkViscosity) const |
Calls GasMixture function to loop through all variations of partial viscosity and calculates the varation in bulk viscosity according to bulk_viscosity_mode. More... | |
Protected Member Functions inherited from FuelCellShop::Equation::EquationBase< dim > | |
EquationBase (FuelCell::SystemManagement &sys_management, boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > data=boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >()) | |
Constructor. More... | |
virtual | ~EquationBase () |
Destructor. More... | |
virtual void | set_parameters (const std::vector< std::string > &name_dvar, const std::vector< double > &value_dvar, ParameterHandler ¶m) |
Set parameters using the parameter file, in order to run parametric/optimization studies. More... | |
virtual void | make_assemblers_generic_constant_data () |
Function used to initialize variable information that will be needed to assemble matrix and residual and that will remain constant throughout the simulation. More... | |
void | select_cell_assemblers (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
This routine is used to select the make_assembly routines that need to be called inside assemble_cell_matrix to compute. More... | |
virtual void | make_internal_flux_couplings () |
This function fills out internal_flux_couplings (DG FEM only) of a derived equation class. More... | |
virtual void | make_component_materialID_value () |
This function fills out component_materialID_value of a derived equation class. More... | |
virtual void | make_component_boundaryID_value () |
This function fills out component_boundaryID_value of a derived equation class. More... | |
virtual void | make_boundary_types () |
This function fills out boundary_types of a derived equation class. More... | |
virtual void | make_multi_boundary_types () |
This function fills out multi_boundary_types of a derived equation class. More... | |
virtual void | make_output_types () |
This function fills out output_types of a derived equation class. More... | |
virtual void | make_multi_output_types () |
This function fills out multi_output_types of a derived equation class. More... | |
virtual void | assemble_cell_Jacobian_matrix (FuelCell::ApplicationCore::MatrixVector &cell_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble the local Jacobian Matrix for Non-Linear problems. More... | |
virtual void | assemble_bdry_Jacobian_matrix (FuelCell::ApplicationCore::MatrixVector &bdry_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local Jacobian boundary matrix for Non-Linear problems. More... | |
virtual void | assemble_cell_linear_matrix (FuelCell::ApplicationCore::MatrixVector &cell_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble the local cell matrix for Linear problems. More... | |
virtual void | assemble_cell_residual_rhs (FuelCell::ApplicationCore::FEVector &cell_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local cell RHS for nonlinear problems. More... | |
virtual void | assemble_cell_linear_rhs (FuelCell::ApplicationCore::FEVector &cell_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local cell RHS for Linear problems. More... | |
virtual void | assemble_bdry_linear_matrix (FuelCell::ApplicationCore::MatrixVector &bdry_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local boundary matrix for linear problems. More... | |
virtual void | assemble_bdry_linear_rhs (FuelCell::ApplicationCore::FEVector &bdry_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer) |
Assemble local boundary RHS for linear problems. More... | |
void | standard_to_block_wise (FullMatrix< double > &target) const |
This function changes the order of dealii::FullMatrix<double> target from standard to block-wise. More... | |
void | standard_to_block_wise (Vector< double > &target) const |
This function changes the order of dealii::Vector<double> target from standard to block-wise. More... | |
void | dealII_to_appframe (FuelCell::ApplicationCore::MatrixVector &dst, const FullMatrix< double > &src, const std::vector< unsigned int > &matrix_block_indices) const |
This function converts the standard ordered structure dealii::FullMatrix<double> src into the block-wise ordered structure FuelCell::ApplicationCore::MatrixVector dst . More... | |
void | dealII_to_appframe (FuelCell::ApplicationCore::FEVector &dst, const Vector< double > &src, const std::vector< unsigned int > &residual_indices) const |
This function converts the standard ordered structure dealii::Vector<double> src into the block-wise ordered structure FuelCell::ApplicationCore::FEVector dst . More... | |
bool | belongs_to_boundary (const unsigned int &tria_boundary_id, const unsigned int ¶m_boundary_id) const |
This function returns true if a boundary indicator of an external face on the triangulation coincides with a boundary indicator defined in the parameters file of a derived equation class. More... | |
void | print_caller_name (const std::string &caller_name) const |
This function is used to print out the name of another function that has been declared in the scope of this class, but not yet been implemented. More... | |
Protected Attributes | |
Boolean constants and form of the drag force in porous media | |
std::vector< bool > | inertia_in_channels |
This object indicates whether the inertia term is ON or OFF in channels. More... | |
std::vector< bool > | shear_stress_in_channels |
This object indicates whether the shear stress term is ON or OFF in channels. More... | |
std::vector< bool > | gravity_in_channels |
This object indicates whether the gravity term is ON or OFF in channels. More... | |
std::vector< bool > | inertia_in_porous_media |
This object indicates whether the inertia term is ON or OFF in porous media. More... | |
std::vector< bool > | shear_stress_in_porous_media |
This object indicates whether the shear stress term is ON or OFF in porous media. More... | |
std::vector< bool > | gravity_in_porous_media |
This object indicates whether the gravity term is ON or OFF in porous media. More... | |
std::vector< std::string > | drag_in_porous_media |
This object indicates which form of the drag term is supposed to be chosen. More... | |
std::vector< bool > | normal_velocity_is_suppressed_weakly |
Sometimes we implement different types of slip boundary conditions on the curved impermeable walls of the domain or the walls that are not perfectly aligned with the coordinate axes. More... | |
Computational constants | |
unsigned int | n_species |
Number of species, . More... | |
bool | coupled_with_fuel_cell_physics |
std::map< unsigned int, std::string > | gas_species_map |
Map of the gas species to their respective equation species number. More... | |
double | R_universal |
Universal gas constant, . More... | |
double | T_mixture |
Temperature of species mixture, . More... | |
double | P_in |
Pressure at inlet, . More... | |
Tensor< 1, dim > | gravity_acceleration |
Gravitational acceleration, . More... | |
SymmetricTensor< 2, dim > | unit |
Unit tensor, . More... | |
std::vector< double > | molar_mass |
Molar mass of pure gas, . More... | |
std::vector< double > | dynamic_viscosity |
Dynamic viscosity of pure gas, . More... | |
std::vector< double > | bulk_viscosity |
Bulk viscosity of pure gas, . More... | |
std::vector< double > | collision_diameter |
Collision diameter of pure gas, . More... | |
std::vector< double > | theta |
Navier slip coefficient of pure gas, . More... | |
std::vector< double > | eta |
Normal velocity suppression coefficient of pure gas, . More... | |
Table< 2, double > | maxwell_stefan_isobaric_diffusion_coefficient |
Each entry of this structure defines a Maxwell-Stefan isobaric diffusion coefficient of gas in gas , . More... | |
std::vector< double > | inlet_outlet_velocity_max |
Maximum inlet-outlet velocity of pure gas, . More... | |
double | inlet_outlet_velocity_mixture_max |
Maximum inlet-outlet velocity of the mixture, . More... | |
bool | use_parabolic_profile |
double | inlet_outlet_velocity_channel_base |
double | inlet_outlet_velocity_channel_roof |
int | inlet_outlet_boundary_ID |
std::string | press_vel_comp_apply_to |
std::string | inlet_outlet_velocity_mixture_equation |
std::vector< double > | maxwell_constant |
These constants are used in the Maxwell slip boundary condition. More... | |
Local CG FEM based assemblers - constant data (generic) | |
std::vector < FEValuesExtractors::Scalar > | density |
Density extractors. More... | |
std::vector < FEValuesExtractors::Vector > | velocity |
Velocity extractors. More... | |
Local CG FEM based assemblers - variable data (previous Newton iteration - cell) | |
Implementation is in the base class. | |
std::vector< std::vector < double > > | density_cell_old |
Density of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | grad_density_cell_old |
Density gradient of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | velocity_cell_old |
Velocity of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < double > > | div_velocity_cell_old |
Velocity divergence of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < SymmetricTensor< 2, dim > > > | grads_velocity_cell_old |
Velocity symmetric gradient of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | mass_flux_cell_old |
Mass flux of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < SymmetricTensor< 2, dim > > > | momentum_flux_cell_old |
Momentum flux of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < double > > | pressure_cell_old |
Partial pressure of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < SymmetricTensor< 2, dim > > > | shear_stress_cell_old |
Shear stress of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < double > > | partial_viscosity_old |
Partial viscosity of mixture of each species, , in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < double > > | bulk_viscosity_old |
Bulk viscosity of mixture of each species, . More... | |
std::vector< std::vector < std::vector< double > > > | paramMatrix_old |
paramMatrix is used for calculating partial viscosity of mixture in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< FullMatrix< double > > | PInv_old |
PInv_old is used for calculating the inverse of P when using the OmegaKG model for calculating partial viscosity in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | drag_cell_old |
Drag force of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | diffusion_cell_old |
Diffusion force of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | gravity_cell_old |
Gravity force of each species in the quadrature points of a cell at a previous Newton iteration. More... | |
Local CG FEM based assemblers - variable data (current Newton iteration - cell) | |
std::vector< std::vector < std::vector< double > > > | phi_density_cell |
Density shape functions. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | grad_phi_density_cell |
Density shape function gradients. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | phi_velocity_cell |
Velocity shape functions. More... | |
std::vector< std::vector < std::vector< double > > > | div_phi_velocity_cell |
Velocity shape function divergences. More... | |
std::vector< std::vector < std::vector< SymmetricTensor < 2, dim > > > > | grads_phi_velocity_cell |
Velocity shape function symmetric gradients. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | delta_mass_flux_cell |
Mass flux perturbations. More... | |
std::vector< std::vector < std::vector< SymmetricTensor < 2, dim > > > > | delta_momentum_flux_cell |
Momentum flux perturbations. More... | |
std::vector< std::vector < std::vector< double > > > | delta_pressure_cell |
Partial pressure perturbations. More... | |
std::vector< std::vector < std::vector< double > > > | delta_partial_viscosity_cell |
Partial dynamic viscosity perturbations. More... | |
std::vector< std::vector < std::vector< double > > > | delta_bulk_viscosity_cell |
Partial bulk viscosity perturbations. More... | |
std::vector< std::vector < std::vector< SymmetricTensor < 2, dim > > > > | delta_shear_stress_cell |
Shear stress perturbations. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | delta_drag_cell |
Drag force perturbations. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | delta_diffusion_cell |
Diffusion force perturbations. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | delta_gravity_cell |
Gravity force perturbations. More... | |
Local CG FEM based assemblers - variable data (previous Newton iteration - boundary) | |
Implementation is in the base class. | |
std::vector< std::vector < double > > | density_bdry_old |
Density of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | grad_density_bdry_old |
Density gradient of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | velocity_bdry_old |
Velocity of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < double > > | div_velocity_bdry_old |
Velocity divergence of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < SymmetricTensor< 2, dim > > > | grads_velocity_bdry_old |
Velocity symmetric gradient of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | mass_flux_bdry_old |
Mass flux of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < SymmetricTensor< 2, dim > > > | momentum_flux_bdry_old |
Momentum flux of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < double > > | pressure_bdry_old |
Partial pressure of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < SymmetricTensor< 2, dim > > > | shear_stress_bdry_old |
Shear stress of each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < double > > | partial_viscosity_bdry_old |
Partial dynamic viscosity, , for each species in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< std::vector < double > > | bulk_viscosity_bdry_old |
Bulk viscosity of mixture of each species, . More... | |
std::vector< std::vector < std::vector< double > > > | paramMatrix_bdry_old |
paramMatrix is used for calculating partial viscosity of mixture in the quadrature points of a boundary at a previous Newton iteration. More... | |
std::vector< FullMatrix< double > > | PInv_bdry_old |
PInv_bdry_old is used for calculating the inverse of P when using the OmegaKG model for calculating partial viscosity in the quadrature points of a boundary at a previous Newton iteration. More... | |
Local CG FEM based assemblers - variable data (current Newton iteration - boundary) | |
std::vector< std::vector < std::vector< double > > > | phi_density_bdry |
Density shape functions. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | grad_phi_density_bdry |
Density shape function gradients. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | phi_velocity_bdry |
Velocity shape functions. More... | |
std::vector< std::vector < std::vector< double > > > | div_phi_velocity_bdry |
Velocity shape function divergences. More... | |
std::vector< std::vector < std::vector< SymmetricTensor < 2, dim > > > > | grads_phi_velocity_bdry |
Velocity shape function symmetric gradients. More... | |
std::vector< std::vector < std::vector< Tensor< 1, dim > > > > | delta_mass_flux_bdry |
Mass flux perturbations. More... | |
std::vector< std::vector < std::vector< SymmetricTensor < 2, dim > > > > | delta_momentum_flux_bdry |
Momentum flux perturbations. More... | |
std::vector< std::vector < std::vector< double > > > | delta_pressure_bdry |
Partial pressure perturbations. More... | |
std::vector< std::vector < std::vector< SymmetricTensor < 2, dim > > > > | delta_shear_stress_bdry |
Shear stress perturbations. More... | |
std::vector< std::vector < std::vector< double > > > | delta_partial_viscosity_bdry |
Partial dyanmic viscosity perturbations. More... | |
std::vector< std::vector < std::vector< double > > > | delta_bulk_viscosity_bdry |
Partial bulk viscosity perturbations. More... | |
Local CG FEM based assemblers - variable data (other - boundary) | |
std::vector< std::vector < std::vector< SymmetricTensor < 2, dim > > > > | n_BY_phi_velocity_bdry_S |
Symmetrized tensor product of normal_vectors [ q ] and phi_velocity_bdry [ s ] [ q ] [ k ] in the quadrature points of a boundary. More... | |
Protected Attributes inherited from FuelCellShop::Equation::EquationBase< dim > | |
unsigned int | dofs_per_cell |
Number of degrees of freedom per cell. More... | |
unsigned int | n_q_points_cell |
Number of quadrature points per cell. More... | |
unsigned int | n_q_points_bdry |
Number of quadrature points per boundary. More... | |
DoFHandler< dim > ::active_cell_iterator | cell |
Currently active DoFHandler<dim> active cell iterator. More... | |
DoFHandler< dim > ::active_face_iterator | bdry |
Currently active DoFHandler<dim> active boundary iterator. More... | |
std::vector< double > | JxW_cell |
Jacobian of mapping by Weight in the quadrature points of a cell. More... | |
std::vector< double > | JxW_bdry |
Jacobian of mapping by Weight in the quadrature points of a boundary. More... | |
std::vector< Point< dim > > | normal_vectors |
Normal vectors in the quadrature points of a boundary. More... | |
std::vector< std::vector < Point< dim > > > | tangential_vectors |
Tangential vectors in the quadrature points of a boundary. More... | |
FuelCell::SystemManagement * | system_management |
Pointer to the external YourApplication<dim>::system_management object. More... | |
couplings_map | internal_cell_couplings |
This object contains the info on how the equations and solution variables of a derived equation class are coupled. More... | |
couplings_map | internal_flux_couplings |
This object contains the info on how the "X" and "Y" of a derived equation class are coupled (DG FEM only). More... | |
component_materialID_value_map | component_materialID_value |
This object reflects the following structure (see FuelCell::InitialAndBoundaryData namespace docs): More... | |
component_boundaryID_value_map | component_boundaryID_value |
This object reflects the following structure (see FuelCell::InitialAndBoundaryData namespace docs): More... | |
std::vector< BoundaryType > | boundary_types |
The list of boundary types of a derived equation class. More... | |
std::vector< std::vector < BoundaryType > > | multi_boundary_types |
The list of multiple boundary types of a derived equation class. More... | |
std::vector< OutputType > | output_types |
The list of output types of a derived equation class. More... | |
std::vector< std::vector < OutputType > > | multi_output_types |
The list of multiple output types of a derived equation class. More... | |
std::string | equation_name |
The name of a derived equation class. More... | |
std::string | name_base_variable |
Const std::string member storing name of the base solution variable corresponding to the equation represented by this class. More... | |
std::vector< unsigned int > | matrix_block_indices |
The system matrix block indices (a derived equation class) drawn from the global structure (a derived equation class + other active equation classes included into the computation). More... | |
std::vector< unsigned int > | residual_indices |
The residual indices (a derived equation class) drawn from the global structure (a derived equation class + other active equation classes included into the computation). More... | |
std::vector< bool > | counter |
This vector contains the collection of internal "counters" used by the derived equation classes. More... | |
EquationFlags | assemble_flags |
This vector contains a collection of internal flags to tell derived equation classes what needs to be re-computed. More... | |
boost::shared_ptr < FuelCell::ApplicationCore::ApplicationData > | data |
Data object for the application data to be passed to the equation classes. More... | |
std::string | solution_vector_name |
The name of the solution vector in FEVectors. More... | |
std::string | residual_vector_name |
The name of the residual vector name in FEVectors. More... | |
Additional Inherited Members | |
Public Attributes inherited from FuelCellShop::Equation::EquationBase< dim > | |
bool | variable_initial_data |
true , if variable initial data is prescribed on a part of the domain. More... | |
bool | variable_boundary_data |
true , if variable Dirichlet boundary conditions are prescribed on a part of the boundary. More... | |
This class implements the multi-component mass transport equations proposed by Kerkhof-Geboers for fluid transport.
These equations are
They take the following form:
ReactionSourceTerms
and SorptionSourceTerms
classes.For any species these equations can be written as
-th mass conservation equation
-th momentum conservation equation
where
-th mass flux
-th momentum flux
-th partial pressure
-th shear stress
-th drag force
-th diffusion force
the inverse of -th, -th Maxwell–Stefan isobaric diffusion coefficient
To be well-posed, these equations are equipped with appropriate boundary conditions. Below we consider several types of boundaries and several types of boundary conditions assigned to those boundaries (here denotes some portion of a boundary of a certain type).
All these boundaries and boundary conditions can be specified in the parameters file as follows:
where
where
where
These equations utilize
Channel
,GasDiffusionLayer
,MicroPorousLayer
,CatalystLayer
as layer classes and
GasMixture
,PureGas
as material classes and have the following output:
We usually do not solve these equations in the membranes of fuel cells which is fully justified under the normal conditions of operation.
The whole problem is solved by linearizing the governing equations at each Newton iteration with subsequent CG FEM discretization in space.
The implementation is performed by means of using only one weak formulation for all scalar governing equations. This approach guarantees the dim-independent programming of the equations at hand.
FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::CompressibleMultiComponentKGEquationsCoupled | ( | FuelCell::SystemManagement & | system_management, |
boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > | data = boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >() |
||
) |
Constructor.
|
virtual |
Destructor.
|
virtual |
Assemble local boundary matrix.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
virtual |
Assemble local boundary residual.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
virtual |
Assemble local cell matrix.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
virtual |
Assemble local cell residual.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
virtual |
Declare parameters.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
inline |
get_gas_species_map()
returns a map relating the different KG species equations to the desired gas species being simulated; e.g 1:oxygen, 2:water vapour, 3:nitrogen means that the first 3 KG equations (Mass Cons., Mom.
Cons. X, Mom. Cons. Y) are for species oxygen, the next three are for water vapour, and the final three are for nitrogen.
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::gas_species_map.
|
inline |
This function tells the application the boundary id to apply velocity profile equation to.
TODO: in the future it might be nice if OpenFCST could detect this based on boundary conditions
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::inlet_outlet_boundary_ID.
|
inline |
This function tells the application the height [cm] of the base of the channel with respect to origin.
|
inline |
This function tells the application the height [cm] of the roof of the channel with respect to origin.
|
inline |
This function returns inlet_outlet_velocity_max
.
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::inlet_outlet_velocity_max.
|
inline |
If get_use_parabolic_profile()
is false then user can specify an equation to use for velocity profile at boundary.
Ex. inlet-outlet
velocity
of
the
mixture
equation
= -625*y*y + 66.25*y - 0.755625
|
inline |
This function returns inlet_outlet_velocity_mixture_max
.
|
inline |
If returns the number of species to use in the parameter file.
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::n_species.
|
inline |
Outputs which component to apply variable boundary equation to (i.e.
Pressure, VelocityX, VelocityY, VelocityZ).
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::press_vel_comp_apply_to.
|
inline |
This function tells the application if the user wishes to use a parabolic profile for a boundary condition.
If false then user must specify equation to use with parameter inlet-outlet
velocity
of
the
mixture
equation
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::use_parabolic_profile.
|
virtual |
Initialize parameters.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protectedvirtual |
This function computes Local CG FEM based assemblers - constant data (boundary) and allocates the memory for.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protectedvirtual |
This function computes.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protectedvirtual |
This function computes Local CG FEM based assemblers - constant data (cell) and allocates the memory for.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protectedvirtual |
This function computes.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protected |
This function computes Computational constants and Local CG FEM based assemblers - constant data (generic).
The template parameter INFO is either typename FuelCell::ApplicationCore::DoFApplication<dim>::CellInfo or typename FuelCell::ApplicationCore::DoFApplication<dim>::FaceInfo.
|
protectedvirtual |
This function fills out internal_cell_couplings
.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protectedvirtual |
This function fills out matrix_block_indices
.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protectedvirtual |
This function fills out residual_indices
.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
virtual |
This function prints out the equations info.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
inlineprotected |
Calls GasMixture function to loop through all variations of partial viscosity and calculates the varation in bulk viscosity according to bulk_viscosity_mode.
All vectors must be sized correctly before using this function.
|
inlineprotected |
This private member function is used to calculate the variation in partial viscosity of the mixture.
Parameters:
prt | - Pointer to the layer we would like to compute the coefficient for |
partialViscosity | - vector of vectors of partial viscosity of mixture; first indice is quadrature point, second is species; [g/cm*s] |
paramMatrix | - vector of vector of vectors used in calculation of partialViscosity; first indice is quadrature point, second is species 1 and second is species 2 |
porosity | - porosity in cell, vector of quadrature points |
deltaDensity | - variation in density; same unit standard as density |
density | - vector of vectors of density; first indice is species and second is quadrature point; [g/cm^3] |
deltaPartialViscosity | - vector returned by reference contains variation of partial viscosity; first indice is species, second is quadrature point, and third is dof |
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::collision_diameter, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::density, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::dynamic_viscosity, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::molar_mass, and FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::T_mixture.
|
inlineprotected |
Inline member function used to update.
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::collision_diameter, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::dynamic_viscosity, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::maxwell_stefan_isobaric_diffusion_coefficient, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::molar_mass, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::n_species, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::P_in, and FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::T_mixture.
|
inlineprotected |
This private member function is used in make_assemblers_cell_variable_data in order to compute the inverse of the product of pressure and molecular diffusion coefficient.
This value is only a function of temperature.
Parameters:
References Units::C_UNIT2, Units::convert(), FuelCellShop::Equation::EquationBase< dim >::n_q_points_cell, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::n_species, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::P_in, and Units::UNIT2.
|
inlineprotected |
This private member function is used to calculate the partial viscosity of the mixture and bulk viscosity, from dynamicViscosity, density, molar_mass.
As well, it returns paramMatrix for possible later calculations (like when calculating the variation). Since this function is used by make_assemblers_cell_variable_data and make_assemblers_bdry_variable_data which have different quadrature points one of the members is the number of quadrature points. NOTE: no need to size vectors paramMatrix, PInv, partialViscosity, and bulkViscosity sizing is done by this function.
Parameters:
prt | - Pointer to the layer we would like to compute the coefficient for |
quadraturePoints | - number of quadrature points in cell/boundary |
porosity | - porosity in cell, vector of quadrature points |
density | - vector of vectors of density; first indice is species and second is quadrature point; [g/cm^3] |
paramMatrix | - used internally; vector of vector of vectors used in calculation of partialViscosity; first indice is quadrature point, second is species 1 and second is species 2 |
PInv | - used internally; this parameter provides the inverse of P for the calculation of OmegaKG viscosity and variation of, |
partialViscosity | - vector of vectors of partial viscosity of mixture; first indice is quadrature point, second is species; [g/cm*s] |
bulkViscosity,: | - vector of vectors of bulk viscosities; first indice is quadrature point, second is species; [g/cm*s] |
References FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::collision_diameter, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::density, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::dynamic_viscosity, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::molar_mass, FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::n_species, and FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::T_mixture.
|
inlineprotected |
Inline member function used to obtain the inverse permeability and the Forchheimer permeability for a given layer.
|
inlineprotected |
Inline function to obtain update the porosity for a given layer class.
|
inlineprotected |
Inline member function used to obtain the tortuosity for a given layer.
|
protected |
Bulk viscosity of pure gas, .
|
protected |
Bulk viscosity of mixture of each species, .
in the quadrature points of a boundary at a previous Newton iteration. NOTE: first index is the quadrature point and second index is the species.
|
protected |
Bulk viscosity of mixture of each species, .
in the quadrature points of a cell at a previous Newton iteration. NOTE: first index is the quadrature point and second index is the species.
|
protected |
Collision diameter of pure gas, .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_delta_partial_viscosities(), FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_gas_properties(), and FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_partial_viscosities().
|
protected |
|
protected |
Partial bulk viscosity perturbations.
delta_viscosity_mix_bdry
[
s
]
[
q
]
[
k
] denotes -th partial bulk viscosity perturbation computed in -th quadrature point of a cell for species .
|
protected |
Partial bulk viscosity perturbations.
delta_bulk_viscosity_cell
[
s
]
[
q
]
[
k
] denotes -th partial bulk velocity perturbation computed in -th quadrature point of a cell for species .
|
protected |
Diffusion force perturbations.
delta_diffusion_cell
[
s
]
[
q
]
[
k
] denotes -th diffusion force perturbation computed in -th quadrature point of a cell for species .
|
protected |
Drag force perturbations.
delta_drag_cell
[
s
]
[
q
]
[
k
] denotes -th drag force perturbation computed in -th quadrature point of a cell for species .
|
protected |
Gravity force perturbations.
delta_gravity_cell
[
s
]
[
q
]
[
k
] denotes -th gravity force perturbation computed in -th quadrature point of a cell for species .
|
protected |
Mass flux perturbations.
delta_mass_flux_bdry
[
s
]
[
q
]
[
k
] denotes -th mass flux perturbation computed in -th quadrature point of a boundary for species .
|
protected |
Mass flux perturbations.
delta_mass_flux_cell
[
s
]
[
q
]
[
k
] denotes -th mass flux perturbation computed in -th quadrature point of a cell for species .
|
protected |
Momentum flux perturbations.
delta_momentum_flux_bdry
[
s
]
[
q
]
[
k
] denotes -th momentum flux perturbation computed in -th quadrature point of a boundary for species .
|
protected |
Momentum flux perturbations.
delta_momentum_flux_cell
[
s
]
[
q
]
[
k
] denotes -th momentum flux perturbation computed in -th quadrature point of a cell for species .
|
protected |
Partial dyanmic viscosity perturbations.
delta_viscosity_mix_bdry
[
s
]
[
q
]
[
k
] denotes -th partial dyanmic viscosity perturbation computed in -th quadrature point of a cell for species .
|
protected |
Partial dynamic viscosity perturbations.
delta_partial_viscosity_cell
[
s
]
[
q
]
[
k
] denotes -th partial dynamic velocity perturbation computed in -th quadrature point of a cell for species .
|
protected |
Partial pressure perturbations.
delta_pressure_bdry
[
s
]
[
q
]
[
k
] denotes -th partial pressure perturbation computed in -th quadrature point of a boundary for species .
|
protected |
Partial pressure perturbations.
delta_pressure_cell
[
s
]
[
q
]
[
k
] denotes -th partial pressure perturbation computed in -th quadrature point of a cell for species .
|
protected |
Shear stress perturbations.
delta_shear_stress_bdry
[
s
]
[
q
]
[
k
] denotes -th shear stress perturbation computed in -th quadrature point of a boundary for species .
|
protected |
Shear stress perturbations.
delta_shear_stress_cell
[
s
]
[
q
]
[
k
] denotes -th shear stress perturbation computed in -th quadrature point of a cell for species .
|
protected |
|
protected |
Density of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Density of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
Diffusion force of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
Velocity shape function divergences.
div_phi_velocity_bdry
[
s
]
[
q
]
[
k
] denotes -th velocity shape function divergence computed in -th quadrature point of a boundary for species .
|
protected |
Velocity shape function divergences.
div_phi_velocity_cell
[
s
]
[
q
]
[
k
] denotes -th velocity shape function divergence computed in -th quadrature point of a cell for species .
|
protected |
Velocity divergence of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Velocity divergence of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
Drag force of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
This object indicates which form of the drag term is supposed to be chosen.
There are four options currently implemented in FCST. These options are
|
protected |
Dynamic viscosity of pure gas, .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_delta_partial_viscosities(), FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_gas_properties(), and FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_partial_viscosities().
|
protected |
Normal velocity suppression coefficient of pure gas, .
|
protected |
Map of the gas species to their respective equation species number.
Possible names for gasses include: oxygen, water, nitrogen. .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::get_gas_species_map().
|
protected |
Density gradient of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Density gradient of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
Density shape function gradients.
grad_phi_density_bdry
[
s
]
[
q
]
[
k
] denotes -th density shape function gradient computed in -th quadrature point of a boundary for species .
|
protected |
Density shape function gradients.
grad_phi_density_cell
[
s
]
[
q
]
[
k
] denotes -th density shape function gradient computed in -th quadrature point of a cell for species .
|
protected |
Velocity shape function symmetric gradients.
grads_phi_velocity_bdry
[
s
]
[
q
]
[
k
] denotes -th velocity shape function symmetric gradient computed in -th quadrature point of a boundary for species .
|
protected |
Velocity shape function symmetric gradients.
grads_phi_velocity_cell
[
s
]
[
q
]
[
k
] denotes -th velocity shape function symmetric gradient computed in -th quadrature point of a cell for species .
|
protected |
Velocity symmetric gradient of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Velocity symmetric gradient of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
Gravitational acceleration, .
|
protected |
Gravity force of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
This object indicates whether the gravity term is ON or OFF in channels.
|
protected |
This object indicates whether the gravity term is ON or OFF in porous media.
|
protected |
This object indicates whether the inertia term is ON or OFF in channels.
|
protected |
This object indicates whether the inertia term is ON or OFF in porous media.
|
protected |
|
protected |
|
protected |
|
protected |
Maximum inlet-outlet velocity of pure gas, .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::get_inlet_outlet_velocity_max().
|
protected |
|
protected |
Maximum inlet-outlet velocity of the mixture, .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::get_inlet_outlet_velocity_mixture_max().
|
protected |
Mass flux of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Mass flux of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
These constants are used in the Maxwell slip boundary condition.
|
protected |
Each entry of this structure defines a Maxwell-Stefan isobaric diffusion coefficient of gas in gas , .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_gas_properties().
|
protected |
Molar mass of pure gas, .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_delta_partial_viscosities(), FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_gas_properties(), and FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_partial_viscosities().
|
protected |
Momentum flux of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Momentum flux of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
Symmetrized tensor product of normal_vectors
[
q
] and
phi_velocity_bdry
[
s
]
[
q
]
[
k
] in the quadrature points of a boundary.
|
protected |
Number of species, .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::get_num_of_species(), FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_gas_properties(), FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_inv_diffusion_coefficients(), and FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_partial_viscosities().
|
protected |
Sometimes we implement different types of slip boundary conditions on the curved impermeable walls of the domain or the walls that are not perfectly aligned with the coordinate axes.
In this case, it might be quite problematic (but still possible) to strictly enforce the normal component of velocity be equal to 0.
The alternative approach is doing that by introducing a penalization term with into the weak formulation of the equations at hand.
This object indicates whether the penalty method is used.
|
protected |
|
protected |
paramMatrix is used for calculating partial viscosity of mixture in the quadrature points of a boundary at a previous Newton iteration.
NOTE: first index is the quadrature point and second and third indices are the species.
|
protected |
paramMatrix is used for calculating partial viscosity of mixture in the quadrature points of a cell at a previous Newton iteration.
NOTE: first index is the quadrature point and second and third indices are the species.
|
protected |
Partial dynamic viscosity, , for each species in the quadrature points of a boundary at a previous Newton iteration.
NOTE: first index is the quadrature point and second index is the species.
|
protected |
Partial viscosity of mixture of each species, , in the quadrature points of a cell at a previous Newton iteration.
NOTE: first index is the quadrature point and second index is the species.
|
protected |
Density shape functions.
phi_density_bdry
[
s
]
[
q
]
[
k
] denotes -th density shape function computed in -th quadrature point of a boundary for species .
|
protected |
Density shape functions.
phi_density_cell
[
s
]
[
q
]
[
k
] denotes -th density shape function computed in -th quadrature point of a cell for species .
|
protected |
Velocity shape functions.
phi_velocity_bdry
[
s
]
[
q
]
[
k
] denotes -th velocity shape function computed in -th quadrature point of a boundary for species .
|
protected |
Velocity shape functions.
phi_velocity_cell
[
s
]
[
q
]
[
k
] denotes -th velocity shape function computed in -th quadrature point of a cell for species .
|
protected |
PInv_bdry_old is used for calculating the inverse of P when using the OmegaKG model for calculating partial viscosity in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
PInv_old is used for calculating the inverse of P when using the OmegaKG model for calculating partial viscosity in the quadrature points of a cell at a previous Newton iteration.
|
protected |
|
protected |
Partial pressure of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Partial pressure of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
Universal gas constant, .
|
protected |
Shear stress of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Shear stress of each species in the quadrature points of a cell at a previous Newton iteration.
|
protected |
This object indicates whether the shear stress term is ON or OFF in channels.
|
protected |
This object indicates whether the shear stress term is ON or OFF in porous media.
|
protected |
Temperature of species mixture, .
Referenced by FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_delta_partial_viscosities(), FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_gas_properties(), and FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >::update_partial_viscosities().
|
protected |
Navier slip coefficient of pure gas, .
|
protected |
Unit tensor, .
|
protected |
|
protected |
Velocity extractors.
|
protected |
Velocity of each species in the quadrature points of a boundary at a previous Newton iteration.
|
protected |
Velocity of each species in the quadrature points of a cell at a previous Newton iteration.