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

nonNewtonianIcoFoam.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     nonNewtonianIcoFoam
00026 
00027 Description
00028     Transient solver for incompressible, laminar flow of non-Newtonian fluids.
00029 
00030 Usage
00031     - nonNewtonianIcoFoam [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 
00053 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00054 
00055 int main(int argc, char *argv[])
00056 {
00057     #include <OpenFOAM/setRootCase.H>
00058 
00059     #include <OpenFOAM/createTime.H>
00060     #include <OpenFOAM/createMeshNoClear.H>
00061     #include "createFields.H"
00062     #include <finiteVolume/initContinuityErrs.H>
00063 
00064     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00065 
00066     Info<< "\nStarting time loop\n" << endl;
00067 
00068     while (runTime.loop())
00069     {
00070         Info<< "Time = " << runTime.timeName() << nl << endl;
00071 
00072         #include <finiteVolume/readPISOControls.H>
00073         #include <finiteVolume/CourantNo.H>
00074 
00075         fluid.correct();
00076 
00077         fvVectorMatrix UEqn
00078         (
00079             fvm::ddt(U)
00080           + fvm::div(phi, U)
00081           - fvm::laplacian(fluid.nu(), U)
00082         );
00083 
00084         solve(UEqn == -fvc::grad(p));
00085 
00086         // --- PISO loop
00087 
00088         for (int corr=0; corr<nCorr; corr++)
00089         {
00090             volScalarField rUA = 1.0/UEqn.A();
00091 
00092             U = rUA*UEqn.H();
00093             phi = (fvc::interpolate(U) & mesh.Sf())
00094                 + fvc::ddtPhiCorr(rUA, U, phi);
00095 
00096             adjustPhi(phi, U, p);
00097 
00098             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00099             {
00100                 fvScalarMatrix pEqn
00101                 (
00102                     fvm::laplacian(rUA, p) == fvc::div(phi)
00103                 );
00104 
00105                 pEqn.setReference(pRefCell, pRefValue);
00106                 pEqn.solve();
00107 
00108                 if (nonOrth == nNonOrthCorr)
00109                 {
00110                     phi -= pEqn.flux();
00111                 }
00112             }
00113 
00114             #include <finiteVolume/continuityErrs.H>
00115 
00116             U -= rUA*fvc::grad(p);
00117             U.correctBoundaryConditions();
00118         }
00119 
00120         runTime.write();
00121 
00122         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
00123             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
00124             << nl << endl;
00125     }
00126 
00127     Info<< "End\n" << endl;
00128 
00129     return 0;
00130 }
00131 
00132 
00133 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines