Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2019-2023 OpenCFD Ltd.
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "streamFunction.H"
#include "surfaceFields.H"
#include "pointFields.H"
#include "emptyPolyPatch.H"
#include "symmetryPlanePolyPatch.H"
#include "symmetryPolyPatch.H"
#include "wedgePolyPatch.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace functionObjects
{
defineTypeNameAndDebug(streamFunction, 0);
addToRunTimeSelectionTable(functionObject, streamFunction, dictionary);
}
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
Foam::tmp<Foam::pointScalarField> Foam::functionObjects::streamFunction::calc
(
const surfaceScalarField& phi
) const
{
Log << " functionObjects::" << type() << " " << name()
<< " calculating stream-function" << endl;
Vector<label> slabNormal((Vector<label>::one - mesh_.geometricD())/2);
const direction slabDir
(
slabNormal
& Vector<label>(Vector<label>::X, Vector<label>::Y, Vector<label>::Z)
);
const pointMesh& pMesh = pointMesh::New(mesh_);
auto tstreamFunction = tmp<pointScalarField>::New
(
(
"streamFunction",
time_.timeName(),
mesh_
),
pMesh,
dimensionedScalar(phi.dimensions(), Zero)
);
pointScalarField& streamFunction = tstreamFunction.ref();
labelList visitedPoint(mesh_.nPoints(), Zero);
label nVisited = 0;
label nVisitedOld = 0;
const faceUList& faces = mesh_.faces();
const pointField& points = mesh_.points();
label nInternalFaces = mesh_.nInternalFaces();
vectorField unitAreas(mesh_.faceAreas());
unitAreas /= mag(unitAreas);
const polyPatchList& patches = mesh_.boundaryMesh();
bool finished = true;
// Find the boundary face with zero flux. Set the stream function
// to zero on that face
bool found = false;
do
{
found = false;
forAll(patches, patchi)
{
const primitivePatch& bouFaces = patches[patchi];
if (!isType<emptyPolyPatch>(patches[patchi]))
{
forAll(bouFaces, facei)
{
if (magSqr(phi.boundaryField()[patchi][facei]) < SMALL)
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
{
const labelList& zeroPoints = bouFaces[facei];
// Zero flux face found
found = true;
forAll(zeroPoints, pointi)
{
if (visitedPoint[zeroPoints[pointi]] == 1)
{
found = false;
break;
}
}
if (found)
{
Log << " Zero face: patch: " << patchi
<< " face: " << facei << endl;
forAll(zeroPoints, pointi)
{
streamFunction[zeroPoints[pointi]] = 0;
visitedPoint[zeroPoints[pointi]] = 1;
nVisited++;
}
break;
}
}
}
}
if (found) break;
}
if (!found)
{
Log << " Zero flux boundary face not found. "
<< "Using cell as a reference."
<< endl;
const cellList& c = mesh_.cells();
{
labelList zeroPoints = c[ci].labels(mesh_.faces());
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
bool found = true;
forAll(zeroPoints, pointi)
{
if (visitedPoint[zeroPoints[pointi]] == 1)
{
found = false;
break;
}
}
if (found)
{
forAll(zeroPoints, pointi)
{
streamFunction[zeroPoints[pointi]] = 0.0;
visitedPoint[zeroPoints[pointi]] = 1;
nVisited++;
}
break;
}
else
{
FatalErrorInFunction
<< "Cannot find initialisation face or a cell."
<< exit(FatalError);
}
}
}
// Loop through all faces. If one of the points on
// the face has the streamfunction value different
// from -1, all points with -1 ont that face have the
// streamfunction value equal to the face flux in
// that point plus the value in the visited point
do
{
finished = true;
for (label facei = nInternalFaces; facei<faces.size(); facei++)
{
const labelList& curBPoints = faces[facei];
bool bPointFound = false;
scalar currentBStream = 0.0;
vector currentBStreamPoint(0, 0, 0);
forAll(curBPoints, pointi)
{
// Check if the point has been visited
if (visitedPoint[curBPoints[pointi]] == 1)
{
// The point has been visited
currentBStream = streamFunction[curBPoints[pointi]];
currentBStreamPoint = points[curBPoints[pointi]];
bPointFound = true;
break;
}
}
if (bPointFound)
{
// Sort out other points on the face
forAll(curBPoints, pointi)
{
// Check if the point has been visited
if (visitedPoint[curBPoints[pointi]] == 0)
{
label patchNo =
mesh_.boundaryMesh().whichPatch(facei);
if
(
!isType<emptyPolyPatch>(patches[patchNo])
&& !isType<symmetryPlanePolyPatch>
(patches[patchNo])
&& !isType<symmetryPolyPatch>(patches[patchNo])
&& !isType<wedgePolyPatch>(patches[patchNo])
)
{
label faceNo =
mesh_.boundaryMesh()[patchNo]
.whichFace(facei);
vector edgeHat =
points[curBPoints[pointi]]
- currentBStreamPoint;
edgeHat.replace(slabDir, 0);
edgeHat.normalise();
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
vector nHat = unitAreas[facei];
if (edgeHat.y() > VSMALL)
{
visitedPoint[curBPoints[pointi]] = 1;
nVisited++;
streamFunction[curBPoints[pointi]] =
currentBStream
+ phi.boundaryField()[patchNo][faceNo]
*sign(nHat.x());
}
else if (edgeHat.y() < -VSMALL)
{
visitedPoint[curBPoints[pointi]] = 1;
nVisited++;
streamFunction[curBPoints[pointi]] =
currentBStream
- phi.boundaryField()[patchNo][faceNo]
*sign(nHat.x());
}
else
{
if (edgeHat.x() > VSMALL)
{
visitedPoint[curBPoints[pointi]] = 1;
nVisited++;
streamFunction[curBPoints[pointi]] =
currentBStream
+ phi.boundaryField()[patchNo][faceNo]
*sign(nHat.y());
}
else if (edgeHat.x() < -VSMALL)
{
visitedPoint[curBPoints[pointi]] = 1;
nVisited++;
streamFunction[curBPoints[pointi]] =
currentBStream
- phi.boundaryField()[patchNo][faceNo]
*sign(nHat.y());
}
}
}
}
}
}
else
{
finished = false;
}
}
for (label facei=0; facei<nInternalFaces; facei++)
{
// Get the list of point labels for the face
const labelList& curPoints = faces[facei];
bool pointFound = false;
scalar currentStream = 0.0;
point currentStreamPoint(0, 0, 0);
forAll(curPoints, pointi)
{
// Check if the point has been visited
if (visitedPoint[curPoints[pointi]] == 1)
{
// The point has been visited
currentStream = streamFunction[curPoints[pointi]];
currentStreamPoint = points[curPoints[pointi]];
pointFound = true;
break;
}
}
if (pointFound)
{
// Sort out other points on the face
forAll(curPoints, pointi)
{
// Check if the point has been visited
if (visitedPoint[curPoints[pointi]] == 0)
{
vector edgeHat =
points[curPoints[pointi]] - currentStreamPoint;
edgeHat.replace(slabDir, 0);
edgeHat.normalise();
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
vector nHat = unitAreas[facei];
if (edgeHat.y() > VSMALL)
{
visitedPoint[curPoints[pointi]] = 1;
nVisited++;
streamFunction[curPoints[pointi]] =
currentStream
+ phi[facei]*sign(nHat.x());
}
else if (edgeHat.y() < -VSMALL)
{
visitedPoint[curPoints[pointi]] = 1;
nVisited++;
streamFunction[curPoints[pointi]] =
currentStream
- phi[facei]*sign(nHat.x());
}
}
}
}
else
{
finished = false;
}
}
if (nVisited == nVisitedOld)
{
// Find new seed. This must be a
// multiply connected domain
Log << " Exhausted a seed, looking for new seed "
<< "(this is correct for multiply connected domains).";
break;
}
else
{
nVisitedOld = nVisited;
}
} while (!finished);
} while (!finished);
// Normalise the stream-function by the 2D mesh thickness
// calculate thickness here to avoid compiler oddness (#2603)
const scalar thickness = vector(slabNormal) & mesh_.bounds().span();
streamFunction /= thickness;
streamFunction.boundaryFieldRef() = 0.0;
return tstreamFunction;
}
bool Foam::functionObjects::streamFunction::calc()
{
const auto* phiPtr = findObject<surfaceScalarField>(fieldName_);
if (phiPtr)
{
const surfaceScalarField& phi = *phiPtr;
return store(resultName_, calc(phi));
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::functionObjects::streamFunction::streamFunction
(
const word& name,
const Time& runTime,
const dictionary& dict
)
:
fieldExpression(name, runTime, dict, "phi")
{
setResultName(typeName, "phi");
const label nD = mesh_.nGeometricD();
if (nD != 2)
{
FatalErrorInFunction
<< "Case is not 2D, stream-function cannot be computed"
<< exit(FatalError);
}
}
// ************************************************************************* //