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-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  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 "dynamicFvMesh.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 "createDynamicFvMesh.H"
61  #include "initContinuityErrs.H"
62  #include "createDyMControls.H"
63  #include "createFields.H"
64  #include "createFieldRefs.H"
65  #include "createAlphaFluxes.H"
66  #include "initCorrectPhi.H"
67  #include "createUfIfPresent.H"
68 
69  turbulence->validate();
70 
71  if (!LTS)
72  {
73  #include "CourantNo.H"
74  #include "setInitialDeltaT.H"
75  }
76 
77  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78  Info<< "\nStarting time loop\n" << endl;
79 
80  while (pimple.run(runTime))
81  {
82  #include "readDyMControls.H"
83 
84  if (LTS)
85  {
86  #include "setRDeltaT.H"
87  }
88  else
89  {
90  #include "CourantNo.H"
91  #include "alphaCourantNo.H"
92  #include "setDeltaT.H"
93  }
94 
95  runTime++;
96 
97  Info<< "Time = " << runTime.timeName() << nl << endl;
98 
99  // --- Pressure-velocity PIMPLE corrector loop
100  while (pimple.loop())
101  {
102  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
103  {
104  // Store divU from the previous mesh so that it can be mapped
105  // and used in correctPhi to ensure the corrected phi has the
106  // same divergence
107  tmp<volScalarField> divU;
108 
109  if
110  (
111  correctPhi
112  && !isType<twoPhaseChangeModels::noPhaseChange>(phaseChange)
113  )
114  {
115  // Construct and register divU for mapping
116  divU = new volScalarField
117  (
118  "divU0",
120  );
121  }
122 
124 
125  mesh.update();
126 
127  if (mesh.changing())
128  {
129  // Do not apply previous time-step mesh compression flux
130  // if the mesh topology changed
131  if (mesh.topoChanging())
132  {
133  talphaPhi1Corr0.clear();
134  }
135 
136  gh = (g & mesh.C()) - ghRef;
137  ghf = (g & mesh.Cf()) - ghRef;
138 
139  MRF.update();
140 
141  if (correctPhi)
142  {
143  #include "correctPhi.H"
144  }
145 
146  mixture.correct();
147 
148  if (checkMeshCourantNo)
149  {
150  #include "meshCourantNo.H"
151  }
152  }
153 
154  divU.clear();
155  }
156 
157  fvModels.correct();
158 
160  (
161  IOobject
162  (
163  "rhoPhi",
164  runTime.timeName(),
165  mesh
166  ),
167  mesh,
169  );
170 
171  #include "alphaControls.H"
172  #include "alphaEqnSubCycle.H"
173 
174  mixture.correct();
175 
176  #include "UEqn.H"
177 
178  // --- Pressure corrector loop
179  while (pimple.correct())
180  {
181  #include "pEqn.H"
182  }
183 
184  if (pimple.turbCorr())
185  {
186  turbulence->correct();
187  }
188  }
189 
190  runTime.write();
191 
192  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
193  << " ClockTime = " << runTime.elapsedClockTime() << " s"
194  << nl << endl;
195  }
196 
197  Info<< "End\n" << endl;
198 
199  return 0;
200 }
201 
202 
203 // ************************************************************************* //
pimpleNoLoopControl & pimple
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:316
engineTime & runTime
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:240
correctPhi
checkMeshCourantNo
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimTime
dynamicFvMesh & mesh
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:188
U
Definition: pEqn.H:72
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 FaceCellWave algorithm to smooth and redistr...
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Creates and initialises the velocity field Uf if required.
tmp< surfaceScalarField > talphaPhi1Corr0