Removed unneeded lib dependency from mdInitialise
[foam-extend-3.2.git] / applications / solvers / engine / turbDyMEngineFoam / turbDyMEngineFoam.C
blob1fd7d362dccc52e01f0505d74c2150cb8ea1736c
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     turbDyMFoam
28 Description
29     Transient solver for incompressible, turbulent flow of Newtonian fluids
30     with moving mesh.
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "singlePhaseTransportModel.H"
36 #include "turbulenceModel.H"
37 #include "dynamicFvMesh.H"
38 #include "engineTime.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 int main(int argc, char *argv[])
45 #   include "setRootCase.H"
47 #   include "createEngineTime.H"
48 #   include "createDynamicFvMesh.H"
49 #   include "initContinuityErrs.H"
50 #   include "createFields.H"
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54     Info<< "\nStarting time loop\n" << endl;
56     while (runTime.run())
57     {
58 #       include "readControls.H"
59 #       include "checkTotalVolume.H"
60 #       include "CourantNo.H"
62         // Make the fluxes absolute
63         fvc::makeAbsolute(phi, U);
65 #       include "setDeltaT.H"
67         runTime++;
69         Info<< "Time = " << runTime.timeName() << nl << endl;
71         bool meshChanged = mesh.update();
72         reduce(meshChanged, orOp<bool>());
74         if (meshChanged)
75         {
76 #           include "checkTotalVolume.H"
77 #           include "correctPhi.H"
78 #           include "CourantNo.H"
79         }
81         // Make the fluxes relative to the mesh motion
82         fvc::makeRelative(phi, U);
84         if (checkMeshCourantNo)
85         {
86 #           include "meshCourantNo.H"
87         }
89 #       include "UEqn.H"
91         // --- PISO loop
92         for (int corr = 0; corr < nCorr; corr++)
93         {
94             rUA = 1.0/UEqn.A();
96             U = rUA*UEqn.H();
97             phi = (fvc::interpolate(U) & mesh.Sf());
98               //+ fvc::ddtPhiCorr(rUA, U, phi);
100             adjustPhi(phi, U, p);
102             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
103             {
104                 fvScalarMatrix pEqn
105                 (
106                     fvm::laplacian(rUA, p) == fvc::div(phi)
107                 );
109                 pEqn.setReference(pRefCell, pRefValue);
111                 if (corr == nCorr - 1 && nonOrth == nNonOrthCorr)
112                 {
113                     pEqn.solve(mesh.solutionDict().solver(p.name() + "Final"));
114                 }
115                 else
116                 {
117                     pEqn.solve(mesh.solutionDict().solver(p.name()));
118                 }
120                 if (nonOrth == nNonOrthCorr)
121                 {
122                     phi -= pEqn.flux();
123                 }
124             }
126 #           include "continuityErrs.H"
128             // Make the fluxes relative to the mesh motion
129             fvc::makeRelative(phi, U);
131             U -= rUA*fvc::grad(p);
132             U.correctBoundaryConditions();
133         }
135         turbulence->correct();
137         runTime.write();
139         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
140             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
141             << nl << endl;
142     }
144     Info<< "End\n" << endl;
146     return(0);
150 // ************************************************************************* //