简介
OpenFoam可以用于处理场类型的数据,但是如果采用“求解器”手脚架(见OpenFOAM编程基本结构)将会引入"setRootCase.H"、“createTime.H”、"createMesh.H"三个代码片段文件,引入了严重的“副作用”,难以将代码模块化。如果同时在多个线程上并行计算往往会出现严重错误。
本文将编写极简的OpenFoam源文件,剔除不必要的“副作用”,为并行计算提供条件
代码
#include "fvCFD.H"
int main()
{
//第一步:构造runTime对象,用于控制时间
Foam::Time runTime(
fileName("."), //可执行程序的根目录
fileName("./pitzDaily")); //案例文件的相对目录
//第二步:构造mesh对象,用于读取网格
const Foam::fvMesh mesh(
Foam::IOobject(
Foam::fvMesh::defaultRegion,
runTime.timeName(),
runTime,
Foam::IOobject::MUST_READ));
//第三步:开始计算
for (int i = 0; i < 20; i++) //以i为自变量步进
{
//初始化自定义场(均匀场)
volScalarField IndexOfCell(
IOobject(
"IndexOfCell",
runTime.timeName(),
mesh,
IOobject::NO_READ, //不从文件中读取,而是直接在代码中生成
IOobject::AUTO_WRITE),
mesh,
dimensionedScalar("tmp", dimless, 0), //初始化为无量纲的0场
"zeroGradient"); //初始化为自然边界条件(零通量边界条件)
//用户自定义的计算,这里给每个单元的数值赋值为单元序号
forAll(IndexOfCell, cellI)
{
IndexOfCell[cellI] = cellI;
}
IndexOfCell.correctBoundaryConditions(); //更新边界处的结果
//输出场
runTime.setTime(i, i); //设置时间
runTime.writeTimeDict(); //输出时间
IndexOfCell.write(); //输出场
}
return 0;
}
代码的核心思想是,对于固定网格问题,我们只需要完成一次runTime
初始化和mesh
初始化,之后在每一步计算中生成我们需要的场,使用各种方法加工这个场,再将结果输出为文件。
本代码使用官方pitzDaily案例,计算完成后,案例文件目录如下
.
├── 0
│ ├── IndexOfCell
│ ├── T
│ ├── U
│ └── uniform
│ └── time
├── 1
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 10
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 11
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 12
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 13
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 14
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 15
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 16
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 17
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 18
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 19
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 2
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 3
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 4
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 5
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 6
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 7
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 8
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── 9
│ ├── IndexOfCell
│ └── uniform
│ └── time
├── constant
│ ├── polyMesh
│ │ ├── boundary
│ │ ├── faces
│ │ ├── neighbour
│ │ ├── owner
│ │ └── points
│ └── transportProperties
└── system
├── blockMeshDict
├── controlDict
├── fvSchemes
└── fvSolution
可以看出我们已经将IndexOfCell场
输出到了对应的时间文件夹,使用paraview实现可视化,如下图
该图可视化地显示了每个单元的标号