This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/incompressible/channelFoam/channelFoam.C

168 lines
4.9 KiB
C
Raw Normal View History

/*---------------------------------------------------------------------------*\
========= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
2016-06-20 15:00:40 +00:00
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
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
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.
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/>.
Application
2010-08-26 14:22:03 +00:00
channelFoam
Description
Incompressible LES solver for flow in a channel.
2017-01-04 13:37:08 +00:00
Consistent formulation without time-step and relaxation dependence by
Jasak and Tukovic.
Author
Hrvoje Jasak, Wikki Ltd. All rights reserved
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
2010-08-26 14:22:03 +00:00
#include "singlePhaseTransportModel.H"
#include "LESModel.H"
#include "IFstream.H"
#include "OFstream.H"
#include "Random.H"
#include "pisoControl.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
pisoControl piso(mesh);
#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())
{
Info<< "Time = " << runTime.timeName() << nl << endl;
#include "CourantNo.H"
sgsModel->correct();
2017-01-04 13:37:08 +00:00
// Time derivative matrix
fvVectorMatrix ddtUEqn(fvm::ddt(U));
// Convection-diffusion matrix
fvVectorMatrix HUEqn
(
fvm::div(phi, U)
+ sgsModel->divDevBeff()
==
flowDirection*gradP
);
if (piso.momentumPredictor())
{
solve(ddtUEqn + HUEqn == -fvc::grad(p));
}
2017-01-04 13:37:08 +00:00
// Prepare clean 1/a_p without time derivative contribution
volScalarField rAU = 1.0/HUEqn.A();
// --- PISO loop
while (piso.correct())
{
2017-01-04 13:37:08 +00:00
// Calculate U from convection-diffusion matrix
U = rAU*HUEqn.H();
// Consistently calculate flux
piso.calcTransientConsistentFlux(phi, U, rAU, ddtUEqn);
adjustPhi(phi, U, p);
while (piso.correctNonOrthogonal())
{
fvScalarMatrix pEqn
(
2017-01-04 13:37:08 +00:00
fvm::laplacian
(
fvc::interpolate(rAU)/piso.aCoeff(),
p,
"laplacian(rAU," + p.name() + ')'
)
==
fvc::div(phi)
);
pEqn.setReference(pRefCell, pRefValue);
pEqn.solve
(
mesh.solutionDict().solver(p.select(piso.finalInnerIter()))
);
if (piso.finalNonOrthogonalIter())
{
phi -= pEqn.flux();
}
}
# include "continuityErrs.H"
2017-01-04 13:37:08 +00:00
// Consistently reconstruct velocity after pressure equation
piso.reconstructTransientVelocity(U, phi, ddtUEqn, rAU, p);
}
// 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 =
2017-01-04 13:37:08 +00:00
(magUbar - magUbarStar)/rAU.weightedAverage(mesh.V());
2017-01-04 13:37:08 +00:00
U += gragPplus*rAU*flowDirection;
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;
}
// ************************************************************************* //