mirror of
https://github.com/OpenFOAM/OpenFOAM-6.git
synced 2025-12-08 06:57:46 +00:00
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:
@ -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"
|
||||
|
||||
@ -8,7 +8,7 @@ application=`getApplication`
|
||||
|
||||
cp -r 0.org 0 > /dev/null 2>&1
|
||||
runApplication blockMesh
|
||||
runApplication $application
|
||||
runApplication $application -writePhi -writep
|
||||
runApplication streamFunction
|
||||
|
||||
# ----------------------------------------------------------------- end-of-file
|
||||
|
||||
@ -28,6 +28,8 @@ gradSchemes
|
||||
divSchemes
|
||||
{
|
||||
default none;
|
||||
div(phi,U) bounded Gauss linear;
|
||||
div(div(phi,U)) Gauss linear;
|
||||
}
|
||||
|
||||
laplacianSchemes
|
||||
|
||||
@ -19,10 +19,20 @@ solvers
|
||||
{
|
||||
Phi
|
||||
{
|
||||
solver PCG;
|
||||
preconditioner DIC;
|
||||
solver GAMG;
|
||||
smoother DIC;
|
||||
cacheAgglomeration on;
|
||||
agglomerator faceAreaPair;
|
||||
nCellsInCoarsestLevel 10;
|
||||
mergeLevels 1;
|
||||
|
||||
tolerance 1e-06;
|
||||
relTol 0;
|
||||
relTol 0.01;
|
||||
}
|
||||
|
||||
p
|
||||
{
|
||||
$Phi;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
@ -8,7 +8,7 @@ application=`getApplication`
|
||||
|
||||
cp -r 0.org 0 > /dev/null 2>&1
|
||||
runApplication blockMesh
|
||||
runApplication $application
|
||||
runApplication $application -writePhi -writep
|
||||
runApplication streamFunction
|
||||
|
||||
# ----------------------------------------------------------------- end-of-file
|
||||
|
||||
@ -28,6 +28,8 @@ gradSchemes
|
||||
divSchemes
|
||||
{
|
||||
default none;
|
||||
div(phi,U) bounded Gauss linear;
|
||||
div(div(phi,U)) Gauss linear;
|
||||
}
|
||||
|
||||
laplacianSchemes
|
||||
|
||||
@ -19,17 +19,26 @@ solvers
|
||||
{
|
||||
Phi
|
||||
{
|
||||
solver PCG;
|
||||
preconditioner DIC;
|
||||
solver GAMG;
|
||||
smoother DIC;
|
||||
cacheAgglomeration on;
|
||||
agglomerator faceAreaPair;
|
||||
nCellsInCoarsestLevel 10;
|
||||
mergeLevels 1;
|
||||
|
||||
tolerance 1e-06;
|
||||
relTol 0;
|
||||
relTol 0.01;
|
||||
}
|
||||
|
||||
p
|
||||
{
|
||||
$Phi;
|
||||
}
|
||||
}
|
||||
|
||||
potentialFlow
|
||||
{
|
||||
nNonOrthogonalCorrectors 0;
|
||||
nNonOrthogonalCorrectors 2;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
|
||||
Reference in New Issue
Block a user