Go to the documentation of this file.
42 namespace functionObjects
57 Log <<
" functionObjects::" <<
type() <<
" " <<
name()
58 <<
" calculating stream-function" <<
endl;
87 label nVisitedOld = 0;
95 unitAreas /=
mag(unitAreas);
113 if (!isType<emptyPolyPatch>(
patches[patchi]))
117 if (
magSqr(
phi.boundaryField()[patchi][facei]) < SMALL)
119 const labelList& zeroPoints = bouFaces[facei];
124 forAll(zeroPoints, pointi)
126 if (visitedPoint[zeroPoints[pointi]] == 1)
135 Log <<
" Zero face: patch: " << patchi
136 <<
" face: " << facei <<
endl;
138 forAll(zeroPoints, pointi)
141 visitedPoint[zeroPoints[pointi]] = 1;
156 Log <<
" Zero flux boundary face not found. "
157 <<
"Using cell as a reference."
168 forAll(zeroPoints, pointi)
170 if (visitedPoint[zeroPoints[pointi]] == 1)
179 forAll(zeroPoints, pointi)
182 visitedPoint[zeroPoints[pointi]] = 1;
191 <<
"Cannot find initialisation face or a cell."
206 for (label facei = nInternalFaces; facei<faces.size(); facei++)
208 const labelList& curBPoints = faces[facei];
209 bool bPointFound =
false;
211 scalar currentBStream = 0.0;
212 vector currentBStreamPoint(0, 0, 0);
214 forAll(curBPoints, pointi)
217 if (visitedPoint[curBPoints[pointi]] == 1)
221 currentBStreamPoint =
points[curBPoints[pointi]];
232 forAll(curBPoints, pointi)
235 if (visitedPoint[curBPoints[pointi]] == 0)
242 !isType<emptyPolyPatch>(
patches[patchNo])
243 && !isType<symmetryPlanePolyPatch>
245 && !isType<symmetryPolyPatch>(
patches[patchNo])
246 && !isType<wedgePolyPatch>(
patches[patchNo])
254 points[curBPoints[pointi]]
255 - currentBStreamPoint;
256 edgeHat.replace(slabDir, 0);
259 vector nHat = unitAreas[facei];
261 if (edgeHat.y() > VSMALL)
263 visitedPoint[curBPoints[pointi]] = 1;
268 +
phi.boundaryField()[patchNo][faceNo]
271 else if (edgeHat.y() < -VSMALL)
273 visitedPoint[curBPoints[pointi]] = 1;
278 -
phi.boundaryField()[patchNo][faceNo]
283 if (edgeHat.x() > VSMALL)
285 visitedPoint[curBPoints[pointi]] = 1;
290 +
phi.boundaryField()[patchNo][faceNo]
293 else if (edgeHat.x() < -VSMALL)
295 visitedPoint[curBPoints[pointi]] = 1;
300 -
phi.boundaryField()[patchNo][faceNo]
314 for (label facei=0; facei<nInternalFaces; facei++)
317 const labelList& curPoints = faces[facei];
319 bool pointFound =
false;
321 scalar currentStream = 0.0;
322 point currentStreamPoint(0, 0, 0);
327 if (visitedPoint[curPoints[pointi]] == 1)
331 currentStreamPoint =
points[curPoints[pointi]];
344 if (visitedPoint[curPoints[pointi]] == 0)
347 points[curPoints[pointi]] - currentStreamPoint;
349 edgeHat.replace(slabDir, 0);
352 vector nHat = unitAreas[facei];
354 if (edgeHat.y() > VSMALL)
356 visitedPoint[curPoints[pointi]] = 1;
363 else if (edgeHat.y() < -VSMALL)
365 visitedPoint[curPoints[pointi]] = 1;
381 if (nVisited == nVisitedOld)
385 Log <<
" Exhausted a seed, looking for new seed "
386 <<
"(this is correct for multiply connected domains).";
392 nVisitedOld = nVisited;
401 return tstreamFunction;
405 bool Foam::functionObjects::streamFunction::calc()
407 if (foundObject<surfaceScalarField>(fieldName_))
412 return store(resultName_, calc(
phi));
437 <<
"Case is not 2D, stream-function cannot be computed"
List< label > labelList
A List of labels.
vectorField pointField
pointField is a vectorField.
virtual const pointField & points() const
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
label nGeometricD() const
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
static word timeName(const scalar t, const int precision=precision_)
const cellList & cells() const
streamFunction(const word &name, const Time &runTime, const dictionary &dict)
const polyBoundaryMesh & boundaryMesh() const
Ostream & endl(Ostream &os)
dimensionedScalar sign(const dimensionedScalar &ds)
const Vector< label > & geometricD() const
label nPoints() const noexcept
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
void setResultName(const word &typeName, const word &defaultArg)
PtrList< polyPatch > polyPatchList
container classes for polyPatch
List< cell > cellList
A List of cells.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label whichPatch(const label faceIndex) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
Vector< scalar > vector
A scalar version of the templated Vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
const boundBox & bounds() const
virtual const faceList & faces() const
#define FatalErrorInFunction
const word & name() const noexcept
UList< face > faceUList
A UList of faces.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual const word & type() const =0
label nInternalFaces() const noexcept
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static tmp< T > New(Args &&... args)
defineTypeNameAndDebug(ObukhovLength, 0)
const polyBoundaryMesh & patches
const dimensionedScalar c
word name(const expressions::valueTypeCode typeCode)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
vector point
Point is a vector.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
const vectorField & faceAreas() const