更新日志
- 2019.01.10
- 修正
incompressible::turbulenceModel::New
实际调用函数的描述
上篇文章中已经介绍了湍流模型框架的顶层设计,本篇讲解求解器中需要调用到的三个类模板。
OpenFOAM 将湍流模拟分为三类:laminar,RAS 及 LES。laminar 表示直接求解 Navier-Stokes 方程,不引入任何湍流模拟,是一种 quasi-DNS 方法。RAS 表示用雷诺平均方法模拟湍流。LES 表示用大涡模拟方法模拟湍流。这三种方法用三个类模板表示:laminarModel
、RASModel
和 LESModel
。
注意
目前只有 IncompressibleTurbulenceModel
和 CompressibleTurbulenceModel
这两种湍流流动是按照上述原则划分。至于 PhaseIncompressibleTurbulenceModel
和 PhaseCompressibleTurbulenceModel
,其代码完成度不高,并没有采取以上原则划分。
IncompressibleTurbulenceModel
在 incompressible
命名空间下特化后得到 incompressible::turbulenceModel
(类型别名),CompressibleTurbulenceModel
在 compressible
命名空间下特化得到 compressible::turbulenceModel
(类型别名)。下面分别分析这两个类的框架(基于OpenFOAM-dev/tree/20181203)。
不可压 turbulenceModel
incompressible::turbulenceModel
被用作模板参数对三个类模板进行特化,同时分别为它们定义类型别名 laminarModel
、RASModel
和 LESModel
(注意类型别名和模板类的名称相同):
1
2
3
4
5
6
7
8
9
|
// $FOAM_SRC/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModel.H
namespace incompressible
{
typedef IncompressibleTurbulenceModel<transportModel> turbulenceModel;
typedef laminarModel<turbulenceModel> laminarModel;
typedef RASModel<turbulenceModel> RASModel;
typedef LESModel<turbulenceModel> LESModel;
// ...
|
这三个均为基于原则设计的类,继承自唯一的模板参数。它们的声明如下:
1
2
3
4
5
6
7
|
// $FOAM_SRC/TurbulenceModels/turbulenceModels/laminar/laminarModel/laminarModel.H
template<class BasicTurbulenceModel>
class laminarModel
:
public BasicTurbulenceModel
{
// ...
|
1
2
3
4
5
6
7
|
// $FOAM_SRC/TurbulenceModels/turbulenceModels/RAS/RASModel/RASModel.H
template<class BasicTurbulenceModel>
class RASModel
:
public BasicTurbulenceModel
{
// ...
|
1
2
3
4
5
6
7
|
// $FOAM_SRC/TurbulenceModels/turbulenceModels/LES/LESModel/LESModel.H
template<class BasicTurbulenceModel>
class LESModel
:
public BasicTurbulenceModel
{
// ...
|
三个类特化时采用相同的模板参数 incompressible::turbulenceModel
,均为 incompressible::turbulenceModel
的派生类,其继承关系如下:
incompressible::turbulenceModel 的继承关系(点击图片放大)
可以通过一个 autoPtr<incompressible::turbulenceModel>
来描述任意一个派生类实例化后的对象。OpenFOAM 通过一个返回 autoPtr<incompressible::turbulenceModel>
的 selector,实现了运行时选择,从字典文件读取关键字并构造对应的湍流模型实例。
incompressible::turbulenceModel
中,incompressible
是命名空间,turbulenceModel
是类型别名,对应 IncompressibleTurbulenceModel<transportModel>
,selector 的对应代码为:
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
|
// $FOAM_SRC/TurbulenceModels/incompressible/IncompressibleTurbulenceModel/IncompressibleTurbulenceModel.C
template<class TransportModel>
Foam::autoPtr<Foam::IncompressibleTurbulenceModel<TransportModel>>
Foam::IncompressibleTurbulenceModel<TransportModel>::New
(
const volVectorField& U,
const surfaceScalarField& phi,
const TransportModel& transport,
const word& propertiesName
)
{
return autoPtr<IncompressibleTurbulenceModel>
(
static_cast<IncompressibleTurbulenceModel*>(
TurbulenceModel
<
geometricOneField,
geometricOneField,
incompressibleTurbulenceModel,
TransportModel
>::New
(
geometricOneField(),
geometricOneField(),
U,
phi,
phi,
transport,
propertiesName
).ptr())
);
}
|
这是一个函数模板,在各求解器的源码中特化并调用,如 pimpleFoam 中的相关代码:
1
2
3
4
5
|
// FOAM_SOLVERS/incompressible/pimpleFoam/createFields.H
autoPtr<incompressible::turbulenceModel> turbulence
(
incompressible::turbulenceModel::New(U, phi, laminarTransport)
);
|
在两相流求解器 interFoam 中也采用了相同的定义,密度是一个 GeometricOneField
。这样的做法存在较为严重的问题,它无法考虑由密度梯度引起的浮力效应(buoyancy effect)。用 interFoam 模拟波浪问题时,如果采用湍流模型,那么在自由面附近会出现较大的湍流粘度,从而导致波浪衰减。若存在波浪翻卷、破碎等现象,波浪衰减问题更为明显。可参考相关文献 。
可压 turbulenceModel
和不可压类似,可压缩版本的 compressible::turbulenceModel
也被作为模板参数对三个类模板进行特化,同时分别定义类型别名:
1
2
3
4
5
6
7
8
9
10
|
// $FOAM_SRC/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModel.H
namespace incompressible
{
typedef ThermalDiffusivity<CompressibleTurbulenceModel<fluidThermo>>
turbulenceModel;
typedef laminarModel<turbulenceModel> laminarModel;
typedef RASModel<turbulenceModel> RASModel;
typedef LESModel<turbulenceModel> LESModel;
// ...
|
这三个类由于采用相同的模板参数,其基类也相同。其继承关系如下:
compressible::turbulenceModel 的继承关系(点击图片放大)
这里也采用了一个返回 autoPtr<compressible::turbulenceModel>
的 selector,实现了运行时选择,从字典文件读取关键字并构造对应的湍流模型实例:
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
|
// $FOAM_SRC/TurbulenceModels/compressible/CompressibleTurbulenceModel/CompressibleTurbulenceModel.C
template<class TransportModel>
Foam::autoPtr<Foam::CompressibleTurbulenceModel<TransportModel>>
Foam::CompressibleTurbulenceModel<TransportModel>::New
(
const volScalarField& rho,
const volVectorField& U,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName
)
{
return autoPtr<CompressibleTurbulenceModel>
(
static_cast<CompressibleTurbulenceModel*>(
TurbulenceModel
<
geometricOneField,
volScalarField,
compressibleTurbulenceModel,
transportModel
>::New
(
geometricOneField(),
rho,
U,
phi,
phi,
transport,
propertiesName
).ptr())
);
}
|
这个 selector 在 combustion、compressible、heatTransfer、lagrangian 和 multiphase 等类型的求解器中特化并调用。如:
1
2
3
4
5
6
7
8
9
10
11
12
|
// $FOAM_SOLVERS/compressible/rhoPimpleFoam/createFields.H
Info<< "Creating turbulence model\n" << endl;
autoPtr<compressible::turbulenceModel> turbulence
(
compressible::turbulenceModel::New
(
rho,
U,
phi,
thermo
)
);
|
又如:
1
2
3
4
5
6
7
8
9
10
11
12
|
// $FOAM_SOLVERS/multiphase/compressibleMultiphaseInterFoam/createFields.H
// Construct compressible turbulence model
autoPtr<compressible::turbulenceModel> turbulence
(
compressible::turbulenceModel::New
(
rho,
U,
mixture.rhoPhi(),
mixture
)
);
|
三个类模板的使用
事实上,laminarModel
、RASModel
和 LESModel
不会单独使用,而是作为模板参数去特化一些更具体的类模板。比如目前有三个类模板用 laminarModel
特化,它们分别是 Stokes
、generalizedNewtonian
、Maxwell
。而用 RASModel
和 LESModel
特化的类模板更多,这两个类的实现更为复杂,会在后续文章中介绍。
派生类
所有具体的流动模型都可以看作是由 laminarModel
、RASModel
或 LESModel
派生出来的类。以 Maxwell
为例,其声明如下:
1
2
3
4
5
6
7
|
// $FOAM_SRC/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.H
template<class BasicTurbulenceModel>
class Maxwell
:
public laminarModel<BasicTurbulenceModel>
{
// ...
|
又比如 Stokes
,其声明如下:
1
2
3
4
5
6
7
|
// $FOAM_SRC/TurbulenceModels/turbulenceModels/laminar/Stokes/Stokes.H
template<class BasicTurbulenceModel>
class Maxwell
:
public linearViscousStress<laminarModel<BasicTurbulenceModel>>
{
// ...
|
这里的 linearViscousStress
也是基于原则设计的类,继承自其模板参数 laminarModel
:
1
2
3
4
5
6
7
|
// $FOAM_SRC/TurbulenceModels/turbulenceModels/linearViscousStress/linearViscousStress.H
template<class BasicTurbulenceModel>
class linearViscousStress
:
public BasicTurbulenceModel
{
// ...
|
laminarModel 的继承关系(点击图片放大)
上图为 laminarModel
的继承关系,其中三个表示具体流动类型的类(Maxwell
、generalizedNewtonian
和 Stokes
)需要进行模板特化才可以使用。
思考
对这些类模板进行模板特化时,模板参数是什么?
派生类的特化
将不可压缩和可压缩流分开考虑,代入不同的模板参数:
思考
transportModelIncompressibleTurbulenceModel
和 fluidThermoCompressibleTurbulenceModel
这两个类的定义在哪?
这两个类通过 makeTurbulenceModelTypes
这个宏完成的定义。
不可压缩流体
对于不可压流,makeTurbulenceModelTypes
宏命令如下:
1
2
3
4
5
6
7
8
9
10
|
// $FOAM_SRC/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C
makeBaseTurbulenceModel
(
geometricOneField,
geometricOneField,
incompressibleTurbulenceModel,
IncompressibleTurbulenceModel,
transportModel
);
|
makeTurbulenceModelTypes
宏的具体定义如下:
1
2
3
4
5
6
7
8
9
10
11
|
// $FOAM_SRC/TurbulenceModels/turbulenceModels/makeTurbulenceModel.H
#define makeTurbulenceModelTypes(Alpha, Rho, baseModel, BaseModel, Transport) \
\
namespace Foam \
{ \
typedef BaseModel<Transport> Transport##BaseModel; \
typedef laminarModel<Transport##BaseModel> \
laminar##Transport##BaseModel; \
typedef RASModel<Transport##BaseModel> RAS##Transport##BaseModel; \
typedef LESModel<Transport##BaseModel> LES##Transport##BaseModel; \
}
|
宏命令展开后,会定义4个类,分别如下:
transportModelIncompressibleTurbulenceModel
laminartransportModelIncompressibleTurbulenceModel
RAStransportModelIncompressibleTurbulenceModel
LEStransportModelIncompressibleTurbulenceModel
可压缩流体
对于可压缩流,makeTurbulenceModelTypes
宏命令如下:
1
2
3
4
5
6
7
8
9
10
|
// $FOAM_SRC/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C
makeBaseTurbulenceModel
(
geometricOneField,
volScalarField,
compressibleTurbulenceModel,
CompressibleTurbulenceModel,
ThermalDiffusivity,
fluidThermo
);
|
makeTurbulenceModelTypes
宏的具体定义如下:
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
|
// $FOAM_SRC/TurbulenceModel/compressible/turbulentFluidThermoModels/makeTurbulenceModel.H
#define makeTurbulenceModelTypes( \
Alpha, Rho, baseModel, BaseModel, TDModel, Transport) \
\
namespace Foam \
{ \
typedef TDModel<BaseModel<Transport>> \
Transport##BaseModel; \
typedef laminarModel<Transport##BaseModel> \
laminar##Transport##BaseModel; \
typedef RASModel<EddyDiffusivity<Transport##BaseModel>> \
RAS##Transport##BaseModel; \
typedef LESModel<EddyDiffusivity<Transport##BaseModel>> \
LES##Transport##BaseModel; \
}
|
宏命令展开后,会定义4个类,分别如下:
fluidThermoCompressibleTurbulenceModel
laminarfluidThermoCompressibleTurbulenceModel
RASfluidThermoCompressibleTurbulenceModel
LESfluidThermoCompressibleTurbulenceModel
用宏命令特化
各具体模型类的模板特化用到了多个宏定义,其中最上层的宏为 makeLaminarModel
、makeRASModel
和 makeLESModel
。这三个宏的功能类似,这里选取 makeLaminarModel
为例介绍。
不可压缩流体
不可压缩流的 makeLaminarModel
宏命令的代码如下:
1
2
3
4
5
6
7
8
9
10
11
12
13
|
// $FOAM_SRC/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C
// -------------------------------------------------------------------------- //
// Laminar models
// -------------------------------------------------------------------------- //
#include "Stokes.H"
makeLaminarModel(Stokes);
#include "generalizedNewtonian.H"
makeLaminarModel(generalizedNewtonian);
#include "Maxwell.H"
makeLaminarModel(Maxwell);
|
makeLaminarModel
宏接受一个参数,宏的定义如下:
1
2
3
4
|
// $FOAM_SRC/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.H
#define makeLaminarModel(Type) \
makeTemplatedTurbulenceModel \
(transportModelIncompressibleTurbulenceModel, laminar, Type)
|
可见它调用了另外一个宏 makeTemplatedTurbulenceModel
,这个宏接受三个参数,定义如下:
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
|
// $FOAM_SRC/TurbulenceModels/turbulenceModels/makeTurbulenceModel.H
#define makeTemplatedTurbulenceModel(BaseModel, SType, Type) \
defineNamedTemplateTypeNameAndDebug \
(Foam::SType##Models::Type<Foam::BaseModel>, 0); \
\
namespace Foam \
{ \
namespace SType##Models \
{ \
typedef Type<BaseModel> Type##SType##BaseModel; \
\
addToRunTimeSelectionTable \
( \
SType##BaseModel, \
Type##SType##BaseModel, \
dictionary \
); \
} \
}
|
以 Stokes
为例,不可压缩版本的 makeLaminarModel(Stokes)
宏展开会生成如下代码:
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
|
defineNamedTemplateTypeNameAndDebug
(Foam::laminarModels::Stokes<Foam::transportModelIncompressibleTurbulenceModel>, 0);
namespace Foam
{
namespace laminarModels
{
typedef Stokes<transportModelIncompressibleTurbulenceModel>
StokeslaminartransportModelIncompressibleTurbulenceModel;
addToRunTimeSelectionTable
(
laminartransportModelIncompressibleTurbulenceModel,
StokeslaminartransportModelIncompressibleTurbulenceModel,
dictionary
);
}
}
|
上面的 defineNamedTemplateTypeNameAndDebug
以及 addToRunTimeSelectionTable
是和运行时选择有关的宏。这样便定义了一个名为 StokeslaminartransportModelIncompressibleTurbulenceModel
的类,并将这个类的相关信息加入到哈希表中。
可压缩流体
可压缩流体的 makeLaminarModel
宏命令的代码和不可压缩流相同,只不过在不同的文件中。具体如下:
1
2
3
4
5
6
7
8
9
10
11
12
13
|
// $FOAM_SRC/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C
// -------------------------------------------------------------------------- //
// Laminar models
// -------------------------------------------------------------------------- //
#include "Stokes.H"
makeLaminarModel(Stokes);
#include "generalizedNewtonian.H"
makeLaminarModel(generalizedNewtonian);
#include "Maxwell.H"
makeLaminarModel(Maxwell);
|
makeLaminarModel
宏的定义如下:
1
2
3
4
|
// $FOAM_SRC/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.H
#define makeLaminarModel(Type) \
makeTemplatedLaminarModel \
(fluidThermoCompressibleTurbulenceModel, laminar, Type)
|
它调用了和不可压版本相同的宏命令 makeTemplatedLaminarModel
,不过传递的第一个参数不一样。
同样以 Stokes
为例,可压缩版本的 makeLaminarModel(Stokes)
宏展开会生成如下代码:
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
|
defineNamedTemplateTypeNameAndDebug
(Foam::laminarModels::Stokes<Foam::fluidThermoCompressibleTurbulenceModel>, 0);
namespace Foam
{
namespace laminarModels
{
typedef Stokes<fluidThermoCompressibleTurbulenceModel>
StokeslaminarfluidThermoCompressibleTurbulenceModel;
addToRunTimeSelectionTable
(
laminarfluidThermoCompressibleTurbulenceModel,
StokeslaminarfluidThermoCompressibleTurbulenceModel,
dictionary
);
}
}
|
关于 laminarModel
、RASModel
和 LESModel
三个类模板的框架设计就先介绍到这里,下篇文章我们将介绍雷诺平均 Navier-Stokes (RANS)方法和大涡模拟方法(LES)在该框架中的实现。
参考资料