Loading...
Searching...
No Matches
streamFunction.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016 OpenFOAM Foundation
9 Copyright (C) 2019-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "streamFunction.H"
30#include "surfaceFields.H"
31#include "pointFields.H"
32#include "emptyPolyPatch.H"
34#include "symmetryPolyPatch.H"
35#include "wedgePolyPatch.H"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
42namespace functionObjects
43{
46}
47}
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
53(
55) const
56{
57 Log << " functionObjects::" << type() << " " << name()
58 << " calculating stream-function" << endl;
59
61 const direction slabDir
62 (
63 slabNormal
65 );
66
67 const pointMesh& pMesh = pointMesh::New(mesh_);
68
69 auto tstreamFunction = tmp<pointScalarField>::New
70 (
72 (
73 "streamFunction",
75 mesh_
76 ),
77 pMesh,
78 dimensionedScalar(phi.dimensions(), Zero)
79 );
80 auto& streamFunction = tstreamFunction.ref();
81
82
83 bitSet visitedPoint(mesh_.nPoints());
84
85 label nVisited = 0;
86 label nVisitedOld = 0;
87
88 const faceUList& faces = mesh_.faces();
89 const pointField& points = mesh_.points();
90
91 const label nInternalFaces = mesh_.nInternalFaces();
92
93 vectorField unitAreas(mesh_.faceAreas());
94 unitAreas.normalise();
95
97
98 bool finished = true;
99
100 // Find the boundary face with zero flux. Set the stream function
101 // to zero on that face
102 bool found = false;
103
104 do
105 {
106 found = false;
107
108 // Check boundary faces first
109 forAll(patches, patchi)
110 {
111 const auto& pp = patches[patchi];
112 const auto& patchPhi = phi.boundaryField()[patchi];
113
114 // Skip empty, symmetry patches etc
115 if
116 (
117 patchPhi.empty()
122 )
123 {
124 continue;
125 }
126
127 forAll(pp, facei)
128 {
129 const auto& f = pp[facei];
130
131 if (magSqr(patchPhi[facei]) < SMALL)
132 {
133 // Zero flux face found
134 found = true;
135
136 for (const label pointi : f)
137 {
138 if (visitedPoint.test(pointi))
139 {
140 found = false;
141 break;
142 }
143 }
144
145 if (found)
146 {
147 Log << " Zero face: patch: " << patchi
148 << " face: " << facei << endl;
149
150 for (const label pointi : f)
151 {
152 visitedPoint.set(pointi);
153 ++nVisited;
154
155 streamFunction[pointi] = 0;
156 }
157
158 break;
159 }
160 }
161 }
162
163 if (found) break;
164 }
165
166 if (!found)
167 {
168 Log << " Zero flux boundary face not found. "
169 << "Using cell as a reference." << endl;
170
171 for (const cell& c : mesh_.cells())
172 {
173 labelList zeroPoints = c.labels(mesh_.faces());
174
175 bool found = true;
176
177 for (const label pointi : zeroPoints)
178 {
179 if (visitedPoint.test(pointi))
180 {
181 found = false;
182 break;
183 }
184 }
185
186 if (found)
187 {
188 for (const label pointi : zeroPoints)
189 {
190 visitedPoint.set(pointi);
191 ++nVisited;
192
193 streamFunction[pointi] = 0;
194 }
195
196 break;
197 }
198 else
199 {
201 << "Cannot find initialisation face or a cell."
202 << exit(FatalError);
203 }
204 }
205 }
206
207 // Loop through all faces. If one of the points on
208 // the face has the streamfunction value different
209 // from -1, all points with -1 ont that face have the
210 // streamfunction value equal to the face flux in
211 // that point plus the value in the visited point
212 do
213 {
214 finished = true;
215
216 scalar currentStreamValue(0);
217 point currentStreamPoint(Zero);
218
219 // Boundary faces first
220 forAll(patches, patchi)
221 {
222 const auto& pp = patches[patchi];
223 const auto& patchPhi = phi.boundaryField()[patchi];
224
225 // Skip empty, symmetry patches etc
226 if
227 (
228 patchPhi.empty()
233 )
234 {
235 continue;
236 }
237
238 forAll(pp, facei)
239 {
240 const auto& f = pp[facei];
241
242 // Check if the point has been visited
243 bool pointFound = false;
244
245 for (const label pointi : f)
246 {
247 if (visitedPoint.test(pointi))
248 {
249 // The point has been visited
250 currentStreamValue = streamFunction[pointi];
251 currentStreamPoint = points[pointi];
252
253 pointFound = true;
254 break;
255 }
256 }
257
258 if (!pointFound)
259 {
260 finished = false;
261 continue;
262 }
263
264
265 // Sort out other points on the face
266 for (const label pointi : f)
267 {
268 // If the point has not yet been visited
269 if (!visitedPoint.test(pointi))
270 {
271 vector edgeHat =
272 (
273 points[pointi] - currentStreamPoint
274 );
275 edgeHat.replace(slabDir, 0);
276 edgeHat.normalise();
277
278 const vector& nHat = unitAreas[facei];
279
280 if (edgeHat.y() > VSMALL)
281 {
282 visitedPoint.set(pointi);
283 ++nVisited;
284
285 streamFunction[pointi] =
286 (
287 currentStreamValue
288 + patchPhi[facei]*sign(nHat.x())
289 );
290 }
291 else if (edgeHat.y() < -VSMALL)
292 {
293 visitedPoint.set(pointi);
294 ++nVisited;
295
296 streamFunction[pointi] =
297 (
298 currentStreamValue
299 - patchPhi[facei]*sign(nHat.x())
300 );
301 }
302 else
303 {
304 if (edgeHat.x() > VSMALL)
305 {
306 visitedPoint.set(pointi);
307 ++nVisited;
308
309 streamFunction[pointi] =
310 (
311 currentStreamValue
312 + patchPhi[facei]*sign(nHat.y())
313 );
314 }
315 else if (edgeHat.x() < -VSMALL)
316 {
317 visitedPoint.set(pointi);
318 ++nVisited;
319
320 streamFunction[pointi] =
321 (
322 currentStreamValue
323 - patchPhi[facei]*sign(nHat.y())
324 );
325 }
326 }
327 }
328 }
329 }
330 }
331
332 // Internal faces next
333 for (label facei = 0; facei < nInternalFaces; ++facei)
334 {
335 const auto& f = faces[facei];
336
337 bool pointFound = false;
338
339 for (const label pointi : f)
340 {
341 // Check if the point has been visited
342 if (visitedPoint.test(pointi))
343 {
344 currentStreamValue = streamFunction[pointi];
345 currentStreamPoint = points[pointi];
346 pointFound = true;
347
348 break;
349 }
350 }
351
352 if (pointFound)
353 {
354 // Sort out other points on the face
355 for (const label pointi : f)
356 {
357 // If the point has not yet been visited
358 if (!visitedPoint.test(pointi))
359 {
360 vector edgeHat =
361 (
362 points[pointi] - currentStreamPoint
363 );
364
365 edgeHat.replace(slabDir, 0);
366 edgeHat.normalise();
367
368 const vector& nHat = unitAreas[facei];
369
370 if (edgeHat.y() > VSMALL)
371 {
372 visitedPoint.set(pointi);
373 ++nVisited;
374
375 streamFunction[pointi] =
376 (
377 currentStreamValue
378 + phi[facei]*sign(nHat.x())
379 );
380 }
381 else if (edgeHat.y() < -VSMALL)
382 {
383 visitedPoint.set(pointi);
384 ++nVisited;
385
386 streamFunction[pointi] =
387 (
388 currentStreamValue
389 - phi[facei]*sign(nHat.x())
390 );
391 }
392 }
393 }
394 }
395 else
396 {
397 finished = false;
398 }
399 }
400
401 if (nVisited == nVisitedOld)
402 {
403 // Find new seed. This must be a
404 // multiply connected domain
405 Log << " Exhausted a seed, looking for new seed "
406 << "(this is correct for multiply connected domains).";
407
408 break;
409 }
410 else
411 {
412 nVisitedOld = nVisited;
413 }
414 } while (!finished);
415 } while (!finished);
416
417 // Normalise the stream-function by the 2D mesh thickness
418 // calculate thickness here to avoid compiler oddness (#2603)
419 const scalar thickness = vector(slabNormal) & mesh_.bounds().span();
420
421 streamFunction /= thickness;
422 streamFunction.boundaryFieldRef() = Zero;
423
424 return tstreamFunction;
425}
426
427
429{
430 const auto* phiPtr = findObject<surfaceScalarField>(fieldName_);
431
432 if (phiPtr)
433 {
434 const surfaceScalarField& phi = *phiPtr;
435
436 return store(resultName_, calc(phi));
437 }
439 return false;
440}
441
442
443// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
444
446(
447 const word& name,
448 const Time& runTime,
449 const dictionary& dict
450)
451:
453{
454 setResultName(typeName, "phi");
455
456 const label nD = mesh_.nGeometricD();
457
458 if (nD != 2)
459 {
461 << "Case is not 2D, stream-function cannot be computed"
462 << exit(FatalError);
463 }
464}
465
466
467// ************************************************************************* //
#define Log
Definition PDRblock.C:28
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
virtual const word & type() const =0
Runtime type information.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
fieldExpression(const word &name, const Time &runTime, const dictionary &dict, const word &fieldName=word::null, const word &resultName=word::null)
Construct from name, Time and dictionary.
virtual bool calc()=0
Calculate the components of the field and return true if successful.
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
const fvMesh & mesh_
Reference to the fvMesh.
Computes the stream function (i.e. https://w.wiki/Ncm).
streamFunction(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
const Time & time_
Reference to the time database.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition polyMesh.C:861
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition polyMesh.C:850
label nInternalFaces() const noexcept
Number of internal faces.
label nPoints() const noexcept
Number of mesh points.
const vectorField & faceAreas() const
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
const pointField & points
const dimensionedScalar c
Speed of light in a vacuum.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool isType(const U &obj)
Check if typeid of the object and Type are identical.
Definition typeInfo.H:112
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
UList< face > faceUList
UList of faces.
Definition faceListFwd.H:43
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.