Skip to content
Snippets Groups Projects
Commit 1002f8d0 authored by henry's avatar henry
Browse files

Changed to use p rather than pd.

parent 918a34dc
No related merge requests found
......@@ -24,10 +24,10 @@
==
fvc::reconstruct
(
(
fvc::interpolate(rho)*(g & mesh.Sf())
+ (
fvc::interpolate(interface.sigmaK())*fvc::snGrad(alpha1)
- ghf*fvc::snGrad(rho)
- fvc::snGrad(pd)
- fvc::snGrad(p)
) * mesh.magSf()
)
);
......
{
# include "continuityErrs.H"
wordList pcorrTypes(pd.boundaryField().types());
wordList pcorrTypes(p.boundaryField().types());
for (label i=0; i<pd.boundaryField().size(); i++)
for (label i=0; i<p.boundaryField().size(); i++)
{
if (pd.boundaryField()[i].fixesValue())
if (p.boundaryField()[i].fixesValue())
{
pcorrTypes[i] = fixedValueFvPatchScalarField::typeName;
}
......@@ -22,7 +22,7 @@
IOobject::NO_WRITE
),
mesh,
dimensionedScalar("pcorr", pd.dimensions(), 0.0),
dimensionedScalar("pcorr", p.dimensions(), 0.0),
pcorrTypes
);
......@@ -37,7 +37,7 @@
fvm::laplacian(rUAf, pcorr) == fvc::div(phi)
);
pcorrEqn.setReference(pdRefCell, pdRefValue);
pcorrEqn.setReference(pRefCell, pRefValue);
pcorrEqn.solve();
if (nonOrth == nNonOrthCorr)
......
Info<< "Reading field pd\n" << endl;
volScalarField pd
Info<< "Reading field p\n" << endl;
volScalarField p
(
IOobject
(
"pd",
"p",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
......@@ -83,45 +83,9 @@
);
Info<< "Calculating field g.h\n" << endl;
volScalarField gh("gh", g & mesh.C());
surfaceScalarField ghf("ghf", g & mesh.Cf());
volScalarField p
(
IOobject
(
"p",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
pd + rho*gh
);
label pdRefCell = 0;
scalar pdRefValue = 0.0;
setRefCell(pd, mesh.solutionDict().subDict("PISO"), pdRefCell, pdRefValue);
label pRefCell = 0;
scalar pRefValue = 0.0;
if (pd.needReference())
{
pRefValue = readScalar
(
mesh.solutionDict().subDict("PISO").lookup("pRefValue")
);
p += dimensionedScalar
(
"p",
p.dimensions(),
pRefValue - getRefCellValue(p, pdRefCell)
);
}
setRefCell(p, mesh.solutionDict().subDict("PISO"), pRefCell, pRefValue);
// Construct interface from alpha1 distribution
......
......@@ -89,18 +89,6 @@ int main(int argc, char *argv[])
#include "continuityErrs.H"
p = pd + rho*gh;
if (pd.needReference())
{
p += dimensionedScalar
(
"p",
p.dimensions(),
pRefValue - getRefCellValue(p, pdRefCell)
);
}
turbulence->correct();
runTime.write();
......
......@@ -12,33 +12,33 @@
phi = phiU +
(
fvc::interpolate(interface.sigmaK())*fvc::snGrad(alpha1)
- ghf*fvc::snGrad(rho)
)*rUAf*mesh.magSf();
fvc::interpolate(interface.sigmaK())*fvc::snGrad(alpha1)*mesh.magSf()
+ fvc::interpolate(rho)*(g & mesh.Sf())
)*rUAf;
adjustPhi(phi, U, pd);
adjustPhi(phi, U, p);
for(int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
fvScalarMatrix pdEqn
fvScalarMatrix pEqn
(
fvm::laplacian(rUAf, pd) == fvc::div(phi)
fvm::laplacian(rUAf, p) == fvc::div(phi)
);
pdEqn.setReference(pdRefCell, pdRefValue);
pEqn.setReference(pRefCell, pRefValue);
if (corr == nCorr-1 && nonOrth == nNonOrthCorr)
{
pdEqn.solve(mesh.solver(pd.name() + "Final"));
pEqn.solve(mesh.solver(p.name() + "Final"));
}
else
{
pdEqn.solve(mesh.solver(pd.name()));
pEqn.solve(mesh.solver(p.name()));
}
if (nonOrth == nNonOrthCorr)
{
phi -= pdEqn.flux();
phi -= pEqn.flux();
}
}
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment