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-2021 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 "dynamicFvMesh.H"
38 #include "fluidThermo.H"
41 #include "pimpleControl.H"
42 #include "pressureReference.H"
43 #include "CorrectPhi.H"
44 #include "fvModels.H"
45 #include "fvConstraints.H"
46 #include "localEulerDdtScheme.H"
47 #include "fvcSmooth.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 int main(int argc, char *argv[])
52 {
53  #include "postProcess.H"
54 
55  #include "setRootCaseLists.H"
56  #include "createTime.H"
57  #include "createDynamicFvMesh.H"
58  #include "createDyMControls.H"
59  #include "initContinuityErrs.H"
60  #include "createFields.H"
61  #include "createFieldRefs.H"
62  #include "createRhoUfIfPresent.H"
63 
64  turbulence->validate();
65 
66  if (!LTS)
67  {
68  #include "compressibleCourantNo.H"
69  #include "setInitialDeltaT.H"
70  }
71 
72  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74  Info<< "\nStarting time loop\n" << endl;
75 
76  while (pimple.run(runTime))
77  {
78  #include "readDyMControls.H"
79 
80  // Store divrhoU from the previous mesh so that it can be mapped
81  // and used in correctPhi to ensure the corrected phi has the
82  // same divergence
83  autoPtr<volScalarField> divrhoU;
84  if (correctPhi)
85  {
86  divrhoU = new volScalarField
87  (
88  "divrhoU",
90  );
91  }
92 
93  if (LTS)
94  {
95  #include "setRDeltaT.H"
96  }
97  else
98  {
99  #include "compressibleCourantNo.H"
100  #include "setDeltaT.H"
101  }
102 
103  runTime++;
104 
105  Info<< "Time = " << runTime.timeName() << nl << endl;
106 
107  // --- Pressure-velocity PIMPLE corrector loop
108  while (pimple.loop())
109  {
110  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
111  {
112  // Store momentum to set rhoUf for introduced faces.
113  autoPtr<volVectorField> rhoU;
114  if (rhoUf.valid())
115  {
116  rhoU = new volVectorField("rhoU", rho*U);
117  }
118 
120 
121  // Do any mesh changes
122  mesh.update();
123 
124  if (mesh.changing())
125  {
126  MRF.update();
127 
128  if (correctPhi)
129  {
130  #include "correctPhi.H"
131  }
132 
133  if (checkMeshCourantNo)
134  {
135  #include "meshCourantNo.H"
136  }
137  }
138  }
139 
140  if
141  (
142  !mesh.steady()
143  && !pimple.simpleRho()
144  && pimple.firstPimpleIter()
145  )
146  {
147  #include "rhoEqn.H"
148  }
149 
150  fvModels.correct();
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  thermophysicalTransport->correct();
165  }
166  }
167 
168  if (!mesh.steady())
169  {
170  rho = thermo.rho();
171  }
172 
173  runTime.write();
174 
175  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
176  << " ClockTime = " << runTime.elapsedClockTime() << " s"
177  << nl << endl;
178  }
179 
180  Info<< "End\n" << endl;
181 
182  return 0;
183 }
184 
185 
186 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
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:316
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:240
correctPhi
checkMeshCourantNo
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
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:188
U
Definition: pEqn.H:72
messageStream Info
Execute application functionObjects to post-process existing results.
autoPtr< surfaceVectorField > rhoUf
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...