interFoam.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  interFoam
26 
27 Description
28  Solver for 2 incompressible, isothermal immiscible fluids using a VOF
29  (volume of fluid) phase-fraction based interface capturing approach,
30  with optional mesh motion and mesh topology changes including adaptive
31  re-meshing.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "fvCFD.H"
36 #include "interfaceCompression.H"
37 #include "CMULES.H"
38 #include "EulerDdtScheme.H"
39 #include "localEulerDdtScheme.H"
40 #include "CrankNicolsonDdtScheme.H"
41 #include "subCycle.H"
43 #include "noPhaseChange.H"
45 #include "pimpleControl.H"
46 #include "pressureReference.H"
47 #include "fvModels.H"
48 #include "fvConstraints.H"
49 #include "CorrectPhi.H"
50 #include "fvcSmooth.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 int main(int argc, char *argv[])
55 {
56  #include "postProcess.H"
57 
58  #include "setRootCaseLists.H"
59  #include "createTime.H"
60  #include "createMesh.H"
61  #include "initContinuityErrs.H"
62  #include "createDyMControls.H"
63  #include "createFields.H"
64  #include "createFieldRefs.H"
65  #include "initCorrectPhi.H"
66  #include "createUfIfPresent.H"
67 
68  if (!LTS)
69  {
70  #include "CourantNo.H"
71  #include "setInitialDeltaT.H"
72  }
73 
74  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75  Info<< "\nStarting time loop\n" << endl;
76 
77  while (pimple.run(runTime))
78  {
79  #include "readDyMControls.H"
80 
81  if (LTS)
82  {
83  #include "setRDeltaT.H"
84  }
85  else
86  {
87  #include "CourantNo.H"
88  #include "alphaCourantNo.H"
89  #include "setDeltaT.H"
90  }
91 
93 
94  // Store divU from the previous mesh so that it can be mapped
95  // and used in correctPhi to ensure the corrected phi has the
96  // same divergence
97  tmp<volScalarField> divU;
98 
99  if
100  (
101  correctPhi
102  && !isType<twoPhaseChangeModels::noPhaseChange>(phaseChange)
103  && mesh.topoChanged()
104  )
105  {
106  // Construct and register divU for correctPhi
107  divU = new volScalarField
108  (
109  "divU0",
111  );
112  }
113 
114  // Update the mesh for topology change, mesh to mesh mapping
115  bool topoChanged = mesh.update();
116 
117  // Do not apply previous time-step mesh compression flux
118  // if the mesh topology changed
119  if (topoChanged)
120  {
121  talphaPhi1Corr0.clear();
122  }
123 
124  runTime++;
125 
126  Info<< "Time = " << runTime.userTimeName() << nl << endl;
127 
128  // --- Pressure-velocity PIMPLE corrector loop
129  while (pimple.loop())
130  {
131  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
132  {
133  if
134  (
135  correctPhi
136  && !isType<twoPhaseChangeModels::noPhaseChange>(phaseChange)
137  && !divU.valid()
138  )
139  {
140  // Construct and register divU for correctPhi
141  divU = new volScalarField
142  (
143  "divU0",
145  );
146  }
147 
148  // Move the mesh
149  mesh.move();
150 
151  if (mesh.changing())
152  {
153  gh = (g & mesh.C()) - ghRef;
154  ghf = (g & mesh.Cf()) - ghRef;
155 
156  MRF.update();
157 
158  if (correctPhi)
159  {
160  #include "correctPhi.H"
161  }
162 
163  mixture.correct();
164 
165  if (checkMeshCourantNo)
166  {
167  #include "meshCourantNo.H"
168  }
169  }
170 
171  divU.clear();
172  }
173 
174  fvModels.correct();
175 
177  (
178  IOobject
179  (
180  "rhoPhi",
181  runTime.timeName(),
182  mesh
183  ),
184  mesh,
186  );
187 
188  #include "alphaControls.H"
189  #include "alphaEqnSubCycle.H"
190 
191  turbulence.correctPhasePhi();
192 
193  mixture.correct();
194 
195  #include "UEqn.H"
196 
197  // --- Pressure corrector loop
198  while (pimple.correct())
199  {
200  #include "pEqn.H"
201  }
202 
203  if (pimple.turbCorr())
204  {
205  turbulence.correct();
206  }
207  }
208 
209  runTime.write();
210 
211  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
212  << " ClockTime = " << runTime.elapsedClockTime() << " s"
213  << nl << endl;
214  }
215 
216  Info<< "End\n" << endl;
217 
218  return 0;
219 }
220 
221 
222 // ************************************************************************* //
pimpleNoLoopControl & pimple
U
Definition: pEqn.H:72
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:353
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
twoPhaseChangeModel & phaseChange
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:255
correctPhi
checkMeshCourantNo
fvMesh & mesh
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet dimTime
phi
Definition: correctPhi.H:3
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
const surfaceScalarField & ghf
static const char nl
Definition: Ostream.H:260
bool LTS
Definition: createRDeltaT.H:1
moveMeshOuterCorrectors
const dimensionSet dimMass
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
Foam::fvModels & fvModels
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:202
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
rhoPhi
Definition: rhoEqn.H:10
const volScalarField & gh
messageStream Info
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedVector & g
Execute application functionObjects to post-process existing results.
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Creates and initialises the velocity field Uf if required.
tmp< surfaceScalarField > talphaPhi1Corr0