initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / solvers / incompressible / channelFoam / channelFoam.C
blob009a30016164a37f65379807f68f3a20b43f605d
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     channelFoam
28 Description
29     Incompressible LES solver for flow in a channel.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
34 #include "singlePhaseTransportModel.H"
35 #include "LESModel.H"
36 #include "IFstream.H"
37 #include "OFstream.H"
38 #include "Random.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 int main(int argc, char *argv[])
44     #include "setRootCase.H"
45     #include "createTime.H"
46     #include "createMesh.H"
47     #include "readTransportProperties.H"
48     #include "createFields.H"
49     #include "initContinuityErrs.H"
50     #include "createGradP.H"
52     Info<< "\nStarting time loop\n" << endl;
54     while (runTime.loop())
55     {
56         Info<< "Time = " << runTime.timeName() << nl << endl;
58         #include "readPISOControls.H"
60         #include "CourantNo.H"
62         sgsModel->correct();
64         fvVectorMatrix UEqn
65         (
66             fvm::ddt(U)
67           + fvm::div(phi, U)
68           + sgsModel->divDevBeff(U)
69          ==
70             flowDirection*gradP
71         );
73         if (momentumPredictor)
74         {
75             solve(UEqn == -fvc::grad(p));
76         }
79         // --- PISO loop
81         volScalarField rUA = 1.0/UEqn.A();
83         for (int corr=0; corr<nCorr; corr++)
84         {
85             U = rUA*UEqn.H();
86             phi = (fvc::interpolate(U) & mesh.Sf())
87                 + fvc::ddtPhiCorr(rUA, U, phi);
89             adjustPhi(phi, U, p);
91             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
92             {
93                 fvScalarMatrix pEqn
94                 (
95                     fvm::laplacian(rUA, p) == fvc::div(phi)
96                 );
98                 pEqn.setReference(pRefCell, pRefValue);
100                 if (corr == nCorr-1 && nonOrth == nNonOrthCorr)
101                 {
102                     pEqn.solve(mesh.solver(p.name() + "Final"));
103                 }
104                 else
105                 {
106                     pEqn.solve(mesh.solver(p.name()));
107                 }
109                 if (nonOrth == nNonOrthCorr)
110                 {
111                     phi -= pEqn.flux();
112                 }
113             }
115             #include "continuityErrs.H"
117             U -= rUA*fvc::grad(p);
118             U.correctBoundaryConditions();
119         }
122         // Correct driving force for a constant mass flow rate
124         // Extract the velocity in the flow direction
125         dimensionedScalar magUbarStar =
126             (flowDirection & U)().weightedAverage(mesh.V());
128         // Calculate the pressure gradient increment needed to
129         // adjust the average flow-rate to the correct value
130         dimensionedScalar gragPplus =
131             (magUbar - magUbarStar)/rUA.weightedAverage(mesh.V());
133         U += flowDirection*rUA*gragPplus;
135         gradP += gragPplus;
137         Info<< "Uncorrected Ubar = " << magUbarStar.value() << tab
138             << "pressure gradient = " << gradP.value() << endl;
140         runTime.write();
142         #include "writeGradP.H"
144         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
145             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
146             << nl << endl;
147     }
149     Info<< "End\n" << endl;
151     return 0;
155 // ************************************************************************* //