initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / solvers / multiphase / interDyMFoam / interDyMFoam.C
blobd353012e48baa57e3b99f07b701e87390f20a848
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 "readPISOControls.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     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61     Info<< "\nStarting time loop\n" << endl;
63     while (runTime.run())
64     {
65         #include "readControls.H"
66         #include "CourantNo.H"
68         // Make the fluxes absolute
69         fvc::makeAbsolute(phi, U);
71         #include "setDeltaT.H"
73         runTime++;
75         Info<< "Time = " << runTime.timeName() << nl << endl;
77         scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
79         // Do any mesh changes
80         mesh.update();
82         if (mesh.changing())
83         {
84             Info<< "Execution time for mesh.update() = "
85                 << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
86                 << " s" << endl;
87         }
89         if (mesh.changing() && correctPhi)
90         {
91             #include "correctPhi.H"
92         }
94         // Make the fluxes relative to the mesh motion
95         fvc::makeRelative(phi, U);
97         if (mesh.changing() && checkMeshCourantNo)
98         {
99             #include "meshCourantNo.H"
100         }
102         twoPhaseProperties.correct();
104         #include "alphaEqnSubCycle.H"
106         #include "UEqn.H"
108         // --- PISO loop
109         for (int corr=0; corr<nCorr; corr++)
110         {
111             #include "pEqn.H"
112         }
114         turbulence->correct();
116         runTime.write();
118         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
119             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
120             << nl << endl;
121     }
123     Info<< "End\n" << endl;
125     return 0;
129 // ************************************************************************* //