Cantera  3.2.0a2
Loading...
Searching...
No Matches
Transport Class Reference

Base class for transport property managers. More...

#include <Transport.h>

Inheritance diagram for Transport:
[legend]

Detailed Description

Base class for transport property managers.

All classes that compute transport properties for a single phase derive from this class. Class Transport is meant to be used as a base class only. It is possible to instantiate it, but its methods throw exceptions if called.

Relationship of the Transport class to the ThermoPhase Class

This section describes how calculations are carried out within the Transport class. The Transport class and derived classes of the the Transport class necessarily use the ThermoPhase class to obtain the list of species and the thermodynamic state of the phase.

No state information is stored within Transport classes. Queries to the underlying ThermoPhase object must be made to obtain the state of the system.

An exception to this however is the state information concerning the the gradients of variables. This information is not stored within the ThermoPhase objects. It may be collected within the Transport objects. In fact, the meaning of const operations within the Transport class refers to calculations which do not change the state of the system nor the state of the first order gradients of the system.

When a const operation is evoked within the Transport class, it is also implicitly assumed that the underlying state within the ThermoPhase object has not changed its values.

Todo:
Provide a general mechanism to store the gradients of state variables within the system.

Definition at line 71 of file Transport.h.

Public Member Functions

 Transport ()=default
 Constructor.
 
 Transport (const Transport &)=delete
 
Transportoperator= (const Transport &)=delete
 
virtual string transportModel () const
 Identifies the model represented by this Transport object.
 
ThermoPhasethermo ()
 Phase object.
 
void checkSpeciesIndex (size_t k) const
 Check that the specified species index is in range.
 
void checkSpeciesArraySize (size_t kk) const
 Check that an array size is at least m_nsp.
 
virtual void getSpeciesFluxes (size_t ndim, const double *const grad_T, size_t ldx, const double *const grad_X, size_t ldf, double *const fluxes)
 Get the species diffusive mass fluxes [kg/m²/s] with respect to the specified solution averaged velocity, given the mole fraction and temperature gradients.
 
virtual void getMolarFluxes (const double *const state1, const double *const state2, const double delta, double *const cfluxes)
 Get the molar fluxes [kmol/m²/s], given the thermodynamic state at two nearby points.
 
virtual void getMassFluxes (const double *state1, const double *state2, double delta, double *mfluxes)
 Get the mass fluxes [kg/m²/s], given the thermodynamic state at two nearby points.
 
virtual void getThermalDiffCoeffs (double *const dt)
 Return a vector of thermal diffusion coefficients [kg/m/s].
 
virtual void getBinaryDiffCoeffs (const size_t ld, double *const d)
 Returns the matrix of binary diffusion coefficients [m²/s].
 
virtual void getMultiDiffCoeffs (const size_t ld, double *const d)
 Return the multicomponent diffusion coefficients [m²/s].
 
virtual void getMixDiffCoeffs (double *const d)
 Return a vector of mixture averaged diffusion coefficients [m²/s].
 
virtual void getMixDiffCoeffsMole (double *const d)
 Returns a vector of mixture averaged diffusion coefficients [m²/s].
 
virtual void getMixDiffCoeffsMass (double *const d)
 Returns a vector of mixture averaged diffusion coefficients [m²/s].
 
virtual void getViscosityPolynomial (size_t i, double *coeffs) const
 Return the polynomial fits to the viscosity of species i.
 
virtual void getConductivityPolynomial (size_t i, double *coeffs) const
 Return the temperature fits of the heat conductivity of species i.
 
virtual void getBinDiffusivityPolynomial (size_t i, size_t j, double *coeffs) const
 Return the polynomial fits to the binary diffusivity of species pair (i, j)
 
virtual void getCollisionIntegralPolynomial (size_t i, size_t j, double *astar_coeffs, double *bstar_coeffs, double *cstar_coeffs) const
 Return the polynomial fits to the collision integral of species pair (i, j)
 
virtual void setViscosityPolynomial (size_t i, double *coeffs)
 Modify the polynomial fits to the viscosity of species i
 
virtual void setConductivityPolynomial (size_t i, double *coeffs)
 Modify the temperature fits of the heat conductivity of species i
 
virtual void setBinDiffusivityPolynomial (size_t i, size_t j, double *coeffs)
 Modify the polynomial fits to the binary diffusivity of species pair (i, j)
 
virtual void setCollisionIntegralPolynomial (size_t i, size_t j, double *astar_coeffs, double *bstar_coeffs, double *cstar_coeffs, bool flag)
 Modify the polynomial fits to the collision integral of species pair (i, j)
 
AnyMap parameters () const
 Return the parameters for a phase definition which are needed to reconstruct an identical object using the newTransport() function.
 
AnyMap fittingErrors () const
 Get error metrics about any functional fits calculated for pure species transport properties.
 
virtual void invalidateCache ()
 Invalidate any cached values which are normally updated only when a change in state is detected.
 
Transport Properties
virtual double viscosity ()
 Get the dynamic viscosity [Pa·s].
 
virtual void getSpeciesViscosities (double *const visc)
 Get the pure species viscosities [Pa·s].
 
virtual double bulkViscosity ()
 The bulk viscosity [Pa·s].
 
virtual double thermalConductivity ()
 Get the mixture thermal conductivity [W/m/K].
 
virtual double electricalConductivity ()
 Get the electrical conductivity [siemens/m].
 
virtual void getMobilities (double *const mobil_e)
 Get the electrical mobilities [m²/V/s].
 
Transport manager construction

These methods are used during construction.

virtual void init (ThermoPhase *thermo, int mode=0)
 Initialize a transport manager.
 
virtual bool CKMode () const
 Boolean indicating the form of the transport properties polynomial fits.
 

Protected Attributes

ThermoPhasem_thermo
 pointer to the object representing the phase
 
size_t m_nsp = 0
 Number of species in the phase.
 
AnyMap m_fittingErrors
 Maximum errors associated with fitting pure species transport properties.
 

Constructor & Destructor Documentation

◆ Transport()

Transport ( )
default

Constructor.

New transport managers should be created using TransportFactory, not by calling the constructor directly.

See also
TransportFactory

◆ ~Transport()

virtual ~Transport ( )
inlinevirtual

Definition at line 83 of file Transport.h.

Member Function Documentation

◆ transportModel()

virtual string transportModel ( ) const
inlinevirtual

Identifies the model represented by this Transport object.

Each derived class should override this method to return a meaningful identifier.

Since
New in Cantera 3.0. The name returned by this method corresponds to the canonical name used in the YAML input format.

Reimplemented in DustyGasTransport, HighPressureGasTransport, ChungHighPressureGasTransport, IonGasTransport, MixTransport, MultiTransport, UnityLewisTransport, and WaterTransport.

Definition at line 93 of file Transport.h.

◆ thermo()

ThermoPhase & thermo ( )
inline

Phase object.

Every transport manager is designed to compute properties for a specific phase of a mixture, which might be a liquid solution, a gas mixture, a surface, etc. This method returns a reference to the object representing the phase itself.

Definition at line 103 of file Transport.h.

◆ checkSpeciesIndex()

void checkSpeciesIndex ( size_t  k) const

Check that the specified species index is in range.

Throws an exception if k is greater than m_nsp.

Definition at line 16 of file Transport.cpp.

◆ checkSpeciesArraySize()

void checkSpeciesArraySize ( size_t  kk) const

Check that an array size is at least m_nsp.

Throws an exception if kk is less than m_nsp. Used before calls which take an array pointer.

Definition at line 23 of file Transport.cpp.

◆ viscosity()

virtual double viscosity ( )
inlinevirtual

Get the dynamic viscosity [Pa·s].

Reimplemented in GasTransport, HighPressureGasTransport, ChungHighPressureGasTransport, IonGasTransport, and WaterTransport.

Definition at line 120 of file Transport.h.

◆ getSpeciesViscosities()

virtual void getSpeciesViscosities ( double *const  visc)
inlinevirtual

Get the pure species viscosities [Pa·s].

Parameters
viscVector of viscosities; length is the number of species

Reimplemented in GasTransport.

Definition at line 129 of file Transport.h.

◆ bulkViscosity()

virtual double bulkViscosity ( )
inlinevirtual

The bulk viscosity [Pa·s].

The bulk viscosity is only non-zero in rare cases. Most transport managers either overload this method to return zero, or do not implement it, in which case an exception is thrown if called.

Reimplemented in WaterTransport.

Definition at line 140 of file Transport.h.

◆ thermalConductivity()

virtual double thermalConductivity ( )
inlinevirtual

Get the mixture thermal conductivity [W/m/K].

Reimplemented in HighPressureGasTransport, ChungHighPressureGasTransport, IonGasTransport, MixTransport, MultiTransport, and WaterTransport.

Definition at line 146 of file Transport.h.

◆ electricalConductivity()

virtual double electricalConductivity ( )
inlinevirtual

Get the electrical conductivity [siemens/m].

Reimplemented in IonGasTransport.

Definition at line 152 of file Transport.h.

◆ getMobilities()

virtual void getMobilities ( double *const  mobil_e)
inlinevirtual

Get the electrical mobilities [m²/V/s].

This function returns the mobilities. In some formulations this is equal to the normal mobility multiplied by Faraday's constant.

Frequently, but not always, the mobility is calculated from the diffusion coefficient using the Einstein relation

\[ \mu^e_k = \frac{F D_k}{R T} \]

Parameters
mobil_eReturns the mobilities of the species in array mobil_e. The array must be dimensioned at least as large as the number of species.

Reimplemented in IonGasTransport, and MixTransport.

Definition at line 173 of file Transport.h.

◆ getSpeciesFluxes()

virtual void getSpeciesFluxes ( size_t  ndim,
const double *const  grad_T,
size_t  ldx,
const double *const  grad_X,
size_t  ldf,
double *const  fluxes 
)
inlinevirtual

Get the species diffusive mass fluxes [kg/m²/s] with respect to the specified solution averaged velocity, given the mole fraction and temperature gradients.

Usually the specified solution average velocity is the mass averaged velocity. This is changed in some subclasses, however.

Parameters
ndimNumber of dimensions in the flux expressions
[in]grad_TGradient of the temperature (length ndim)
ldxLeading dimension of the grad_X array (usually equal to the number of species)
[in]grad_XGradients of the mole fractions; flattened matrix such that \( dX_k/dx_n = \tt{ grad\_X[n*ldx+k]} \) is the gradient of species k in dimension n. Length is ldx * ndim.
ldfLeading dimension of the fluxes array (usually equal to the number of species)
[out]fluxesThe diffusive mass fluxes; flattened matrix such that \( j_{kn} = \tt{ fluxes[n*ldf+k]} \) is the flux of species k in dimension n. Length is ldf * ndim.

Reimplemented in MixTransport, and MultiTransport.

Definition at line 199 of file Transport.h.

◆ getMolarFluxes()

virtual void getMolarFluxes ( const double *const  state1,
const double *const  state2,
const double  delta,
double *const  cfluxes 
)
inlinevirtual

Get the molar fluxes [kmol/m²/s], given the thermodynamic state at two nearby points.

Parameters
[in]state1Array of temperature, density, and mass fractions for state 1.
[in]state2Array of temperature, density, and mass fractions for state 2.
[in]deltaDistance [m] from state 1 to state 2.
[out]cfluxesArray containing the diffusive molar fluxes of species from state1 to state2; Length is number of species.

Reimplemented in DustyGasTransport, and MultiTransport.

Definition at line 217 of file Transport.h.

◆ getMassFluxes()

virtual void getMassFluxes ( const double *  state1,
const double *  state2,
double  delta,
double *  mfluxes 
)
inlinevirtual

Get the mass fluxes [kg/m²/s], given the thermodynamic state at two nearby points.

Parameters
[in]state1Array of temperature, density, and mass fractions for state 1.
[in]state2Array of temperature, density, and mass fractions for state 2.
[in]deltaDistance [m] from state 1 to state 2.
[out]mfluxesArray containing the diffusive mass fluxes of species from state1 to state2; length is number of species.

Reimplemented in MultiTransport.

Definition at line 235 of file Transport.h.

◆ getThermalDiffCoeffs()

virtual void getThermalDiffCoeffs ( double *const  dt)
inlinevirtual

Return a vector of thermal diffusion coefficients [kg/m/s].

The thermal diffusion coefficient \( D^T_k \) is defined so that the diffusive mass flux of species k induced by the local temperature gradient is given by:

\[ \mathbf{j}_k = -D^T_k \frac{\nabla T}{T}. \]

The thermal diffusion coefficient can be either positive or negative.

Parameters
dtOn return, dt will contain the species thermal diffusion coefficients. Dimension dt at least as large as the number of species.

Reimplemented in MixTransport, MultiTransport, and UnityLewisTransport.

Definition at line 257 of file Transport.h.

◆ getBinaryDiffCoeffs()

virtual void getBinaryDiffCoeffs ( const size_t  ld,
double *const  d 
)
inlinevirtual

Returns the matrix of binary diffusion coefficients [m²/s].

Parameters
[in]ldLeading dimension of the flattened array d used to store the diffusion coefficient matrix; usually equal to the number of species.
[out]dDiffusion coefficient matrix stored in column-major (Fortran) order, such that \( \mathcal{D}_{ij} = \tt{d[ld*j + i]} \); must be at least ld times the number of species in length.
See also
GasTransport::fitDiffCoeffs()

Reimplemented in GasTransport, and HighPressureGasTransportBase.

Definition at line 272 of file Transport.h.

◆ getMultiDiffCoeffs()

virtual void getMultiDiffCoeffs ( const size_t  ld,
double *const  d 
)
inlinevirtual

Return the multicomponent diffusion coefficients [m²/s].

If the transport manager implements a multicomponent diffusion model, then this method returns the array of multicomponent diffusion coefficients. Otherwise it throws an exception.

Parameters
[in]ldLeading dimension of the flattened array d used to store the diffusion coefficient matrix; usually equal to the number of species.
[out]dDiffusion coefficient matrix stored in column-major (Fortran) order, such that \( D_{ij} = \tt{d[ld*j + i]} \) is the diffusion coefficient for species i due to concentration gradients in species j; must be at least ld times the number of species in length.

Reimplemented in DustyGasTransport, and MultiTransport.

Definition at line 292 of file Transport.h.

◆ getMixDiffCoeffs()

virtual void getMixDiffCoeffs ( double *const  d)
inlinevirtual

Return a vector of mixture averaged diffusion coefficients [m²/s].

Mixture-averaged diffusion coefficients [m^2/s]. If the transport manager implements a mixture-averaged diffusion model, then this method returns the array of mixture-averaged diffusion coefficients. Otherwise it throws an exception.

Parameters
dReturn vector of mixture averaged diffusion coefficients; length is the number of species.

Reimplemented in GasTransport, HighPressureGasTransportBase, IonGasTransport, and UnityLewisTransport.

Definition at line 307 of file Transport.h.

◆ getMixDiffCoeffsMole()

virtual void getMixDiffCoeffsMole ( double *const  d)
inlinevirtual

Returns a vector of mixture averaged diffusion coefficients [m²/s].

Reimplemented in GasTransport, HighPressureGasTransportBase, and UnityLewisTransport.

Definition at line 313 of file Transport.h.

◆ getMixDiffCoeffsMass()

virtual void getMixDiffCoeffsMass ( double *const  d)
inlinevirtual

Returns a vector of mixture averaged diffusion coefficients [m²/s].

Reimplemented in GasTransport, HighPressureGasTransportBase, and UnityLewisTransport.

Definition at line 319 of file Transport.h.

◆ getViscosityPolynomial()

virtual void getViscosityPolynomial ( size_t  i,
double *  coeffs 
) const
inlinevirtual

Return the polynomial fits to the viscosity of species i.

Reimplemented in GasTransport.

Definition at line 325 of file Transport.h.

◆ getConductivityPolynomial()

virtual void getConductivityPolynomial ( size_t  i,
double *  coeffs 
) const
inlinevirtual

Return the temperature fits of the heat conductivity of species i.

Reimplemented in GasTransport.

Definition at line 331 of file Transport.h.

◆ getBinDiffusivityPolynomial()

virtual void getBinDiffusivityPolynomial ( size_t  i,
size_t  j,
double *  coeffs 
) const
inlinevirtual

Return the polynomial fits to the binary diffusivity of species pair (i, j)

Reimplemented in GasTransport.

Definition at line 337 of file Transport.h.

◆ getCollisionIntegralPolynomial()

virtual void getCollisionIntegralPolynomial ( size_t  i,
size_t  j,
double *  astar_coeffs,
double *  bstar_coeffs,
double *  cstar_coeffs 
) const
inlinevirtual

Return the polynomial fits to the collision integral of species pair (i, j)

Reimplemented in GasTransport.

Definition at line 343 of file Transport.h.

◆ setViscosityPolynomial()

virtual void setViscosityPolynomial ( size_t  i,
double *  coeffs 
)
inlinevirtual

Modify the polynomial fits to the viscosity of species i

Reimplemented in GasTransport.

Definition at line 352 of file Transport.h.

◆ setConductivityPolynomial()

virtual void setConductivityPolynomial ( size_t  i,
double *  coeffs 
)
inlinevirtual

Modify the temperature fits of the heat conductivity of species i

Reimplemented in GasTransport.

Definition at line 358 of file Transport.h.

◆ setBinDiffusivityPolynomial()

virtual void setBinDiffusivityPolynomial ( size_t  i,
size_t  j,
double *  coeffs 
)
inlinevirtual

Modify the polynomial fits to the binary diffusivity of species pair (i, j)

Reimplemented in GasTransport.

Definition at line 364 of file Transport.h.

◆ setCollisionIntegralPolynomial()

virtual void setCollisionIntegralPolynomial ( size_t  i,
size_t  j,
double *  astar_coeffs,
double *  bstar_coeffs,
double *  cstar_coeffs,
bool  flag 
)
inlinevirtual

Modify the polynomial fits to the collision integral of species pair (i, j)

Reimplemented in GasTransport.

Definition at line 370 of file Transport.h.

◆ parameters()

AnyMap parameters ( ) const

Return the parameters for a phase definition which are needed to reconstruct an identical object using the newTransport() function.

This excludes the individual species transport properties, which are handled separately.

Definition at line 30 of file Transport.cpp.

◆ fittingErrors()

AnyMap fittingErrors ( ) const
inline

Get error metrics about any functional fits calculated for pure species transport properties.

See GasTransport::fitDiffCoeffs and GasTransport::fitProperties.

Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.
Since
New in Cantera 3.1.

Definition at line 392 of file Transport.h.

◆ init()

virtual void init ( ThermoPhase thermo,
int  mode = 0 
)
inlinevirtual

Initialize a transport manager.

This routine sets up a transport manager. It calculates the collision integrals and populates species-dependent data structures.

Parameters
thermoPointer to the ThermoPhase object
modeChemkin compatible mode or not. This alters the specification of the collision integrals. defaults to no.

Reimplemented in IonGasTransport, GasTransport, HighPressureGasTransport, ChungHighPressureGasTransport, MixTransport, MultiTransport, and WaterTransport.

Definition at line 408 of file Transport.h.

◆ CKMode()

virtual bool CKMode ( ) const
inlinevirtual

Boolean indicating the form of the transport properties polynomial fits.

Returns true if the Chemkin form is used.

Reimplemented in GasTransport.

Definition at line 412 of file Transport.h.

◆ invalidateCache()

virtual void invalidateCache ( )
inlinevirtual

Invalidate any cached values which are normally updated only when a change in state is detected.

Since
New in Cantera 3.1.

Reimplemented in GasTransport, and MultiTransport.

Definition at line 422 of file Transport.h.

Member Data Documentation

◆ m_thermo

ThermoPhase* m_thermo
protected

pointer to the object representing the phase

Definition at line 426 of file Transport.h.

◆ m_nsp

size_t m_nsp = 0
protected

Number of species in the phase.

Definition at line 429 of file Transport.h.

◆ m_fittingErrors

AnyMap m_fittingErrors
protected

Maximum errors associated with fitting pure species transport properties.

Definition at line 432 of file Transport.h.


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