mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
potentialFoam: Solve for velocity potential named Phi rather than using the pressure field for this purpose
The Phi field is read if available otherwise created automatically with boundary conditions obtained automatically from the pressure field if available (with optional name) otherwise inferred from the velocity field. Phi Laplacian scheme and solver specification are required. See tutorials for examples.
This commit is contained in:
@ -1,21 +1,4 @@
|
||||
Info<< "Reading field p\n" << endl;
|
||||
volScalarField p
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"p",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::MUST_READ,
|
||||
IOobject::NO_WRITE
|
||||
),
|
||||
mesh
|
||||
);
|
||||
|
||||
p = dimensionedScalar("zero", p.dimensions(), 0.0);
|
||||
|
||||
|
||||
Info<< "Reading field U\n" << endl;
|
||||
Info<< "Reading velocity field U\n" << endl;
|
||||
volVectorField U
|
||||
(
|
||||
IOobject
|
||||
@ -51,12 +34,65 @@
|
||||
}
|
||||
|
||||
|
||||
label pRefCell = 0;
|
||||
scalar pRefValue = 0.0;
|
||||
// Default name for the pressure field
|
||||
word pName("p");
|
||||
|
||||
// Update name of the pressure field from the command-line option
|
||||
args.optionReadIfPresent("pName", pName);
|
||||
|
||||
// Infer the pressure BCs from the velocity BCs
|
||||
wordList pBCTypes
|
||||
(
|
||||
U.boundaryField().size(),
|
||||
fixedValueFvPatchScalarField::typeName
|
||||
);
|
||||
|
||||
forAll(U.boundaryField(), patchi)
|
||||
{
|
||||
if (U.boundaryField()[patchi].fixesValue())
|
||||
{
|
||||
pBCTypes[patchi] = zeroGradientFvPatchScalarField::typeName;
|
||||
}
|
||||
}
|
||||
|
||||
Info<< "Constructing pressure field " << pName << nl << endl;
|
||||
volScalarField p
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
pName,
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::READ_IF_PRESENT,
|
||||
IOobject::NO_WRITE
|
||||
),
|
||||
mesh,
|
||||
dimensionedScalar(pName, sqr(dimVelocity), 0),
|
||||
pBCTypes
|
||||
);
|
||||
|
||||
Info<< "Constructing velocity potential field Phi\n" << endl;
|
||||
volScalarField Phi
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"Phi",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::READ_IF_PRESENT,
|
||||
IOobject::NO_WRITE
|
||||
),
|
||||
mesh,
|
||||
dimensionedScalar("Phi", dimLength*dimVelocity, 0),
|
||||
p.boundaryField().types()
|
||||
);
|
||||
|
||||
label PhiRefCell = 0;
|
||||
scalar PhiRefValue = 0;
|
||||
setRefCell
|
||||
(
|
||||
p,
|
||||
potentialFlow,
|
||||
pRefCell,
|
||||
pRefValue
|
||||
Phi,
|
||||
potentialFlow.dict(),
|
||||
PhiRefCell,
|
||||
PhiRefValue
|
||||
);
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -25,29 +25,48 @@ Application
|
||||
potentialFoam
|
||||
|
||||
Description
|
||||
Simple potential flow solver which can be used to generate starting fields
|
||||
for full Navier-Stokes codes.
|
||||
Potential flow solver which solves for the velocity potential
|
||||
from which the flux-field is obtained and velocity field by reconstructing
|
||||
the flux.
|
||||
|
||||
This application is particularly useful to generate starting fields for
|
||||
Navier-Stokes codes.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "fvCFD.H"
|
||||
#include "fvIOoptionList.H"
|
||||
#include "pisoControl.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
int main(int argc, char *argv[])
|
||||
{
|
||||
argList::addBoolOption("writep", "write the final pressure field");
|
||||
argList::addOption
|
||||
(
|
||||
"pName",
|
||||
"pName",
|
||||
"Name of the pressure field"
|
||||
);
|
||||
|
||||
argList::addBoolOption
|
||||
(
|
||||
"initialiseUBCs",
|
||||
"initialise U boundary conditions"
|
||||
"Initialise U boundary conditions"
|
||||
);
|
||||
|
||||
argList::addBoolOption
|
||||
(
|
||||
"writePhi",
|
||||
"Write the velocity potential field"
|
||||
);
|
||||
|
||||
#include "setRootCase.H"
|
||||
#include "createTime.H"
|
||||
#include "createMesh.H"
|
||||
#include "readControls.H"
|
||||
|
||||
pisoControl potentialFlow(mesh, "potentialFlow");
|
||||
|
||||
#include "createFields.H"
|
||||
#include "createFvOptions.H"
|
||||
|
||||
@ -63,55 +82,47 @@ int main(int argc, char *argv[])
|
||||
|
||||
adjustPhi(phi, U, p);
|
||||
|
||||
|
||||
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
|
||||
// Non-orthogonal velocity potential corrector loop
|
||||
while (potentialFlow.correctNonOrthogonal())
|
||||
{
|
||||
fvScalarMatrix pEqn
|
||||
fvScalarMatrix PhiEqn
|
||||
(
|
||||
fvm::laplacian
|
||||
(
|
||||
dimensionedScalar
|
||||
(
|
||||
"1",
|
||||
dimTime/p.dimensions()*dimensionSet(0, 2, -2, 0, 0),
|
||||
1
|
||||
),
|
||||
p
|
||||
)
|
||||
fvm::laplacian(dimensionedScalar("1", dimless, 1), Phi)
|
||||
==
|
||||
fvc::div(phi)
|
||||
);
|
||||
|
||||
pEqn.setReference(pRefCell, pRefValue);
|
||||
pEqn.solve();
|
||||
PhiEqn.setReference(PhiRefCell, PhiRefValue);
|
||||
PhiEqn.solve();
|
||||
|
||||
if (nonOrth == nNonOrthCorr)
|
||||
if (potentialFlow.finalNonOrthogonalIter())
|
||||
{
|
||||
phi -= pEqn.flux();
|
||||
phi -= PhiEqn.flux();
|
||||
}
|
||||
}
|
||||
|
||||
fvOptions.makeAbsolute(phi);
|
||||
|
||||
Info<< "continuity error = "
|
||||
Info<< "Continuity error = "
|
||||
<< mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
|
||||
<< endl;
|
||||
|
||||
U = fvc::reconstruct(phi);
|
||||
U.correctBoundaryConditions();
|
||||
|
||||
Info<< "Interpolated U error = "
|
||||
Info<< "Interpolated velocity error = "
|
||||
<< (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
|
||||
/sum(mesh.magSf())).value()
|
||||
<< endl;
|
||||
|
||||
// Force the write
|
||||
// Write U and phi
|
||||
U.write();
|
||||
phi.write();
|
||||
|
||||
if (args.optionFound("writep"))
|
||||
// Optionally write Phi
|
||||
if (args.optionFound("writePhi"))
|
||||
{
|
||||
p.write();
|
||||
Phi.write();
|
||||
}
|
||||
|
||||
runTime.functionObjects().end();
|
||||
|
||||
Reference in New Issue
Block a user