potentialFoam: Added new method to estimate the static pressure field from the velocity

Uses a form of the Euler equation in which only variation along the streamlines is considered
This commit is contained in:
Henry
2015-02-19 19:05:17 +00:00
parent f43abd9f3b
commit e593fef659
7 changed files with 87 additions and 11 deletions

View File

@ -35,8 +35,8 @@ Description
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "fvIOoptionList.H"
#include "pisoControl.H"
#include "fvIOoptionList.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -61,6 +61,12 @@ int main(int argc, char *argv[])
"Write the velocity potential field"
);
argList::addBoolOption
(
"writep",
"Calculate and write the pressure field"
);
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
@ -125,6 +131,53 @@ int main(int argc, char *argv[])
Phi.write();
}
// Calculate the pressure field
if (args.optionFound("writep"))
{
Info<< nl << "Calculating approximate pressure field" << endl;
label pRefCell = 0;
scalar pRefValue = 0.0;
setRefCell
(
p,
potentialFlow.dict(),
pRefCell,
pRefValue
);
// Calculate the flow-direction filter tensor
volScalarField magSqrU(magSqr(U));
volSymmTensorField F(sqr(U)/(magSqrU + SMALL*average(magSqrU)));
// Calculate the divergence of the flow-direction filtered div(U*U)
// Filtering with the flow-direction generates a more reasonable
// pressure distribution in regions of high velocity gradient in the
// direction of the flow
volScalarField divDivUU
(
fvc::div
(
F & fvc::div(phi, U),
"div(div(phi,U))"
)
);
// Solve a Poisson equation for the approximate pressure
while (potentialFlow.correctNonOrthogonal())
{
fvScalarMatrix pEqn
(
fvm::laplacian(p) + divDivUU
);
pEqn.setReference(pRefCell, pRefValue);
pEqn.solve();
}
p.write();
}
runTime.functionObjects().end();
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"