All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
multiphaseInterFoam.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  multiphaseInterFoam
26 
27 Description
28  Solver for n incompressible fluids which captures the interfaces and
29  includes surface-tension and contact-angle effects for each phase, with
30  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 "multiphaseMixture.H"
39 #include "pimpleControl.H"
40 #include "pressureReference.H"
41 #include "fvModels.H"
42 #include "fvConstraints.H"
43 #include "CorrectPhi.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  #include "postProcess.H"
50 
51  #include "setRootCaseLists.H"
52  #include "createTime.H"
53  #include "createMesh.H"
54  #include "initContinuityErrs.H"
55  #include "createDyMControls.H"
56  #include "createFields.H"
57  #include "initCorrectPhi.H"
58  #include "createUfIfPresent.H"
59 
60  turbulence->validate();
61 
62  #include "CourantNo.H"
63  #include "setInitialDeltaT.H"
64 
65  const surfaceScalarField& rhoPhi(mixture.rhoPhi());
66 
67  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69  Info<< "\nStarting time loop\n" << endl;
70 
71  while (pimple.run(runTime))
72  {
73  #include "readDyMControls.H"
74  #include "CourantNo.H"
75  #include "alphaCourantNo.H"
76  #include "setDeltaT.H"
77 
79 
80  // Update the mesh for topology change, mesh to mesh mapping
81  mesh.update();
82 
83  runTime++;
84 
85  Info<< "Time = " << runTime.userTimeName() << nl << endl;
86 
87  // --- Pressure-velocity PIMPLE corrector loop
88  while (pimple.loop())
89  {
90  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
91  {
92  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
93 
94  // Move the mesh
95  mesh.move();
96 
97  if (mesh.changing())
98  {
99  Info<< "Execution time for mesh.update() = "
100  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
101  << " s" << endl;
102 
103  gh = (g & mesh.C()) - ghRef;
104  ghf = (g & mesh.Cf()) - ghRef;
105 
106  MRF.update();
107 
108  if (correctPhi)
109  {
110  #include "correctPhi.H"
111  }
112 
113  mixture.correct();
114 
115  if (checkMeshCourantNo)
116  {
117  #include "meshCourantNo.H"
118  }
119  }
120  }
121 
122  fvModels.correct();
123 
124  mixture.solve();
125  rho = mixture.rho();
126 
127  #include "UEqn.H"
128 
129  // --- Pressure corrector loop
130  while (pimple.correct())
131  {
132  #include "pEqn.H"
133  }
134 
135  if (pimple.turbCorr())
136  {
137  turbulence->correct();
138  }
139  }
140 
141  runTime.write();
142 
143  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
144  << " ClockTime = " << runTime.elapsedClockTime() << " s"
145  << nl << endl;
146  }
147 
148  Info<< "End\n" << endl;
149 
150  return 0;
151 }
152 
153 
154 // ************************************************************************* //
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 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
moveMeshOuterCorrectors
Foam::fvModels & fvModels
Calculates and outputs the mean and maximum Courant Numbers.
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.
Creates and initialises the velocity field Uf if required.