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-2019 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"
40 #include "pimpleControl.H"
41 #include "pressureControl.H"
42 #include "CorrectPhi.H"
43 #include "fvOptions.H"
44 #include "localEulerDdtScheme.H"
45 #include "fvcSmooth.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
51  #include "postProcess.H"
52 
53  #include "setRootCaseLists.H"
54  #include "createTime.H"
55  #include "createDynamicFvMesh.H"
56  #include "createDyMControls.H"
57  #include "initContinuityErrs.H"
58  #include "createFields.H"
59  #include "createFieldRefs.H"
60  #include "createRhoUfIfPresent.H"
61 
62  turbulence->validate();
63 
64  if (!LTS)
65  {
66  #include "compressibleCourantNo.H"
67  #include "setInitialDeltaT.H"
68  }
69 
70  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72  Info<< "\nStarting time loop\n" << endl;
73 
74  while (runTime.run())
75  {
76  #include "readDyMControls.H"
77 
78  // Store divrhoU from the previous mesh so that it can be mapped
79  // and used in correctPhi to ensure the corrected phi has the
80  // same divergence
81  autoPtr<volScalarField> divrhoU;
82  if (correctPhi)
83  {
84  divrhoU = new volScalarField
85  (
86  "divrhoU",
88  );
89  }
90 
91  if (LTS)
92  {
93  #include "setRDeltaT.H"
94  }
95  else
96  {
97  #include "compressibleCourantNo.H"
98  #include "setDeltaT.H"
99  }
100 
101  runTime++;
102 
103  Info<< "Time = " << runTime.timeName() << nl << endl;
104 
105  // --- Pressure-velocity PIMPLE corrector loop
106  while (pimple.loop())
107  {
108  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
109  {
110  // Store momentum to set rhoUf for introduced faces.
111  autoPtr<volVectorField> rhoU;
112  if (rhoUf.valid())
113  {
114  rhoU = new volVectorField("rhoU", rho*U);
115  }
116 
117  // Do any mesh changes
118  mesh.update();
119 
120  if (mesh.changing())
121  {
122  MRF.update();
123 
124  if (correctPhi)
125  {
126  // Calculate absolute flux
127  // from the mapped surface velocity
128  phi = mesh.Sf() & rhoUf();
129 
130  #include "correctPhi.H"
131 
132  // Make the fluxes relative to the mesh-motion
134  }
135 
136  if (checkMeshCourantNo)
137  {
138  #include "meshCourantNo.H"
139  }
140  }
141  }
142 
143  if (pimple.firstPimpleIter() && !pimple.simpleRho())
144  {
145  #include "rhoEqn.H"
146  }
147 
148  #include "UEqn.H"
149  #include "EEqn.H"
150 
151  // --- Pressure corrector loop
152  while (pimple.correct())
153  {
154  if (pimple.consistent())
155  {
156  #include "pcEqn.H"
157  }
158  else
159  {
160  #include "pEqn.H"
161  }
162  }
163 
164  if (pimple.turbCorr())
165  {
166  turbulence->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
surfaceScalarField & phi
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:256
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::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
correctPhi
checkMeshCourantNo
rhoReactionThermo & thermo
Definition: createFields.H:28
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:265
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...