Removed unneeded lib dependency from mdInitialise
[foam-extend-3.2.git] / applications / solvers / multiphase / interDyMFoam / interDyMFoam.C
blob74f38f63243c5c7ac2cb665f541da8813ab69efd
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     interDyMFoam
28 Description
29     Solver for 2 incompressible, isothermal immiscible fluids using a VOF
30     (volume of fluid) phase-fraction based interface capturing approach,
31     with optional mesh motion and mesh topology changes including adaptive
32     re-meshing.
34 \*---------------------------------------------------------------------------*/
36 #include "fvCFD.H"
37 #include "dynamicFvMesh.H"
38 #include "MULES.H"
39 #include "subCycle.H"
40 #include "interfaceProperties.H"
41 #include "twoPhaseMixture.H"
42 #include "turbulenceModel.H"
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 int main(int argc, char *argv[])
48 #   include "setRootCase.H"
49 #   include "createTime.H"
50 #   include "createDynamicFvMesh.H"
51 #   include "readGravitationalAcceleration.H"
52 #   include "readPIMPLEControls.H"
53 #   include "initContinuityErrs.H"
54 #   include "createFields.H"
55 #   include "readTimeControls.H"
56 #   include "correctPhi.H"
57 #   include "CourantNo.H"
58 #   include "setInitialDeltaT.H"
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62     Info<< "\nStarting time loop\n" << endl;
64     while (runTime.run())
65     {
66 #       include "readControls.H"
67 #       include "CourantNo.H"
69         // Make the fluxes absolute
70         fvc::makeAbsolute(phi, U);
72 #       include "setDeltaT.H"
74         runTime++;
76         Info<< "Time = " << runTime.timeName() << nl << endl;
78         bool meshChanged = mesh.update();
79         reduce(meshChanged, orOp<bool>());
81 #       include "volContinuity.H"
83         volScalarField gh("gh", g & mesh.C());
84         surfaceScalarField ghf("ghf", g & mesh.Cf());
86         if (correctPhi && meshChanged)
87         {
88 #           include "correctPhi.H"
89         }
91         // Make the fluxes relative to the mesh motion
92         fvc::makeRelative(phi, U);
94         if (checkMeshCourantNo)
95         {
96 #           include "meshCourantNo.H"
97         }
99         // Pressure-velocity corrector
100         int oCorr = 0;
101         do
102         {
103             twoPhaseProperties.correct();
105 #           include "alphaEqnSubCycle.H"
107 #           include "UEqn.H"
109             // --- PISO loop
110             for (int corr = 0; corr < nCorr; corr++)
111             {
112                     #           include "pEqn.H"
113             }
115             p = pd + rho*gh;
117             if (pd.needReference())
118             {
119                 p += dimensionedScalar
120                 (
121                     "p",
122                     p.dimensions(),
123                     pRefValue - getRefCellValue(p, pdRefCell)
124                 );
125             }
127             turbulence->correct();
128         } while (++oCorr < nOuterCorr);
130         runTime.write();
132         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
133             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
134             << nl << endl;
135     }
137     Info<< "End\n" << endl;
139     return 0;
143 // ************************************************************************* //