overRhoPimpleDyMFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  overRhoPimpleDyMFoam
29 
30 Group
31  grpCompressibleSolvers grpMovingMeshSolvers
32 
33 Description
34  Transient solver for laminar or turbulent flow of compressible fluids
35  for HVAC and similar applications.
36 
37  Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
38  pseudo-transient simulations.
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #include "fvCFD.H"
43 #include "dynamicFvMesh.H"
44 #include "fluidThermo.H"
46 #include "pimpleControl.H"
47 #include "pressureControl.H"
48 #include "CorrectPhi.H"
49 #include "fvOptions.H"
50 #include "localEulerDdtScheme.H"
51 #include "fvcSmooth.H"
52 #include "cellCellStencilObject.H"
53 #include "localMin.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 int main(int argc, char *argv[])
58 {
59  argList::addNote
60  (
61  "Transient solver for compressible turbulent flow.\n"
62  "With optional mesh motion and mesh topology changes."
63  );
64 
65  #include "setRootCaseLists.H"
66  #include "createTime.H"
67  #include "createDynamicFvMesh.H"
68  #include "createDyMControls.H"
69  #include "createRDeltaT.H"
70  #include "initContinuityErrs.H"
71  #include "createFields.H"
72  #include "createMRF.H"
73  #include "createFvOptions.H"
74  #include "createRhoUfIfPresent.H"
75  #include "createControls.H"
76 
77  turbulence->validate();
78 
79  if (!LTS)
80  {
81  #include "compressibleCourantNo.H"
82  #include "setInitialDeltaT.H"
83  }
84 
85  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87  Info<< "\nStarting time loop\n" << endl;
88 
89  while (runTime.run())
90  {
91  #include "readDyMControls.H"
92 
93  // Store divrhoU from the previous mesh so that it can be mapped
94  // and used in correctPhi to ensure the corrected phi has the
95  // same divergence
96  autoPtr<volScalarField> divrhoU;
97  if (correctPhi)
98  {
99  divrhoU.reset
100  (
101  new volScalarField
102  (
103  "divrhoU",
105  )
106  );
107  }
108 
109  if (LTS)
110  {
111  #include "setRDeltaT.H"
112  }
113  else
114  {
115  #include "compressibleCourantNo.H"
116  #include "setDeltaT.H"
117  }
118 
119  ++runTime;
120 
121  Info<< "Time = " << runTime.timeName() << nl << endl;
122 
123  // --- Pressure-velocity PIMPLE corrector loop
124  while (pimple.loop())
125  {
126  if (pimple.firstIter() || moveMeshOuterCorrectors)
127  {
128  // Do any mesh changes
129  mesh.update();
130 
131  if (mesh.changing())
132  {
133  MRF.update();
134 
135  #include "setCellMask.H"
136  #include "setInterpolatedCells.H"
137  #include "correctRhoPhiFaceMask.H"
138 
139  if (correctPhi)
140  {
141  // Corrects flux on separated regions
142  #include "correctPhi.H"
143  }
144 
145  // Make the fluxes relative to the mesh-motion
147 
148  if (checkMeshCourantNo)
149  {
150  #include "meshCourantNo.H"
151  }
152  }
153  }
154 
155  if (pimple.firstIter() && !pimple.SIMPLErho())
156  {
157  #include "rhoEqn.H"
158  }
159 
160  #include "UEqn.H"
161  #include "EEqn.H"
162 
163  // --- Pressure corrector loop
164  while (pimple.correct())
165  {
166  #include "pEqn.H"
167  }
168 
169  if (pimple.turbCorr())
170  {
171  turbulence->correct();
172  }
173  }
174 
175  rho = thermo.rho();
176 
177  runTime.write();
178 
179  runTime.printExecutionTime(Info);
180  }
181 
182  Info<< "End\n" << endl;
183 
184  return 0;
185 }
186 
187 
188 // ************************************************************************* //
Sets blocked cells mask field.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
correctPhi
checkMeshCourantNo
Sets blocked cells mask field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
psiReactionThermo & thermo
Definition: createFields.H:28
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
IOMRFZoneList & MRF
bool LTS
Definition: createRDeltaT.H:1
pimpleControl & pimple
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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", 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
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:183
U
Definition: pEqn.H:72
messageStream Info
Information stream (stdout output on master, null elsewhere)
Correct phi on new faces C-I faces.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:70
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...