Fixing indentation in applications/solvers/solidMechanics
[foam-extend-3.2.git] / applications / solvers / multiphase / interDyMFoam / interDyMFoam.C
blob4e3d84eba640b621870980eb523175322820c05d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Application
25     interDyMFoam
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.
33 \*---------------------------------------------------------------------------*/
35 #include "fvCFD.H"
36 #include "dynamicFvMesh.H"
37 #include "MULES.H"
38 #include "subCycle.H"
39 #include "interfaceProperties.H"
40 #include "twoPhaseMixture.H"
41 #include "turbulenceModel.H"
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 int main(int argc, char *argv[])
47 #   include "setRootCase.H"
48 #   include "createTime.H"
49 #   include "createDynamicFvMesh.H"
50 #   include "readGravitationalAcceleration.H"
51 #   include "readPIMPLEControls.H"
52 #   include "initContinuityErrs.H"
53 #   include "createFields.H"
54 #   include "readTimeControls.H"
55 #   include "correctPhi.H"
56 #   include "CourantNo.H"
57 #   include "setInitialDeltaT.H"
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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         bool meshChanged = mesh.update();
78         reduce(meshChanged, orOp<bool>());
80 #       include "volContinuity.H"
82         volScalarField gh("gh", g & mesh.C());
83         surfaceScalarField ghf("ghf", g & mesh.Cf());
85         if (correctPhi && meshChanged)
86         {
87 #           include "correctPhi.H"
88         }
90         // Make the fluxes relative to the mesh motion
91         fvc::makeRelative(phi, U);
93         if (checkMeshCourantNo)
94         {
95 #           include "meshCourantNo.H"
96         }
98         // Pressure-velocity corrector
99         int oCorr = 0;
100         do
101         {
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             p = pd + rho*gh;
116             if (pd.needReference())
117             {
118                 p += dimensionedScalar
119                 (
120                     "p",
121                     p.dimensions(),
122                     pRefValue - getRefCellValue(p, pdRefCell)
123                 );
124             }
126             turbulence->correct();
127         } while (++oCorr < nOuterCorr);
129         runTime.write();
131         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
132             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
133             << nl << endl;
134     }
136     Info<< "End\n" << endl;
138     return 0;
142 // ************************************************************************* //