1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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
29 Incompressible LES solver for flow in a channel.
31 \*---------------------------------------------------------------------------*/
34 #include "singlePhaseTransportModel.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())
56 Info<< "Time = " << runTime.timeName() << nl << endl;
58 #include "readPISOControls.H"
60 #include "CourantNo.H"
68 + sgsModel->divDevBeff(U)
73 if (momentumPredictor)
75 solve(UEqn == -fvc::grad(p));
81 volScalarField rUA = 1.0/UEqn.A();
83 for (int corr=0; corr<nCorr; corr++)
86 phi = (fvc::interpolate(U) & mesh.Sf())
87 + fvc::ddtPhiCorr(rUA, U, phi);
91 for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
95 fvm::laplacian(rUA, p) == fvc::div(phi)
98 pEqn.setReference(pRefCell, pRefValue);
100 if (corr == nCorr-1 && nonOrth == nNonOrthCorr)
102 pEqn.solve(mesh.solver(p.name() + "Final"));
106 pEqn.solve(mesh.solver(p.name()));
109 if (nonOrth == nNonOrthCorr)
115 #include "continuityErrs.H"
117 U -= rUA*fvc::grad(p);
118 U.correctBoundaryConditions();
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;
137 Info<< "Uncorrected Ubar = " << magUbarStar.value() << tab
138 << "pressure gradient = " << gradP.value() << endl;
142 #include "writeGradP.H"
144 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
145 << " ClockTime = " << runTime.elapsedClockTime() << " s"
149 Info<< "End\n" << endl;
155 // ************************************************************************* //