1 2 3 4 5 6 7 8 9 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 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 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 160 161 162 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 255 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 348 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 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491
|
/*---------------------------------------------------------------------------* \
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
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/>.
Application
streamFunction
Description
Calculates and writes the stream function of velocity field U at each
time.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "pointFields.H"
#include "emptyPolyPatch.H"
#include "symmetryPlanePolyPatch.H"
#include "symmetryPolyPatch.H"
#include "wedgePolyPatch.H"
#include "OSspecific.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
timeSelector::addOptions();
#include "addRegionOption.H"
#include "setRootCase.H"
#include "createTime.H"
instantList timeDirs = timeSelector::select0(runTime, args);
#include "createNamedMesh.H"
label nD = mesh.nGeometricD();
if (nD != 2)
{
FatalErrorInFunction
<< "Case is not 2D, stream-function cannot be computed"
<< exit(FatalError);
}
Vector<label> slabNormal((Vector<label>::one - mesh.geometricD())/2);
const direction slabDir
(
slabNormal
& Vector<label>(Vector<label>::X, Vector<label>::Y, Vector<label>::Z)
);
scalar thickness = vector(slabNormal) & mesh.bounds().span();
const pointMesh& pMesh = pointMesh::New(mesh);
forAll(timeDirs, timeI)
{
runTime.setTime(timeDirs[timeI], timeI);
Info<< nl << "Time: " << runTime.timeName() << endl;
IOobject phiHeader
(
"phi",
runTime.timeName(),
mesh,
IOobject::NO_READ
);
if (phiHeader.headerOk())
{
mesh.readUpdate();
Info<< nl << "Reading field phi" << endl;
surfaceScalarField phi
(
IOobject
(
"phi",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
),
mesh
);
pointScalarField streamFunction
(
IOobject
(
"streamFunction",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
pMesh,
dimensionedScalar("zero", phi.dimensions(), 0.0)
);
labelList visitedPoint(mesh.nPoints());
forAll(visitedPoint, pointi)
{
visitedPoint[pointi] = 0;
}
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
)
{
const labelList& zeroPoints = bouFaces[facei];
// Zero flux face found
found = true;
forAll(zeroPoints, pointi)
{
if (visitedPoint[zeroPoints[pointi]] == 1)
{
found = false;
break;
}
}
if (found)
{
Info<< "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)
{
Info<< "zero flux boundary face not found. "
<< "Using cell as a reference."
<< endl;
const cellList& c = mesh.cells();
forAll(c, cI)
{
labelList zeroPoints = c[cI].labels(mesh.faces());
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."
<< abort(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 /= mag(edgeHat);
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 /= mag(edgeHat);
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;
}
}
Info<< ".";
if (nVisited == nVisitedOld)
{
// Find new seed. This must be a
// multiply connected domain
Info<< nl << "Exhausted a seed. Looking for new seed "
<< "(this is correct for multiply connected "
<< "domains).";
break;
}
else
{
nVisitedOld = nVisited;
}
} while (!finished);
Info<< endl;
} while (!finished);
// Normalise the stream-function by the 2D mesh thickness
streamFunction /= thickness;
streamFunction.boundaryFieldRef() = 0.0;
streamFunction.write();
}
else
{
WarningInFunction
<< "Flux field does not exist."
<< " Stream function not calculated" << endl;
}
}
Info<< "\nEnd\n" << endl;
return 0;
}
// ************************************************************************* //
|