Merge branch 'upstream/OpenFOAM' into master
[freefoam.git] / applications / solvers / compressible / sonicLiquidFoam / sonicLiquidFoam.C
blob87575743580e0806783703f73845754ea836e00e
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     sonicLiquidFoam
28 Description
29     Transient solver for trans-sonic/supersonic, laminar flow of a
30     compressible liquid.
32 Usage
33     - sonicLiquidFoam [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>
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 int main(int argc, char *argv[])
58     #include <OpenFOAM/setRootCase.H>
59     #include <OpenFOAM/createTime.H>
60     #include <OpenFOAM/createMesh.H>
61     #include "readThermodynamicProperties.H"
62     #include "readTransportProperties.H"
63     #include "createFields.H"
64     #include <finiteVolume/initContinuityErrs.H>
66     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68     Info<< "\nStarting time loop\n" << endl;
70     while (runTime.loop())
71     {
72         Info<< "Time = " << runTime.timeName() << nl << endl;
74         #include <finiteVolume/readPISOControls.H>
75         #include <finiteVolume/compressibleCourantNo.H>
77         #include <finiteVolume/rhoEqn.H>
79         fvVectorMatrix UEqn
80         (
81             fvm::ddt(rho, U)
82           + fvm::div(phi, U)
83           - fvm::laplacian(mu, U)
84         );
86         solve(UEqn == -fvc::grad(p));
89         // --- PISO loop
91         for (int corr=0; corr<nCorr; corr++)
92         {
93             volScalarField rUA = 1.0/UEqn.A();
94             U = rUA*UEqn.H();
96             surfaceScalarField phid
97             (
98                 "phid",
99                 psi
100                *(
101                     (fvc::interpolate(U) & mesh.Sf())
102                   + fvc::ddtPhiCorr(rUA, rho, U, phi)
103                 )
104             );
106             phi = (rhoO/psi)*phid;
108             fvScalarMatrix pEqn
109             (
110                 fvm::ddt(psi, p)
111               + fvc::div(phi)
112               + fvm::div(phid, p)
113               - fvm::laplacian(rho*rUA, p)
114             );
116             pEqn.solve();
118             phi += pEqn.flux();
120             #include "compressibleContinuityErrs.H"
122             U -= rUA*fvc::grad(p);
123             U.correctBoundaryConditions();
124         }
126         rho = rhoO + psi*p;
128         runTime.write();
130         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
131             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
132             << nl << endl;
133     }
135     Info<< "End\n" << endl;
137     return 0;
141 // ************************ vim: set sw=4 sts=4 et: ************************ //