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

channelFoam.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     channelFoam
00026 
00027 Description
00028     Incompressible LES solver for flow in a channel.
00029 
00030 Usage
00031     - channelFoam [OPTION]
00032 
00033     @param -case <dir> \n
00034     Specify the case directory
00035 
00036     @param -parallel \n
00037     Run the case in parallel
00038 
00039     @param -help \n
00040     Display short usage message
00041 
00042     @param -doc \n
00043     Display Doxygen documentation page
00044 
00045     @param -srcDoc \n
00046     Display source code
00047 
00048 \*---------------------------------------------------------------------------*/
00049 
00050 #include <finiteVolume/fvCFD.H>
00051 #include <incompressibleTransportModels/singlePhaseTransportModel.H>
00052 #include <incompressibleLESModels/LESModel.H>
00053 #include <OpenFOAM/IFstream.H>
00054 #include <OpenFOAM/OFstream.H>
00055 #include <OpenFOAM/Random.H>
00056 
00057 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00058 
00059 int main(int argc, char *argv[])
00060 {
00061     #include <OpenFOAM/setRootCase.H>
00062     #include <OpenFOAM/createTime.H>
00063     #include <OpenFOAM/createMesh.H>
00064     #include "readTransportProperties.H"
00065     #include "createFields.H"
00066     #include <finiteVolume/initContinuityErrs.H>
00067     #include "createGradP.H"
00068 
00069     Info<< "\nStarting time loop\n" << endl;
00070 
00071     while (runTime.loop())
00072     {
00073         Info<< "Time = " << runTime.timeName() << nl << endl;
00074 
00075         #include <finiteVolume/readPISOControls.H>
00076 
00077         #include <finiteVolume/CourantNo.H>
00078 
00079         sgsModel->correct();
00080 
00081         fvVectorMatrix UEqn
00082         (
00083             fvm::ddt(U)
00084           + fvm::div(phi, U)
00085           + sgsModel->divDevBeff(U)
00086          ==
00087             flowDirection*gradP
00088         );
00089 
00090         if (momentumPredictor)
00091         {
00092             solve(UEqn == -fvc::grad(p));
00093         }
00094 
00095 
00096         // --- PISO loop
00097 
00098         volScalarField rUA = 1.0/UEqn.A();
00099 
00100         for (int corr=0; corr<nCorr; corr++)
00101         {
00102             U = rUA*UEqn.H();
00103             phi = (fvc::interpolate(U) & mesh.Sf())
00104                 + fvc::ddtPhiCorr(rUA, U, phi);
00105 
00106             adjustPhi(phi, U, p);
00107 
00108             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00109             {
00110                 fvScalarMatrix pEqn
00111                 (
00112                     fvm::laplacian(rUA, p) == fvc::div(phi)
00113                 );
00114 
00115                 pEqn.setReference(pRefCell, pRefValue);
00116 
00117                 if (corr == nCorr-1 && nonOrth == nNonOrthCorr)
00118                 {
00119                     pEqn.solve(mesh.solver(p.name() + "Final"));
00120                 }
00121                 else
00122                 {
00123                     pEqn.solve(mesh.solver(p.name()));
00124                 }
00125 
00126                 if (nonOrth == nNonOrthCorr)
00127                 {
00128                     phi -= pEqn.flux();
00129                 }
00130             }
00131 
00132             #include <finiteVolume/continuityErrs.H>
00133 
00134             U -= rUA*fvc::grad(p);
00135             U.correctBoundaryConditions();
00136         }
00137 
00138 
00139         // Correct driving force for a constant mass flow rate
00140 
00141         // Extract the velocity in the flow direction
00142         dimensionedScalar magUbarStar =
00143             (flowDirection & U)().weightedAverage(mesh.V());
00144 
00145         // Calculate the pressure gradient increment needed to
00146         // adjust the average flow-rate to the correct value
00147         dimensionedScalar gragPplus =
00148             (magUbar - magUbarStar)/rUA.weightedAverage(mesh.V());
00149 
00150         U += flowDirection*rUA*gragPplus;
00151 
00152         gradP += gragPplus;
00153 
00154         Info<< "Uncorrected Ubar = " << magUbarStar.value() << tab
00155             << "pressure gradient = " << gradP.value() << endl;
00156 
00157         runTime.write();
00158 
00159         #include "writeGradP.H"
00160 
00161         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
00162             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
00163             << nl << endl;
00164     }
00165 
00166     Info<< "End\n" << endl;
00167 
00168     return 0;
00169 }
00170 
00171 
00172 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines