OpenFOAMV9中的isoSurface类可以用来提取等值面。
该类的实例化方式为:
sampledSurfaces::isoSurface isosurf = sampledSurfaces::isoSurface(
"isoSurface",
mesh,
isoSurfaceDict);
"isoSurface"
是一个自定义的名字(一般取对象名),mesh
是所研究问题的网格,isoSurfaceDict
是一个数据字典,该数据字典的内容如下
/*--------------------------------*- C++ -*----------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Version: 9
\\/ M anipulation |
-------------------------------------------------------------------------------
Description
Writes out iso-surface files with interpolated field data in VTK format.
\*---------------------------------------------------------------------------*/
FoamFile
{
format ascii;
class dictionary;
location "system";
object isoSurfaceDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
type isoSurface;
isoField p; //需要提取等值线的场,这里是压力场
isoValue 10; //等值线的数值,这里三10Pa
filter full;
interpolate yes;
// ************************************************************************* //
有了isosurf
这个对象之后,可以用以下一行代码完成等值面提取的工作
isosurf.sample(p);
提取完成等值面后,可以获得等值面的面单元向量和等值面的顶点坐标
//等值面的面单元向量
faceList faces = isosurf.faces();
//等值面的顶点坐标
pointField points = isosurf.points();
points描述了该等值面是由哪些点构成的,提供了这些点的坐标信息。faces描述了这些点的连接关系。如果要访问第faceI
个单元面的面积,可以使用以下代码:
mag(faces[faceI].area(points))
为了方便可视化,还可以将得到的等值面输出为VTK文件
vtkSurfaceWriter vtkWriter = vtkSurfaceWriter(IOstream::streamFormat::ASCII);
vtkWriter.write("postProcess",
"someContours",
points,
faces);
以教程案例pitzDaily为例子,将其复制到自己的文件夹。使用simpleFoam求解器完成求解,结果如下:
在simpleFoam求解器源代码的基础上添加有关等值面的内容,改进后的求解器main函数如下
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see .
Application
isoSimpleFoam
Description
Steady-state solver for incompressible, turbulent flow, using the SIMPLE
algorithm.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "singlePhaseTransportModel.H"
#include "kinematicMomentumTransportModel.H"
#include "simpleControl.H"
#include "pressureReference.H"
#include "fvModels.H"
#include "fvConstraints.H"
//有关等值面的头文件
#include "sampledIsoSurface.H"
#include "vtkSurfaceWriter.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "postProcess.H"
#include "setRootCaseLists.H"
#include "createTime.H"
#include "createMesh.H"
#include "createControl.H"
#include "createFields.H"
#include "initContinuityErrs.H"
//读取等值线字典
dictionary isoSurfaceDict = IOdictionary(IOobject(
"isoSurfaceDict",
mesh.time().system(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE));
turbulence->validate();
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info << "\nStarting time loop\n"
<< endl;
while (simple.loop(runTime))
{
Info << "Time = " << runTime.timeName() << nl << endl;
fvModels.correct();
// --- Pressure-velocity SIMPLE corrector
{
#include "UEqn.H"
#include "pEqn.H"
}
laminarTransport.correct();
turbulence->correct();
//修改等值线字典中的isoValue数值
isoSurfaceDict.set("isoValue", 7);
//实例化isoSurface对象
sampledSurfaces::isoSurface isosurf = sampledSurfaces::isoSurface(
"isoSurface",
mesh,
isoSurfaceDict);
//提取等值面
isosurf.sample(p);
//等值面的面单元向量
faceList faces = isosurf.faces();
//等值面的顶点坐标
pointField points = isosurf.points();
//计算等值面的面积
scalar area = 0;
forAll(faces, faceI)
{
area += mag(faces[faceI].area(points));
}
Info << "面积:" << area << endl;
//将等值面输出为VTK文件
vtkSurfaceWriter vtkWriter = vtkSurfaceWriter(IOstream::streamFormat::ASCII);
vtkWriter.write("postProcess",
"someContours",
points,
faces);
runTime.write();
Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info << "End\n"
<< endl;
return 0;
}
// ************************************************************************* //
编译的过程见OpenFOAM用户手册,这里仅介绍基本步骤。
新建Make文件夹,在Make文件夹中新建files文件,写入编译信息
isoSimpleFoam.C
EXE = $(FOAM_USER_APPBIN)/isoSimpleFoam
在Make文件夹中新建options文件,写入依赖信息
EXE_INC = \
-I$(LIB_SRC)/MomentumTransportModels/momentumTransportModels/lnInclude \
-I$(LIB_SRC)/MomentumTransportModels/incompressible/lnInclude \
-I$(LIB_SRC)/transportModels/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude \
-I$(LIB_SRC)/sampling/lnInclude \
-I$(LIB_SRC)/surfMesh/lnInclude
EXE_LIBS = \
-lmomentumTransportModels \
-lincompressibleMomentumTransportModels \
-ltransportModels \
-lfiniteVolume \
-lmeshTools \
-lfvModels \
-lfvConstraints \
-lsampling\
-lsurfMesh
使用wmake
命令完成编译,得到isoSimpleFoam
求解器。该求解器将计算pitzDaily案例中压力为7Pa的等值面面积,并输出为VTK文件。
计算完成后,使用paraview读取postProcess
文件夹下的someContours.vtk
文件,可得到下图
这就是提取得到的等值面,其面积为9.043e-05