chtMultiRegionSimpleFoam.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  chtMultiRegionSimpleFoam
29 
30 Group
31  grpHeatTransferSolvers
32 
33 Description
34  Steady-state solver for buoyant, turbulent fluid flow and solid heat
35  conduction with conjugate heat transfer between solid and fluid regions.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "fvCFD.H"
40 #include "rhoThermo.H"
43 #include "regionProperties.H"
44 #include "solidThermo.H"
45 #include "radiationModel.H"
46 #include "fvOptions.H"
47 #include "coordinateSystem.H"
48 #include "loopControl.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 int main(int argc, char *argv[])
53 {
54  argList::addNote
55  (
56  "Steady-state solver for buoyant, turbulent fluid flow and solid heat"
57  " conduction with conjugate heat transfer"
58  " between solid and fluid regions."
59  );
60 
61  #define NO_CONTROL
62  #define CREATE_MESH createMeshesPostProcess.H
63  #include "postProcess.H"
64 
65  #include "setRootCaseLists.H"
66  #include "createTime.H"
67  #include "createMeshes.H"
68  #include "createFields.H"
69  #include "createCoupledRegions.H"
70  #include "initContinuityErrs.H"
71 
72  while (runTime.loop())
73  {
74  Info<< "Time = " << runTime.timeName() << nl << endl;
75 
77  {
78  Info<< "\nSolving for fluid region "
79  << fluidRegions[i].name() << endl;
80  #include "setRegionFluidFields.H"
82  #include "solveFluid.H"
83  }
84 
86  {
87  #include "setRegionSolidFields.H"
89  #include "solveSolid.H"
90  }
91 
92 
93  if (coupled)
94  {
95  Info<< "\nSolving energy coupled regions" << endl;
96  fvMatrixAssemblyPtr->solve();
97  #include "correctThermos.H"
98 
100  {
101  #include "setRegionFluidFields.H"
103  if (!frozenFlow)
104  {
105  #include "pEqn.H"
106  turb.correct();
107  }
108  }
109 
110  fvMatrixAssemblyPtr->clear();
111  }
112 
113  // Additional loops for energy solution only
114  {
115  loopControl looping(runTime, "SIMPLE", "energyCoupling");
116 
117  while (looping.loop())
118  {
119  Info<< nl << looping << nl;
120 
121  forAll(fluidRegions, i)
122  {
123  Info<< "\nSolving for fluid region "
124  << fluidRegions[i].name() << endl;
125  #include "setRegionFluidFields.H"
127  frozenFlow = true;
128  #include "solveFluid.H"
129  }
130 
131  forAll(solidRegions, i)
132  {
133  Info<< "\nSolving for solid region "
134  << solidRegions[i].name() << endl;
135  #include "setRegionSolidFields.H"
137  #include "solveSolid.H"
138  }
139 
140  if (coupled)
141  {
142  Info<< "\nSolving energy coupled regions.. " << endl;
143  fvMatrixAssemblyPtr->solve();
144  #include "correctThermos.H"
145 
146  fvMatrixAssemblyPtr->clear();
147  }
148  }
149  }
150 
151  runTime.write();
152 
153  runTime.printExecutionTime(Info);
154  }
155 
156  Info<< "End\n" << endl;
157 
158  return 0;
159 }
160 
161 
162 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
fvOptions.H
regionProperties.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
setRegionSolidFields.H
coordinateSystem.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
solidThermo.H
rhoThermo.H
setRootCaseLists.H
readSolidMultiRegionSIMPLEControls.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
createCoupledRegions.H
readFluidMultiRegionSIMPLEControls.H
postProcess.H
Execute application functionObjects to post-process existing results.
fixedGradientFvPatchFields.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
loopControl.H
fluidRegions
PtrList< fvMesh > fluidRegions(fluidNames.size())
solidRegions
PtrList< fvMesh > solidRegions(solidNames.size())
correctThermos.H
createTime.H
createMeshes.H
fvCFD.H
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
frozenFlow
bool frozenFlow
Definition: setRegionFluidFields.H:34
radiationModel.H
turbulentFluidThermoModel.H
fvMatrixAssemblyPtr
autoPtr< fvMatrix< scalar > > fvMatrixAssemblyPtr
Definition: createCoupledRegions.H:5