1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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
29 Transient solver for trans-sonic/supersonic, laminar flow of a
33 - sonicLiquidFoam [OPTION]
35 @param -case \<dir\> \n
36 Specify the case directory
39 Run the case in parallel
42 Display short usage message
45 Display Doxygen documentation page
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())
72 Info<< "Time = " << runTime.timeName() << nl << endl;
74 #include <finiteVolume/readPISOControls.H>
75 #include <finiteVolume/compressibleCourantNo.H>
77 #include <finiteVolume/rhoEqn.H>
83 - fvm::laplacian(mu, U)
86 solve(UEqn == -fvc::grad(p));
91 for (int corr=0; corr<nCorr; corr++)
93 volScalarField rUA = 1.0/UEqn.A();
96 surfaceScalarField phid
101 (fvc::interpolate(U) & mesh.Sf())
102 + fvc::ddtPhiCorr(rUA, rho, U, phi)
106 phi = (rhoO/psi)*phid;
113 - fvm::laplacian(rho*rUA, p)
120 #include "compressibleContinuityErrs.H"
122 U -= rUA*fvc::grad(p);
123 U.correctBoundaryConditions();
130 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
131 << " ClockTime = " << runTime.elapsedClockTime() << " s"
135 Info<< "End\n" << endl;
141 // ************************ vim: set sw=4 sts=4 et: ************************ //