buoyantFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  buoyantFoam
26 
27 Description
28  Solver for steady or transient buoyant, turbulent flow of compressible
29  fluids for ventilation and heat-transfer, with optional mesh motion and mesh
30  topology changes.
31 
32  Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
33  pseudo-transient simulations.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
38 #include "fluidThermo.H"
41 #include "pimpleControl.H"
42 #include "pressureReference.H"
44 #include "CorrectPhi.H"
45 #include "fvModels.H"
46 #include "fvConstraints.H"
47 #include "localEulerDdtScheme.H"
48 #include "fvcSmooth.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 int main(int argc, char *argv[])
53 {
54  #include "postProcess.H"
55 
56  #include "setRootCaseLists.H"
57  #include "createTime.H"
58  #include "createMesh.H"
59  #include "createDyMControls.H"
60  #include "initContinuityErrs.H"
61  #include "createFields.H"
62  #include "createFieldRefs.H"
63  #include "createRhoUfIfPresent.H"
64 
65  turbulence->validate();
66 
67  if (!LTS)
68  {
69  #include "compressibleCourantNo.H"
70  #include "setInitialDeltaT.H"
71  }
72 
73  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75  Info<< "\nStarting time loop\n" << endl;
76 
77  while (pimple.run(runTime))
78  {
79  #include "readDyMControls.H"
80 
81  // Store divrhoU from the previous mesh so that it can be mapped
82  // and used in correctPhi to ensure the corrected phi has the
83  // same divergence
84  autoPtr<volScalarField> divrhoU;
85  if (correctPhi)
86  {
87  divrhoU = new volScalarField
88  (
89  "divrhoU",
91  );
92  }
93 
94  if (LTS)
95  {
96  #include "setRDeltaT.H"
97  }
98  else
99  {
100  #include "compressibleCourantNo.H"
101  #include "setDeltaT.H"
102  }
103 
105 
106  // Store momentum to set rhoUf for introduced faces.
107  autoPtr<volVectorField> rhoU;
108  if (rhoUf.valid())
109  {
110  rhoU = new volVectorField("rhoU", rho*U);
111  }
112 
113  // Update the mesh for topology change, mesh to mesh mapping
114  mesh.update();
115 
116  runTime++;
117 
118  Info<< "Time = " << runTime.userTimeName() << nl << endl;
119 
120  // --- Pressure-velocity PIMPLE corrector loop
121  while (pimple.loop())
122  {
123  if (!pimple.flow())
124  {
125  if (pimple.models())
126  {
127  fvModels.correct();
128  }
129 
130  if (pimple.thermophysics())
131  {
132  #include "EEqn.H"
133  }
134  }
135  else
136  {
137  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
138  {
139  // Move the mesh
140  mesh.move();
141 
142  if (mesh.changing())
143  {
144  gh = (g & mesh.C()) - ghRef;
145  ghf = (g & mesh.Cf()) - ghRef;
146 
147  MRF.update();
148 
149  if (correctPhi)
150  {
151  #include "correctPhi.H"
152  }
153 
154  if (checkMeshCourantNo)
155  {
156  #include "meshCourantNo.H"
157  }
158  }
159  }
160 
161  if (pimple.firstPimpleIter() && !pimple.simpleRho())
162  {
163  #include "rhoEqn.H"
164  }
165 
166  if (pimple.models())
167  {
168  fvModels.correct();
169  }
170 
171  #include "UEqn.H"
172 
173  if (pimple.thermophysics())
174  {
175  #include "EEqn.H"
176  }
177 
178  // --- Pressure corrector loop
179  while (pimple.correct())
180  {
181  #include "pEqn.H"
182  }
183 
184  if (pimple.turbCorr())
185  {
186  turbulence->correct();
187  thermophysicalTransport->correct();
188  }
189  }
190  }
191 
192  rho = thermo.rho();
193 
194  runTime.write();
195 
196  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
197  << " ClockTime = " << runTime.elapsedClockTime() << " s"
198  << nl << endl;
199  }
200 
201  Info<< "End\n" << endl;
202 
203  return 0;
204 }
205 
206 
207 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
U
Definition: pEqn.H:72
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:353
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:255
correctPhi
checkMeshCourantNo
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Creates and initialises the velocity field rhoUf if required.
phi
Definition: correctPhi.H:3
Info<< "Reading field U\"<< 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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
const surfaceScalarField & ghf
static const char nl
Definition: Ostream.H:260
bool LTS
Definition: createRDeltaT.H:1
fluidReactionThermophysicalTransportModel & thermophysicalTransport
moveMeshOuterCorrectors
Foam::fvModels & fvModels
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
const volScalarField & gh
messageStream Info
const dimensionedVector & g
Execute application functionObjects to post-process existing results.
autoPtr< surfaceVectorField > rhoUf
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...