All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
cavitatingFoam.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  cavitatingFoam
26 
27 Description
28  Transient cavitation code based on the homogeneous equilibrium model
29  from which the compressibility of the liquid/vapour "mixture" is obtained,
30  with optional mesh motion and mesh topology changes.
31 
32  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "dynamicFvMesh.H"
41 #include "CorrectPhi.H"
42 #include "pimpleControl.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 int main(int argc, char *argv[])
47 {
48  #include "postProcess.H"
49 
50  #include "setRootCaseLists.H"
51  #include "createTime.H"
52  #include "createDynamicFvMesh.H"
53  #include "createControls.H"
54  #include "createFields.H"
55  #include "createUfIfPresent.H"
56  #include "createPcorrTypes.H"
57  #include "CourantNo.H"
58  #include "setInitialDeltaT.H"
59 
60  turbulence->validate();
61 
62  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64  Info<< "\nStarting time loop\n" << endl;
65 
66  while (pimple.run(runTime))
67  {
68  #include "readControls.H"
69 
70  {
71  #include "CourantNo.H"
72  #include "setDeltaT.H"
73 
74  runTime++;
75 
76  Info<< "Time = " << runTime.timeName() << nl << endl;
77 
78  // Do any mesh changes
79  mesh.update();
80 
81  if (mesh.changing() && correctPhi)
82  {
83  #include "correctPhi.H"
84  }
85  }
86 
87  // --- Pressure-velocity PIMPLE corrector loop
88  while (pimple.loop())
89  {
90  #include "rhoEqn.H"
91  #include "alphavPsi.H"
92  #include "UEqn.H"
93 
94  // --- Pressure corrector loop
95  while (pimple.correct())
96  {
97  #include "pEqn.H"
98  }
99 
100  if (pimple.turbCorr())
101  {
102  turbulence->correct();
103  }
104  }
105 
106  runTime.write();
107 
108  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
109  << " ClockTime = " << runTime.elapsedClockTime() << " s"
110  << nl << endl;
111  }
112 
113  Info<< "End\n" << endl;
114 
115  return 0;
116 }
117 
118 
119 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
correctPhi
dynamicFvMesh & 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
messageStream Info
Execute application functionObjects to post-process existing results.
Creates and initialises the velocity field Uf if required.