43 int main(
int argc,
char *argv[])
45 timeSelector::addOptions();
60 <<
"Case is not 2D, stream-function cannot be computed" 71 scalar thickness =
vector(slabNormal) &
mesh.bounds().span();
77 runTime.setTime(timeDirs[timeI], timeI);
79 Info<<
nl <<
"Time: " << runTime.timeName() <<
endl;
89 if (phiHeader.headerOk())
123 forAll(visitedPoint, pointi)
125 visitedPoint[pointi] = 0;
128 label nVisitedOld = 0;
133 label nInternalFaces =
mesh.nInternalFaces();
136 unitAreas /=
mag(unitAreas);
140 bool finished =
true;
154 if (!isType<emptyPolyPatch>(patches[patchi]))
164 const labelList& zeroPoints = bouFaces[facei];
169 forAll(zeroPoints, pointi)
171 if (visitedPoint[zeroPoints[pointi]] == 1)
180 Info<<
"Zero face: patch: " << patchi
181 <<
" face: " << facei <<
endl;
183 forAll(zeroPoints, pointi)
185 streamFunction[zeroPoints[pointi]] = 0;
186 visitedPoint[zeroPoints[pointi]] = 1;
201 Info<<
"zero flux boundary face not found. " 202 <<
"Using cell as a reference." 213 forAll(zeroPoints, pointi)
215 if (visitedPoint[zeroPoints[pointi]] == 1)
224 forAll(zeroPoints, pointi)
226 streamFunction[zeroPoints[pointi]] = 0.0;
227 visitedPoint[zeroPoints[pointi]] = 1;
236 <<
"Cannot find initialisation face or a cell." 253 label facei = nInternalFaces;
258 const labelList& curBPoints = faces[facei];
259 bool bPointFound =
false;
261 scalar currentBStream = 0.0;
262 vector currentBStreamPoint(0, 0, 0);
264 forAll(curBPoints, pointi)
267 if (visitedPoint[curBPoints[pointi]] == 1)
271 streamFunction[curBPoints[pointi]];
272 currentBStreamPoint =
273 points[curBPoints[pointi]];
284 forAll(curBPoints, pointi)
287 if (visitedPoint[curBPoints[pointi]] == 0)
290 mesh.boundaryMesh().whichPatch(facei);
294 !isType<emptyPolyPatch>
296 && !isType<symmetryPlanePolyPatch>
298 && !isType<symmetryPolyPatch>
300 && !isType<wedgePolyPatch>
305 mesh.boundaryMesh()[patchNo]
309 points[curBPoints[pointi]]
310 - currentBStreamPoint;
312 edgeHat /=
mag(edgeHat);
314 vector nHat = unitAreas[facei];
316 if (edgeHat.y() > VSMALL)
318 visitedPoint[curBPoints[pointi]] =
322 streamFunction[curBPoints[pointi]]
325 +
phi.boundaryField()
329 else if (edgeHat.y() < -VSMALL)
331 visitedPoint[curBPoints[pointi]] =
335 streamFunction[curBPoints[pointi]]
338 -
phi.boundaryField()
344 if (edgeHat.x() > VSMALL)
347 [curBPoints[pointi]] = 1;
351 [curBPoints[pointi]] =
353 +
phi.boundaryField()
357 else if (edgeHat.x() < -VSMALL)
360 [curBPoints[pointi]] = 1;
364 [curBPoints[pointi]] =
366 -
phi.boundaryField()
381 for (
label facei=0; facei<nInternalFaces; facei++)
384 const labelList& curPoints = faces[facei];
386 bool pointFound =
false;
388 scalar currentStream = 0.0;
389 point currentStreamPoint(0, 0, 0);
394 if (visitedPoint[curPoints[pointi]] == 1)
398 streamFunction[curPoints[pointi]];
400 points[curPoints[pointi]];
413 if (visitedPoint[curPoints[pointi]] == 0)
416 points[curPoints[pointi]]
417 - currentStreamPoint;
420 edgeHat /=
mag(edgeHat);
422 vector nHat = unitAreas[facei];
424 if (edgeHat.y() > VSMALL)
426 visitedPoint[curPoints[pointi]] = 1;
429 streamFunction[curPoints[pointi]] =
433 else if (edgeHat.y() < -VSMALL)
435 visitedPoint[curPoints[pointi]] = 1;
438 streamFunction[curPoints[pointi]] =
453 if (nVisited == nVisitedOld)
457 Info<<
nl <<
"Exhausted a seed. Looking for new seed " 458 <<
"(this is correct for multiply connected " 465 nVisitedOld = nVisited;
473 streamFunction /= thickness;
474 streamFunction.boundaryFieldRef() = 0.0;
475 streamFunction.
write();
480 <<
"Flux field does not exist." 481 <<
" Stream function not calculated" <<
endl;
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
PtrList< polyPatch > polyPatchList
container classes for polyPatch
dimensionedScalar sign(const dimensionedScalar &ds)
List< instant > instantList
List of instants.
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Ostream & endl(Ostream &os)
Add newline and flush stream.
volVectorField vectorField(fieldObject, mesh)
static const Vector< label > one
Vector< scalar > vector
A scalar version of the templated Vector.
void replace(const direction, const Cmpt &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
vectorField pointField
pointField is a vectorField.
List< label > labelList
A List of labels.
errorManip< error > abort(error &err)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
PrimitivePatch< face, SubList, const pointField & > primitivePatch
Addressing for a faceList slice.
vector point
Point is a vector.
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual Ostream & write(const token &)=0
Write next token to stream.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Foam::argList args(argc, argv)
List< cell > cellList
list of cells