Cosmetic improvements.
[OpenFOAM-1.6.x.git] / applications / solvers / compressible / sonicDyMFoam / sonicDyMFoam.C
blob9551fcadd8a1f8ace50ffba04559739cb5620f07
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     sonicDyMFoam
28 Description
29     Transient solver for trans-sonic/supersonic, laminar or turbulent flow
30     of a compressible gas with mesh motion..
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "basicPsiThermo.H"
36 #include "turbulenceModel.H"
37 #include "motionSolver.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43     #include "setRootCase.H"
44     #include "createTime.H"
45     #include "createMesh.H"
46     #include "createFields.H"
47     #include "initContinuityErrs.H"
49     autoPtr<Foam::motionSolver> motionPtr = motionSolver::New(mesh);
51     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53     Info<< "\nStarting time loop\n" << endl;
55     while (runTime.loop())
56     {
57         Info<< "Time = " << runTime.timeName() << nl << endl;
59         #include "readPISOControls.H"
60         #include "compressibleCourantNo.H"
62         mesh.movePoints(motionPtr->newPoints());
64         #include "rhoEqn.H"
66         fvVectorMatrix UEqn
67         (
68             fvm::ddt(rho, U)
69           + fvm::div(phi, U)
70           + turbulence->divDevRhoReff(U)
71         );
73         solve(UEqn == -fvc::grad(p));
75         #include "eEqn.H"
78         // --- PISO loop
80         for (int corr=0; corr<nCorr; corr++)
81         {
82             U = UEqn.H()/UEqn.A();
84             surfaceScalarField phid
85             (
86                 "phid",
87                 fvc::interpolate(psi)
88                *(
89                     (fvc::interpolate(U) & mesh.Sf()) - fvc::meshPhi(rho, U)
90                 )
91             );
93             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
94             {
95                 fvScalarMatrix pEqn
96                 (
97                     fvm::ddt(psi, p)
98                   + fvm::div(phid, p)
99                   - fvm::laplacian(rho/UEqn.A(), p)
100                 );
102                 pEqn.solve();
104                 phi = pEqn.flux();
105             }
107             #include "compressibleContinuityErrs.H"
109             U -= fvc::grad(p)/UEqn.A();
110             U.correctBoundaryConditions();
111         }
113         DpDt = fvc::DDt
114         (
115             surfaceScalarField
116             (
117                 "phiU",
118                 phi/fvc::interpolate(rho) + fvc::meshPhi(rho, U)
119             ),
120             p
121         );
123         turbulence->correct();
125         rho = psi*p;
127         runTime.write();
129         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
130             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
131             << nl << endl;
132     }
134     Info<< "End\n" << endl;
136     return 0;
140 // ************************************************************************* //