吾生有涯 学海无涯
析模有界 知识无界

OpenFOAM编程案例|10 标量方程求解

本案例演示利用OpenFOAM创建一个标量输运方程求解器。

本案例要求解的方程为:

没有瞬态项和源项。

1 文件框架

利用foamNewApp快速创建文件框架。

run
foamNewApp demo10 && cd demo10

执行完毕后的文件框架如下图所示。

本案例比较简单,只需要一个源文件demo10.C即可。

2 源代码

源文件demo10.C的内容如下。

#include "fvCFD.H"
 
int main(int argc, char *argv[])
{
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
 
    Info << "读取物理场beta" << endl;
    volScalarField beta
    (
        IOobject
        (
            "beta",
            runTime.timeName(),
            mesh,
            IOobject::MUST_READ,
            IOobject::AUTO_WRITE
        ),
        mesh
    )
;
 
    Info << "读取物理场U" << endl;
    volVectorField U
    (
        IOobject
        (
            "U",
            runTime.timeName(),
            mesh,
            IOobject::MUST_READ,
            IOobject::AUTO_WRITE
        ),
        mesh
    )
;
 
    // 读取字典文件
    Info << "读取字典文件transportProperties" << endl;
 
    IOdictionary transportProperties
    (
        IOobject
        (
            "transportProperties",
            runTime.constant(), //- 此文件在constant文件夹下面
            mesh,
            IOobject::MUST_READ_IF_MODIFIED,
            IOobject::NO_WRITE
        )
    )
;
 
    Info << "从字典文件中读取扩散系数gamma" << endl;
    dimensionedScalar gamma("gamma", dimViscosity, transportProperties);
 
    Info << "读取/计算面通量场phi" << endl;
    surfaceScalarField phi
    (
        IOobject
        (
            "phi",
            runTime.timeName(),
            mesh,
            IOobject::READ_IF_PRESENT, //若phi文件存在,则读取,否则采用速度与面积相乘计算得到
            IOobject::AUTO_WRITE       // 保存该变量
        ),
        fvc::interpolate(U) & mesh.Sf() //速度向量与面积向量点积
    )
;
 
    //求解控制方程
    solve(fvm::div(phi, beta) - fvm::laplacian(gamma, beta));
 
    // 创建一个标量场用于存储结果,也可以直接用beta
    volScalarField result
    (
        IOobject
        (
            "result",
            runTime.timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::AUTO_WRITE
        ),
        beta     //将beta的值拷贝给result
    )
;
 
    result.write();
 
    Info << nl << "计算完毕" << endl;
    runTime.printExecutionTime(Info);
    return 0;
}

编译并测试运行,如下图所示。

3 测试计算

案例文件需要根据求解器的需求来准备。上面准备的求解器中需要求解物理量beta。

  • constant/transportProperties文件中需要指定扩散系数gamma
FoamFile
{
    version     2.0;
    format      ascii;
    class       dictionary;
    location    "constant";
    object      transportProperties;
}
// * * * * * * * * * * * * * * * * * //
gamma            0.01;

由于在读取字典文件时指定了量纲,因此在字典文件中无需重复指定量纲。

  • system/fvSolution文件中需要指定beta的求解方法
FoamFile
{
    version     2.0;
    format      ascii;
    class       dictionary;
    location    "system";
    object      fvSolution;
}
// * * * * * * * * * * * * * * * * * //
solvers
{
    beta
    {
        solver        GAMG;
        smoother    DILUGaussSeidel;
        tolerance    1e-6;
        relTol        0;
    }
}
  • system/fvSchemes文件需要指定控制方程中对流项、梯度项以及laplacian项的离散格式
FoamFile
{
    version     2.0;
    format      ascii;
    class       dictionary;
    location    "system";
    object      fvSchemes;
}
// * * * * * * * * * * * * * * //
 
ddtSchemes
{
    default none;
}
// 梯度项
gradSchemes
{
    default none;
    grad(beta) Gauss linear;
}
//散度项
divSchemes
{
    default none;
    div(phi,beta) bounded Gauss upwind;
}
// laplacian项
laplacianSchemes
{
    default none;
    laplacian(gamma,beta) Gauss linear corrected;
}
// 通量插值
interpolationSchemes
{
    default none;
    interpolate(U) linear;
}
 
snGradSchemes
{
    default none;
}
 
  • system/controlDict按常规设置
FoamFile
{
    version     2.0;
    format      ascii;
    class       dictionary;
    location    "system";
    object      controlDict;
}
// * * * * * * * * * * * * * * * * //
startFrom       startTime;
startTime       0;
stopAt          endTime;
endTime         1;
deltaT          1;
writeControl    timeStep;
writeInterval   1;
purgeWrite      0;
writeFormat     ascii;
writePrecision  6;
writeCompression off;
timeFormat      general;
timePrecision   6;
runTimeModifiable true;
  • 0/beta文件指定物理场beta的边界条件与初始条件
FoamFile
{
    version     2.0;
    format      ascii;
    class       volScalarField;
    object      beta;
}
// 注意beta是标量
dimensions      [0 0 0 0 0 0 0];
internalField   uniform 0;
boundaryField
{
    left
    {
        type            fixedValue;
        value           uniform 1;
    }
 
    lower
    {
        type            fixedValue;
        value           uniform 0;
    }
 
    "(upper|right)"
    {
        type zeroGradient;
    }
 
    frontAndBack
    {
        type            empty;
    }
}
  • 0/U文件
FoamFile
{
    version     2.0;
    format      ascii;
    class       volVectorField;
    object      U;
}
// * * * * * * * * * * * * * //
 
dimensions      [0 1 -1 0 0 0 0];
 
internalField   uniform (1 1 0);
 
boundaryField
{
    "(lower|upper|left|right)"
    {
        type            fixedValue;
        value           $internalField;
    }
 
    frontAndBack
    {
        type            empty;
    }
}

在案例根路径下执行下面的命令进行计算。

blockMesh
demo10

很快就计算完了。

计算结果如下图所示。


(完毕)

本篇文章来源于微信公众号: CFD之道

赞(0) 打赏
版权声明:未经允许,请勿随意用于商业用途。
文章名称:《OpenFOAM编程案例|10 标量方程求解》
文章链接:https://www.topcfd.cn/18617/
本站资源仅供个人学习交流,请于下载后24小时内删除,不允许用于商业用途,否则法律问题自行承担。
分享到

说两句 抢沙发

评论前必须登录!

 

觉得文章有用就打赏一下文章作者吧

非常感谢你的打赏,我们将继续给力更多优质内容,让我们一起创建更加美好的网络世界!

支付宝扫一扫

微信扫一扫

登录

找回密码

注册