FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

sonicLiquidFoam.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Application
00025     sonicLiquidFoam
00026 
00027 Description
00028     Transient solver for trans-sonic/supersonic, laminar flow of a
00029     compressible liquid.
00030 
00031 Usage
00032     - sonicLiquidFoam [OPTION]
00033 
00034     @param -case <dir> \n
00035     Specify the case directory
00036 
00037     @param -parallel \n
00038     Run the case in parallel
00039 
00040     @param -help \n
00041     Display short usage message
00042 
00043     @param -doc \n
00044     Display Doxygen documentation page
00045 
00046     @param -srcDoc \n
00047     Display source code
00048 
00049 \*---------------------------------------------------------------------------*/
00050 
00051 #include <finiteVolume/fvCFD.H>
00052 
00053 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00054 
00055 int main(int argc, char *argv[])
00056 {
00057     #include <OpenFOAM/setRootCase.H>
00058     #include <OpenFOAM/createTime.H>
00059     #include <OpenFOAM/createMesh.H>
00060     #include "readThermodynamicProperties.H"
00061     #include "readTransportProperties.H"
00062     #include "createFields.H"
00063     #include <finiteVolume/initContinuityErrs.H>
00064 
00065     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00066 
00067     Info<< "\nStarting time loop\n" << endl;
00068 
00069     while (runTime.loop())
00070     {
00071         Info<< "Time = " << runTime.timeName() << nl << endl;
00072 
00073         #include <finiteVolume/readPISOControls.H>
00074         #include <finiteVolume/compressibleCourantNo.H>
00075 
00076         #include <finiteVolume/rhoEqn.H>
00077 
00078         fvVectorMatrix UEqn
00079         (
00080             fvm::ddt(rho, U)
00081           + fvm::div(phi, U)
00082           - fvm::laplacian(mu, U)
00083         );
00084 
00085         solve(UEqn == -fvc::grad(p));
00086 
00087 
00088         // --- PISO loop
00089 
00090         for (int corr=0; corr<nCorr; corr++)
00091         {
00092             volScalarField rUA = 1.0/UEqn.A();
00093             U = rUA*UEqn.H();
00094 
00095             surfaceScalarField phid
00096             (
00097                 "phid",
00098                 psi
00099                *(
00100                     (fvc::interpolate(U) & mesh.Sf())
00101                   + fvc::ddtPhiCorr(rUA, rho, U, phi)
00102                 )
00103             );
00104 
00105             phi = (rhoO/psi)*phid;
00106 
00107             fvScalarMatrix pEqn
00108             (
00109                 fvm::ddt(psi, p)
00110               + fvc::div(phi)
00111               + fvm::div(phid, p)
00112               - fvm::laplacian(rho*rUA, p)
00113             );
00114 
00115             pEqn.solve();
00116 
00117             phi += pEqn.flux();
00118 
00119             #include "compressibleContinuityErrs.H"
00120 
00121             U -= rUA*fvc::grad(p);
00122             U.correctBoundaryConditions();
00123         }
00124 
00125         rho = rhoO + psi*p;
00126 
00127         runTime.write();
00128 
00129         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
00130             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
00131             << nl << endl;
00132     }
00133 
00134     Info<< "End\n" << endl;
00135 
00136     return 0;
00137 }
00138 
00139 
00140 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines