Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / solvers / engine / turbDyMEngineFoam / turbDyMEngineFoam.C
blob9edd1497660c7d637daff8757c918dfbf41118aa
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     turbDyMFoam
27 Description
28     Transient solver for incompressible, turbulent flow of Newtonian fluids
29     with moving mesh.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
34 #include "singlePhaseTransportModel.H"
35 #include "turbulenceModel.H"
36 #include "dynamicFvMesh.H"
37 #include "engineTime.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
44 #   include "setRootCase.H"
46 #   include "createEngineTime.H"
47 #   include "createDynamicFvMesh.H"
48 #   include "initContinuityErrs.H"
49 #   include "createFields.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53     Info<< "\nStarting time loop\n" << endl;
55     while (runTime.run())
56     {
57 #       include "readControls.H"
58 #       include "checkTotalVolume.H"
59 #       include "CourantNo.H"
61         // Make the fluxes absolute
62         fvc::makeAbsolute(phi, U);
64 #       include "setDeltaT.H"
66         runTime++;
68         Info<< "Time = " << runTime.timeName() << nl << endl;
70         bool meshChanged = mesh.update();
71         reduce(meshChanged, orOp<bool>());
73         if (meshChanged)
74         {
75 #           include "checkTotalVolume.H"
76 #           include "correctPhi.H"
77 #           include "CourantNo.H"
78         }
80         // Make the fluxes relative to the mesh motion
81         fvc::makeRelative(phi, U);
83         if (checkMeshCourantNo)
84         {
85 #           include "meshCourantNo.H"
86         }
88 #       include "UEqn.H"
90         // --- PISO loop
91         for (int corr = 0; corr < nCorr; corr++)
92         {
93             rUA = 1.0/UEqn.A();
95             U = rUA*UEqn.H();
96             phi = (fvc::interpolate(U) & mesh.Sf());
97               //+ fvc::ddtPhiCorr(rUA, U, phi);
99             adjustPhi(phi, U, p);
101             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
102             {
103                 fvScalarMatrix pEqn
104                 (
105                     fvm::laplacian(rUA, p) == fvc::div(phi)
106                 );
108                 pEqn.setReference(pRefCell, pRefValue);
110                 if (corr == nCorr - 1 && nonOrth == nNonOrthCorr)
111                 {
112                     pEqn.solve(mesh.solutionDict().solver(p.name() + "Final"));
113                 }
114                 else
115                 {
116                     pEqn.solve(mesh.solutionDict().solver(p.name()));
117                 }
119                 if (nonOrth == nNonOrthCorr)
120                 {
121                     phi -= pEqn.flux();
122                 }
123             }
125 #           include "continuityErrs.H"
127             // Make the fluxes relative to the mesh motion
128             fvc::makeRelative(phi, U);
130             U -= rUA*fvc::grad(p);
131             U.correctBoundaryConditions();
132         }
134         turbulence->correct();
136         runTime.write();
138         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
139             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
140             << nl << endl;
141     }
143     Info<< "End\n" << endl;
145     return(0);
149 // ************************************************************************* //