Removed unnecessary return statement
[foam-extend-3.2.git] / applications / solvers / incompressible / channelFoam / channelFoam.C
blobba2b7703869726df6c958283f511a4d0bdf328c4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Application
25     channelFoam
27 Description
28     Incompressible LES solver for flow in a channel.
30 \*---------------------------------------------------------------------------*/
32 #include "fvCFD.H"
33 #include "singlePhaseTransportModel.H"
34 #include "LESModel.H"
35 #include "IFstream.H"
36 #include "OFstream.H"
37 #include "Random.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43     #include "setRootCase.H"
44     #include "createTime.H"
45     #include "createMesh.H"
46     #include "readTransportProperties.H"
47     #include "createFields.H"
48     #include "initContinuityErrs.H"
49     #include "createGradP.H"
51     Info<< "\nStarting time loop\n" << endl;
53     while (runTime.loop())
54     {
55         Info<< "Time = " << runTime.timeName() << nl << endl;
57         #include "readPISOControls.H"
59         #include "CourantNo.H"
61         sgsModel->correct();
63         fvVectorMatrix UEqn
64         (
65             fvm::ddt(U)
66           + fvm::div(phi, U)
67           + sgsModel->divDevBeff(U)
68          ==
69             flowDirection*gradP
70         );
72         if (momentumPredictor)
73         {
74             solve(UEqn == -fvc::grad(p));
75         }
78         // --- PISO loop
80         volScalarField rUA = 1.0/UEqn.A();
82         for (int corr = 0; corr < nCorr; corr++)
83         {
84             U = rUA*UEqn.H();
85             phi = (fvc::interpolate(U) & mesh.Sf())
86                 + fvc::ddtPhiCorr(rUA, U, phi);
88             adjustPhi(phi, U, p);
90             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
91             {
92                 fvScalarMatrix pEqn
93                 (
94                     fvm::laplacian(rUA, p) == fvc::div(phi)
95                 );
97                 pEqn.setReference(pRefCell, pRefValue);
99                 if (corr == nCorr-1 && nonOrth == nNonOrthCorr)
100                 {
101                     pEqn.solve(mesh.solutionDict().solver(p.name() + "Final"));
102                 }
103                 else
104                 {
105                     pEqn.solve(mesh.solutionDict().solver(p.name()));
106                 }
108                 if (nonOrth == nNonOrthCorr)
109                 {
110                     phi -= pEqn.flux();
111                 }
112             }
114             #include "continuityErrs.H"
116             U -= rUA*fvc::grad(p);
117             U.correctBoundaryConditions();
118         }
121         // Correct driving force for a constant mass flow rate
123         // Extract the velocity in the flow direction
124         dimensionedScalar magUbarStar =
125             (flowDirection & U)().weightedAverage(mesh.V());
127         // Calculate the pressure gradient increment needed to
128         // adjust the average flow-rate to the correct value
129         dimensionedScalar gragPplus =
130             (magUbar - magUbarStar)/rUA.weightedAverage(mesh.V());
132         U += flowDirection*rUA*gragPplus;
134         gradP += gragPplus;
136         Info<< "Uncorrected Ubar = " << magUbarStar.value() << tab
137             << "pressure gradient = " << gradP.value() << endl;
139         runTime.write();
141         #include "writeGradP.H"
143         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
144             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
145             << nl << endl;
146     }
148     Info<< "End\n" << endl;
150     return 0;
154 // ************************************************************************* //