Backported PDRFoam, rhoReactingFoam, sonicTurbDyMEngineFoam, turbDyMEngineFoam (vanil...
[foam-extend-4.0.git] / applications / solvers / engine / turbDyMEngineFoam / turbDyMEngineFoam.C
blob1a8d4646378bf9879864b8fc04213ce09f840f92
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"
38 #include "pisoControl.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 int main(int argc, char *argv[])
45 #   include "setRootCase.H"
47 #   include "createEngineTime.H"
48 #   include "createDynamicFvMesh.H"
50     pisoControl piso(mesh);
52 #   include "initContinuityErrs.H"
53 #   include "createFields.H"
54 #   include "createControls.H"
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58     Info<< "\nStarting time loop\n" << endl;
60     while (runTime.run())
61     {
62 #       include "readControls.H"
63 #       include "checkTotalVolume.H"
64 #       include "CourantNo.H"
66         // Make the fluxes absolute
67         fvc::makeAbsolute(phi, U);
69 #       include "setDeltaT.H"
71         runTime++;
73         Info<< "Time = " << runTime.timeName() << nl << endl;
75         bool meshChanged = mesh.update();
76         reduce(meshChanged, orOp<bool>());
78         if (meshChanged)
79         {
80 #           include "checkTotalVolume.H"
81 #           include "correctPhi.H"
82 #           include "CourantNo.H"
83         }
85         // Make the fluxes relative to the mesh motion
86         fvc::makeRelative(phi, U);
88         if (checkMeshCourantNo)
89         {
90 #           include "meshCourantNo.H"
91         }
93 #       include "UEqn.H"
95         // --- PISO loop
96         while (piso.correct())
97         {
98             rUA = 1.0/UEqn.A();
100             U = rUA*UEqn.H();
101             phi = (fvc::interpolate(U) & mesh.Sf());
102               //+ fvc::ddtPhiCorr(rUA, U, phi);
104             adjustPhi(phi, U, p);
106             while (piso.correctNonOrthogonal())
107             {
108                 fvScalarMatrix pEqn
109                 (
110                     fvm::laplacian(rUA, p) == fvc::div(phi)
111                 );
113                 pEqn.setReference(pRefCell, pRefValue);
115                 pEqn.solve
116                 (
117                     mesh.solutionDict().solver(p.select(piso.finalInnerIter()))
118                 );
120                 if (piso.finalNonOrthogonalIter())
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 // ************************************************************************* //