🐺

【OpenFOAM】手法を変えて処理を追加④ : 継承クラスに処理を記述する

2020/12/10に公開

オープンCAEアドベントカレンダー 9日目のエントリーです。

はじめに

『【OpenFOAM】手法を変えて処理を追加③ : 関数やクラスとして処理を追加するの続きです。OpenFOAMの計算に手法を変えて処理を追加していきます。
 前回は関数やクラスとして処理を記述しました。
今回はそのクラスによる方法の発展型です。継承クラスに処理を記述してみましょう。

手法一覧

この記事は様々な手法を紹介するシリーズの一つであり、全体としては以下の方法で行う予定です。実行確認済みのファイルはgithub.com/inabower/getScoreTutorialsで公開しています。

  • 例題を用意する ()(キーワード:pimpleFoam
  • 既存機能でやるとしたら ()(キーワード:functionObjects
  1. 計算ファイルに処理を書き込む ()(キーワード:codedFunctionObject
  2. ソルバーに処理を追加する ()
  3. 関数やクラスとして処理を追加する()
  4. 継承クラスに処理を追加 ★今ここ
  5. 設定項目を実行時に読み込む(キーワード:IOdictionary
  6. 継承クラスを実行時に選択する(キーワード:runTimeSelectionTable
  7. クラスを別でコンパイルする(キーワード:wmake libso
  8. 継承クラスの一つとしてcodedを追加(キーワード:codedTemplete

動作環境

  • OS : Ubuntu 20.04 LTS (WSL2 on Windows 10)
  • OpenFOAM : ESI版 v2006 (from source)
  • gcc : 9.3.0

継承クラスの構成

実際のコードはgithub.com/inabower/getScoreTutorials/004.inheritanceに公開しています。

前回の【OpenFOAM】手法を変えて処理を追加③ : 関数やクラスとして処理を追加するの後半で作成したscoreClassPimpleFoam(GitHub)と比較する形で中身をを見ていきます。

ディレクトリ構造

今回のディレクトリ構造は以下のようになっています。

.
├── Make
│   ├── files
│   └── options
├── UEqn.H
├── correctPhi.H
├── cpCorrectPhi.H
├── createFields.H
├── getScore
│   ├── getScore.C
│   ├── getScore.H
│   ├── sumInletOutletScore.C #NEW!!
│   └── sumInletOutletScore.H #NEW!!
├── pEqn.H
└── pimpleFoam.C

前回に比べて、sumInletOutletScore.(C|H)というファイルが増えました。この中での役割はgetScore基底クラスsumInletOutletScore継承クラスとなります。これにより、sumInletOutletScoregetScoreの機能を使用したり、一部を上書きしたりすることができます。
 前回のgetScoreクラスには、メンバ変数score_を計算するためのcalculate()関数と、そのscore_を呼び出すためのvalue()関数が備わっていました。今回はこのうちの『calculate()関数についてはは継承クラスで記述する』という仕組みを作ってみます。

基底クラスgetScore

基底クラスとして以下のようなgetScoreを作成します。

getScore/getScore.H

#ifndef getScore_H
#define getScore_H

#include "fvMesh.H"

namespace Foam
{
class getScore
{
protected:
    // Protected data
        const fvMesh& mesh_;
        scalar score_;
public:
    // Constructors
        getScore(const fvMesh& mesh);
    //- Destructor
        // - ~getScore()
        virtual ~getScore()
        {}
    // Member Functions
        // - void calculate();
        virtual void calculate();

        const scalar& value() const
        {
            return score_;
        }
};
}

#endif

getScore/getScore.C

#include "getScore.H"
#include "surfaceFields.H"

Foam::getScore::getScore
(
    const fvMesh& mesh
)
:
    mesh_(mesh)
{
    Info << "Set constructor" << endl;
}

void Foam::getScore::calculate()
{
    // ここに書かれた内容は継承先で上書きされる。(ここを継承先から呼び出すこともできる)
    Info << "Here will be overwritten by inheritance function" << endl;
}

前回と異なる点はvirtualです。~getScore()void calculate();の頭にvirtualがついています。デストラクタと継承先で上書きしたい関数については関数の宣言の際に先頭にvirtualとつけることになります。
 このように基底クラスとしてgetScoreを作成したのちに、今度は継承クラスを作成します。

継承クラス sumInletOutletScore

継承クラスも基底クラスと同様に"*.H"ファイルで宣言を行い、"*.C"ファイルで定義を行うという構成になります。それぞれの中身は以下のようになっています。
 なお継承クラスは機能ごとに作成することが多いため、その機能を示すような名前にすることが多いです。今回は出入り口の流量の総和を求めるという機能であるため、こんな名前にしてあります。

getScore/sumInletOutletScore.H

#ifndef sumInletOutletScore_H
#define sumInletOutletScore_H

#include "getScore.H"

namespace Foam
{
class sumInletOutletScore
: public getScore
{
protected:
    // Protected data
public:
    // Constructors
        sumInletOutletScore(const fvMesh& mesh);
    //- Destructor
        ~sumInletOutletScore()
        {}
    // Member Functions
        virtual void calculate();
};
}

#endif

この中ではsumInletOutletScoreクラスを作成しています。そのためコンストラクタとデストラクタを以下のように宣言しています。

public:
    // Constructors
        sumInletOutletScore(const fvMesh& mesh);
    //- Destructor
        ~sumInletOutletScore()
        {}

継承を扱う上に特徴的な部分を説明します。まず以下のように基底クラスの宣言を読み込んでいます。

#include "getScore.H"

また以下の部分で継承を示しています。: public getScoregetScoreを継承する事を示しています。

class sumInletOutletScore
: public getScore
{

また上書きするcalculate()関数をvirtual付きで宣言しています。これにより、基底クラスのcalculate()関数をこれで上書きすることができます。

    // Member Functions
        virtual void calculate();

次に定義を見ていきましょう。

getScore/sumInletOutletScore.C

#include "sumInletOutletScore.H"
#include "surfaceFields.H"

Foam::sumInletOutletScore::sumInletOutletScore
(
    const fvMesh& mesh
)
:
    getScore(mesh)
{
    Info << "Set constructor" << endl;
}

void Foam::sumInletOutletScore::calculate()
{

    Info << "Running inheritance function" << endl;

    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];
        if (!(patch.name()=="inlet" || patch.name()=="outlet1" || patch.name()=="outlet2")) continue;
        Info << "    sum(" << patch.name() << ") of phi = " << gSum(pphi) << endl;
        score_ += gSum(pphi);
    }
}

ここではコンストラクタとcalculate()関数を定義しています。

継承クラスのコンストラクタは以下のようになっています。このうち:getScore(mesh)部分では基底クラスのコンストラクタを呼び出しています。こうすることでgetScore(mesh)で呼び出した場合のgetScoreクラスとしての機能を備えることができます。

Foam::sumInletOutletScore::sumInletOutletScore
(
    const fvMesh& mesh
)
:
    getScore(mesh)
{
    Info << "Set constructor" << endl;
}

また以下のようにすることでsumInletOutletScoreクラスとしてのcalculate()関数を定義しています。

void Foam::sumInletOutletScore::calculate()
{
    // iroiro
}

main()側で継承クラスを呼び出す

main()関数では以下のように呼び出しています。前回は#include "getScore.H"としていたところをこのように変更します。getScoreを使用する予定がない場合はここで#includeする必要はありません。

#include "sumInletOutletScore.H"
  • 110行目

sumInletOutletScoreクラスとしてmyScoreという名前のオブジェクトを作成します。

  sumInletOutletScore myScore(mesh);
  • 174行目

最後にmyScoreオブジェクトのcalculate()value()を実行します。ここはgetScoreクラスの時と同じです。ここで実行されるcalculate()sumInletOutletScoreクラスとして再定義されて上書きされたものになります。

        myScore.calculate();
        Info << nl << "score: " << runTime.value() << tab << myScore.value() << nl << endl;

コンパイル

コンパイルのためのMake/の設定は設定は以下のようになります。

Make/files

ここではgetScore.CsumInletOutletScore.Cの両方をコンパイルします。

getScore/getScore.C
getScore/sumInletOutletScore.C
pimpleFoam.C

EXE = $(FOAM_USER_APPBIN)/scoreInheritancePimpleFoam

Make/options

ここはgetScoreのみの時と同じです。もし継承クラスを別のディレクトリに作成した場合にはEXE_INCに追加する必要があります。

EXE_INC = \
    -IgetScore \
    -I$(LIB_SRC)/finiteVolume/lnInclude \
    -I$(LIB_SRC)/meshTools/lnInclude \
    -I$(LIB_SRC)/sampling/lnInclude \
    -I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
    -I$(LIB_SRC)/TurbulenceModels/incompressible/lnInclude \
    -I$(LIB_SRC)/transportModels \
    -I$(LIB_SRC)/transportModels/incompressible/singlePhaseTransportModel \
    -I$(LIB_SRC)/dynamicMesh/lnInclude \
    -I$(LIB_SRC)/dynamicFvMesh/lnInclude

EXE_LIBS = \
    -lfiniteVolume \
    -lfvOptions \
    -lmeshTools \
    -lsampling \
    -lturbulenceModels \
    -lincompressibleTurbulenceModels \
    -lincompressibleTransportModels \
    -ldynamicMesh \
    -ldynamicFvMesh \
    -ltopoChangerFvMesh \
    -latmosphericModels

動作確認

以上のコードをコンパイルして実行すると前回と同様にscoreに関する出力がされると思います。

最後に

今回はgetScoreクラスのcalculate()関数について継承クラスで作成するということを行いました。これにより、例えば別々のcalculate()を行うクラスを複数作成したい場合などには、共通部分をgetScoreクラスに作っておくことで、コードがシンプルになり可読性が向上します。
 次回はコードではなく計算ケースの方で設定を追加する方法を見ていきます。

Discussion