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/simpleSRFFoam/simpleSRFFoam.C
Henrik Rusche b858c51ff2 Merge Request #30: Time consistent incompressible solvers. Author: Hrvoje Jasak. Merge: Henrik Rusche
https://sourceforge.net/p/foam-extend/foam-extend-3.2/merge-requests/30/

Conflicts:
	applications/solvers/incompressible/simpleFoam/pEqn.H
	applications/solvers/incompressible/simpleSRFFoam/simpleSRFFoam.C
	tutorials/incompressible/pimpleDyMFoam/movingCone/system/fvSchemes
	tutorials/incompressible/simpleFoam/motorBike/system/fvSchemes
2016-05-24 13:42:50 +02:00

132 lines
3.9 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 3.2
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
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
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
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
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Application
simpleSRFFoam
Description
Steady-state solver for incompressible, turbulent flow of Newtonian
fluids with single rotating frame.
Consistent formulation without time-step and relaxation dependence by Jasak
Author
Hrvoje Jasak, Wikki Ltd. All rights reserved
\*---------------------------------------------------------------------------*/
#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));
p.boundaryField().updateCoeffs();
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);
}
// ************************************************************************* //