rhoPimplecFoam: New compressible solver featuring PIMPLEC

This commit is contained in:
Henry
2012-05-03 14:14:40 +01:00
parent 664fb30912
commit 26c9d330e7
21 changed files with 1211 additions and 0 deletions

View File

@ -0,0 +1,3 @@
rhoPimplecFoam.C
EXE = $(FOAM_APPBIN)/rhoPimplecFoam

View File

@ -0,0 +1,15 @@
EXE_INC = \
-I.. \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/turbulenceModels/compressible/turbulenceModel \
-I$(LIB_SRC)/finiteVolume/cfdTools \
-I$(LIB_SRC)/finiteVolume/lnInclude
EXE_LIBS = \
-lbasicThermophysicalModels \
-lspecie \
-lcompressibleTurbulenceModel \
-lcompressibleRASModels \
-lcompressibleLESModels \
-lfiniteVolume \
-lmeshTools

View File

@ -0,0 +1,117 @@
rho = thermo.rho();
rho = max(rho, rhoMin);
rho = min(rho, rhoMax);
rho.relax();
volScalarField rAU(1.0/UEqn().A());
volScalarField rAtU(1.0/(1.0/rAU - UEqn().H1()));
volVectorField HbyA("HbyA", U);
HbyA = rAU*UEqn().H();
if (pimple.nCorrPIMPLE() <= 1)
{
UEqn.clear();
}
if (pimple.transonic())
{
surfaceScalarField phid
(
"phid",
fvc::interpolate(psi)
*(
(fvc::interpolate(HbyA) & mesh.Sf())
+ fvc::ddtPhiCorr(rAU, rho, U, phi)
)
);
surfaceScalarField phic
(
"phic",
fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
);
HbyA -= (rAU - rAtU)*fvc::grad(p);
volScalarField Dp("Dp", rho*rAtU);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix pEqn
(
fvm::ddt(psi, p)
+ fvm::div(phid, p)
+ fvc::div(phic)
- fvm::laplacian(Dp, p)
);
// Relax the pressure equation to maintain diagonal dominance
pEqn.relax();
pEqn.solve();
if (pimple.finalNonOrthogonalIter())
{
phi == phic + pEqn.flux();
}
}
}
else
{
surfaceScalarField phiHbyA
(
"phiHbyA",
fvc::interpolate(rho)
*(
(fvc::interpolate(HbyA) & mesh.Sf())
+ fvc::ddtPhiCorr(rAU, rho, U, phi)
)
);
phiHbyA += fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf();
HbyA -= (rAU - rAtU)*fvc::grad(p);
volScalarField Dp("Dp", rho*rAtU);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix pEqn
(
fvm::ddt(psi, p)
+ fvc::div(phiHbyA)
- fvm::laplacian(Dp, p)
);
pEqn.solve();
if (pimple.finalNonOrthogonalIter())
{
phi = phiHbyA + pEqn.flux();
}
}
}
#include "rhoEqn.H"
#include "compressibleContinuityErrs.H"
// Explicitly relax pressure for momentum corrector
p.relax();
U = HbyA - rAtU*fvc::grad(p);
U.correctBoundaryConditions();
K = 0.5*magSqr(U);
dpdt = fvc::ddt(p);
// Recalculate density from the relaxed pressure
rho = thermo.rho();
rho = max(rho, rhoMin);
rho = min(rho, rhoMax);
if (!pimple.transonic())
{
rho.relax();
}
Info<< "rho max/min : " << max(rho).value() << " " << min(rho).value() << endl;

View File

@ -0,0 +1,105 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM 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.
OpenFOAM 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 OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
rhoPimplecFoam
Description
Transient solver for laminar or turbulent flow of compressible fluids
for HVAC and similar applications.
Uses the flexible PIMPLEC (PISOC-SIMPLEC) solution for time-resolved and
pseudo-transient simulations.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "basicPsiThermo.H"
#include "turbulenceModel.H"
#include "bound.H"
#include "pimpleControl.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
pimpleControl pimple(mesh);
#include "createFields.H"
#include "initContinuityErrs.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (runTime.run())
{
#include "readTimeControls.H"
#include "compressibleCourantNo.H"
#include "setDeltaT.H"
runTime++;
Info<< "Time = " << runTime.timeName() << nl << endl;
if (pimple.nCorrPIMPLE() <= 1)
{
#include "rhoEqn.H"
}
// --- Pressure-velocity PIMPLE corrector loop
while (pimple.loop())
{
#include "UEqn.H"
#include "hEqn.H"
// --- Pressure corrector loop
while (pimple.correct())
{
#include "pEqn.H"
}
if (pimple.turbCorr())
{
turbulence->correct();
}
}
runTime.write();
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //