rhoPimpleFoam.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  rhoPimpleFoam
26 
27 Description
28  Transient solver for turbulent flow of compressible fluids for HVAC and
29  similar applications, with optional mesh motion and mesh topology changes.
30 
31  Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
32  pseudo-transient simulations.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "fluidThermo.H"
40 #include "pimpleControl.H"
41 #include "pressureReference.H"
42 #include "CorrectPhi.H"
43 #include "fvModels.H"
44 #include "fvConstraints.H"
45 #include "localEulerDdtScheme.H"
46 #include "fvcSmooth.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 int main(int argc, char *argv[])
51 {
52  #include "postProcess.H"
53 
54  #include "setRootCaseLists.H"
55  #include "createTime.H"
56  #include "createMesh.H"
57  #include "createDyMControls.H"
58  #include "initContinuityErrs.H"
59  #include "createFields.H"
60  #include "createFieldRefs.H"
61  #include "createRhoUfIfPresent.H"
62 
63  turbulence->validate();
64 
65  if (!LTS)
66  {
67  #include "compressibleCourantNo.H"
68  #include "setInitialDeltaT.H"
69  }
70 
71  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73  Info<< "\nStarting time loop\n" << endl;
74 
75  while (pimple.run(runTime))
76  {
77  #include "readDyMControls.H"
78 
79  // Store divrhoU from the previous mesh so that it can be mapped
80  // and used in correctPhi to ensure the corrected phi has the
81  // same divergence
82  autoPtr<volScalarField> divrhoU;
83  if (correctPhi)
84  {
85  divrhoU = new volScalarField
86  (
87  "divrhoU",
89  );
90  }
91 
92  if (LTS)
93  {
94  #include "setRDeltaT.H"
95  }
96  else
97  {
98  #include "compressibleCourantNo.H"
99  #include "setDeltaT.H"
100  }
101 
103 
104  // Store momentum to set rhoUf for introduced faces.
105  autoPtr<volVectorField> rhoU;
106  if (rhoUf.valid())
107  {
108  rhoU = new volVectorField("rhoU", rho*U);
109  }
110 
111  // Update the mesh for topology change, mesh to mesh mapping
112  mesh.update();
113 
114  runTime++;
115 
116  Info<< "Time = " << runTime.userTimeName() << nl << endl;
117 
118  // --- Pressure-velocity PIMPLE corrector loop
119  while (pimple.loop())
120  {
121  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
122  {
123  // Move the mesh
124  mesh.move();
125 
126  if (mesh.changing())
127  {
128  MRF.update();
129 
130  if (correctPhi)
131  {
132  #include "correctPhi.H"
133  }
134 
135  if (checkMeshCourantNo)
136  {
137  #include "meshCourantNo.H"
138  }
139  }
140  }
141 
142  if
143  (
144  !mesh.schemes().steady()
145  && !pimple.simpleRho()
146  && pimple.firstPimpleIter()
147  )
148  {
149  #include "rhoEqn.H"
150  }
151 
152  fvModels.correct();
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  if (!mesh.schemes().steady())
171  {
172  rho = thermo.rho();
173  }
174 
175  runTime.write();
176 
177  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
178  << " ClockTime = " << runTime.elapsedClockTime() << " s"
179  << nl << endl;
180  }
181 
182  Info<< "End\n" << endl;
183 
184  return 0;
185 }
186 
187 
188 // ************************************************************************* //
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
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
messageStream Info
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...