2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
2013-12-11 16:09:41 +00:00
|
|
|
\\ / F ield | foam-extend: Open Source CFD
|
2015-05-17 13:32:07 +00:00
|
|
|
\\ / O peration | Version: 3.2
|
|
|
|
\\ / A nd | Web: http://www.foam-extend.org
|
|
|
|
\\/ M anipulation | For copyright notice see file Copyright
|
2010-05-12 13:27:55 +00:00
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
2013-12-11 16:09:41 +00:00
|
|
|
This file is part of foam-extend.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
2010-05-12 13:27:55 +00:00
|
|
|
under the terms of the GNU General Public License as published by the
|
2013-12-11 16:09:41 +00:00
|
|
|
Free Software Foundation, either version 3 of the License, or (at your
|
2010-05-12 13:27:55 +00:00
|
|
|
option) any later version.
|
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is distributed in the hope that it will be useful, but
|
|
|
|
WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
|
|
|
|
General Public License for more details.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
2013-12-11 16:09:41 +00:00
|
|
|
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Application
|
2010-08-26 14:22:03 +00:00
|
|
|
channelFoam
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Description
|
|
|
|
Incompressible LES solver for flow in a channel.
|
2016-04-05 11:53:56 +00:00
|
|
|
Consistent formulation without time-step and relaxation dependence by Jasak
|
|
|
|
|
|
|
|
Author
|
|
|
|
Hrvoje Jasak, Wikki Ltd. All rights reserved
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
2010-08-26 14:22:03 +00:00
|
|
|
#include "singlePhaseTransportModel.H"
|
|
|
|
#include "LESModel.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
#include "IFstream.H"
|
|
|
|
#include "OFstream.H"
|
|
|
|
#include "Random.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
|
|
|
#include "setRootCase.H"
|
|
|
|
#include "createTime.H"
|
|
|
|
#include "createMesh.H"
|
|
|
|
#include "readTransportProperties.H"
|
|
|
|
#include "createFields.H"
|
|
|
|
#include "initContinuityErrs.H"
|
|
|
|
#include "createGradP.H"
|
|
|
|
|
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
while (runTime.loop())
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
|
|
|
|
#include "readPISOControls.H"
|
|
|
|
|
|
|
|
#include "CourantNo.H"
|
|
|
|
|
|
|
|
sgsModel->correct();
|
|
|
|
|
2016-04-05 11:53:56 +00:00
|
|
|
// Convection-diffusion matrix
|
|
|
|
fvVectorMatrix HUEqn
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2016-04-05 11:53:56 +00:00
|
|
|
fvm::div(phi, U)
|
2010-05-12 13:27:55 +00:00
|
|
|
+ sgsModel->divDevBeff(U)
|
|
|
|
==
|
|
|
|
flowDirection*gradP
|
|
|
|
);
|
|
|
|
|
2016-04-05 11:53:56 +00:00
|
|
|
// Time derivative matrix
|
|
|
|
fvVectorMatrix ddtUEqn(fvm::ddt(U));
|
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
if (momentumPredictor)
|
|
|
|
{
|
2016-04-05 11:53:56 +00:00
|
|
|
solve(ddtUEqn + HUEqn == -fvc::grad(p));
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
2016-04-05 11:53:56 +00:00
|
|
|
// Prepare clean Ap without time derivative contribution
|
|
|
|
// HJ, 26/Oct/2015
|
|
|
|
volScalarField aU = HUEqn.A();
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// --- PISO loop
|
|
|
|
|
2011-10-12 10:41:39 +00:00
|
|
|
for (int corr = 0; corr < nCorr; corr++)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2016-04-05 11:53:56 +00:00
|
|
|
U = HUEqn.H()/aU;
|
|
|
|
phi = (fvc::interpolate(U) & mesh.Sf());
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
adjustPhi(phi, U, p);
|
|
|
|
|
2016-04-05 11:53:56 +00:00
|
|
|
for (int nonOrth = 0; nonOrth <= nNonOrthCorr; nonOrth++)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
fvScalarMatrix pEqn
|
|
|
|
(
|
2016-04-05 11:53:56 +00:00
|
|
|
fvm::laplacian(1/aU, p) == fvc::div(phi)
|
2010-05-12 13:27:55 +00:00
|
|
|
);
|
|
|
|
|
|
|
|
pEqn.setReference(pRefCell, pRefValue);
|
|
|
|
|
|
|
|
if (corr == nCorr-1 && nonOrth == nNonOrthCorr)
|
|
|
|
{
|
2011-08-14 16:39:59 +00:00
|
|
|
pEqn.solve(mesh.solutionDict().solver(p.name() + "Final"));
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
2011-08-14 16:39:59 +00:00
|
|
|
pEqn.solve(mesh.solutionDict().solver(p.name()));
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
if (nonOrth == nNonOrthCorr)
|
|
|
|
{
|
|
|
|
phi -= pEqn.flux();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2016-04-05 11:53:56 +00:00
|
|
|
# include "continuityErrs.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2016-04-05 11:53:56 +00:00
|
|
|
// Note: cannot call H(U) here because the velocity is not complete
|
|
|
|
// HJ, 22/Jan/2016
|
|
|
|
U = 1.0/(aU + ddtUEqn.A())*
|
|
|
|
(
|
|
|
|
U*aU - fvc::grad(p) + ddtUEqn.H()
|
|
|
|
);
|
2010-05-12 13:27:55 +00:00
|
|
|
U.correctBoundaryConditions();
|
|
|
|
}
|
|
|
|
|
|
|
|
// Correct driving force for a constant mass flow rate
|
|
|
|
|
|
|
|
// Extract the velocity in the flow direction
|
|
|
|
dimensionedScalar magUbarStar =
|
|
|
|
(flowDirection & U)().weightedAverage(mesh.V());
|
|
|
|
|
|
|
|
// Calculate the pressure gradient increment needed to
|
|
|
|
// adjust the average flow-rate to the correct value
|
|
|
|
dimensionedScalar gragPplus =
|
2016-04-05 11:53:56 +00:00
|
|
|
(magUbar - magUbarStar)*aU.weightedAverage(mesh.V());
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2016-04-05 11:53:56 +00:00
|
|
|
U += gragPplus/aU*flowDirection;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
gradP += gragPplus;
|
|
|
|
|
|
|
|
Info<< "Uncorrected Ubar = " << magUbarStar.value() << tab
|
|
|
|
<< "pressure gradient = " << gradP.value() << endl;
|
|
|
|
|
|
|
|
runTime.write();
|
|
|
|
|
|
|
|
#include "writeGradP.H"
|
|
|
|
|
|
|
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
|
|
<< nl << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
return 0;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|