buoyantPimpleFoam.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-2020 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  buoyantPimpleFoam
26 
27 Description
28  Transient solver for buoyant, turbulent flow of compressible fluids for
29  ventilation and heat-transfer, with optional mesh motion and
30  mesh 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 "dynamicFvMesh.H"
39 #include "fluidThermo.H"
42 #include "radiationModel.H"
43 #include "pimpleControl.H"
44 #include "pressureControl.H"
45 #include "CorrectPhi.H"
46 #include "fvOptions.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 "createDynamicFvMesh.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 
104  runTime++;
105 
106  Info<< "Time = " << runTime.timeName() << nl << endl;
107 
108  // --- Pressure-velocity PIMPLE corrector loop
109  while (pimple.loop())
110  {
111  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
112  {
113  // Store momentum to set rhoUf for introduced faces.
114  autoPtr<volVectorField> rhoU;
115  if (rhoUf.valid())
116  {
117  rhoU = new volVectorField("rhoU", rho*U);
118  }
119 
120  // Do any mesh changes
121  mesh.update();
122 
123  if (mesh.changing())
124  {
125  gh = (g & mesh.C()) - ghRef;
126  ghf = (g & mesh.Cf()) - ghRef;
127 
128  MRF.update();
129 
130  if (correctPhi)
131  {
132  // Calculate absolute flux
133  // from the mapped surface velocity
134  phi = mesh.Sf() & rhoUf();
135 
136  #include "correctPhi.H"
137 
138  // Make the fluxes relative to the mesh-motion
140  }
141 
142  if (checkMeshCourantNo)
143  {
144  #include "meshCourantNo.H"
145  }
146  }
147  }
148 
149  if (pimple.firstPimpleIter() && !pimple.simpleRho())
150  {
151  #include "rhoEqn.H"
152  }
153 
154  #include "UEqn.H"
155  #include "EEqn.H"
156 
157  // --- Pressure corrector loop
158  while (pimple.correct())
159  {
160  #include "pEqn.H"
161  }
162 
163  if (pimple.turbCorr())
164  {
165  turbulence->correct();
166  thermophysicalTransport->correct();
167  }
168  }
169 
170  rho = thermo.rho();
171 
172  runTime.write();
173 
174  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
175  << " ClockTime = " << runTime.elapsedClockTime() << " s"
176  << nl << endl;
177  }
178 
179  Info<< "End\n" << endl;
180 
181  return 0;
182 }
183 
184 
185 // ************************************************************************* //
pimpleNoLoopControl & pimple
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
correctPhi
checkMeshCourantNo
rhoReactionThermo & thermo
Definition: createFields.H:28
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
phi
Definition: pEqn.H:104
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
rhoReactionThermophysicalTransportModel & thermophysicalTransport
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
moveMeshOuterCorrectors
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:188
U
Definition: pEqn.H:72
const volScalarField & gh
messageStream Info
const dimensionedVector & g
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
autoPtr< surfaceVectorField > rhoUf
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...