こんにちは(@t_kun_kamakiri)
今回はOpenFOAMで以下のような配管内の水の流れを解析し、境界と任意の断面での流量を計算します。
境界の「inlet」「outlet1」「outlet2」での流量はOpenFOAMのfunction objectの機能を使うことで簡単に出力することができます。
【OpenFOAM】境界面での流量計算(volFieldValue)
※本記事のモデル作成方法はこちらをご参考ください。
しかし、c1~c5の断面での流量を出力するにはどうしたらいいのかなと思って自分なりに試してみた内容を書いています。
OpenFOAMで任意の断面での流量を計算する。
【使用環境】
OpenFOAM-v2012(WSL2)
流量を出力する方針
流量を計算する方法はパッと考えると以下の方法があるでしょう。
- 流量計測したい断面にあらかじめ面かセル領域を用意しておく
- cellZoneSetを使ってfunction objectを使う
- codedFunctionObjectという機能で直接コードを埋め込む
①の方法はあらかじめ計算する前に用意する必要があると考えて今回は試しませんでした。
今回は解析計算が終わった後でも出力できる方法として②と③を試したのでご紹介します。
③が本記事の内容
とはいえ、任意の断面での流量を直接出力するのはわからなかったので任意の領域を通過する平均流速を求めて面積を掛けて流量$Q=vA$とする方針にしました。
※任意の面を通過する流量とすると、入ってきて出ていく量の合計を計算するのことになって結局0になってしまいました。なので、いったん平均流速を求めて面積を掛けて流量を求めればいいやと思った次第です。
③codedFunctionObject機能で直接コードを埋め込む
③は直接コードを書いて出力する方法です。
②の方法のようにtopoSetとかわざわざ作ったりするのが面倒だと感じたり、既存機能だけではやりたいことができない場合に役に立ちます。
codedFunctionObjectはsystem/controlDictのfunction objectでC++コードを直接埋め込んで計算実行と同時にコンパイルしてくれる機能です。
OpenFOAM用のC++コードに慣れないと難しく感じるかもしれません。
以下に参考記事を載せておきます。
今回はc1領域の平均流速、平均圧力、平均流量を計算するコードだけを追加します。
system/controlDictのfunctionsに以下を追加します。
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 |
codevolFieldValue { type coded; name codevolFieldValue; libs (utilityFunctionObjects); codeExecute #{ Info << "============= codevolFieldValue ===============" << endl; volVectorField C = mesh().C();//cell infomation scalarField V = mesh().V();//cell volume // intitial label Vi = 0; scalar Volsum = 0.0; vector UVsum(0.0, 0.0, 0.0); scalar pVsum = 0.0; scalar phisum = 0.0; //c1での半径と断面積 scalar radius = 0.02; //[m] scalar Area = radius*radius*3.14; // U,p, phiを定義 const volVectorField& U = mesh().lookupObject<volVectorField>("U"); const volScalarField& p = mesh().lookupObject<volScalarField>("p"); const surfaceScalarField& phi = mesh().lookupObject<surfaceScalarField>("phi"); //https://caefn.com/openfoam/ja-cellsource //https://caefn.com/openfoam/ja-vector-operations forAll(C, i) { //x,y,x座標を定義 1つ目の引数がcellのid、2つ目の引数がx,y,z座標 scalar x = C[i][0]; scalar y = C[i][1]; scalar z = C[i][2]; //x,y,x座標を定義 1つ目の引数がcellのid、2つ目の引数がx,y,z座標 if ( (z >=0.05) && (z <=0.06) ) //zが0.05~0.06の領域のみ(つまりc1領域だけ) { if( sqrt(pow(x, 2) + pow(y, 2)) <= 0.02) { /* idが取得できたらU,p,phiにアクセスすると値を取得できる */ vector pos(x, y, z); label cellNo = mesh().findCell(pos); Vi += 1; Volsum += V[i]; UVsum += U[cellNo]*V[i]; pVsum += p[cellNo]*V[i]; phisum += phi[cellNo]*V[i]; } } } // 出力 Info << "volFieldValue codevolFieldValue write:" << endl; Info << "volAverage(codevolFieldValue) of U =" << UVsum/Volsum << endl; Info << "volAverage(codevolFieldValue) of p =" << pVsum/Volsum << endl; Info << "volAverage(codevolFieldValue) of phi = U*A = " << mag(UVsum/Volsum*Area) << endl; Info << "Vi " << Vi << endl; #}; } |
流れは以下のようになっています。
- cellとvolumeを定義
- c1領域の座標のみのcellの値を取得するif文
- 座標を取得
- 座標からidを取得
- idからU.p.phiの値を取得
- 最後に領域内の平均値を計算
計算を実行します。
1 |
simpleFoam > log.simpleFoam & |
計算が終わった後にcontrolDictのfunction objectだけを実行したい場合は、
1 |
simpleFoam -postProcess |
とします。
lob.simpleFoamのログファイルに、
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 |
(省略) Time = 2000 smoothSolver: Solving for Ux, Initial residual = 0.0192337, Final residual = 0.000478075, No Iterations 2 smoothSolver: Solving for Uy, Initial residual = 0.00328875, Final residual = 7.75909e-05, No Iterations 2 smoothSolver: Solving for Uz, Initial residual = 0.00138017, Final residual = 3.70736e-05, No Iterations 2 GAMG: Solving for p, Initial residual = 0.000104941, Final residual = 2.85289e-06, No Iterations 2 time step continuity errors : sum local = 0.00324322, global = 4.34941e-06, cumulative = 0.189192 ExecutionTime = 390.66 s ClockTime = 430 s continuityError continuityError1 write: local = 0.00324322 global = 4.34941e-06 cumulative = 0.00447026 surfaceFieldValue inlet write: sum(inlet) of phi = -0.000125382 surfaceFieldValue outlet1 write: sum(outlet1) of phi = 0.000111621 surfaceFieldValue outlet2 write: sum(outlet2) of phi = 1.37635e-05 volFieldValue volFieldValue1 write: volAverage(c1) of U = (-6.71016e-06 -4.74256e-06 0.100544) volAverage(c1) of p = 18.1052 volFieldValue volFieldValue2 write: volAverage(c2) of U = (3.25019e-07 -6.6332e-05 2.52346) volAverage(c2) of p = 7.04559 volFieldValue volFieldValue3 write: volAverage(c3) of U = (-1.4252e-05 0.0107654 2.24645) volAverage(c3) of p = 2.30779 volFieldValue volFieldValue4 write: volAverage(c4) of U = (-0.00569703 -0.00613122 0.0919273) volAverage(c4) of p = -0.289308 volFieldValue volFieldValue5 write: volAverage(c5) of U = (2.92012e-07 -0.276999 1.76482e-06) volAverage(c5) of p = 1.98803 ============= codevolFieldValue =============== volFieldValue codevolFieldValue write: volAverage(codevolFieldValue) of U =(-6.71016e-06 -4.74256e-06 0.100544) volAverage(codevolFieldValue) of p =18.1052 volAverage(codevolFieldValue) of phi = U*A = 0.000126284 Vi 1216 End |
注目すべきは以下の部分です。
同じc1領域を計算したものです。
- topoSetDictで作った領域からvolFieldValueで計算したのが上
- codedFunctionObjectで計算したのが下
1 2 3 4 5 6 7 8 9 10 11 |
volFieldValue volFieldValue1 write: volAverage(c1) of U = (-6.71016e-06 -4.74256e-06 0.100544) volAverage(c1) of p = 18.1052 (省略) ============= codevolFieldValue =============== volFieldValue codevolFieldValue write: volAverage(codevolFieldValue) of U =(-6.71016e-06 -4.74256e-06 0.100544) volAverage(codevolFieldValue) of p =18.1052 volAverage(codevolFieldValue) of phi = U*A = 0.000126284 Vi 1216 End |
平均流速、平均圧力ともに同じ値になっていますね。
また、topoSetで作ったセルの数が1216でしたが、コードで計算したものも同じ1216です。
さらにコードの中で平均流速に面積を掛けて流量を計算して$Q=vA=0.000126284$と出しています。
結局「controlDict」は色々書き込んで以下のようになりました。
system/controlDict
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 |
/*--------------------------------*- C++ -*----------------------------------*\ | ========= | | | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | | \\ / O peration | Version: v2012 | | \\ / A nd | Website: www.openfoam.com | | \\/ M anipulation | | \*---------------------------------------------------------------------------*/ FoamFile { version 2.0; format ascii; class dictionary; location "system"; object controlDict; } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // application simpleFoam; startFrom latestTime; startTime 0; stopAt endTime; endTime 2000; deltaT 1; writeControl timeStep; writeInterval 100; purgeWrite 0; writeFormat ascii; writePrecision 6; writeCompression off; timeFormat general; timePrecision 6; runTimeModifiable true; functions { // #includeFunc streamlines // #include sampleDict //#includeFunc solverInfo continuityError1 { type continuityError; libs ( fieldFunctionObjects ); // Mandatory entries (unmodifiable)// Optional entries (runtime modifiable) phi phi; // Optional (inherited) entries writePrecision 8; writeToFile yes; useUserTime yes; region region0; enabled yes; log yes; timeStart 0; timeEnd 10000; executeControl timeStep; executeInterval 1; writeControl timeStep; writeInterval 50; } inlet { type surfaceFieldValue; libs ( "libfieldFunctionObjects.so" ); writeControl timeStep; log yes; writeFields no; regionType patch; name inlet; operation sum; fields ( phi ); // Optional (inherited) entries writePrecision 8; writeToFile yes; useUserTime yes; region region0; enabled yes; timeStart 0; timeEnd 10000; executeControl timeStep; executeInterval 1; writeInterval 100; } outlet1 { $inlet; name outlet1; } outlet2 { $inlet; name outlet2; } residuals { type solverInfo; //libs("libutilityFunctionObjects.so"); libs (utilityFunctionObjects); fields (".*" ); } //Add volFieldValue1 { // Mandatory entries (unmodifiable) type volFieldValue; libs (fieldFunctionObjects); // Mandatory entries (runtime modifiable) fields (U p); operation volAverage; regionType cellZone; name c1; // Optional entries (runtime modifiable) postOperation none; weightField alpha1; // Optional (inherited) entries writeFields false; scalingFactor 1.0; writePrecision 8; writeToFile true; useUserTime true; region region0; enabled true; log true; timeStart 0; timeEnd 10000; executeControl timeStep; executeInterval 1; writeControl timeStep; writeInterval 100; } volFieldValue2 { $volFieldValue1 name c2; } volFieldValue3 { $volFieldValue1 name c3; } volFieldValue4 { $volFieldValue1 name c4; } volFieldValue5 { $volFieldValue1 name c5; } //https://www.openfoam.com/documentation/guides/latest/api/classFoam_1_1functionObjects_1_1codedFunctionObject.html codevolFieldValue { type coded; name codevolFieldValue; libs (utilityFunctionObjects); codeExecute #{ Info << "============= codevolFieldValue ===============" << endl; volVectorField C = mesh().C();//cell infomation scalarField V = mesh().V();//cell volume // intitial label Vi = 0; scalar Volsum = 0.0; vector UVsum(0.0, 0.0, 0.0); scalar pVsum = 0.0; scalar phisum = 0.0; //c1での半径と断面積 scalar radius = 0.02; //[m] scalar Area = (radius*radius*3.14); // U,p, phiを定義 const volVectorField& U = mesh().lookupObject<volVectorField>("U"); const volScalarField& p = mesh().lookupObject<volScalarField>("p"); const surfaceScalarField& phi = mesh().lookupObject<surfaceScalarField>("phi"); //https://caefn.com/openfoam/ja-cellsource //https://caefn.com/openfoam/ja-vector-operations forAll(C, i) { //x,y,x座標を定義 1つ目の引数がcellのid、2つ目の引数がx,y,z座標 scalar x = C[i][0]; scalar y = C[i][1]; scalar z = C[i][2]; //x,y,x座標を定義 1つ目の引数がcellのid、2つ目の引数がx,y,z座標 if ( (z >=0.05) && (z <=0.06) ) //zが0.05~0.06の領域のみ(つまりc1領域だけ) { if( sqrt(pow(x, 2) + pow(y, 2)) <= 0.02) { /* idが取得できたらU,p,phiにアクセスすると値を取得できる */ vector pos(x, y, z); label cellNo = mesh().findCell(pos); Vi += 1; Volsum += V[i]; UVsum += U[cellNo]*V[i]; pVsum += p[cellNo]*V[i]; phisum += phi[cellNo]*V[i]; } } } // 出力 Info << "volFieldValue codevolFieldValue write:" << endl; Info << "volAverage(codevolFieldValue) of U =" << UVsum/Volsum << endl; Info << "volAverage(codevolFieldValue) of p =" << pVsum/Volsum << endl; Info << "volAverage(codevolFieldValue) of phi = U*A = " << mag(UVsum/Volsum*Area) << endl; Info << "Vi " << Vi << endl; #}; } } // ************************************************************************* // |
ファイルにまとめてincludeした方がいいかもしれないですね。
このままでは少し煩雑・・・
codedFUnctionObjectを別ファイルから読み込む
controlDictに直接コードを書くとファイルが長くなって見にくいのでcodedFunctionObjectは別ファイルに保存して読み込むようにします。
チュートリアルで同じようなことをやっている例題がないか探します。
以下のコマンドで「coded」という文字を使っているファイルをチュートリアルから探します。
1 |
find $FOAM_TUTORIALS -type f | xargs grep -ln "coded" |
候補が出てきました。
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 |
/opt/OpenFOAM/OpenFOAM-v2012/tutorials/multiphase/overInterDyMFoam/floatingBody/background/system/controlDict /opt/OpenFOAM/OpenFOAM-v2012/tutorials/multiphase/overInterDyMFoam/twoSimpleRotors/system/controlDict /opt/OpenFOAM/OpenFOAM-v2012/tutorials/multiphase/icoReactingMultiPhaseInterFoam/poolEvaporation/0.orig/T /opt/OpenFOAM/OpenFOAM-v2012/tutorials/multiphase/reactingTwoPhaseEulerFoam/laminar/bubbleColumnEvaporating/system/continuityFunctions /opt/OpenFOAM/OpenFOAM-v2012/tutorials/heatTransfer/chtMultiRegionFoam/externalCoupledHeater/externalSolver /opt/OpenFOAM/OpenFOAM-v2012/tutorials/heatTransfer/chtMultiRegionSimpleFoam/externalCoupledHeater/externalSolver /opt/OpenFOAM/OpenFOAM-v2012/tutorials/combustion/XiDyMFoam/oscillatingCylinder/0/pointDisplacementy /opt/OpenFOAM/OpenFOAM-v2012/tutorials/combustion/XiEngineFoam/kivaTest/system/controlDict /opt/OpenFOAM/OpenFOAM-v2012/tutorials/incompressible/simpleFoam/pipeCyclic/0.orig/U /opt/OpenFOAM/OpenFOAM-v2012/tutorials/mesh/moveDynamicMesh/twistingColumn/constant/dynamicMeshDict /opt/OpenFOAM/OpenFOAM-v2012/tutorials/compressible/rhoPimpleFoam/RAS/externalCoupledSquareBendLiq/externalSolver /opt/OpenFOAM/OpenFOAM-v2012/tutorials/compressible/rhoPorousSimpleFoam/angledDuct/implicit/system/sampling /opt/OpenFOAM/OpenFOAM-v2012/tutorials/compressible/rhoPimpleAdiabaticFoam/rutlandVortex2D/system/preProcess /opt/OpenFOAM/OpenFOAM-v2012/tutorials/compressible/overRhoPimpleDyMFoam/twoSimpleRotors/system/controlDict /opt/OpenFOAM/OpenFOAM-v2012/tutorials/basic/potentialFoam/cylinder/system/controlDict /opt/OpenFOAM/OpenFOAM-v2012/tutorials/resources/geometry/nacaAirfoil/nacaAirfoil.inp |
おそらく最後から4番目がやりたいチュートリアルかと思いますので「codedvolField」というファイル名で「system」フォルダにコピーします。
1 |
cp -r /opt/OpenFOAM/OpenFOAM-v2012/tutorials/compressible/rhoPimpleAdiabaticFoam/rutlandVortex2D/system/preProcess system/codedvolField |
※こちらのチュートリアルにあるAllrunを見ると実行コマンドのヒントになります。
では、先ほどのコードをコピーします。
system/codedvolField
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 |
/*--------------------------------*- C++ -*----------------------------------*\ | ========= | | | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | | \\ / O peration | Version: v2012 | | \\ / A nd | Website: www.openfoam.com | | \\/ M anipulation | | \*---------------------------------------------------------------------------*/ FoamFile { version 2.0; format ascii; class dictionary; location "system"; object controlDict; } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // functions { //https://www.openfoam.com/documentation/guides/latest/api/classFoam_1_1functionObjects_1_1codedFunctionObject.html codevolFieldValue { type coded; name codevolFieldValue; libs (utilityFunctionObjects); codeExecute #{ Info << "============= codevolFieldValue ===============" << endl; volVectorField C = mesh().C();//cell infomation scalarField V = mesh().V();//cell volume // intitial label Vi = 0; scalar Volsum = 0.0; vector UVsum(0.0, 0.0, 0.0); scalar pVsum = 0.0; scalar phisum = 0.0; //c1での半径と断面積 scalar radius = 0.02; //[m] scalar Area = (radius*radius*3.14); // U,p, phiを定義 const volVectorField& U = mesh().lookupObject<volVectorField>("U"); const volScalarField& p = mesh().lookupObject<volScalarField>("p"); const surfaceScalarField& phi = mesh().lookupObject<surfaceScalarField>("phi"); //https://caefn.com/openfoam/ja-cellsource //https://caefn.com/openfoam/ja-vector-operations forAll(C, i) { //x,y,x座標を定義 1つ目の引数がcellのid、2つ目の引数がx,y,z座標 scalar x = C[i][0]; scalar y = C[i][1]; scalar z = C[i][2]; //x,y,x座標を定義 1つ目の引数がcellのid、2つ目の引数がx,y,z座標 if ( (z >=0.05) && (z <=0.06) ) //zが0.05~0.06の領域のみ(つまりc1領域だけ) { if( sqrt(pow(x, 2) + pow(y, 2)) <= 0.02) { /* idが取得できたらU,p,phiにアクセスすると値を取得できる */ vector pos(x, y, z); label cellNo = mesh().findCell(pos); Vi += 1; Volsum += V[i]; UVsum += U[cellNo]*V[i]; pVsum += p[cellNo]*V[i]; phisum += phi[cellNo]*V[i]; } } } // 出力 Info << "volFieldValue codevolFieldValue write:" << endl; Info << "volAverage(codevolFieldValue) of U =" << UVsum/Volsum << endl; Info << "volAverage(codevolFieldValue) of p =" << pVsum/Volsum << endl; Info << "volAverage(codevolFieldValue) of phi = U*A = " << mag(UVsum/Volsum*Area) << endl; Info << "Vi " << Vi << endl; #}; } } |
これで計算を実行すればOKです。
計算終了後に出力したい場合は以下のコマンドでも出力できます。
1 |
simpleFoam -postProcess -dict system/codedvolField -latestTime |
出力設定をいじりたい場合には便利ですね。
codedFunctionObject機能(境界面の場を取得)
codedFunctionObjectは直接C++のコードを書くことができる機能なので任意のセル情報以外にも、名前を指定した境界面の情報も取得できます。
補足で少し解説しておきますのでご参考ください。
まず、実行すると「dynamicCode/codevolFieldValue」に以下のファイルが生成されコンパイルが開始されます。
1 2 3 4 5 6 7 8 9 10 11 12 |
├── Make │ ├── SHA1Digest │ ├── files │ ├── linux64Gcc63DPInt32Opt │ │ ├── functionObjectTemplate.C.dep │ │ ├── functionObjectTemplate.o │ │ ├── options │ │ ├── sourceFiles │ │ └── variables │ └── options ├── functionObjectTemplate.C ├── functionObjectTemplate.H |
ポイントは、
- functionObjectTemplate.Hに変数や関数の定義
- functionObjectTemplate.Cに計算処理内容
です。
大まかな情報を得るのには主にヘッダーファイル「functionObjectTemplate.H」を見れば良いです。
何を引数に入れて何が返ってくるのかさえ押さえておけばだいたい使えます。OpenFOAMには異なる型を引数に入れても動作するように関数をオーバーロード(上書き)して同じ名前の関数でも引数違いの関数を複数定義しています。
細かい計算内容を知りたい場合は拡張子.Cの方も見ていけば良いです。
上記のコードは特定の座標に対応するセル中心の物理量の場(流速ベクトル、圧力)を取得する方法でしたが、他にも色々な情報を取得できます。
以下、境界面の値の取得例を出しておきます。
まず注目すべきは、
functionObjectTemplate.H
1 2 |
//- Cast reference of objectRegistry to fvMesh const fvMesh& mesh() const; |
ここです。
- mesh().V():体積ボリューム
- mesh().Sf():セル表面ベクトル
- mesh().magSf():セル表面の大きさ
- mesh().C():セル中心の座標情報
- mesh().Cf():セル表面中心
など色々とセル情報の取得ができます。
さらに以下のに書けば名前が付いた境界面の値も取得できます。
1 2 3 4 5 6 7 8 9 |
Info << "============= Boundary U value ===============" << endl; label inletvelID = mesh().boundaryMesh().findPatchID("inlet"); Info << "inletvelID = " << inletvelID << endl; const fvPatchVectorField& UinletList = U.boundaryField()[inletvelID]; forAll(UinletList, fi) { Info << "U[" << fi << "] =" << UinletList[fi] << endl; } |
fvMeshが継承しているpolyMeshクラスで用意されている関数を調べるとわかります。
mesh().boundaryMesh()はpolyMeshクラスの中のメンバ関数boundaryMeshで以下のように、polyBoundaryMeshを返します。
const polyBoundaryMesh & boundaryMesh()const
まとめ
今回は流量を出力する方法をまとめました。
- 任意の断面の流量計算(topoSetとcodedFuction)
topoSetとcodedFuctionの結果は同じになりましたので、慣れればcodedFuctionが汎用性が効きますね。
参考書
有限体積法の勉強にはこの本で間違いないです。
対流項目が中心差分で安定しないことを理論的に評価しておりとても参考になります。
OpenFOAMの基本的な設定を学ぶことができる書籍として以下の参考書がお勧めです。
改訂新版 OpenFOAMの歩き方 (技術の泉シリーズ(NextPublishing))
本書はESI版ではなくFoundation版でありOpenFOAMv2.xで書かれたやや古い内容です。その点に注意してESI版のお使いのversionとの違いを丁寧に調べることができる方であれば、十分参考になる書籍です。
OpenFOAM全般の設定などまとめられた書籍として以下のものがあります。
こちらもFoundation版なので、ESI版のお使いのversionとの違いを丁寧に調べることができる方であれば、十分参考になる書籍です。
こちらのオープンCAE学会から技術書典にOpenFOAMでメッシュ作成【PDF版】が出ています。難しい内容ではありますが、OpenFOAMのメッシュに関する内容が1000円で学べます。