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
|
|
|
|
rhopSonicFoam
|
|
|
|
|
|
|
|
Description
|
|
|
|
Pressure-density-based compressible flow solver.
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
|
|
|
#include "weighted.H"
|
|
|
|
#include "gaussConvectionScheme.H"
|
|
|
|
#include "multivariateGaussConvectionScheme.H"
|
|
|
|
#include "MUSCL.H"
|
|
|
|
#include "LimitedScheme.H"
|
|
|
|
#include "boundaryTypes.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
|
|
|
|
|
|
|
# include "setRootCase.H"
|
|
|
|
# include "createTime.H"
|
|
|
|
# include "createMesh.H"
|
|
|
|
# include "readThermodynamicProperties.H"
|
|
|
|
# include "createFields.H"
|
2016-05-06 21:31:48 +00:00
|
|
|
# include "createTimeControls.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
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.value() << nl << endl;
|
|
|
|
|
|
|
|
# include "readPISOControls.H"
|
|
|
|
scalar HbyAblend = readScalar(piso.lookup("HbyAblend"));
|
|
|
|
|
|
|
|
# include "readTimeControls.H"
|
|
|
|
|
|
|
|
scalar CoNum = max
|
|
|
|
(
|
|
|
|
mesh.surfaceInterpolation::deltaCoeffs()
|
|
|
|
*mag(phiv)/mesh.magSf()
|
|
|
|
).value()*runTime.deltaT().value();
|
|
|
|
|
|
|
|
Info<< "Max Courant Number = " << CoNum << endl;
|
|
|
|
|
|
|
|
# include "setDeltaT.H"
|
|
|
|
|
2011-10-12 10:41:39 +00:00
|
|
|
for (int outerCorr = 0; outerCorr < nOuterCorr; outerCorr++)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
magRhoU = mag(rhoU);
|
|
|
|
H = (rhoE + p)/rho;
|
|
|
|
|
|
|
|
fv::multivariateGaussConvectionScheme<scalar> mvConvection
|
|
|
|
(
|
|
|
|
mesh,
|
|
|
|
fields,
|
|
|
|
phiv,
|
2011-08-14 16:39:59 +00:00
|
|
|
mesh.schemesDict().divScheme("div(phiv,rhoUH)")
|
2010-05-12 13:27:55 +00:00
|
|
|
);
|
|
|
|
|
|
|
|
solve
|
|
|
|
(
|
|
|
|
fvm::ddt(rho)
|
|
|
|
+ mvConvection.fvmDiv(phiv, rho)
|
|
|
|
);
|
|
|
|
|
|
|
|
surfaceScalarField rhoUWeights =
|
|
|
|
mvConvection.interpolationScheme()()(magRhoU)()
|
|
|
|
.weights(magRhoU);
|
|
|
|
|
|
|
|
weighted<vector> rhoUScheme(rhoUWeights);
|
|
|
|
|
|
|
|
fvVectorMatrix rhoUEqn
|
|
|
|
(
|
|
|
|
fvm::ddt(rhoU)
|
|
|
|
+ fv::gaussConvectionScheme<vector>(mesh, phiv, rhoUScheme)
|
|
|
|
.fvmDiv(phiv, rhoU)
|
|
|
|
);
|
|
|
|
|
|
|
|
solve(rhoUEqn == -fvc::grad(p));
|
|
|
|
|
|
|
|
solve
|
|
|
|
(
|
|
|
|
fvm::ddt(rhoE)
|
|
|
|
+ mvConvection.fvmDiv(phiv, rhoE)
|
|
|
|
==
|
|
|
|
- mvConvection.fvcDiv(phiv, p)
|
|
|
|
);
|
|
|
|
|
|
|
|
T = (rhoE - 0.5*rho*magSqr(rhoU/rho))/Cv/rho;
|
|
|
|
psi = 1.0/(R*T);
|
|
|
|
p = rho/psi;
|
|
|
|
|
2011-10-12 10:41:39 +00:00
|
|
|
for (int corr = 0; corr < nCorr; corr++)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
volScalarField rrhoUA = 1.0/rhoUEqn.A();
|
|
|
|
surfaceScalarField rrhoUAf("rrhoUAf", fvc::interpolate(rrhoUA));
|
|
|
|
volVectorField HbyA = rrhoUA*rhoUEqn.H();
|
|
|
|
|
|
|
|
surfaceScalarField HbyAWeights =
|
|
|
|
HbyAblend*mesh.weights()
|
|
|
|
+ (1.0 - HbyAblend)*
|
|
|
|
LimitedScheme
|
|
|
|
<vector, MUSCLLimiter<NVDTVD>, limitFuncs::magSqr>
|
|
|
|
(mesh, phi, IStringStream("HbyA")()).weights(HbyA);
|
|
|
|
|
|
|
|
phi =
|
|
|
|
(
|
|
|
|
surfaceInterpolationScheme<vector>::interpolate
|
|
|
|
(HbyA, HbyAWeights) & mesh.Sf()
|
|
|
|
)
|
|
|
|
+ HbyAblend*fvc::ddtPhiCorr(rrhoUA, rho, rhoU, phi);
|
|
|
|
|
|
|
|
surfaceScalarField phiGradp =
|
|
|
|
rrhoUAf*mesh.magSf()*fvc::snGrad(p);
|
|
|
|
|
|
|
|
phi -= phiGradp;
|
|
|
|
|
|
|
|
# include "resetPhiPatches.H"
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
surfaceScalarField rhof =
|
2010-05-12 13:27:55 +00:00
|
|
|
mvConvection.interpolationScheme()()(rho)()
|
|
|
|
.interpolate(rho);
|
|
|
|
|
|
|
|
phiv = phi/rhof;
|
|
|
|
|
|
|
|
fvScalarMatrix pEqn
|
|
|
|
(
|
|
|
|
fvm::ddt(psi, p)
|
|
|
|
+ mvConvection.fvcDiv(phiv, rho)
|
|
|
|
+ fvc::div(phiGradp)
|
|
|
|
- fvm::laplacian(rrhoUAf, p)
|
|
|
|
);
|
|
|
|
|
|
|
|
pEqn.solve();
|
|
|
|
|
|
|
|
phi += phiGradp + pEqn.flux();
|
|
|
|
rho = psi*p;
|
2010-08-26 14:22:03 +00:00
|
|
|
rhof =
|
2010-05-12 13:27:55 +00:00
|
|
|
mvConvection.interpolationScheme()()(rho)()
|
|
|
|
.interpolate(rho);
|
|
|
|
phiv = phi/rhof;
|
|
|
|
|
|
|
|
rhoU = HbyA - rrhoUA*fvc::grad(p);
|
|
|
|
rhoU.correctBoundaryConditions();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
U = rhoU/rho;
|
|
|
|
|
|
|
|
runTime.write();
|
|
|
|
|
|
|
|
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
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|