6
7

Delete article

Deleted articles cannot be recovered.

Draft of this article would be also deleted.

Are you sure you want to delete this article?

More than 5 years have passed since last update.

OpenFOAMAdvent Calendar 2017

Day 23

OpenFOAMカスタマイズ入門:レイノルズ平均型乱流モデルのカスタマイズに挑戦

Last updated at Posted at 2017-12-16

OpenFOAMカスタマイズ入門:レイノルズ平均型乱流モデルのカスタマイズに挑戦

version.20171221-2011

0. 表記方法

端末から実行するコマンドは,次のように表記する。

command

ソースファイル等に記載する内容は,次のように表記する。

    void main { return 0; }

1. 概要

OpenFOAMのライブラリを改造する手順を学ぶ。RANS型乱流モデルのkOmegaSSTモデルに,旋回/曲率影響を考慮する機能を追加する。オープンCAEシンポジウム2017で中山勝之氏が発表したカスタマイズモデルを,ハンズオン形式で作成します。
http://www.opencae.or.jp/activity/symposium/opencae_symposium2017/

入門とはいいながら,手数は多くなりそうです。コード構成の改修が継続的に進められている乱流モデルのコードを確認します。

2. 環境

推奨環境:DEXCS2017 for OpenFOAM
http://mogura7.zenno.info/~et/wordpress/ocse/?p=2848

ハンズオンに必要な準備:OpenFOAM v1706がインストールされていること。OpenFOAMのソースコードがコンパイルできる環境であること。

対象となるOpenFOAMのバージョン:OpenFOAM v1706(v1612でも可)(OpenFOAM 4, 5 については,違いを説明する予定)

ハンズオンで想定する受講者:OpenFOAMの使用経験がある。何らかのプログラミング言語を使える。Linuxが使える。OpenFOAMのコードを読もうとトライしたことがある。講習内容が難しかったり,講師の説明が不十分だったとしても,それを許容できる広い心を備えている。

想定する受講者(聴講のみ):OpenFOAMに関心がある。講習内容が難しかったり,講師の説明が不十分だったとしても,それを許容できる広い心を備えている。

3. 基本方針

kOmegaSSTクラスをベースにする。これを継承して,必要な機能だけを追加する。

まずは,カスタムライブラリの作成手順を確認する。既存クラスとまったく同じものを,別名として作成し,コンパイルと実行が可能なことを確認する。

その後,kOmegaSST派生クラスを参考にしながら,補正モデルを作成する。

kOmegaSSTクラスを継承する方法は,kOmegaSSTLMクラスなどを参考にする。

kOmegaSSTRCH内での式の書き方などは,kOmegaSSTBaseクラスを参考にする。

4. 基礎式など

乱流モデルの呼称は,NASA Langley Research Center, Turbulence Modeling Resource を基準とする。
https://turbmodels.larc.nasa.gov/sst.html

OpenFOAM付属のkOmegaSSTモデルは,Menter SST Two-Equation Model from 2003 (SST-2003) である。

Image.png

Image [1].png

今回は,このモデルに簡略化した旋回/曲率補正モデルである Menter SST Two-Equation Model with Hellsten's Simplified Rotation/Curvature Correction (SST-RC-Hellsten) を作成する。このモデルでは,上記のomega式のdestruction項 Image [2].pngImage [3].png に置き換える。ここで,補正係数F4は下記の式から求める。

Image [4].png
Image [5].png
Image [6].png
Image [7].png
Image [8].png
Image [9].png
Image [10].png

5. オリジナルのコードを使ってコンパイルの練習

まず手始めに,OpenFOAMに標準で搭載されている kOmegaSST モデルをコピーして,名前だけを変更したモデル kOmegaSSTTest を作成してみる。コードの構成や,コンパイル方法を確認する。

5.1 そもそも,オリジナルのRANS乱流モデルはどこにある?

OpenFOAM-v1706/src/TurbulenceModels/
    ├── compressible/        //LIB = $(FOAM_LIBBIN)/libcompressibleTurbulenceModels
    ├── incompressible/        //LIB = $(FOAM_LIBBIN)/libincompressibleTurbulenceModels
    ├── phaseCompressible/
    ├── phaseIncompressible/
    ├── schemes/        //LIB = $(FOAM_LIBBIN)/libturbulenceModelSchemes
    └── turbulenceModels/        //LIB = $(FOAM_LIBBIN)/libturbulenceModels

kOmegaSST関係のモデルは下記のディレクトリに収納されている。
OpenFOAM-v1706/src/TurbulenceModels/turbulenceModels/

kOmegaSSTの基準となるクラス(kOmegaSSTBase)
Base/kOmegaSST/

kOmegaSSTクラス
RAS/kOmegaSST/

kOmegaSSTクラスからの派生クラス例
RAS/kOmegaSSTLM/

5.2 そもそも,オリジナルのRANS乱流モデルはどこからコンパイルされる?

圧縮性用乱流モデル(turbulentTransportModel)の場合,実体の作成は下記で行なわれる。
//OpenFOAM-v1706/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C

非圧縮性用乱流モデル(turbulentFluidThermoModel)の場合,実体の作成は下記で行なわれる。
//OpenFOAM-v1706/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C

圧縮性,非圧縮性の場合のどちらも,下記の様なマクロを使って乱流モデルが登録されている。

#include "kOmegaSST.H"
makeRASModel(kOmegaSST);

マクロは,圧縮性と非圧縮性とで異なる。圧縮性の場合には,下記に存在する。

makeBaseTurbulenceModelなどのマクロ

// /OpenFOAM-v1706 /src/TurbulenceModels/turbulenceModels/makeTurbulenceModel.H

makeRASModelなどのマクロ

//OpenFOAM-v1706 /src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.H

C++のマクロに関するメモ
        #define  マクロ名 実体(複数行可,\で繋ぐ)
        マクロ内の ## は,前後の文字列を連結する演算子

        typedef と #difine とでは,順番が逆になる

OpenFOAM 4系とOpenFOAM v1706系とで,マクロの構成が異なる。詳細については,本資料の付録に説明がある。

5.3 ソースコードディレクトリの作成

mkdir -p $WM_PROJECT_USER_DIR/src/myTurbulenceModels/kOmegaSSTTest

5.4 元となるコードのコピー

ソースコードディレクトリ $WM_PROJECT_USER_DIR/src/myTurbulenceModels/kOmegaSSTTest へ移動する。

cd $WM_PROJECT_USER_DIR/src/myTurbulenceModels/kOmegaSSTTest

オリジナルコードを kOmegaSSTTest ディレクトリへコピーする。

cp $WM_PROJECT_DIR/src/TurbulenceModels/turbulenceModels/RAS/kOmegaSST/kOmegaSST.* ./

ソースコードディレクトリ $WM_PROJECT_USER_DIR/src/myTurbulenceModels へ移動する。

cd $WM_PROJECT_USER_DIR/src/myTurbulenceModels/

オリジナルコードを myTurbulenceModels ディレクトリへコピーする。

cp $WM_PROJECT_DIR/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C ./

cp $WM_PROJECT_DIR/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.H ./

cp $WM_PROJECT_DIR/src/TurbulenceModels/turbulenceModels/makeTurbulenceModel.H ./

cp -r $WM_PROJECT_DIR/src/TurbulenceModels/incompressible/Make/ ./

5.5 ファイル名の変更

cd /$WM_PROJECT_USER_DIR/src/myTurbulenceModels/kOmegaSSTTest/

mv kOmegaSST.H kOmegaSSTTest.H

mv kOmegaSST.C kOmegaSSTTest.C

5.6 クラス名の変更

kOmegaSSTTest.H および kOmegaSSTTest.C ファイルにて,kOmegaSSTクラスを表す文字列”kOmegaSST”を"kOmegaSSTTest"に書き換える。

注意:kOmegaSSTBase の部分は,変更してはいけない。置換時に注意が必要である。

sed -i -e ' /kOmegaSSTBase/!s/kOmegaSST/kOmegaSSTTest/g' kOmegaSSTTest.H

sed -i -e ' s/kOmegaSSTTest_H/kOmegaSSTTEST_H/g' kOmegaSSTTest.H

sed -i -e ' /kOmegaSSTBase/!s/kOmegaSST/kOmegaSSTTest/g' kOmegaSSTTest.C

5.7 乱流モデル作成用makeTurbulenceModels.C の作成

$WM_PROJECT_USER_DIR/src/myTurbulenceModels に,makeTurbulenceModels.C というファイルを作成する。

touch $WM_PROJECT_USER_DIR/src/myTurbulenceModels/makeTurbulenceModels.C

makeTurbulenceModel.H のマクロをコピペする。不要部分を削除し,下記の内容とする。

#define makeTurbulenceModelTypes(Alpha, Rho, baseModel, BaseModel, Transport)  \
                                                                               \
    namespace Foam                                                             \
    {                                                                          \
        typedef BaseModel<Transport> Transport##BaseModel;                     \
        typedef RASModel<Transport##BaseModel> RAS##Transport##BaseModel;      \
    }

#define makeBaseTurbulenceModel(Alpha, Rho, baseModel, BaseModel, Transport)   \
                                                                               \
    namespace Foam                                                             \
    {                                                                          \
        typedef TurbulenceModel                                                \
        <                                                                      \
            Alpha,                                                             \
            Rho,                                                               \
            baseModel,                                                         \
            Transport                                                          \
        > Transport##baseModel;                                                \
                                                                               \
        defineTemplateRunTimeSelectionTable                                    \
        (                                                                      \
            Transport##baseModel,                                              \
            dictionary                                                         \
        );                                                                     \
                                                                               \
                                                                               \
        defineNamedTemplateTypeNameAndDebug(RAS##Transport##BaseModel, 0);     \
                                                                               \
        defineTemplateRunTimeSelectionTable                                    \
        (RAS##Transport##BaseModel, dictionary);                               \
                                                                               \
                                                                               \
    }

#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                                                     \
            );                                                                 \
        }                                                                      \
    }

makeTurbulenceModels.Cに追記していく。

turbulentTransportModels.HのmakeTurbulenceModelTypesとmakeRASModel(Type) ,turbulentTransportModels.CのmakeBaseTurbulenceModelをコピペする。
最後に下記を書いて,kOmegaSSTTestモデルの作成と登録を実施する。

#include "kOmegaSSTTest.H"
makeRASModel(kOmegaSSTTest);

なお,makeTurbulenceModels.Cの冒頭では,他のクラスが利用するいくつかの定義を読み込んでおく必要がある。下記の内容をファイルの冒頭に記載する。

#include "IncompressibleTurbulenceModel.H"
#include "incompressible/transportModel/transportModel.H"
#include "addToRunTimeSelectionTable.H"

#include "laminarModel.H"
#include "RASModel.H"

#include "fvOptions.H"

5.8 Makeディレクトリの修正

filesファイルの内容を下記に修正する。

makeTurbulenceModels.C
LIB = $(FOAM_USER_LIBBIN)/libmyTurbulenceModels

optionsファイルの内容を下記に修正する。

EXE_INC = \
    -I$(LIB_SRC)/finiteVolume/lnInclude \
    -I$(LIB_SRC)/meshTools/lnInclude \
    -I$(LIB_SRC)/fvOptions/lnInclude \
    -I$(LIB_SRC)/transportModels \
    -I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
    -I$(LIB_SRC)/TurbulenceModels/incompressible/lnInclude

LIB_LIBS = \
    -lincompressibleTurbulenceModels \
    -lturbulenceModels \
    -lfvOptions \
    -lfiniteVolume \
    -lmeshTools

5.9 コンパイル

プログラムがコンパイルでき,動作することを確認する。

cd $WM_PROJECT_USER_DIR/src/myTurbulenceModels/

wmake

5.10 テスト用例題ケースの作成

チュートリアルディレクトリから,simpleFoam/pitzDaily/例題を作業用ディレクトリにコピーする。

cp -r $WM_PROJECT_DIR/tutorials/incompressible/simpleFoam/pitzDaily/ $WM_PROJECT_USER_DIR/run/

例題ディレクトリにAllcleanスクリプトファイルを作成する。ファイルの内容は下記とする。ファイルに実行可能権限を与える。

#!/bin/sh
cd ${0%/*} || exit 1    # Run from this directory

# Source tutorial clean functions
. $WM_PROJECT_DIR/bin/tools/CleanFunctions

cleanCase

例題ディレクトリにAllrunスクリプトファイルを作成する。ファイルの内容は下記とする。ファイルに実行可能権限を与える。

#!/bin/sh
cd ${0%/*} || exit 1    # Run from this directory

# Source tutorial run functions
. $WM_PROJECT_DIR/bin/tools/RunFunctions

application=`getApplication`

runApplication blockMesh
runApplication $application

カスタマイズした乱流モデルを利用可能とするために,system/controlDictにライブラリを追加する。controlDict の最下部に下記を追記する。

libs
(
    "libmyTurbulenceModels.so"
);

5.11 テスト用例題ケースの実行

選択できる乱流モデル一覧を表示させるために,constant/turbulencePropertiesファイル中のRASModel を kEpsilon から BANANA に変更する。

sed -i -e s/"(RASModel[ \t]*) kEpsilon;"/"\1 BANANA;"/g constant/turbulenceProperties

./Allrun スクリプトを実行する。log.simpleFoamファイルを開き,選択できる乱流モデル一覧の部分に,kOmegaSSTTestが表示されることを確認する。

kOmegaSSTTestモデルを使って実行する。先ほどと同様にして,BANANAモデルをkOmegaSSTTestモデルに修正する。

sed -i -e s/"(RASModel[ \t]*) BANANA;"/"\1 kOmegaSSTTest;"/g constant/turbulenceProperties

./Allrun スクリプトを実行する。log.simpleFoamファイルを開き,エラーなく実行できていることを確認する。

6. RCHモデルの作成

コンパイル方法などが分かったので,ようやく,カスタマイズした乱流モデルの作成に挑戦する。今回作成するモデルは,kOmegaSSTをベースとして,一部の項を修正するものである。そのため,標準のkOmegaSSTクラスを継承し,異なる部分のみのコードを書くことにする。

標準の乱流モデルの中には,kOmegaSSTクラスを継承して作成されているモデル(例えば,kOmegaSSTLM など)が存在する。そこで,kOmegaSSTLM のソースコードをコピーして,必要な変更を加えることにする。

cd $WM_PROJECT_USER_DIR/src/myTurbulenceModels/

cp -p -r $WM_PROJECT_DIR/src/TurbulenceModels/turbulenceModels/RAS/kOmegaSSTLM/ kOmegaSSTRCH/

6.1 ファイル名の変更

cd $WM_PROJECT_USER_DIR/src/myTurbulenceModels/kOmegaSSTRCH/

mv kOmegaSSTLM.H kOmegaSSTRCH.H

mv kOmegaSSTLM.C kOmegaSSTRCH.C

6.2 クラス名の変更

kOmegaSSTRCH.H および kOmegaSSTRCH.C ファイルにて,kOmegaSSTLMクラスを表す文字列”kOmegaSSTLM”を"kOmegaSSTRCH"に書き換える。

sed -i 's/kOmegaSSTLM/kOmegaSSTRCH/g' kOmegaSSTRCH.H

sed -i 's/kOmegaSSTLM/kOmegaSSTRCH/g' kOmegaSSTRCH.C

6.3 コードの修正:宣言ファイル kOmegaSSTRCH.H

元になるkOmegaSSTLMクラスのソースファイルには,不要な部分が多い。116行目付近からの protected: の内容を削除し,下記の宣言だけを加える。

追加するもの:ヘッダファイル

  • dimensionedScalar cRC_;
  • volScalarField S_;
  • volScalarField W_;
  • volScalarField Ri_;
  • tmp<volScalarField> F4() const;

さらに,170行付近の // Member Functions では,下記だけを残すようにして,他は削除する。

  • virtual bool read();
  • virtual void correct();

6.4 コードの修正:定義ファイル kOmegaSSTRCH.C

先ほどの宣言ファイルに記載した関数 F4( ),read( ), correct( ) の定義を記述する。さらに,変数の初期化を行なうようにコンストラクタを修正する。

追加するもの:定義ファイル

  • F4( )関数
  • コンストラクタでの初期化
  • read( )関数
  • correct( )関数

6.4.1 Private Member Functions:関数 F4( )

kOmegaSSTRCH.Cファイルの中で,F1関数の定義部分を下記に書き換える。その他の Private Member Functions は削除する。

// * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //

template<class BasicTurbulenceModel>
tmp<volScalarField>
kOmegaSSTRCH<BasicTurbulenceModel>::kOmegaSSTRCH::F4() const
{
    return 1 / (1 + cRC_*Ri_);
}

6.4.2 コンストラクタ

コンストラクトの初期化部分では,まず始めに,継承元のkOmegaSSTクラスの初期化が行なわれる。この中で,引数にtypeを与えるように書き加える。これによって,係数ファイルの名前がここで与えるtype(kOmegaSSTRCH)をベースに作成されることになる。

さらに,cRC_, S_, W_, Ri_ の各変数を初期化部に書き入れる。

まとめると,Constructors の初期化部は下記となる。

    kOmegaSST<BasicTurbulenceModel>
    (
        alpha,
        rho,
        U,
        alphaRhoPhi,
        phi,
        transport,
        propertiesName,
        type
    ),
    cRC_
    (
        dimensioned<scalar>::lookupOrAddToDict
        (
            "cRC",
            this->coeffDict_,
            1.4
        )
    ),
    S_
    (
        IOobject
        (
            "S",
            this->runTime_.timeName(),
            this->mesh_,
            IOobject::NO_READ, //READ_IF_PRESENT
            IOobject::AUTO_WRITE
        ),
        this->mesh_,
        dimensionedScalar("S",dimless/dimTime,0.)
    ),
    W_
    (
        IOobject
        (
            "W",
            this->runTime_.timeName(),
            this->mesh_,
            IOobject::NO_READ, //READ_IF_PRESENT,
            IOobject::AUTO_WRITE
        ),
        this->mesh_,
        dimensionedScalar("W",dimless/dimTime,0.)
    ),
    Ri_
    (
        IOobject
        (
            "Ri", 
            this->runTime_.timeName(),
            this->mesh_,
            IOobject::NO_READ, //READ_IF_PRESENT,
            IOobject::AUTO_WRITE
        ),
        this->mesh_,  //scalar(0.)
        dimensionedScalar("Ri",dimless,scalar(0.0)) 
    )

さらに,Constructors の実行部分( {}の内部)に,係数を出力するコードを記載する。

    if (type == typeName)
    {
        this->printCoeffs(type);
    }

6.4.3 Member Functions:関数 read( ), correct( )

read関数とcorrect関数以外の関数の定義を削除する。

設定ファイルからcRcの値を読み込むために,read関数の内容を次のように書き換える。

    if (kOmegaSST<BasicTurbulenceModel>::read())
    {
        cRC_.readIfPresent(this->coeffDict());

        return true;
    }
    else
    {
        return false;
    }

F4関数を組み込んだomega式を解くために,correct関数を変更する。

correct関数で解くomegaの式は,kOmegaSSTの基本の式とほとんど同じである。そのため,すでにコードが作成されているものをコピーすることにする。そのために,src/TurbulenceModels/turbulenceModels/Base/kOmegaSST/kOmegaSSTBase.C ファイルを開く。このファイルのcorrect関数の中身をコピーして,kOmegaSSTRCH.C の定義として貼付ける。

fv::optionを定義している行の下辺りに,kOmegaSSTBaseで定義されているkやomegaを使うためのLocal referencesを定義する。

    volScalarField& k_ = this->k_;
    volScalarField& omega_ = this->omega_;
    const dimensionedScalar& alphaOmega2_ = this->alphaOmega2_;

correct関数をBasicEddyViscosityModelではなく,BasicTurbulenceModelのものを使うように書き換える。

    BasicTurbulenceModel::correct();

G の定義の下辺りに,W_, S_, Ri_ を求める式を記載する。

    W_ = sqrt(2*magSqr(skew(tgradU())));
    S_ = sqrt(S2);
    Ri_ = W_/(S_+dimensionedScalar ("VSMALL", dimless/dimTime, VSMALL))
         *(W_/(S_+dimensionedScalar ("VSMALL", dimless/dimTime, VSMALL))-1);

フィールド値としてF23を定義している下に,volScalarFieldとしてF4を定義し,初期値をF4関数から与える。

    volScalarField F4(this->F4());

omegaEqnの右辺第3項に,F4をかける。

          - fvm::Sp(F4*alpha()*rho()*beta*omega_(), omega_)    // F4 is included.

DomegaEff, GbyNu 関数は,kOmegaSSTBaseで定義されている。これを利用するためにthisポインタを追加する。

同じくomegaEqnの右辺において,kOmegaSSTBaseで定義されている関数QsasおよびomegaSoureceを利用するためのthis->を追加する。

          + this->Qsas(S2(), gamma, beta)
          + this->omegaSource()

kEqnは,Baseと同じである。ただし,DkEff, Pk, epsilonByk, kSource といった関数については,thisポインタを追加する。

さらに,correctNut関数にもthisポインタを追加する。

6.5 コードの修正: makeTurbulenceModel.C

makeTurbulenceModel.C に下記を追加して,ライブラリに新しいモデルを登録する。

#include "kOmegaSSTRCH.H"
makeRASModel(kOmegaSSTRCH);

6.6 コンパイル

プログラムがコンパイルでき,動作することを確認する。

cd $WM_PROJECT_USER_DIR/src/myTurbulenceModels/

wmake

これによって,先に作成した kOmegaSSTTest と kOmegaSSTRCH の2つのモデルが,RASモデルとして RCturbulenceModelsライブラリに含まれることになる。

例題で,どちらのモデルも使用できることを確認する。

7. 付録

7.1 OpenFOAM バージョン リリース 乱流モデルの遷移

Foundation系

26th July 2017
OpenFOAM 5.0 Released

13th October 2016
OpenFOAM 4.1 Released

17 Sep 2016
DEXCS2016
Build: 4.x-be7fba6cff9b

28th June 2016
OpenFOAM 4.0 Released

2016-11-05 23:30
kOmegaSSTBase: Corrected read() to re-read the base-class settings
https://bugs.openfoam.org/view.php?id=2318
Resolved in OpenFOAM-dev by commit c3fdc191c2e648ba210902bf848c02614af728a1
Resolved in OpenFOAM-4.x by commit c4b8d5e443e91ddec6766d992d830bcbaba60af6

16 Feb 2015
turbulenceModels/RAS/kOmegaSSTSAS/kOmegaSSTSAS: Added the k-omega-SST-SAS model
OpenFOAM-3.0.x
https://github.com/OpenFOAM/OpenFOAM-3.0.x/commit/f99884de0052105fb9b18a96185139ccc1160de6
Qsas追加

ESI-OpenCFD Inc. 系

Jul 28, 2016
https://develop.openfoam.com/Development/OpenFOAM-plus/commit/02dd85167e1ebb56fdc435e169c3d0a5645b641e
TurbulenceModels::kOmegaSST.*: Updated source-terms and associated functions to use volScalarField::Internal
This is more efficient, avoids divide-by-0 when evaluating unnecessary
boundary values and avoids unnecessary communications when running in parallel.

7.2 乱流モデル作成用マクロ等のバージョンによる相違

OpenFOAM v1706

OpenFOAM v1706 では,次のマクロが使われている。

  • makeTurbulenceModelTypes(OpenFOAM 4には存在しない。makeBaseTurbulenceModelに含まれる。)
  • makeBaseTurbulenceModel
  • makeTemplatedTurbulenceModel
  • makeRASModel

これらに関するソースコードは次の通りとなる。

マクロ(makeTurbulenceModelTypes, makeBaseTurbulenceModel,および makeTemplatedTurbulenceModel)の定義:$WM_PROJECT_DIR/src/TurbulenceModels/turbulenceModels/makeTurbulenceModel.H

makeTurbulenceModelTypesマクロの呼びだし と makeRASModel(Type)マクロの定義:$WM_PROJECT_DIR/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.H

makeBaseTurbulenceModelマクロの呼びだし 各乱流モデルの呼びだし:$WM_PROJECT_DIR/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C

OpenFOAM 4

OpenFOAM 4 では,次のマクロが使われている。

  • makeBaseTurbulenceModel
  • makeTemplatedTurbulenceModel
  • makeRASModel

これらに関するソースコードは次の通りとなる。

マクロ(makeBaseTurbulenceModel,および makeTemplatedTurbulenceModel)の定義:$WM_PROJECT_DIR/src/TurbulenceModels/turbulenceModels/makeTurbulenceModel.H

makeBaseTurbulenceModelマクロの呼びだし と makeRASModel(Type)マクロの定義と 各乱流モデルの呼びだし[of4では,makeBaseTurbulenceModelの中で,makeTurbulenceModelTypesと同様の作業(typedef RASModel )が行なわれている。]:$WM_PROJECT_DIR/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C

7.3 チェック

constant/turbulencePropertiesファイルに下記を追加して,係数を変更できるか試す。

    kOmegaSSTRCHCoeffs 
    {
        cRC             1.2; // for test Only
    }
6
7
0

Register as a new user and use Qiita more conveniently

  1. You get articles that match your needs
  2. You can efficiently read back useful information
  3. You can use dark theme
What you can do with signing up
6
7

Delete article

Deleted articles cannot be recovered.

Draft of this article would be also deleted.

Are you sure you want to delete this article?