Loading...
Searching...
No Matches
chtMultiRegionFoam.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2019,2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 chtMultiRegionFoam
29
30Group
31 grpHeatTransferSolvers
32
33Description
34 Transient solver for buoyant, turbulent fluid flow and solid heat
35 conduction with conjugate heat transfer between solid and fluid regions.
36
37 It handles secondary fluid or solid circuits which can be coupled
38 thermally with the main fluid region. i.e radiators, etc.
39
40\*---------------------------------------------------------------------------*/
41
42#include "fvCFD.H"
44#include "rhoReactionThermo.H"
45#include "CombustionModel.H"
47#include "regionProperties.H"
48#include "compressibleCourantNo.H"
49#include "solidRegionDiffNo.H"
50#include "solidThermo.H"
51#include "radiationModel.H"
52#include "fvOptions.H"
53#include "coordinateSystem.H"
54#include "loopControl.H"
55#include "pressureControl.H"
56
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60int main(int argc, char *argv[])
61{
62 argList::addNote
63 (
64 "Transient solver for buoyant, turbulent fluid flow and solid heat"
65 " conduction with conjugate heat transfer"
66 " between solid and fluid regions."
67 );
68
69 #define NO_CONTROL
70 #define CREATE_MESH createMeshesPostProcess.H
71 #include "postProcess.H"
72
73 #include "addCheckCaseOptions.H"
74 #include "setRootCaseLists.H"
75 #include "createTime.H"
76 #include "createMeshes.H"
77 #include "createFields.H"
78 #include "initContinuityErrs.H"
79 #include "createTimeControls.H"
81 #include "compressibleMultiRegionCourantNo.H"
82 #include "solidRegionDiffusionNo.H"
84
85 #include "createCoupledRegions.H"
86
87 while (runTime.run())
88 {
89 #include "readTimeControls.H"
91 #include "readPIMPLEControls.H"
92
93 #include "compressibleMultiRegionCourantNo.H"
94 #include "solidRegionDiffusionNo.H"
95 #include "setMultiRegionDeltaT.H"
96
97 ++runTime;
98
99 Info<< "Time = " << runTime.timeName() << nl << endl;
100
101 if (nOuterCorr != 1)
102 {
104 {
105 #include "storeOldFluidFields.H"
106 }
107 }
108
109 // --- PIMPLE loop
110 for (int oCorr=0; oCorr<nOuterCorr; ++oCorr)
111 {
112 const bool finalIter = (oCorr == nOuterCorr-1);
113
115 {
116 fvMesh& mesh = fluidRegions[i];
117
118 #include "readFluidMultiRegionPIMPLEControls.H"
119 #include "setRegionFluidFields.H"
120 #include "solveFluid.H"
121 }
122
124 {
125 fvMesh& mesh = solidRegions[i];
126
128 #include "setRegionSolidFields.H"
129 #include "solveSolid.H"
130 }
131
132 if (coupled)
133 {
134 Info<< "\nSolving energy coupled regions " << endl;
135 fvMatrixAssemblyPtr->solve();
136 #include "correctThermos.H"
137
139 {
140 fvMesh& mesh = fluidRegions[i];
141
142 #include "readFluidMultiRegionPIMPLEControls.H"
143 #include "setRegionFluidFields.H"
144 if (!frozenFlow)
145 {
146 Info<< "\nSolving for fluid region "
147 << fluidRegions[i].name() << endl;
148 // --- PISO loop
149 for (int corr=0; corr<nCorr; corr++)
150 {
151 #include "pEqn.H"
152 }
153 turbulence.correct();
154 }
155
156 rho = thermo.rho();
157 Info<< "Min/max T:" << min(thermo.T()).value() << ' '
158 << max(thermo.T()).value() << endl;
159 }
160
161 fvMatrixAssemblyPtr->clear();
162 }
163
164 // Additional loops for energy solution only
165 if (!oCorr && nOuterCorr > 1)
166 {
167 loopControl looping(runTime, pimple, "energyCoupling");
168
169 while (looping.loop())
170 {
171 Info<< nl << looping << nl;
172
174 {
175 fvMesh& mesh = fluidRegions[i];
176
177 Info<< "\nSolving for fluid region "
178 << fluidRegions[i].name() << endl;
179 #include "readFluidMultiRegionPIMPLEControls.H"
180 #include "setRegionFluidFields.H"
181 frozenFlow = true;
182 #include "solveFluid.H"
183 }
184
186 {
187 fvMesh& mesh = solidRegions[i];
188
189 Info<< "\nSolving for solid region "
190 << solidRegions[i].name() << endl;
192 #include "setRegionSolidFields.H"
193 #include "solveSolid.H"
194 }
195
196 if (coupled)
197 {
198 Info<< "\nSolving energy coupled regions " << endl;
199 fvMatrixAssemblyPtr->solve();
200 #include "correctThermos.H"
201
203 {
204 #include "setRegionFluidFields.H"
205 rho = thermo.rho();
206 }
207
208 fvMatrixAssemblyPtr->clear();
209 }
210 }
211 }
212 }
213
214 runTime.write();
215
216 runTime.printExecutionTime(Info);
217 }
218
219 Info<< "End\n" << endl;
220
221 return 0;
222}
223
224
225// ************************************************************************* //
Required Classes.
PtrList< fvMesh > fluidRegions(fluidNames.size())
bool frozenFlow
pimpleControl & pimple
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
autoPtr< fvMatrix< scalar > > fvMatrixAssemblyPtr
bool coupled
dynamicFvMesh & mesh
engineTime & runTime
PtrList< fvMesh > solidRegions(solidNames.size())
Read the control parameters used by setDeltaT.
compressible::turbulenceModel & turbulence
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Execute application functionObjects to post-process existing results.
const int nOuterCorr
Read the control parameters used in the solid.
Read the control parameters used by setDeltaT.
Set the initial timestep for the CHT MultiRegion solver.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299