All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
pimpleFoam.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  pimpleFoam
26 
27 Description
28  Transient solver for incompressible, turbulent flow of Newtonian fluids,
29  with optional mesh motion and mesh topology changes.
30 
31  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "fvCFD.H"
36 #include "viscosityModel.H"
38 #include "pimpleControl.H"
39 #include "pressureReference.H"
40 #include "CorrectPhi.H"
41 #include "fvModels.H"
42 #include "fvConstraints.H"
43 #include "localEulerDdtScheme.H"
44 #include "fvcSmooth.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 int main(int argc, char *argv[])
49 {
50  #include "postProcess.H"
51 
52  #include "setRootCaseLists.H"
53  #include "createTime.H"
54  #include "createMesh.H"
55  #include "initContinuityErrs.H"
56  #include "createDyMControls.H"
57  #include "createFields.H"
58  #include "createUfIfPresent.H"
59 
60  turbulence->validate();
61 
62  if (!LTS)
63  {
64  #include "CourantNo.H"
65  #include "setInitialDeltaT.H"
66  }
67 
68  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70  Info<< "\nStarting time loop\n" << endl;
71 
72  while (pimple.run(runTime))
73  {
74  #include "readDyMControls.H"
75 
76  if (LTS)
77  {
78  #include "setRDeltaT.H"
79  }
80  else
81  {
82  #include "CourantNo.H"
83  #include "setDeltaT.H"
84  }
85 
87 
88  // Update the mesh for topology change, mesh to mesh mapping
89  mesh.update();
90 
91  runTime++;
92 
93  Info<< "Time = " << runTime.userTimeName() << nl << endl;
94 
95  // --- Pressure-velocity PIMPLE corrector loop
96  while (pimple.loop())
97  {
98  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
99  {
100  // Move the mesh
101  mesh.move();
102 
103  if (mesh.changing())
104  {
105  MRF.update();
106 
107  if (correctPhi)
108  {
109  #include "correctPhi.H"
110  }
111 
112  if (checkMeshCourantNo)
113  {
114  #include "meshCourantNo.H"
115  }
116  }
117  }
118 
119  fvModels.correct();
120 
121  #include "UEqn.H"
122 
123  // --- Pressure corrector loop
124  while (pimple.correct())
125  {
126  #include "pEqn.H"
127  }
128 
129  if (pimple.turbCorr())
130  {
131  viscosity->correct();
132  turbulence->correct();
133  }
134  }
135 
136  runTime.write();
137 
138  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
139  << " ClockTime = " << runTime.elapsedClockTime() << " s"
140  << nl << endl;
141  }
142 
143  Info<< "End\n" << endl;
144 
145  return 0;
146 }
147 
148 
149 // ************************************************************************* //
pimpleNoLoopControl & pimple
IOMRFZoneList & MRF
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
fvMesh & mesh
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
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
Foam::fvModels & fvModels
Calculates and outputs the mean and maximum Courant Numbers.
messageStream Info
Execute application functionObjects to post-process existing results.
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...
Creates and initialises the velocity field Uf if required.