initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / solvers / incompressible / nonNewtonianIcoFoam / nonNewtonianIcoFoam.C
blobfb7ff1c8904c263cabb3a8dfcdb39266cbefb45d
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     nonNewtonianIcoFoam
28 Description
29     Transient solver for incompressible, laminar flow of non-Newtonian fluids.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
34 #include "singlePhaseTransportModel.H"
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 int main(int argc, char *argv[])
40     #include "setRootCase.H"
42     #include "createTime.H"
43     #include "createMeshNoClear.H"
44     #include "createFields.H"
45     #include "initContinuityErrs.H"
47     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49     Info<< "\nStarting time loop\n" << endl;
51     while (runTime.loop())
52     {
53         Info<< "Time = " << runTime.timeName() << nl << endl;
55         #include "readPISOControls.H"
56         #include "CourantNo.H"
58         fluid.correct();
60         fvVectorMatrix UEqn
61         (
62             fvm::ddt(U)
63           + fvm::div(phi, U)
64           - fvm::laplacian(fluid.nu(), U)
65         );
67         solve(UEqn == -fvc::grad(p));
69         // --- PISO loop
71         for (int corr=0; corr<nCorr; corr++)
72         {
73             volScalarField rUA = 1.0/UEqn.A();
75             U = rUA*UEqn.H();
76             phi = (fvc::interpolate(U) & mesh.Sf())
77                 + fvc::ddtPhiCorr(rUA, U, phi);
79             adjustPhi(phi, U, p);
81             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
82             {
83                 fvScalarMatrix pEqn
84                 (
85                     fvm::laplacian(rUA, p) == fvc::div(phi)
86                 );
88                 pEqn.setReference(pRefCell, pRefValue);
89                 pEqn.solve();
91                 if (nonOrth == nNonOrthCorr)
92                 {
93                     phi -= pEqn.flux();
94                 }
95             }
97             #include "continuityErrs.H"
99             U -= rUA*fvc::grad(p);
100             U.correctBoundaryConditions();
101         }
103         runTime.write();
105         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
106             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
107             << nl << endl;
108     }
110     Info<< "End\n" << endl;
112     return 0;
116 // ************************************************************************* //