2013-09-17 15:04:22 +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
|
2013-09-17 15:04:22 +00:00
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
2013-12-11 16:09:41 +00:00
|
|
|
This file is part of foam-extend.
|
2013-09-17 15:04:22 +00:00
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
2013-09-17 15:04:22 +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
|
2013-09-17 15:04:22 +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.
|
2013-09-17 15:04:22 +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/>.
|
2013-09-17 15:04:22 +00:00
|
|
|
|
|
|
|
Application
|
|
|
|
simpleSRFFoam
|
|
|
|
|
|
|
|
Description
|
2016-04-08 13:19:19 +00:00
|
|
|
Steady-state solver for incompressible, turbulent flow of Newtonian
|
2013-09-17 15:04:22 +00:00
|
|
|
fluids with single rotating frame.
|
2016-04-08 13:19:19 +00:00
|
|
|
Consistent formulation without time-step and relaxation dependence by Jasak
|
|
|
|
|
|
|
|
Author
|
|
|
|
Hrvoje Jasak, Wikki Ltd. All rights reserved
|
2013-09-17 15:04:22 +00:00
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
|
|
|
#include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
|
|
|
|
#include "incompressible/RAS/RASModel/RASModel.H"
|
|
|
|
#include "SRFModel.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
|
|
|
|
|
|
|
# include "setRootCase.H"
|
|
|
|
|
|
|
|
# include "createTime.H"
|
|
|
|
# include "createMesh.H"
|
|
|
|
# include "createFields.H"
|
|
|
|
# include "initContinuityErrs.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
|
|
|
|
|
|
|
for (runTime++; !runTime.end(); runTime++)
|
|
|
|
{
|
|
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
|
|
|
|
# include "readSIMPLEControls.H"
|
|
|
|
|
|
|
|
p.storePrevIter();
|
|
|
|
|
|
|
|
// Pressure-velocity SIMPLE corrector
|
|
|
|
{
|
|
|
|
// Momentum predictor
|
|
|
|
tmp<fvVectorMatrix> UrelEqn
|
|
|
|
(
|
|
|
|
fvm::div(phi, Urel)
|
|
|
|
+ turbulence->divDevReff(Urel)
|
|
|
|
+ SRF->Su()
|
|
|
|
);
|
|
|
|
|
|
|
|
UrelEqn().relax();
|
|
|
|
|
|
|
|
solve(UrelEqn() == -fvc::grad(p));
|
|
|
|
|
2016-05-24 10:55:36 +00:00
|
|
|
p.boundaryField().updateCoeffs();
|
2013-09-17 15:04:22 +00:00
|
|
|
volScalarField AUrel = UrelEqn().A();
|
|
|
|
Urel = UrelEqn().H()/AUrel;
|
|
|
|
UrelEqn.clear();
|
|
|
|
phi = fvc::interpolate(Urel) & mesh.Sf();
|
|
|
|
adjustPhi(phi, Urel, p);
|
|
|
|
|
|
|
|
// Non-orthogonal pressure corrector loop
|
|
|
|
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
|
|
|
|
{
|
|
|
|
fvScalarMatrix pEqn
|
|
|
|
(
|
|
|
|
fvm::laplacian(1.0/AUrel, p) == fvc::div(phi)
|
|
|
|
);
|
|
|
|
|
|
|
|
pEqn.setReference(pRefCell, pRefValue);
|
|
|
|
pEqn.solve();
|
|
|
|
|
|
|
|
if (nonOrth == nNonOrthCorr)
|
|
|
|
{
|
|
|
|
phi -= pEqn.flux();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
# include "continuityErrs.H"
|
|
|
|
|
|
|
|
// Explicitly relax pressure for momentum corrector
|
|
|
|
p.relax();
|
|
|
|
|
|
|
|
// Momentum corrector
|
|
|
|
Urel -= fvc::grad(p)/AUrel;
|
|
|
|
Urel.correctBoundaryConditions();
|
|
|
|
}
|
|
|
|
|
|
|
|
turbulence->correct();
|
|
|
|
|
|
|
|
// Recalculate Uabs
|
|
|
|
Uabs = Urel + SRF->U();
|
|
|
|
|
|
|
|
runTime.write();
|
|
|
|
|
|
|
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
|
|
<< nl << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
|
|
|
|
return(0);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|