【OpenFOAM】手法を変えて処理を追加⑤ : 設定項目を実行時に読み込む (IOdictionary)
オープンCAEアドベントカレンダー 10日目のエントリーです。
はじめに
『【OpenFOAM】手法を変えて処理を追加④ : 継承クラスに処理を記述する』の続きです。OpenFOAMの計算に手法を変えて処理を追加していきます。
前回は継承クラスに処理を記述してみました。色々な処理を行いたい場合には継承を使って類似したクラスをシンプルに作成することができました。ただその処理の内容を微修正したい場合などには毎回コンパイルをする必要があります。
今回はOpenFOAMのIOdictrionary
というクラスも継承することにより、計算ケースで実行時に読み込むファイルの内容を処理に反映させてみたいと思います。
手法一覧
この記事は様々な手法を紹介するシリーズの一つであり、全体としては以下の方法で行う予定です。実行確認済みのファイルはgithub.com/inabower/getScoreTutorialsで公開しています。
- 計算ファイルに処理を書き込む (①)(キーワード:
codedFunctionObject
) - ソルバーに処理を追加する (②)
- 関数やクラスとして処理を追加する(③)
- 継承クラスに処理を追加(④)
- 設定項目を実行時に読み込む ★今ここ (キーワード:
IOdictionary
) - 継承クラスを実行時に選択する(キーワード:
runTimeSelectionTable
) - クラスを別でコンパイルする(キーワード:
wmake libso
) - 継承クラスの一つとしてcodedを追加(キーワード:
codedTemplete
)
動作環境
- OS : Ubuntu 20.04 LTS (WSL2 on Windows 10)
- OpenFOAM : ESI版 v2006 (from source)
- gcc : 9.3.0
実行時に設定を読み込む
今回は実行時のケースに設定を書き込みたい例を考えます。例えば前回までに記述した処理には以下のような分岐がありました。
if (!(patch.name()=="inlet" || patch.name()=="outlet1" || patch.name()=="outlet2")) continue;
これは、『境界の名前がinlet
かoutlet1
かoutlet2
か出ない場合はループを抜ける』という意味ですが、この境界の名前については実行するメッシュによって異なることが一般的です。前回までの作成方法では、この部分をメッシュによって変更してコンパイルする必要がありました。しかし、この境界の名前の情報を実行時に指定してあげることができれば毎回コンパイルを行う必要はなくなります。
今回はこの『境界の名前』を記述した『scoreDict』をcase/constant/
に配置して読み込む方法を考えます。
そのために、前回基底クラスとして作成したgetScore
クラスをOpenFOAMの標準ライブラリであるIOdictionary
クラスから継承するという方法を考えてみます。これにより、getScore
クラスおよびその継承クラスはIOdictionary
クラスとしての機能を使えるようになります。
IOdictionary
IOdictionary
はOpenFOAMのソルバー実行時にケースの中のファイルを読み込む関数です。jsonに似た以下のようなFoam::dictionary
形式のファイルを読み込むことができます。例えばconstant/sampleDict
という以下のようなファイルがあるとします。
- constant/sampleDict
key1 10;
key2 "Hello world";
dict3
{
key4 100;
key5 ( 0 0 0 );
dict6
{
key7 1000;
}
}
この状態でソルバーの中では以下のようにオブジェクトを作成します。
IOdictionary sampleDictObj
(
IOobject
(
"sampleDict",
mesh.time().constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::AUTO_WRITE
)
);
こうすることで、例えば以下のようにconstant/sampleDict
の中身を呼び出すことができます・
- ソルバーの中
scalar val1 = sampleDictObj.get<scalar>("key1"); // 10.0
word val2 = sampleDictObj.get<word>("key2"); // "Hello world"
label val4 = sampleDictObj.subDict("dict3").get<label>("key4"); // 1000
vector val4 = sampleDictObj.subDict("dict3").get<vector>("key5"); // vector(0.0, 0.0, 0.0)
dictionary& dict7 = sampleDictObj.subDict("dict3").subDict("dict7");
scalar val7 = dict7.get<scalar>("key7"); // 1000.0
さらに詳しく動作を知りたい場合は、以下のopenfoam.com/documentationでコードの構成などを確認することができます。
実装してみた
完成済みのコードは(github.com/inabower/getScoreTutorials/tree/main/005.IOdictionary)に公開しています。前回作成した継承クラスをアレンジする形で作成しました。
ディレクトリ構成
ソルバーとケースを含んだディレクトリ構成は以下のようになっています。
├── plainCase
│ ├── Allclean
│ ├── Allrun
│ ├── 0.orig
│ ├── constant
│ │ ├── scoreDict
│ │ ├── transportProperties
│ │ └── turbulenceProperties
│ └── system
│ ├── blockMeshDict
│ ├── controlDict
│ ├── fvSchemes
│ └── fvSolution
└── scoreDictPimpleFoam
├── Make
│ ├── files
│ └── options
├── UEqn.H
├── correctPhi.H
├── cpCorrectPhi.H
├── createFields.H
├── getScore
│ ├── getScore.C
│ ├── getScore.H
│ ├── sumInletOutletScore.C
│ └── sumInletOutletScore.H
├── pEqn.H
└── pimpleFoam.C
読み込むファイルの中身
読み込むためのscoreDict
は以下のようになっています。
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "constant";
object scoreDict;
}
targetPatches (inlet outlet1 outlet2);
message Hello_world;
この中のtargetPatches
が今回扱いたい境界の名前とします。ケースを実行する際にここを読み込むこととします。その下のmessage
は一例として、この文章を出力の際に表示させるようにしてみます。FoamFile
はOpenFOAMのファイルであることを示すヘッダです。中身は重要ではないので他のファイルからのコピペでOKです。
IOdictionaryの継承
まず基底クラスであるgetScore
でIOdictionary
を継承してみます。
getScore.H
まずgetScore.H
は以下のようにしてみました。
#ifndef getScore_H
#define getScore_H
#include "fvMesh.H"
#include "IOdictionary.H"
namespace Foam
{
class getScore
:
public IOdictionary
{
protected:
// Protected data
const fvMesh& mesh_;
scalar score_;
public:
// Constructors
getScore(const fvMesh& mesh);
getScore(const fvMesh& mesh, const IOdictionary& dict);
//- Destructor
virtual ~getScore()
{}
// Member Functions
virtual void calculate();
const scalar& value() const
{
return score_;
}
};
}
#endif
このうち以下の部分でIOdictionsaryを継承していることが分かります。
class getScore
:
public IOdictionary
また、以下のようにコンストラクタを二つ作ってみました。これにより二通りの方法でこのクラスを呼び出せるようになります。詳しくは次の定義のところで説明します。
getScore.C
クラスの定義を行うgetScore.C
は以下のようになっています。
#include "getScore.H"
#include "surfaceFields.H"
Foam::getScore::getScore
(
const fvMesh& mesh
)
:
IOdictionary
(
IOobject
(
"scoreDict",
mesh.time().constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::AUTO_WRITE
)
),
mesh_(mesh)
{
Info << "Set getScore constructor from only mesh" << endl;
}
Foam::getScore::getScore
(
const fvMesh& mesh,
const IOdictionary& dict
)
:
IOdictionary(dict),
mesh_(mesh)
{
Info << "Set getScore constructor from mesh and dict" << endl;
}
void Foam::getScore::calculate()
{
Info << "Here will be overwritten by inheritance function" << endl;
}
まず以下のようにコンストラクタを二つ作成しています。
Foam::getScore::getScore
(
const fvMesh& mesh
)
:
IOdictionary
(
IOobject
(
"scoreDict",
mesh.time().constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::AUTO_WRITE
)
),
mesh_(mesh)
{
Foam::getScore::getScore
(
const fvMesh& mesh,
const IOdictionary& dict
)
:
IOdictionary(dict),
mesh_(mesh)
{
前者はconst fvMesh& mesh
、後者はconst fvMesh& mesh
とconst IOdictionary& dict
を引数としています。これは継承元のIOdictionary
のコンストラクタをどのように呼び出すかを、以下の★部分を変更することで選択できるようになっています。
Foam::getScore::getScore(?)
:
IOdictionary(★),
mesh_(mesh)
{
前者の場合はこの引数を以下のように引数として与えられたconst fvMesh& mesh
から紐づけて、『その中のconstantディレクトリ(mesh.time().constant()
)の"scoreDict"
という名前のファイルを指定』しています。ちなみにIOobject
はファイルの場所や名前や読み込み方などを示すクラスです。
Foam::getScore::getScore
(
const fvMesh& mesh
)
:
IOdictionary
(
IOobject
(
"scoreDict",
mesh.time().constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::AUTO_WRITE
)
),
後者の場合は引数として与えられたIOdictionary
クラスをそのまま引数としています。つまりこの場合は、この時点ではどのファイルを読み込むかは指定されておらず、ソルバー側でgetScoreのコンストラクタが呼び出されるときに指定することになります。
Foam::getScore::getScore
(
const fvMesh& mesh,
const IOdictionary& dict
)
:
IOdictionary(dict),
これでgetScore
はIOdictionary
を継承しました。それでは次にここから更に継承したsumInletOutletScore
を見てみましょう。
sumInletOutletScore.H
クラスの宣言を行うsumInletOutletScore.H
は以下のようになっています。ここでは継承元のgetScoreに合わせてコンストラクタを一つ増やしています。
#ifndef sumInletOutletScore_H
#define sumInletOutletScore_H
#include "getScore.H"
namespace Foam
{
class sumInletOutletScore
: public getScore
{
protected:
// Protected data
public:
// Constructors
sumInletOutletScore(const fvMesh& mesh);
sumInletOutletScore(const fvMesh& mesh, const IOdictionary& dict);
//- Destructor
~sumInletOutletScore()
{}
// Member Functions
virtual void calculate();
};
}
#endif
sumInletOutletScore.C
クラスの定義を行うsumInletOutletScore.C
は以下のようになっています。
#include "sumInletOutletScore.H"
#include "surfaceFields.H"
Foam::sumInletOutletScore::sumInletOutletScore
(
const fvMesh& mesh
)
:
getScore(mesh)
{
Info << "Set sumInletOutletScore constructor from only mesh" << endl;
}
Foam::sumInletOutletScore::sumInletOutletScore
(
const fvMesh& mesh,
const IOdictionary& dict
)
:
getScore(mesh, dict)
{
Info << "Set sumInletOutletScore constructor from mesh and dictionary" << endl;
}
void Foam::sumInletOutletScore::calculate()
{
Info << "Running inheritance function with dictionary" << endl;
const wordList targetPatches = get<wordList>("targetPatches");
const surfaceScalarField& phi = mesh_.lookupObject<surfaceScalarField>("phi");
score_ = 0.0;
forAll(phi.boundaryField(), patchI)
{
const fvPatch& patch = phi.boundaryField()[patchI].patch();
const scalarField& pphi = phi.boundaryField()[patchI];
bool isTargetPatch = false;
forAll(targetPatches, i)
{
if (patch.name()==targetPatches[i]) isTargetPatch = true;
}
if (!isTargetPatch) continue;
// ちなみにPythonだと if not patch.name() in targetPatches: continue
Info << " sum(" << patch.name() << ") of phi = " << gSum(pphi) << endl;
score_ += gSum(pphi);
}
const word message = get<word>("message");
Info << "message in scoreDict is " << message << endl;
}
コンストラクタ二つの後のcalculateの中身に注目します。以下のようにget<type>
という関数を直接使用できていることが分かります。これはsumInletOutletScore
クラスが継承しているgetScore
クラスが継承しているIOdictionary
クラスの持つ関数です。コンストラクタで読み込んだ"scoreDict"
の内容を読み込んでいます。これにより、境界の名前をtargetPatches
、任意の出力をmessage
として読み込み、使用しています。
const wordList targetPatches = get<wordList>("targetPatches");
...
const word message = get<word>("message");
クラスの呼び出し
こうして作成したクラスは以下のように呼び出しています。
- scorePimpleFoam/pimpleFoam.C : 110行目
sumInletOutletScore myScore
(
mesh,
IOdictionary
(
IOobject
(
"scoreDict",
mesh.time().constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::AUTO_WRITE
)
)
);
ここでは呼び出しの時に辞書ファイルの場所と名前を指定するような後者のコンストラクタを使用しています。前者のコンストラクタを使用する場合は以下のように記述します。なお今回の内容だとどちらもファイルを指定していることになります。
sumInletOutletScore myScore(mesh);
コンパイルと実行
wmakeとケースの実行を以前と同じように行います。例えば以下のようにscoreDictに記述したmessageがちゃんと出力されているかを確認してもよいかもしれません。
$ cd scoreDictPimpleFoam
$ wmake
$ cd ../plainCase
$ foamRunTutorials
$ cat log.scoreDictPimpleFoam | grep -e "message in scoreDict is" > log.message
message in scoreDict is hello_world
message in scoreDict is hello_world
message in scoreDict is hello_world
message in scoreDict is hello_world
message in scoreDict is hello_world
最後に
今回はOpenFOAMの標準ライブラリであるIOdictionary
を継承することで、ケース側で指定した情報を読み込むことができました。
次回は『どの継承クラスを使うか』という方法をケース側で指定する方法について紹介します。
Discussion