Loading...
Searching...
No Matches
overBuoyantPimpleDyMFoam.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) 2019-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Application
27 overBuoyantPimpleDymFoam
28
29Group
30 grpHeatTransferSolvers
31
32Description
33 Transient solver for buoyant, turbulent flow of compressible fluids for
34 ventilation and heat-transfer with overset feature
35
36 Turbulence is modelled using a run-time selectable compressible RAS or
37 LES model.
38
39\*---------------------------------------------------------------------------*/
40
41#include "fvCFD.H"
42#include "dynamicFvMesh.H"
43#include "rhoThermo.H"
45#include "radiationModel.H"
46#include "fvOptions.H"
47#include "pimpleControl.H"
48#include "pressureControl.H"
49
50#include "CorrectPhi.H"
52#include "localMin.H"
53#include "oversetAdjustPhi.H"
54
55// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57int main(int argc, char *argv[])
58{
59 argList::addNote
60 (
61 "Transient solver for buoyant, turbulent fluid flow"
62 " of compressible fluids, including radiation."
63 );
64
65 #include "postProcess.H"
66
67 #include "addCheckCaseOptions.H"
68 #include "setRootCaseLists.H"
69 #include "createTime.H"
70 #include "createDynamicFvMesh.H"
71 #include "createDyMControls.H"
72 #include "createFields.H"
73 #include "createFieldRefs.H"
74 #include "initContinuityErrs.H"
75
76 #include "createRhoUfIfPresent.H"
77 #include "createControls.H"
78
79 #include "compressibleCourantNo.H"
80 #include "setInitialDeltaT.H"
81
82 turbulence->validate();
83
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85
86 Info<< "\nStarting time loop\n" << endl;
87
88 while (runTime.run())
89 {
90 #include "readDyMControls.H"
91
92 #include "compressibleCourantNo.H"
93 #include "setDeltaT.H"
94
95 // Store divrhoU from the previous mesh so that it can be mapped
96 // and used in correctPhi to ensure the corrected phi has the
97 // same divergence
98 autoPtr<volScalarField> divrhoU;
99 if (correctPhi)
100 {
101 divrhoU.reset
102 (
103 new volScalarField
104 (
105 "divrhoU",
106 fvc::div(fvc::absolute(phi, rho, U))
107 )
108 );
109 }
110
111
112 ++runTime;
113
114 Info<< "Time = " << runTime.timeName() << nl << endl;
115
116 // --- Pressure-velocity PIMPLE corrector loop
117 while (pimple.loop())
118 {
119 if (pimple.firstIter() || moveMeshOuterCorrectors)
120 {
121 // Do any mesh changes
122 mesh.update();
123
124 if (mesh.changing())
125 {
126 MRF.update();
127
128 #include "setCellMask.H"
129 #include "setInterpolatedCells.H"
130 #include "correctRhoPhiFaceMask.H"
131
132 if (correctPhi)
133 {
134 #include "correctPhi.H"
135 }
136
137 // Make the fluxes relative to the mesh-motion
138 fvc::makeRelative(phi, rho, U);
139 }
140
142 {
143 #include "meshCourantNo.H"
144 }
145 }
146
147 if (pimple.firstIter())
148 {
149 #include "rhoEqn.H"
150 }
151
152 #include "UEqn.H"
153 #include "EEqn.H"
154
155 // --- Pressure corrector loop
156 while (pimple.correct())
157 {
158 #include "pEqn.H"
159 }
160
161 if (pimple.turbCorr())
162 {
163 turbulence->correct();
164 }
165 }
166
167 rho = thermo.rho();
168
169 runTime.write();
170
171 runTime.printExecutionTime(Info);
172 }
173
174 Info<< "End\n" << endl;
175
176 return 0;
177}
178
179
180// ************************************************************************* //
Required Classes.
IOMRFZoneList & MRF
pimpleControl & pimple
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
U
Definition pEqn.H:72
Correct phi on new faces C-I faces.
dynamicFvMesh & mesh
engineTime & runTime
Creates and initialises the velocity field rhoUf if required.
compressible::turbulenceModel & turbulence
Calculates and outputs the mean and maximum Courant Numbers.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity.
Execute application functionObjects to post-process existing results.
checkMeshCourantNo
correctPhi
moveMeshOuterCorrectors
Sets blocked cells mask field.
Sets blocked cells mask field.