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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Application
28  chtMultiRegionFoam
29 
30 Group
31  grpHeatTransferSolvers
32 
33 Description
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 
60 int 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 "setRootCaseLists.H"
74  #include "createTime.H"
75  #include "createMeshes.H"
76  #include "createFields.H"
77  #include "initContinuityErrs.H"
78  #include "createTimeControls.H"
79  #include "readSolidTimeControls.H"
80  #include "compressibleMultiRegionCourantNo.H"
81  #include "solidRegionDiffusionNo.H"
83 
84  #include "createCoupledRegions.H"
85 
86  while (runTime.run())
87  {
88  #include "readTimeControls.H"
89  #include "readSolidTimeControls.H"
90  #include "readPIMPLEControls.H"
91 
92  #include "compressibleMultiRegionCourantNo.H"
93  #include "solidRegionDiffusionNo.H"
94  #include "setMultiRegionDeltaT.H"
95 
96  ++runTime;
97 
98  Info<< "Time = " << runTime.timeName() << nl << endl;
99 
100  if (nOuterCorr != 1)
101  {
102  forAll(fluidRegions, i)
103  {
104  #include "storeOldFluidFields.H"
105  }
106  }
107 
108  // --- PIMPLE loop
109  for (int oCorr=0; oCorr<nOuterCorr; ++oCorr)
110  {
111  const bool finalIter = (oCorr == nOuterCorr-1);
112 
113  forAll(fluidRegions, i)
114  {
115  #include "setRegionFluidFields.H"
116  #include "readFluidMultiRegionPIMPLEControls.H"
117  #include "solveFluid.H"
118  }
119 
120  forAll(solidRegions, i)
121  {
122  #include "setRegionSolidFields.H"
124  #include "solveSolid.H"
125  }
126 
127  if (coupled)
128  {
129  Info<< "\nSolving energy coupled regions " << endl;
130  fvMatrixAssemblyPtr->solve();
131  #include "correctThermos.H"
132 
133  forAll(fluidRegions, i)
134  {
135  #include "setRegionFluidFields.H"
136  #include "readFluidMultiRegionPIMPLEControls.H"
137  if (!frozenFlow)
138  {
139  Info<< "\nSolving for fluid region "
140  << fluidRegions[i].name() << endl;
141  // --- PISO loop
142  for (int corr=0; corr<nCorr; corr++)
143  {
144  #include "pEqn.H"
145  }
146  turbulence.correct();
147  }
148 
149  rho = thermo.rho();
150  Info<< "Min/max T:" << min(thermo.T()).value() << ' '
151  << max(thermo.T()).value() << endl;
152  }
153 
154  fvMatrixAssemblyPtr->clear();
155  }
156 
157  // Additional loops for energy solution only
158  if (!oCorr && nOuterCorr > 1)
159  {
160  loopControl looping(runTime, pimple, "energyCoupling");
161 
162  while (looping.loop())
163  {
164  Info<< nl << looping << nl;
165 
166  forAll(fluidRegions, i)
167  {
168  Info<< "\nSolving for fluid region "
169  << fluidRegions[i].name() << endl;
170  #include "setRegionFluidFields.H"
171  #include "readFluidMultiRegionPIMPLEControls.H"
172  frozenFlow = true;
173  #include "solveFluid.H"
174  }
175 
176  forAll(solidRegions, i)
177  {
178  Info<< "\nSolving for solid region "
179  << solidRegions[i].name() << endl;
180  #include "setRegionSolidFields.H"
182  #include "solveSolid.H"
183  }
184 
185  if (coupled)
186  {
187  Info<< "\nSolving energy coupled regions " << endl;
188  fvMatrixAssemblyPtr->solve();
189  #include "correctThermos.H"
190 
191  forAll(fluidRegions, i)
192  {
193  #include "setRegionFluidFields.H"
194  rho = thermo.rho();
195  }
196 
197  fvMatrixAssemblyPtr->clear();
198  }
199  }
200  }
201  }
202 
203  runTime.write();
204 
205  runTime.printExecutionTime(Info);
206  }
207 
208  Info<< "End\n" << endl;
209 
210  return 0;
211 }
212 
213 
214 // ************************************************************************* //
readSolidMultiRegionPIMPLEControls.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
setInitialMultiRegionDeltaT.H
Set the initial timestep for the CHT MultiRegion solver.
nCorr
const int nCorr
Definition: readFluidMultiRegionPIMPLEControls.H:3
fvOptions.H
regionProperties.H
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
setRegionSolidFields.H
rho
rho
Definition: readInitialConditions.H:88
coordinateSystem.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
solidThermo.H
rhoReactionThermo.H
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
createCoupledRegions.H
CombustionModel.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
solidRegionDiffNo.H
Calculates and outputs the mean and maximum Diffusion Numbers for the solid regions.
readSolidTimeControls.H
Read the control parameters used in the solid.
postProcess.H
Execute application functionObjects to post-process existing results.
nOuterCorr
const int nOuterCorr
Definition: readPIMPLEControls.H:7
fixedGradientFvPatchFields.H
pressureControl.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
loopControl.H
fluidRegions
PtrList< fvMesh > fluidRegions(fluidNames.size())
createTimeControls.H
Read the control parameters used by setDeltaT.
solidRegions
PtrList< fvMesh > solidRegions(solidNames.size())
readTimeControls.H
Read the control parameters used by setDeltaT.
correctThermos.H
createTime.H
createMeshes.H
fvCFD.H
readPIMPLEControls.H
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
frozenFlow
bool frozenFlow
Definition: setRegionFluidFields.H:34
radiationModel.H
turbulentFluidThermoModel.H
fvMatrixAssemblyPtr
autoPtr< fvMatrix< scalar > > fvMatrixAssemblyPtr
Definition: createCoupledRegions.H:5