Transferred copyright to the OpenFOAM Foundation
[OpenFOAM-2.0.x.git] / applications / solvers / incompressible / nonNewtonianIcoFoam / nonNewtonianIcoFoam.C
blob29964a7d4d0d672a0651557cff4c9a63e9be79ca
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Application
25     nonNewtonianIcoFoam
27 Description
28     Transient solver for incompressible, laminar flow of non-Newtonian fluids.
30 \*---------------------------------------------------------------------------*/
32 #include "fvCFD.H"
33 #include "singlePhaseTransportModel.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 int main(int argc, char *argv[])
39     #include "setRootCase.H"
41     #include "createTime.H"
42     #include "createMeshNoClear.H"
43     #include "createFields.H"
44     #include "initContinuityErrs.H"
46     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48     Info<< "\nStarting time loop\n" << endl;
50     while (runTime.loop())
51     {
52         Info<< "Time = " << runTime.timeName() << nl << endl;
54         #include "readPISOControls.H"
55         #include "CourantNo.H"
57         fluid.correct();
59         fvVectorMatrix UEqn
60         (
61             fvm::ddt(U)
62           + fvm::div(phi, U)
63           - fvm::laplacian(fluid.nu(), U)
64         );
66         solve(UEqn == -fvc::grad(p));
68         // --- PISO loop
70         for (int corr=0; corr<nCorr; corr++)
71         {
72             volScalarField rAU(1.0/UEqn.A());
74             U = rAU*UEqn.H();
75             phi = (fvc::interpolate(U) & mesh.Sf())
76                 + fvc::ddtPhiCorr(rAU, U, phi);
78             adjustPhi(phi, U, p);
80             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
81             {
82                 fvScalarMatrix pEqn
83                 (
84                     fvm::laplacian(rAU, p) == fvc::div(phi)
85                 );
87                 pEqn.setReference(pRefCell, pRefValue);
88                 pEqn.solve();
90                 if (nonOrth == nNonOrthCorr)
91                 {
92                     phi -= pEqn.flux();
93                 }
94             }
96             #include "continuityErrs.H"
98             U -= rAU*fvc::grad(p);
99             U.correctBoundaryConditions();
100         }
102         runTime.write();
104         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
105             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
106             << nl << endl;
107     }
109     Info<< "End\n" << endl;
111     return 0;
115 // ************************************************************************* //