OpenFOAM 热物理库深度解析

本帖描述 OpenFOAM 热物理库的框架结构(以OF-7为基础,同时也尝试追踪自OpenFOAM开创以来,整个热物理库的架构变迁过程),希望对你有用。

OF-7

UML 类图

下图是热物理库部分类的 UML 类图。欲获取更佳阅读体验,可以右键在新标签页打开该图。

OF-7 UML

求解器中创建对象

solver 中的 createFields.H 中:

autoPtr<psiReactionThermo> pThermo(psiReactionThermo::New(mesh));
psiReactionThermo& thermo = pThermo();

创建 thermo,在 reactingFoam 自带算例中,创建的是 hePsiThermo 的对象,但是这里 thermo 是基类 psiReactionThermo 的引用,即只能使用 psiReactionThermo 这一继承链的方法。

basicSpecieMixture& composition = thermo.composition();

创建 composition,是指系统中的混合物。可以调用继承链中 mixture 这一分支的方法。

燃烧类的求解器中,一般都没有对单个组分进行操作,这些操作都在热物理库里边进行。

multiComponentMixture 类中的私有数据 speciesData_ 的类型是 PtrList<sutherlandTransport<species::thermo<janafThermo<perfectGas<specie>>,sensibleEnthalpy>>>,它是组分的 list,可以用它获取组分的所有信息。solver 中可以通过 multiComponentMixturesolvercomposition 的类型是 basicSpecieMixture,需要强制转换成 multiComponentMixture 才行 ) 的 speciesData() 或者 getLocalThermo(speciei) 这两个函数来调用。如:

dynamic_cast<const multiComponentMixture<gasHThermoPhysics>&>(composition).speciesData()[speciei].rho(p,T)

也可以通过 multiComponentMixture 的派生类 reactingMixture 来调用。

dynamic_cast<const reactingMixture<gasHThermoPhysics>&>(composition).speciesData()[speciei].rho(p,T)

也可以通过派生类的派生类 SpecieMixture 来调用。

dynamic_cast<const SpecieMixture<reactingMixture<gasHThermoPhysics>>&>(composition).speciesData()[speciei].rho(p,T)

上述操作说明可以使用组分列表调用单个组分的函数(即 sutherlandTransport<species::thermo<janafThermo<perfectGas<specie>>,sensibleEnthalpy>> 继承链中的任何函数),如 speciesData_[i].rho(p, T),这里的 rho(p,T) 来自于 perfectGas 类。

也可以在求解器中重新构造组分列表:

autoPtr<psiReactionThermo> pThermo(psiReactionThermo::New(mesh));
psiReactionThermo& thermo = pThermo();

basicSpecieMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();

PtrList<gasHThermoPhysics> specieData(Y.size());
forAll(specieData, i)
{
//对于每个组分,通过 thermo 调用 Lib 中的 speciesData_ 来构造 solver 中的 specieData
specieData.set
(
i,
new gasHThermoPhysics
(
dynamic_cast<const reactingMixture<gasHThermoPhysics>&>
(thermo).speciesData()[i]
)
);
}

// 每个网格构造一个 mixture,这个 mixture 相当于 multiComponentMixture 中的私有数据 mixture_。
// 它是混合物等效成的一个组分,和单个组分同一类型,即 `gasHThermoPhysics`。
forAll(T, celli)
{
gasHThermoPhysics mixture(0.*specieData[0]);
forAll(Y, speciei)
{
mixture += Y[speciei][celli]*specieData[speciei];
}
}

gasHThermoPhysicssutherlandTransport<species::thermo<janafThermo<perfectGas<specie>>,sensibleEnthalpy>> 的别名,它是单个组分的类型。
事实上,multiComponentMixture 类中的函数 cellMixture 函数可以用于获取 mixture_

template<class ThermoType>
const ThermoType& Foam::multiComponentMixture<ThermoType>::cellMixture
(
const label celli
) const
{
mixture_ = Y_[0][celli]*speciesData_[0];

for (label n=1; n<Y_.size(); n++)
{
mixture_ += Y_[n][celli]*speciesData_[n];
}

return mixture_;
}

使用方法参考 hePsiThermo.C

const scalarField& hCells = this->he_;
const scalarField& pCells = this->p_;

scalarField& TCells = this->T_.primitiveFieldRef();
scalarField& psiCells = this->psi_.primitiveFieldRef();
scalarField& muCells = this->mu_.primitiveFieldRef();
scalarField& alphaCells = this->alpha_.primitiveFieldRef();

forAll(TCells, celli)
{
const typename MixtureType::thermoType& mixture_ =
this->cellMixture(celli);

TCells[celli] = mixture_.THE
(
hCells[celli],
pCells[celli],
TCells[celli]
);

psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);

muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
}

RTS 选择热物理模型

从字典文件的设置开始

下面介绍我们是如何从 thermophysicalProperties 字典文件中的下列设置构建对象的:

thermoType
{
type hePsiThermo;
mixture reactingMixture;
transport sutherland;
thermo janaf;
energy sensibleEnthalpy;
equationOfState perfectGas;
specie specie;
}

solver 中创建模型:

autoPtr<psiReactionThermo> pThermo(psiReactionThermo::New(mesh));

这个函数定义于 psiReactionThermo

Foam::autoPtr<Foam::psiReactionThermo> Foam::psiReactionThermo::New
(
const fvMesh& mesh,
const word& phaseName
)
{
return basicThermo::New<psiReactionThermo>(mesh, phaseName);
}

它调用了 basicThermo 中的 New 函数:

template<class Thermo>
Foam::autoPtr<Thermo> Foam::basicThermo::New(const fvMesh& mesh,const word& phaseName)
{
// IOdictionary thermoDict---thermophysicalProperties
typename Thermo::fvMeshConstructorTable::iterator cstrIter =
lookupThermo<Thermo, typename Thermo::fvMeshConstructorTable>
(
thermoDict,
Thermo::fvMeshConstructorTablePtr_
);
//fvMeshConstructorTablePtr_
//这个是在 runTimeSelectionTables.H 的 declareRunTimeSelectionTable 函数中定义的。
return autoPtr<Thermo>(cstrIter()(mesh, phaseName));
}

调用了 2 参数的 lookupThermo:
basicThermoTemplates.C 中的 lookupThermo 函数:

template<class Thermo, class Table>
typename Table::iterator Foam::basicThermo::lookupThermo
(
const dictionary& thermoDict,
Table* tablePtr //即 fvMeshConstructorTablePtr_,它是哈希表add类型的一个指针,typeName和模型都存放在这里
)
{
const word thermoTypeName
(
word(thermoTypeDict.lookup("type")) + '<'
+ word(thermoTypeDict.lookup("mixture")) + '<'
+ word(thermoTypeDict.lookup("transport")) + '<'
+ word(thermoTypeDict.lookup("thermo")) + '<'
+ word(thermoTypeDict.lookup("equationOfState")) + '<'
+ word(thermoTypeDict.lookup("specie")) + ">>,"
+ word(thermoTypeDict.lookup("energy")) + ">>>"
);
//通过字典文件读入的数据,构造 thermoTypeName: //hePsiThermo<reactingMixture<sutherland<janaf<perfectGas<specie>>,sensibleEnthalpy>>>
//注意这里并不是类名,而是省略版

return lookupThermo<Thermo, Table>
(
thermoTypeDict,
tablePtr,
nCmpt,
cmptNames,
thermoTypeName
);
}

最后返回值是调用 5 个参数的 lookupThermo 函数:

template<class Thermo, class Table>
typename Table::iterator Foam::basicThermo::lookupThermo
(
const dictionary& thermoTypeDict,
Table* tablePtr,
const int nCmpt,
const char* cmptNames[],
const word& thermoTypeName
)
{
// Lookup the thermo package
typename Table::iterator cstrIter = tablePtr->find(thermoTypeName);
return cstrIter;
}
//上述步骤是根据字典文件在备选库中选择我们需要的组合,那么备选库是怎么形成的呢?

第一个 lookupThermo 函数的作用是,将 thermophysicalProperties 字典文件中的设置翻译成:

hePsiThermo<reactingMixture<sutherland<janaf<perfectGas<specie>>,sensibleEnthalpy>>>

但该字符串并不是最终的 class,比如 janaf,实际上类名是 janafThermo

第二个 lookupThermo 函数的作用就是,根据该字符串,在 RTS 库中找到我们选择的类。
那么上面提到的省略版字符串,怎么和 RTS 库中的组合匹配的呢?下面来看 RTS 库的生成过程。

RTS 库的生成

按照Giskard 博客提到的 RTS 流程:

  1. 基类类体里调用 TypeNamedeclareRunTimeSelectionTable 两个函数,类体外面调用 defineTypeNameAndDebugdefineRunTimeSelectionTable 两个函数。
  2. 基类中需要一个静态 New 函数作为 selector
  3. 派生类类体中需要调用 TypeName 函数,类体外调用 defineRunTimeSelectionTableaddToRunTimeSelectionTable 两个宏函数。

我们发现:makeReactionThermo.H 中有两个重要函数:defineThermoPhysicsReactionThermomakeThermoPhysicsReactionThermos

下面一一介绍:

defineThermoPhysicsReactionThermo

其中的 defineThermoPhysicsReactionThermo, 这个函数通过一组形参列表得到 CThermo##Mixture##ThermoPhys,然后再转换成一串简化的字符串。这里的形参列表就是实例化所需要提供的各种类,对应(psiReactionThermo,hePsiThermo,reactingMixture,sutherlandTransport)

#define defineThermoPhysicsReactionThermo(BaseReactionThermo,CThermo,Mixture,ThermoPhys) \
typedef CThermo//hePsiThermo \
< \
BaseReactionThermo,//psiReactionThermo \
SpecieMixture \
< \
Mixture//reactingMixture \
< \
ThermoPhys//sutherlandTransport \
> \
> \
> CThermo##Mixture##ThermoPhys; \
\
defineTemplateTypeNameAndDebugWithName \
( \
CThermo##Mixture##ThermoPhys, \
(#CThermo"<" + Mixture<ThermoPhys>::typeName() + ">").c_str(), \
0 \
)
//作用:把尖括号表示的类变成字符串 CThermo##Mixture##ThermoPhys,然后再变成简化的尖括号。

如:

hePsiThermo
<
psiReactionThermo,
SpecieMixture
<
reactingMixture
<
sutherlandTransport
<
species::thermo
<
janafThermo<perfectGas<specie>>,sensibleEnthalpy
>
>
>
>
>

变成 RTS 库中保存的字符串:
hePsiThermo<reactingMixture<sutherland<janaf<perfectGas<specie>>,sensibleEnthalpy>>>
隐藏了一些信息,如 psiReactionThermoSpecieMixture

具体过程:
typeName() 依次调用 reactingMixturesutherlandTransport, species::thermo, janafThermo, sensibleEnthalpy, perfectGas, specietypeName()。返回值依次是:

"reactingMixture<" + ThermoType::typeName() + '>'
sutherland<" + Thermo::typeName() + '>'
Thermo::typeName() + ',' + Type<thermo<Thermo, Type>>::typeName()
"janaf<" + EquationOfState::typeName() + '>'
"sensibleEnthalpy"
"perfectGas<" + word(Specie::typeName_()) + '>'

最底层的 specie 类只有一个 ClassName("specie");
我们找到一个文件className.H 中:

#define ClassName(TypeNameString)                             \
ClassNameNoDebug(TypeNameString); \
static int debug

#define ClassNameNoDebug(TypeNameString) \
static const char* typeName_() { return TypeNameString; } \
static const ::Foam::word typeName

所以这里的意义就在于用简化版的字符串表示完整的类型。而简化版的字符串则是读取 thermophysicalProperties 之后稍作处理得到的。(在前面的2 参数的 lookupThermo 中)

makeThermoPhysicsReactionThermos

定义如下:

#define makeThermoPhysicsReactionThermos(BaseThermo,BaseReactionThermo,CThermo,Mixture,ThermoPhys)
// 定义于 makeReactionThermo.H
defineThermoPhysicsReactionThermo(BaseReactionThermo,CThermo,Mixture,ThermoPhys);

addThermoPhysicsThermo(basicThermo, CThermo##Mixture##ThermoPhys);
addThermoPhysicsThermo(fluidThermo, CThermo##Mixture##ThermoPhys);
addThermoPhysicsThermo(BaseThermo, CThermo##Mixture##ThermoPhys);
addThermoPhysicsThermo(BaseReactionThermo, CThermo##Mixture##ThermoPhys)

调用了上边的 defineThermoPhysicsReactionThermo 用于创建 RTS 中的 typeName,调用了 addThermoPhysicsThermo(位于 makeThermo.H), 用于写入到 RTS 哈希表。
makeThermoPhysicsReactionThermos 函数的形参:(BaseThermo,BaseReactionThermo,CThermo,Mixture,ThermoPhys),这是实例化需要提供的所有信息。
比如其中的一种组合:makeThermoPhysicsReactionThermos(psiThermo,psiReactionThermo,hePsiThermo,reactingMixture,gasHThermoPhysics)(所有的组合在 psiReactionThermos.C 中列出)。

这里调用 4 个 addThermoPhysicsThermo 函数是为何?假定先分析第一个 addThermoPhysicsThermo 函数:

#define addThermoPhysicsThermo(BaseThermo_,CThermoMixtureThermoPhys_) //定义于makeThermo.H
addToRunTimeSelectionTable
(
BaseThermo_,//基类名
CThermoMixtureThermoPhys_,//派生类名
fvMesh
);

BaseThermo_basicThermoCThermoMixtureThermoPhys_hePsiThermo##reactingMixture##gasHThermoPhysics

这里的 addToRunTimeSelectionTable 函数定义于 addToRunTimeSelectionTable.H 中:

#define addToRunTimeSelectionTable(baseType,thisType,argNames)
baseType::add##argNames##ConstructorToTable<thisType>
add##thisType##argNames##ConstructorTo##baseType##Table_

addfvMeshConstructorToTable 这个类是在 runTimeSelectionTables.H 文件中的 declareRunTimeSelectionTable 宏函数中定义的。

所以 addThermoPhysicsThermo 中创建了一个 baseType::addfvMeshConstructorToTable<thisType> 类的对象 add##thisType##fvMesh##ConstructorTo##baseType##Table_,而该对象构造时会往哈希表中插入键和值,即派生类的名字和 New 函数,此函数返回指向派生类对象的基类指针。
这里的基类 baseTypebasicThermo,派生类 thisTypehePsiThermo##reactingMixture##gasHThermoPhysics,这里是省略版。

OF-8

OF-8 去除了 reactingMixture 这个类!
它的成员 speciesComposition_ 移至 multiComponentMixture 中。
另外 multiComponentMixture 还新增了 readSpeciesDatareadSpeciesComposition 这两个函数,用于初始化 specieThermos_speciesComposition_

之前(OF-7)化学反应的信息是保存在 reactingMixture 中,然后复制引用到 StandardChemistryModel中的 reactions_
现在(OF-8)是直接保存在 StandardChemistryModel 中的 reactions_

化学反应的信息在 StandardChemistryModel 中,以前(OF-7)是

const PtrList<Reaction<ThermoType>>& reactions_;

以前(OF-7)是下面这样,调用的是复制构造函数,因为所有的化学反应的信息都包含在 reactingMixture 中!

reactions_
(
dynamic_cast<const reactingMixture<ThermoType>&>(this->thermo())
)

现在(OF-8)是

const ReactionList<ThermoType> reactions_;

现在在 StandardChemistryModel 中是这样构造:

reactions_
(
dynamic_cast<const multiComponentMixture<ThermoType>&>
(
this->thermo()
).species(),
specieThermo_,
this->mesh(),
*this
)

调用的是这个构造函数:

//- Construct from thermo list, objectRegistry and dictionary
ReactionList
(
const speciesTable& species,
const PtrList<ThermoType>& speciesThermo,
const objectRegistry& ob,
const dictionary& dict
);basicThermo

下图是热物理库部分类的 UML 类图。欲获取更佳阅读体验,可以右键在新标签页打开该图。

OF-8 UML

OF-10

OF-10 UML

OF-1.1

OF-1.1 UML

reactingFoam

autoPtr<hCombustionThermo> thermo
(
hCombustionThermo::New(mesh)
);

combustionMixture& composition = thermo->composition();

CASE/constant/thermodynamicProperties

thermoType hMixtureThermo<reactingMixture>;
typedef sutherlandTransport<specieThermo<janafThermo<perfectGas> > >
reactionThermo;

OF-1.2

OF-1.2 UML
typedef sutherlandTransport<specieThermo<janafThermo<perfectGas> > >
reactionThermo;

reactingFoam:

autoPtr<hCombustionThermo> thermo
(
hCombustionThermo::New(mesh)
);

combustionMixture& composition = thermo->composition();
PtrList<volScalarField>& Y = composition.Y();

CASE/constant/thermophysicalProperties:

thermoType hMixtureThermo<reactingMixture>;

OF-1.3

OF-1.3 UML
typedef sutherlandTransport<specieThermo<janafThermo<perfectGas> > >
reactionThermo;

reactingFoam:

autoPtr<hCombustionThermo> thermo
(
hCombustionThermo::New(mesh)
);

combustionMixture& composition = thermo->composition();
PtrList<volScalarField>& Y = composition.Y();

CASE/constant/thermophysicalProperties:

thermoType hMixtureThermo<reactingMixture>;

OF-1.4

OF-1.4 UML
typedef sutherlandTransport<specieThermo<janafThermo<perfectGas> > >
reactionThermo;

reactingFoam:

autoPtr<hCombustionThermo> thermo
(
hCombustionThermo::New(mesh)
);

combustionMixture& composition = thermo->composition();
PtrList<volScalarField>& Y = composition.Y();

CASE/constant/thermophysicalProperties:

thermoType hMixtureThermo<reactingMixture>;

OF-1.5

OF-1.5 UML
typedef sutherlandTransport<specieThermo<janafThermo<perfectGas> > >
reactionThermo;

reactingFoam:

autoPtr<hCombustionThermo> thermo
(
hCombustionThermo::New(mesh)
);

combustionMixture& composition = thermo->composition();
PtrList<volScalarField>& Y = composition.Y();

CASE/constant/thermophysicalProperties:

thermoType hMixtureThermo<reactingMixture>;

OF1.2到OF1.5继承关系没有发生任何变化。

OF-1.6

OF-1.6 UML
typedef sutherlandTransport<specieThermo<janafThermo<perfectGas> > >
gasThermoPhysics;

reactingFoam:

autoPtr<psiChemistryModel> pChemistry
(
psiChemistryModel::New(mesh)
);
psiChemistryModel& chemistry = pChemistry();

hCombustionThermo& thermo = chemistry.thermo();

basicMultiComponentMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();

CASE/constant/thermophysicalProperties:

thermoType      hPsiMixtureThermo<reactingMixture<gasThermoPhysics>>;

OF-1.7.0

OF-1.7.0 UML
typedef sutherlandTransport<specieThermo<janafThermo<perfectGas> > >
gasThermoPhysics;

reactingFoam:

autoPtr<psiChemistryModel> pChemistry
(
psiChemistryModel::New(mesh)
);
psiChemistryModel& chemistry = pChemistry();

hsCombustionThermo& thermo = chemistry.thermo();

basicMultiComponentMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();

CASE/constant/thermophysicalProperties:

thermoType hsPsiMixtureThermo<reactingMixture<gasThermoPhysics>>;

OF-2.0.x

OF-2.0.x UML

reactingMixture 的模板是 gasThermoPhysics,所以UML图中仍旧写的是 sutherlandTransport

typedef sutherlandTransport<specieThermo<janafThermo<perfectGas> > > gasThermoPhysics;

reactingFoam:

autoPtr<psiChemistryModel> pChemistry
(
psiChemistryModel::New(mesh)
);
psiChemistryModel& chemistry = pChemistry();

hsCombustionThermo& thermo = chemistry.thermo();
basicMultiComponentMixture& composition = thermo.composition();

CASE/constant/thermophysicalProperties:

thermoType hsPsiMixtureThermo<reactingMixture<gasThermoPhysics>>;

OF-2.3.x

OF-2.3.x UML

相比于2.0.x,把显焓 hs 抽象成了一个类 sensibleEnthalpy。

reactingFoam:

autoPtr<combustionModels::psiCombustionModel> reaction
(
combustionModels::psiCombustionModel::New(mesh)
);

psiReactionThermo& thermo = reaction->thermo();
thermo.validate(args.executable(), "h", "e");

basicMultiComponentMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();

CASE/constant/thermophysicalProperties:

thermoType
{
type hePsiThermo;
mixture reactingMixture;
transport sutherland;
thermo janaf;
energy sensibleEnthalpy;
equationOfState perfectGas;
specie specie;
}

OF-4.x

OF-4.x UML

reactingFoam:

autoPtr<combustionModels::psiCombustionModel> reaction
(
combustionModels::psiCombustionModel::New(mesh)
);

psiReactionThermo& thermo = reaction->thermo();
thermo.validate(args.executable(), "h", "e");

basicMultiComponentMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();

CASE/constant/thermophysicalProperties:

thermoType
{
type hePsiThermo;
mixture reactingMixture;
transport sutherland;
thermo janaf;
energy sensibleEnthalpy;
equationOfState perfectGas;
specie specie;
}
文章作者: Yan Zhang
文章链接: https://openfoam.top/thermodynamicLIB/
版权声明: 本博客所有文章除特别声明外,均采用 CC BY-NC-SA 4.0 许可协议。转载请注明来自 OpenFOAM 成长之路
您的肯定会给我更大动力~