Cantera  2.1.2
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
TransportFactory Class Reference

Factory class for creating new instances of classes derived from Transport. More...

#include <TransportFactory.h>

Inheritance diagram for TransportFactory:
[legend]
Collaboration diagram for TransportFactory:
[legend]

Public Member Functions

virtual void deleteFactory ()
 Deletes the statically allocated factory instance. More...
 
virtual LTPspeciesnewLTP (const XML_Node &trNode, const std::string &name, TransportPropertyType tp_ind, thermo_t *thermo)
 Make one of several transport models, and return a base class pointer to it. More...
 
virtual LiquidTranInteractionnewLTI (const XML_Node &trNode, TransportPropertyType tp_ind, LiquidTransportParams &trParam)
 Factory function for the construction of new LiquidTranInteraction objects, which are transport models. More...
 
virtual TransportnewTransport (const std::string &model, thermo_t *thermo, int log_level=0, int ndim=1)
 Build a new transport manager using a transport manager that may not be the same as in the phase description and return a base class pointer to it. More...
 
virtual TransportnewTransport (thermo_t *thermo, int log_level=0)
 Build a new transport manager using the default transport manager in the phase description and return a base class pointer to it. More...
 
virtual void initTransport (Transport *tr, thermo_t *thermo, int mode=0, int log_level=0)
 Initialize an existing transport manager. More...
 
virtual void initLiquidTransport (Transport *tr, thermo_t *thermo, int log_level=0)
 Initialize an existing transport manager for liquid phase. More...
 
- Public Member Functions inherited from FactoryBase
virtual ~FactoryBase ()
 destructor More...
 

Static Public Member Functions

static TransportFactoryfactory ()
 Return a pointer to a TransportFactory instance. More...
 
static std::string modelName (int model)
 Get the name of the transport model corresponding to the specified constant. More...
 
- Static Public Member Functions inherited from FactoryBase
static void deleteFactories ()
 static function that deletes all factories in the internal registry maintained in a static variable More...
 

Private Member Functions

virtual void initSolidTransport (Transport *tr, thermo_t *thermo, int log_level=0)
 Initialize an existing transport manager for solid phase. More...
 
 TransportFactory ()
 The constructor is private; use static method factory() to get a pointer to a factory instance. More...
 
void getTransportData (const std::vector< const XML_Node * > &xspecies, XML_Node &log, const std::vector< std::string > &names, GasTransportParams &tr)
 Read the transport database. More...
 
void getLiquidSpeciesTransportData (const std::vector< const XML_Node * > &db, XML_Node &log, const std::vector< std::string > &names, LiquidTransportParams &tr)
 Read transport property data from a file for a list of species that comprise the phase. More...
 
void getLiquidInteractionsTransportData (const XML_Node &phaseTran_db, XML_Node &log, const std::vector< std::string > &names, LiquidTransportParams &tr)
 Read transport property data from a file for interactions between species. More...
 
void getSolidTransportData (const XML_Node &transportNode, XML_Node &log, const std::string phaseName, SolidTransportData &tr)
 Read transport property data from a file for a solid phase. More...
 
void fitProperties (GasTransportParams &tr, MMCollisionInt &integrals, std::ostream &logfile)
 Generate polynomial fits to the viscosity, conductivity, and the binary diffusion coefficients. More...
 
void fitCollisionIntegrals (std::ostream &logfile, GasTransportParams &tr, MMCollisionInt &integrals)
 Generate polynomial fits to collision integrals. More...
 
void setupMM (std::ostream &flog, const std::vector< const XML_Node * > &transport_database, thermo_t *thermo, int mode, int log_level, GasTransportParams &tr)
 Prepare to build a new kinetic-theory-based transport manager for low-density gases. More...
 
void setupLiquidTransport (std::ostream &flog, thermo_t *thermo, int log_level, LiquidTransportParams &trParam)
 Prepare to build a new transport manager for liquids assuming that viscosity transport data is provided in Arrhenius form. More...
 
void setupSolidTransport (std::ostream &flog, thermo_t *thermo, int log_level, SolidTransportData &trParam)
 Prepare to build a new transport manager for solids. More...
 
void getBinDiffCorrection (doublereal t, const GasTransportParams &tr, MMCollisionInt &integrals, size_t k, size_t j, doublereal xk, doublereal xj, doublereal &fkj, doublereal &fjk)
 Second-order correction to the binary diffusion coefficients. More...
 
void makePolarCorrections (size_t i, size_t j, const GasTransportParams &tr, doublereal &f_eps, doublereal &f_sigma)
 Corrections for polar-nonpolar binary diffusion coefficients. More...
 

Private Attributes

bool m_verbose
 Boolean indicating whether to turn on verbose printing. More...
 
std::map< std::string, int > m_models
 Mapping between between the string name for a transport model and the integer name. More...
 
std::map< int, std::string > m_modelNames
 Inverse mapping of transport models, from integer constant to string. More...
 
std::map< std::string,
TransportPropertyType
m_tranPropMap
 Mapping between between the string name for a transport property and the integer name. More...
 
std::map< std::string,
LTPTemperatureDependenceType
m_LTRmodelMap
 Mapping between between the string name for a species-specific transport property model and the integer name. More...
 
std::map< std::string,
LiquidTranMixingModel
m_LTImodelMap
 Mapping between between the string name for a liquid mixture transport property model and the integer name. More...
 

Static Private Attributes

static TransportFactorys_factory = 0
 Static instance of the factor -> This is the only instance of this object allowed. More...
 
static mutex_t transport_mutex
 Static instance of the mutex used to ensure the proper reading of the transport database. More...
 

Additional Inherited Members

- Protected Member Functions inherited from FactoryBase
 FactoryBase ()
 Constructor. More...
 

Detailed Description

Factory class for creating new instances of classes derived from Transport.

Creates 'transport managers', which are classes derived from class Transport that provide transport properties. TransportFactory handles all initialization required, including evaluation of collision integrals and generating polynomial fits. Transport managers can also be created in other ways.

Definition at line 38 of file TransportFactory.h.

Constructor & Destructor Documentation

TransportFactory ( )
private

The constructor is private; use static method factory() to get a pointer to a factory instance.

The default constructor for this class sets up m_models[], a mapping between the string name for a transport model and the integer name.

Definition at line 166 of file TransportFactory.cpp.

References TransportFactory::m_LTImodelMap, TransportFactory::m_LTRmodelMap, TransportFactory::m_modelNames, TransportFactory::m_models, and TransportFactory::m_tranPropMap.

Referenced by TransportFactory::factory().

Member Function Documentation

static TransportFactory* factory ( )
inlinestatic

Return a pointer to a TransportFactory instance.

TransportFactory is implemented as a 'singleton', which means that at most one instance may be created. The constructor is private. When a TransportFactory instance is required, call static method factory() to return a pointer to the TransportFactory instance.

Definition at line 53 of file TransportFactory.h.

References TransportFactory::s_factory, TransportFactory::transport_mutex, and TransportFactory::TransportFactory().

Referenced by TransportFactory::modelName(), Cantera::newDefaultTransportMgr(), and Cantera::newTransportMgr().

void deleteFactory ( )
virtual

Deletes the statically allocated factory instance.

Implements FactoryBase.

Definition at line 216 of file TransportFactory.cpp.

References TransportFactory::s_factory, and TransportFactory::transport_mutex.

std::string modelName ( int  model)
static

Get the name of the transport model corresponding to the specified constant.

Parameters
modelInteger representing the model name

Definition at line 225 of file TransportFactory.cpp.

References TransportFactory::factory(), and TransportFactory::m_modelNames.

LTPspecies * newLTP ( const XML_Node trNode,
const std::string &  name,
TransportPropertyType  tp_ind,
thermo_t thermo 
)
virtual

Make one of several transport models, and return a base class pointer to it.

This method operates at the level of a single transport property as a function of temperature and possibly composition. It's a factory for LTPspecies classes.

Parameters
trNodeXML node
namereference to the name
tp_indTransportPropertyType class
thermoPointer to the ThermoPhase class

Definition at line 236 of file TransportFactory.cpp.

References Cantera::lowercase(), and TransportFactory::m_LTRmodelMap.

Referenced by TransportFactory::getLiquidSpeciesTransportData(), and TransportFactory::getSolidTransportData().

LiquidTranInteraction * newLTI ( const XML_Node trNode,
TransportPropertyType  tp_ind,
LiquidTransportParams trParam 
)
virtual

Factory function for the construction of new LiquidTranInteraction objects, which are transport models.

This method operates at the level of a single mixture transport property. Individual species transport properties are addressed by the LTPspecies returned by newLTP.

Parameters
trNodeXML_Node containing the information for the interaction
tp_indTransportPropertyType object
trParamreference to the LiquidTransportParams object

Definition at line 261 of file TransportFactory.cpp.

References LiquidTranInteraction::init(), TransportFactory::m_LTImodelMap, and TransportParams::thermo.

Referenced by TransportFactory::getLiquidInteractionsTransportData().

Transport * newTransport ( const std::string &  model,
thermo_t thermo,
int  log_level = 0,
int  ndim = 1 
)
virtual

Build a new transport manager using a transport manager that may not be the same as in the phase description and return a base class pointer to it.

Parameters
modelString name for the transport manager
thermoThermoPhase object
log_levellog level
ndimNumber of dimensions for fluxes

Definition at line 314 of file TransportFactory.cpp.

References DustyGasTransport::initialize(), TransportFactory::initLiquidTransport(), TransportFactory::initSolidTransport(), TransportFactory::initTransport(), TransportFactory::m_models, Phase::restoreState(), Phase::saveState(), and Transport::setThermo().

Referenced by Cantera::newDefaultTransportMgr(), TransportFactory::newTransport(), and Cantera::newTransportMgr().

Transport * newTransport ( thermo_t thermo,
int  log_level = 0 
)
virtual

Build a new transport manager using the default transport manager in the phase description and return a base class pointer to it.

Parameters
thermoThermoPhase object
log_levellog level

Definition at line 388 of file TransportFactory.cpp.

References XML_Node::attrib(), XML_Node::child(), XML_Node::hasChild(), TransportFactory::newTransport(), and Phase::xml().

void initTransport ( Transport tr,
thermo_t thermo,
int  mode = 0,
int  log_level = 0 
)
virtual

Initialize an existing transport manager.

This routine sets up an existing gas-phase transport manager. It calculates the collision integrals and calls the initGas() function to populate the species-dependent data structure.

Parameters
trPointer to the Transport manager
thermoPointer to the ThermoPhase object
modeChemkin compatible mode or not. This alters the specification of the collision integrals. defaults to no.
log_levelDefaults to zero, no logging

In DEBUG_MODE, this routine will create the file transport_log.xml and write informative information to it.

Definition at line 614 of file TransportFactory.cpp.

References Transport::initGas(), TransportFactory::m_verbose, TransportFactory::setupMM(), ThermoPhase::speciesData(), TransportFactory::transport_mutex, and TransportParams::xml.

Referenced by TransportFactory::newTransport().

void initLiquidTransport ( Transport tr,
thermo_t thermo,
int  log_level = 0 
)
virtual

Initialize an existing transport manager for liquid phase.

This routine sets up an existing liquid-phase transport manager. It is similar to initTransport except that it uses the LiquidTransportParams class and calls setupLiquidTransport().

Parameters
trPointer to the Transport manager
thermoPointer to the ThermoPhase object
log_levelDefaults to zero, no logging

In DEBUG_MODE, this routine will create the file transport_log.xml and write informative information to it.

Definition at line 649 of file TransportFactory.cpp.

References Transport::initLiquid(), TransportFactory::m_verbose, TransportFactory::setupLiquidTransport(), and TransportParams::xml.

Referenced by TransportFactory::newTransport().

void initSolidTransport ( Transport tr,
thermo_t thermo,
int  log_level = 0 
)
privatevirtual

Initialize an existing transport manager for solid phase.

This routine sets up an existing solid-phase transport manager. It is similar to initTransport except that it uses the SolidTransportData class and calls setupSolidTransport().

Parameters
trPointer to the Transport manager
thermoPointer to the ThermoPhase object
log_levelDefaults to zero, no logging

In DEBUG_MODE, this routine will create the file transport_log.xml and write informative information to it.

Definition at line 677 of file TransportFactory.cpp.

References Transport::initSolid(), TransportFactory::m_verbose, TransportFactory::setupSolidTransport(), and TransportParams::xml.

Referenced by TransportFactory::newTransport().

void getTransportData ( const std::vector< const XML_Node * > &  xspecies,
XML_Node log,
const std::vector< std::string > &  names,
GasTransportParams tr 
)
private

Read the transport database.

Read transport property data from a file for a list of species. Given the name of a file containing transport property parameters and a list of species names, this method returns an instance of TransportParams containing the transport data for these species read from the file.

Parameters
xspeciesVector of pointers to species XML_Node databases.
logreference to an XML_Node that will contain the log (unused)
namesvector of species names that must be filled in with valid transport parameters
trOutput object containing the transport parameters for the species listed in names (in the order of their listing in names).

Definition at line 779 of file TransportFactory.cpp.

References GasTransportParams::alpha, Cantera::Boltzmann, XML_Node::child(), GasTransportParams::crot, GasTransportParams::dipole, GasTransportParams::eps, ctml::getByTitle(), ctml::getFloat(), GasTransportParams::polar, GasTransportParams::sigma, Cantera::SqrtTen, XML_Node::value(), and GasTransportParams::zrot.

Referenced by TransportFactory::setupMM().

void getLiquidSpeciesTransportData ( const std::vector< const XML_Node * > &  db,
XML_Node log,
const std::vector< std::string > &  names,
LiquidTransportParams tr 
)
private

Read transport property data from a file for a list of species that comprise the phase.

Given a vector of pointers to species XML data bases and a list of species names, this method constructs the LiquidTransport Params object containing the transport data for these species.

It is an error to not find a "transport" XML element within each of the species XML elements listed in the names vector.

Parameters
dbReference to a vector of XML_Node pointers containing the species XML nodes.
logReference to an XML log file. (currently unused)
namesVector of names of species. On output, tr will contain transport data for each of of these names in the order determined by this vector.
trReference to the LiquidTransportParams object that will contain the results.

Definition at line 863 of file TransportFactory.cpp.

References XML_Node::child(), LiquidTransportData::electCond, XML_Node::hasChild(), LiquidTransportData::hydroRadius, LiquidTransportData::ionConductivity, LiquidTransportParams::LTData, TransportFactory::m_tranPropMap, LiquidTransportData::mobilityRatio, XML_Node::name(), XML_Node::nChildren(), TransportFactory::newLTP(), TransportParams::nsp_, CanteraError::save(), LiquidTransportData::selfDiffusion, LiquidTransportData::speciesDiffusivity, Phase::speciesIndex(), LiquidTransportData::speciesName, LiquidTransportData::thermalCond, TransportParams::thermo, and LiquidTransportData::viscosity.

Referenced by TransportFactory::setupLiquidTransport().

void getLiquidInteractionsTransportData ( const XML_Node phaseTran_db,
XML_Node log,
const std::vector< std::string > &  names,
LiquidTransportParams tr 
)
private

Read transport property data from a file for interactions between species.

Given the XML_Node database for transport interactions defined within the current phase and a list of species names within the phase, this method returns an instance of TransportParams containing the transport data for these species read from the file.

This routine reads interaction parameters between species within the phase.

Parameters
phaseTran_dbReference to the transport XML field for the phase
logReference to an XML log file. (currently unused)
namesVector of names of species. On output, tr will contain transport data for each of of these names in the order determined by this vector.
trReference to the LiquidTransportParams object that will contain the results.

Definition at line 1003 of file TransportFactory.cpp.

References XML_Node::attrib(), XML_Node::child(), LiquidTransportParams::electCond, XML_Node::hasChild(), LiquidTransportParams::hydroRadius, LiquidTransportParams::ionConductivity, TransportFactory::m_tranPropMap, LiquidTransportParams::mobilityRatio, XML_Node::name(), XML_Node::nChildren(), TransportFactory::newLTI(), TransportParams::nsp_, LiquidTransportParams::selfDiffusion, LiquidTransportParams::speciesDiffusivity, Phase::speciesIndex(), LiquidTransportParams::thermalCond, TransportParams::thermo, Cantera::VB_MASSAVG, Cantera::VB_MOLEAVG, TransportParams::velocityBasis_, LiquidTransportParams::viscosity, and CanteraError::what().

Referenced by TransportFactory::setupLiquidTransport().

void getSolidTransportData ( const XML_Node transportNode,
XML_Node log,
const std::string  phaseName,
SolidTransportData tr 
)
private

Read transport property data from a file for a solid phase.

Given a phase XML data base, this method constructs the SolidTransportData object containing the transport data for the phase.

Parameters
transportNodeReference to XML_Node containing the phase.
logReference to an XML log file. (currently unused)
phaseNamename of the corresponding phase
trReference to the SolidTransportData object that will contain the results.

Definition at line 1116 of file TransportFactory.cpp.

References XML_Node::child(), SolidTransportData::defectActivity, SolidTransportData::defectDiffusivity, SolidTransportData::electConductivity, SolidTransportData::ionConductivity, TransportFactory::m_tranPropMap, XML_Node::name(), XML_Node::nChildren(), TransportFactory::newLTP(), Cantera::showErrors(), SolidTransportData::thermalConductivity, and TransportParams::thermo.

Referenced by TransportFactory::setupSolidTransport().

void fitProperties ( GasTransportParams tr,
MMCollisionInt integrals,
std::ostream &  logfile 
)
private

Generate polynomial fits to the viscosity, conductivity, and the binary diffusion coefficients.

If CK_mode, then the fits are of the form

\[ \log(\eta(i)) = \sum_{n = 0}^3 a_n(i) (\log T)^n \]

and

\[ \log(D(i,j)) = \sum_{n = 0}^3 a_n(i,j) (\log T)^n \]

Otherwise the fits are of the form

\[ \eta(i)/sqrt(k_BT) = \sum_{n = 0}^4 a_n(i) (\log T)^n \]

and

\[ D(i,j)/sqrt(k_BT)) = \sum_{n = 0}^4 a_n(i,j) (\log T)^n \]

Parameters
trReference to the GasTransportParams object that will contain the results.
logfileReference to an ostream that will contain log information when in DEBUG_MODE
integralsinterpolator for the collision integrals

Definition at line 1173 of file TransportFactory.cpp.

References Cantera::Avogadro, Cantera::Boltzmann, GasTransportParams::condcoeffs, GasTransportParams::crot, DATA_PTR, GasTransportParams::delta, GasTransportParams::diam, GasTransportParams::diffcoeffs, GasTransportParams::eps, GasTransportParams::epsilon, Cantera::FiveSixteenths, Cantera::GasConstant, TransportFactory::getBinDiffCorrection(), TransportParams::log_level, TransportFactory::m_verbose, TransportParams::mode_, TransportParams::mw, TransportParams::nsp_, Cantera::Pi, Cantera::poly3(), Cantera::poly4(), Cantera::polyfit(), GasTransportParams::reducedMass, Phase::setTemperature(), GasTransportParams::sigma, Phase::speciesName(), TransportParams::thermo, TransportParams::tmax, TransportParams::tmin, GasTransportParams::visccoeffs, TransportParams::xml, and GasTransportParams::zrot.

Referenced by TransportFactory::setupMM().

void fitCollisionIntegrals ( std::ostream &  logfile,
GasTransportParams tr,
MMCollisionInt integrals 
)
private

Generate polynomial fits to collision integrals.

Parameters
logfileReference to an ostream that will contain log information when in DEBUG_MODE
trReference to the GasTransportParams object that will contain the results.
integralsinterpolator for the collision integrals

Definition at line 711 of file TransportFactory.cpp.

References GasTransportParams::astar_poly, GasTransportParams::bstar_poly, COLL_INT_POLY_DEGREE, GasTransportParams::cstar_poly, DATA_PTR, GasTransportParams::delta, GasTransportParams::fitlist, TransportParams::log_level, TransportFactory::m_verbose, TransportParams::mode_, TransportParams::nsp_, GasTransportParams::omega22_poly, GasTransportParams::poly, and TransportParams::xml.

Referenced by TransportFactory::setupMM().

void setupMM ( std::ostream &  flog,
const std::vector< const XML_Node * > &  transport_database,
thermo_t thermo,
int  mode,
int  log_level,
GasTransportParams tr 
)
private

Prepare to build a new kinetic-theory-based transport manager for low-density gases.

This class fills up the GastransportParams structure for the current phase

Uses polynomial fits to Monchick & Mason collision integrals. store then in tr

Parameters
flogReference to the ostream for writing log info
transport_databaseReference to a vector of pointers containing the transport database for each species
thermoPointer to the ThermoPhase object
modeMode -> Either it's CK_Mode, chemkin compatibility mode, or it is not We usually run with chemkin compatibility mode turned off.
log_levellog level
trGasTransportParams structure to be filled up with information

Definition at line 407 of file TransportFactory.cpp.

References GasTransportParams::alpha, Cantera::Avogadro, Cantera::Boltzmann, GasTransportParams::crot, GasTransportParams::delta, GasTransportParams::diam, GasTransportParams::dipole, GasTransportParams::eps, GasTransportParams::epsilon, TransportFactory::fitCollisionIntegrals(), TransportFactory::fitProperties(), TransportFactory::getTransportData(), MMCollisionInt::init(), TransportParams::log_level, TransportFactory::m_verbose, TransportFactory::makePolarCorrections(), ThermoPhase::maxTemp(), ThermoPhase::minTemp(), TransportParams::mode_, Phase::molecularWeights(), TransportParams::mw, TransportParams::nsp_, Phase::nSpecies(), GasTransportParams::polar, GasTransportParams::poly, GasTransportParams::reducedMass, DenseMatrix::resize(), GasTransportParams::sigma, Phase::speciesNames(), TransportParams::thermo, TransportParams::tmax, TransportParams::tmin, TransportParams::xml, and GasTransportParams::zrot.

Referenced by TransportFactory::initTransport().

void setupLiquidTransport ( std::ostream &  flog,
thermo_t thermo,
int  log_level,
LiquidTransportParams trParam 
)
private
void setupSolidTransport ( std::ostream &  flog,
thermo_t thermo,
int  log_level,
SolidTransportData trParam 
)
private

Prepare to build a new transport manager for solids.

Parameters
flogReference to the ostream for writing log info
thermoPointer to the ThermoPhase object
log_levellog level
trParamSolidTransportData structure to be filled up with information

Definition at line 577 of file TransportFactory.cpp.

References XML_Node::child(), TransportFactory::getSolidTransportData(), XML_Node::hasChild(), TransportParams::log_level, ThermoPhase::maxTemp(), ThermoPhase::minTemp(), Phase::molecularWeights(), TransportParams::mw, Phase::name(), TransportParams::nsp_, Phase::nSpecies(), TransportParams::thermo, TransportParams::tmax, TransportParams::tmin, and Phase::xml().

Referenced by TransportFactory::initSolidTransport().

void getBinDiffCorrection ( doublereal  t,
const GasTransportParams tr,
MMCollisionInt integrals,
size_t  k,
size_t  j,
doublereal  xk,
doublereal  xj,
doublereal &  fkj,
doublereal &  fjk 
)
private

Second-order correction to the binary diffusion coefficients.

Calculate second-order corrections to binary diffusion coefficient pair (dkj, djk). At first order, the binary diffusion coefficients are independent of composition, and d(k,j) = d(j,k). But at second order, there is a weak dependence on composition, with the result that d(k,j) != d(j,k). This method computes the multiplier by which the first-order binary diffusion coefficient should be multiplied to produce the value correct to second order. The expressions here are taken from Marerro and Mason, J. Phys. Chem. Ref. Data, vol. 1, p. 3 (1972).

Parameters
tTemperature (K)
trTransport parameters
integralsinterpolator for the collision integrals
kindex of first species
jindex of second species
xkMole fraction of species k
xjMole fraction of species j
fkjmultiplier for d(k,j)
fjkmultiplier for d(j,k)
Note
This method is not used currently.

Definition at line 74 of file TransportFactory.cpp.

References Cantera::Boltzmann, GasTransportParams::delta, GasTransportParams::eps, TransportParams::mw, and GasTransportParams::sigma.

Referenced by TransportFactory::fitProperties().

void makePolarCorrections ( size_t  i,
size_t  j,
const GasTransportParams tr,
doublereal &  f_eps,
doublereal &  f_sigma 
)
private

Corrections for polar-nonpolar binary diffusion coefficients.

Calculate corrections to the well depth parameter and the diameter for use in computing the binary diffusion coefficient of polar-nonpolar pairs. For more information about this correction, see Dixon-Lewis, Proc. Royal Society (1968).

Parameters
iSpecies one - this is a bimolecular correction routine
jspecies two - this is a bimolecular correction routine
trDatabase of species properties read in from the input xml file.
f_epsMultiplicative correction factor to be applied to epsilon(i,j)
f_sigmaMultiplicative correction factor to be applied to diam(i,j)

Definition at line 139 of file TransportFactory.cpp.

References GasTransportParams::alpha, GasTransportParams::dipole, GasTransportParams::eps, GasTransportParams::polar, and GasTransportParams::sigma.

Referenced by TransportFactory::setupMM().

Member Data Documentation

TransportFactory * s_factory = 0
staticprivate

Static instance of the factor -> This is the only instance of this object allowed.

Definition at line 169 of file TransportFactory.h.

Referenced by TransportFactory::deleteFactory(), and TransportFactory::factory().

mutex_t transport_mutex
staticprivate

Static instance of the mutex used to ensure the proper reading of the transport database.

Definition at line 172 of file TransportFactory.h.

Referenced by TransportFactory::deleteFactory(), TransportFactory::factory(), and TransportFactory::initTransport().

bool m_verbose
private
std::map<std::string, int> m_models
private

Mapping between between the string name for a transport model and the integer name.

Definition at line 382 of file TransportFactory.h.

Referenced by TransportFactory::newTransport(), and TransportFactory::TransportFactory().

std::map<int, std::string> m_modelNames
private

Inverse mapping of transport models, from integer constant to string.

Definition at line 385 of file TransportFactory.h.

Referenced by TransportFactory::modelName(), and TransportFactory::TransportFactory().

std::map<std::string, TransportPropertyType> m_tranPropMap
private

Mapping between between the string name for a transport property and the integer name.

Definition at line 389 of file TransportFactory.h.

Referenced by TransportFactory::getLiquidInteractionsTransportData(), TransportFactory::getLiquidSpeciesTransportData(), TransportFactory::getSolidTransportData(), and TransportFactory::TransportFactory().

std::map<std::string, LTPTemperatureDependenceType> m_LTRmodelMap
private

Mapping between between the string name for a species-specific transport property model and the integer name.

Definition at line 393 of file TransportFactory.h.

Referenced by TransportFactory::newLTP(), and TransportFactory::TransportFactory().

std::map<std::string, LiquidTranMixingModel> m_LTImodelMap
private

Mapping between between the string name for a liquid mixture transport property model and the integer name.

Definition at line 397 of file TransportFactory.h.

Referenced by TransportFactory::newLTI(), and TransportFactory::TransportFactory().


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