Go to the documentation of this file.
93 int main(
int argc,
char *argv[])
99 "Name of the pressure field"
102 argList::addBoolOption
105 "Initialise U boundary conditions"
108 argList::addBoolOption
111 "Write the final velocity potential field"
114 argList::addBoolOption
117 "Calculate and write the Euler pressure field"
120 argList::addBoolOption
122 "withFunctionObjects",
123 "execute functionObjects"
133 #include "createMRF.H"
137 Info<<
nl <<
"Calculating potential flow" <<
endl;
141 runTime.functionObjects().start();
146 while (potentialFlow.correctNonOrthogonal())
156 PhiEqn.setReference(PhiRefCell, PhiRefValue);
159 if (potentialFlow.finalNonOrthogonalIter())
161 phi -= PhiEqn.flux();
167 Info<<
"Continuity error = "
172 U.correctBoundaryConditions();
174 Info<<
"Interpolated velocity error = "
192 Info<<
nl <<
"Calculating approximate pressure field" <<
endl;
199 potentialFlow.dict(),
222 while (potentialFlow.correctNonOrthogonal())
236 runTime.functionObjects().end();
238 Info<<
"ExecutionTime = " << runTime.elapsedCpuTime() <<
" s"
239 <<
" ClockTime = " << runTime.elapsedClockTime() <<
" s"
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Specialization of the pimpleControl class for PISO control.
const dictionary & dict() const
Return the solution dictionary.
tmp< surfaceScalarField > interpolate(const RhoType &rho)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Ostream & endl(Ostream &os)
Add newline and flush stream.
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
void setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
bool finalNonOrthogonalIter() const
Helper function to identify final non-orthogonal iteration.
const surfaceVectorField & Sf() const
Return cell face area vectors.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
bool optionFound(const word &opt) const
Return true if the named option is found.
int main(int argc, char *argv[])
Generic GeometricField class.
Foam::argList args(argc, argv)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
dimensioned< scalar > magSqr(const dimensioned< Type > &)