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-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  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 "pressureControl.H"
43 #include "CorrectPhi.H"
44 #include "fvOptions.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 "createDynamicFvMesh.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 
102  runTime++;
103 
104  Info<< "Time = " << runTime.timeName() << nl << endl;
105 
106  // --- Pressure-velocity PIMPLE corrector loop
107  while (pimple.loop())
108  {
109  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
110  {
111  // Store momentum to set rhoUf for introduced faces.
112  autoPtr<volVectorField> rhoU;
113  if (rhoUf.valid())
114  {
115  rhoU = new volVectorField("rhoU", rho*U);
116  }
117 
118  // Do any mesh changes
119  mesh.update();
120 
121  if (mesh.changing())
122  {
123  MRF.update();
124 
125  if (correctPhi)
126  {
127  // Calculate absolute flux
128  // from the mapped surface velocity
129  phi = mesh.Sf() & rhoUf();
130 
131  #include "correctPhi.H"
132 
133  // Make the fluxes relative to the mesh-motion
135  }
136 
137  if (checkMeshCourantNo)
138  {
139  #include "meshCourantNo.H"
140  }
141  }
142  }
143 
144  if
145  (
146  !mesh.steady()
147  && !pimple.simpleRho()
148  && pimple.firstPimpleIter()
149  )
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  if (!mesh.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 // ************************************************************************* //
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
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
messageStream Info
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...