Files
OpenFOAM-12/applications/solvers/basic/potentialFoam/createFields.H
Henry Weller ed7e703040 Time::timeName(): no longer needed, calls replaced by name()
The timeName() function simply returns the dimensionedScalar::name() which holds
the user-time name of the current time and now that timeName() is no longer
virtual the dimensionedScalar::name() can be called directly.  The timeName()
function implementation is maintained for backward-compatibility.
2022-11-30 15:53:51 +00:00

117 lines
2.1 KiB
C++

Info<< "Reading velocity field U\n" << endl;
volVectorField U
(
IOobject
(
"U",
runTime.name(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
U = Zero;
surfaceScalarField phi
(
IOobject
(
"phi",
runTime.name(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
fvc::flux(U)
);
if (args.optionFound("initialiseUBCs"))
{
U.correctBoundaryConditions();
phi = fvc::flux(U);
}
// 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
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.name(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::NO_WRITE
),
mesh,
dimensionedScalar(pName, sqr(dimVelocity), 0),
pBCTypes
);
// Infer the velocity potential BCs from the pressure
wordList PhiBCTypes
(
p.boundaryField().size(),
zeroGradientFvPatchScalarField::typeName
);
forAll(p.boundaryField(), patchi)
{
if (p.boundaryField()[patchi].fixesValue())
{
PhiBCTypes[patchi] = fixedValueFvPatchScalarField::typeName;
}
}
Info<< "Constructing velocity potential field Phi\n" << endl;
volScalarField Phi
(
IOobject
(
"Phi",
runTime.name(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::NO_WRITE
),
mesh,
dimensionedScalar(dimLength*dimVelocity, 0),
PhiBCTypes
);
label PhiRefCell = 0;
scalar PhiRefValue = 0;
setRefCell
(
Phi,
potentialFlow.dict(),
PhiRefCell,
PhiRefValue
);
mesh.schemes().setFluxRequired(Phi.name());
#include "createMRF.H"