ENH: Added Usage section to all solvers
[freefoam.git] / applications / solvers / compressible / sonicTurbFoam / sonicTurbFoam.C
blobcec173df2cf16cb9632618d03e8c941f1aec8473
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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     sonicTurbFoam
28 Description
29     Transient solver for trans-sonic/supersonic, turbulent flow of a 
30     compressible gas.
32 Usage
33     - sonicTurbFoam [OPTION]
35     @param -case \<dir\> \n
36     Specify the case directory
38     @param -parallel \n
39     Run the case in parallel
41     @param -help \n
42     Display short usage message
44     @param -doc \n
45     Display Doxygen documentation page
47     @param -srcDoc \n
48     Display source code
50 \*---------------------------------------------------------------------------*/
52 #include <finiteVolume/fvCFD.H>
53 #include <basicThermophysicalModels/basicThermo.H>
54 #include <compressibleRASModels/RASModel.H>
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 int main(int argc, char *argv[])
61 #   include <OpenFOAM/setRootCase.H>
62 #   include <OpenFOAM/createTime.H>
63 #   include <OpenFOAM/createMesh.H>
64 #   include "createFields.H"
65 #   include <finiteVolume/initContinuityErrs.H>
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70     Info<< "\nStarting time loop\n" << endl;
72     for (runTime++; !runTime.end(); runTime++)
73     {
74         Info<< "Time = " << runTime.timeName() << nl << endl;
76 #       include <finiteVolume/readPISOControls.H>
77 #       include <finiteVolume/compressibleCourantNo.H>
79 #       include <finiteVolume/rhoEqn.H>
81         fvVectorMatrix UEqn
82         (
83             fvm::ddt(rho, U)
84           + fvm::div(phi, U)
85           + turbulence->divDevRhoReff(U)
86         );
88         solve(UEqn == -fvc::grad(p));
90 #       include "hEqn.H"
93         // --- PISO loop
95         for (int corr=0; corr<nCorr; corr++)
96         {
97             volScalarField rUA = 1.0/UEqn.A();
98             U = rUA*UEqn.H();
100             surfaceScalarField phid
101             (
102                 "phid",
103                 fvc::interpolate(thermo->psi())
104                *(
105                    (fvc::interpolate(U) & mesh.Sf())
106                  + fvc::ddtPhiCorr(rUA, rho, U, phi)
107                )
108             );
110             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
111             {
112                 fvScalarMatrix pEqn
113                 (
114                     fvm::ddt(psi, p)
115                   + fvm::div(phid, p)
116                   - fvm::laplacian(rho*rUA, p)
117                 );
119                 pEqn.solve();
121                 if (nonOrth == nNonOrthCorr)
122                 {
123                     phi = pEqn.flux();
124                 }
125             }
127 #           include "compressibleContinuityErrs.H"
129             U -= rUA*fvc::grad(p);
130             U.correctBoundaryConditions();
131         }
133         DpDt = 
134             fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p);
136         turbulence->correct();
138         rho = psi*p;
140         runTime.write();
142         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
143             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
144             << nl << endl;
145     }
147     Info<< "End\n" << endl;
149     return(0);
153 // ************************ vim: set sw=4 sts=4 et: ************************ //