FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

streamFunction.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------* \
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Application
00025     streamFunction
00026 
00027 Description
00028     Calculates and writes the stream function of velocity field U at each time
00029 
00030 Usage
00031 
00032     - streamFunction [OPTIONS]
00033 
00034     @param -noZero \n
00035     Ignore timestep 0.
00036 
00037     @param -constant \n
00038     Include the constant directory.
00039 
00040     @param -time <time>\n
00041     Apply only to specific time.
00042 
00043     @param -latestTime \n
00044     Only apply to latest time step.
00045 
00046     @param -case <dir>\n
00047     Case directory.
00048 
00049     @param -parallel \n
00050     Run in parallel.
00051 
00052     @param -help \n
00053     Display help message.
00054 
00055     @param -doc \n
00056     Display Doxygen API documentation page for this application.
00057 
00058     @param -srcDoc \n
00059     Display Doxygen source documentation page for this application.
00060 
00061 \*---------------------------------------------------------------------------*/
00062 
00063 #include <finiteVolume/fvCFD.H>
00064 #include <OpenFOAM/pointFields.H>
00065 #include <OpenFOAM/emptyPolyPatch.H>
00066 #include <OpenFOAM/symmetryPolyPatch.H>
00067 #include <OpenFOAM/wedgePolyPatch.H>
00068 #include <OpenFOAM/OSspecific.H>
00069 
00070 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00071 //  Main program:
00072 
00073 int main(int argc, char *argv[])
00074 {
00075     timeSelector::addOptions();
00076 
00077 #   include <OpenFOAM/setRootCase.H>
00078 #   include <OpenFOAM/createTime.H>
00079 
00080     instantList timeDirs = timeSelector::select0(runTime, args);
00081 
00082 #   include <OpenFOAM/createMeshNoClear.H>
00083 
00084     pointMesh pMesh(mesh);
00085 
00086     forAll(timeDirs, timeI)
00087     {
00088         runTime.setTime(timeDirs[timeI], timeI);
00089 
00090         Info<< nl << "Time: " << runTime.timeName() << endl;
00091 
00092         IOobject phiHeader
00093         (
00094             "phi",
00095             runTime.timeName(),
00096             mesh,
00097             IOobject::NO_READ
00098         );
00099 
00100         if (phiHeader.headerOk())
00101         {
00102             mesh.readUpdate();
00103 
00104             Info<< nl << "Reading field phi" << endl;
00105 
00106             surfaceScalarField phi
00107             (
00108                 IOobject
00109                 (
00110                     "phi",
00111                     runTime.timeName(),
00112                     mesh,
00113                     IOobject::MUST_READ,
00114                     IOobject::NO_WRITE
00115                 ),
00116                 mesh
00117             );
00118 
00119             pointScalarField streamFunction
00120             (
00121                 IOobject
00122                 (
00123                     "streamFunction",
00124                     runTime.timeName(),
00125                     mesh,
00126                     IOobject::NO_READ,
00127                     IOobject::NO_WRITE
00128                 ),
00129                 pMesh,
00130                 dimensionedScalar("zero", phi.dimensions(), 0.0)
00131             );
00132 
00133             labelList visitedPoint(mesh.nPoints());
00134             forAll (visitedPoint, pointI)
00135             {
00136                 visitedPoint[pointI] = 0;
00137             }
00138             label nVisited = 0;
00139             label nVisitedOld = 0;
00140 
00141             const unallocFaceList& faces = mesh.faces();
00142             const pointField& points = mesh.points();
00143 
00144             label nInternalFaces = mesh.nInternalFaces();
00145 
00146             vectorField unitAreas = mesh.faceAreas();
00147             unitAreas /= mag(unitAreas);
00148 
00149             const polyPatchList& patches = mesh.boundaryMesh();
00150 
00151             bool finished = true;
00152 
00153             // Find the boundary face with zero flux. set the stream function
00154             // to zero on that face
00155             bool found = false;
00156 
00157             do
00158             {
00159                 found = false;
00160 
00161                 forAll (patches, patchI)
00162                 {
00163                     const primitivePatch& bouFaces = patches[patchI];
00164 
00165                     if (!isType<emptyPolyPatch>(patches[patchI]))
00166                     {
00167                         forAll (bouFaces, faceI)
00168                         {
00169                             if
00170                             (
00171                                 magSqr(phi.boundaryField()[patchI][faceI])
00172                               < SMALL
00173                             )
00174                             {
00175                                 const labelList& zeroPoints = bouFaces[faceI];
00176 
00177                                 // Zero flux face found
00178                                 found = true;
00179 
00180                                 forAll (zeroPoints, pointI)
00181                                 {
00182                                     if (visitedPoint[zeroPoints[pointI]] == 1)
00183                                     {
00184                                         found = false;
00185                                         break;
00186                                     }
00187                                 }
00188 
00189                                 if (found)
00190                                 {
00191                                     Info<< "Zero face: patch: " << patchI
00192                                         << "    face: " << faceI << endl;
00193 
00194                                     forAll (zeroPoints, pointI)
00195                                     {
00196                                         streamFunction[zeroPoints[pointI]] = 0;
00197                                         visitedPoint[zeroPoints[pointI]] = 1;
00198                                         nVisited++;
00199                                     }
00200 
00201                                     break;
00202                                 }
00203                             }
00204                         }
00205                     }
00206 
00207                     if (found) break;
00208                 }
00209 
00210                 if (!found)
00211                 {
00212                     Info<< "zero flux boundary face not found. "
00213                         << "Using cell as a reference."
00214                         << endl;
00215 
00216                     const cellList& c = mesh.cells();
00217 
00218                     forAll (c, cI)
00219                     {
00220                         labelList zeroPoints = c[cI].labels(mesh.faces());
00221 
00222                         bool found = true;
00223 
00224                         forAll (zeroPoints, pointI)
00225                         {
00226                             if (visitedPoint[zeroPoints[pointI]] == 1)
00227                             {
00228                                 found = false;
00229                                 break;
00230                             }
00231                         }
00232 
00233                         if (found)
00234                         {
00235                             forAll (zeroPoints, pointI)
00236                             {
00237                                 streamFunction[zeroPoints[pointI]] = 0.0;
00238                                 visitedPoint[zeroPoints[pointI]] = 1;
00239                                 nVisited++;
00240                             }
00241 
00242                             break;
00243                         }
00244                         else
00245                         {
00246                             FatalErrorIn(args.executable())
00247                                 << "Cannot find initialisation face or a cell."
00248                                 << abort(FatalError);
00249                         }
00250                     }
00251                 }
00252 
00253                 // Loop through all faces. If one of the points on
00254                 // the face has the streamfunction value different
00255                 // from -1, all points with -1 ont that face have the
00256                 // streamfunction value equal to the face flux in
00257                 // that point plus the value in the visited point
00258                 do
00259                 {
00260                      finished = true;
00261 
00262                      for
00263                      (
00264                          label faceI = nInternalFaces;
00265                          faceI<faces.size();
00266                          faceI++
00267                      )
00268                      {
00269                          const labelList& curBPoints = faces[faceI];
00270                          bool bPointFound = false;
00271 
00272                          scalar currentBStream = 0.0;
00273                          vector currentBStreamPoint(0, 0, 0);
00274 
00275                          forAll (curBPoints, pointI)
00276                          {
00277                              // Check if the point has been visited
00278                              if (visitedPoint[curBPoints[pointI]] == 1)
00279                              {
00280                                  // The point has been visited
00281                                  currentBStream =
00282                                      streamFunction[curBPoints[pointI]];
00283                                  currentBStreamPoint =
00284                                      points[curBPoints[pointI]];
00285 
00286                                  bPointFound = true;
00287 
00288                                  break;
00289                              }
00290                          }
00291 
00292                          if (bPointFound)
00293                          {
00294                              // Sort out other points on the face
00295                              forAll (curBPoints, pointI)
00296                              {
00297                                  // Check if the point has been visited
00298                                  if (visitedPoint[curBPoints[pointI]] == 0)
00299                                  {
00300                                      label patchNo =
00301                                          mesh.boundaryMesh().whichPatch(faceI);
00302 
00303                                      if
00304                                      (
00305                                         !isType<emptyPolyPatch>
00306                                          (patches[patchNo])
00307                                      && !isType<symmetryPolyPatch>
00308                                          (patches[patchNo])
00309                                      && !isType<wedgePolyPatch>
00310                                          (patches[patchNo])
00311                                      )
00312                                      {
00313                                          label faceNo =
00314                                              mesh.boundaryMesh()[patchNo]
00315                                              .whichFace(faceI);
00316 
00317                                          vector edgeHat =
00318                                              points[curBPoints[pointI]]
00319                                              - currentBStreamPoint;
00320                                          edgeHat.replace(vector::Z, 0);
00321                                          edgeHat /= mag(edgeHat);
00322 
00323                                          vector nHat = unitAreas[faceI];
00324 
00325                                          if (edgeHat.y() > VSMALL)
00326                                          {
00327                                              visitedPoint[curBPoints[pointI]] =
00328                                                  1;
00329                                              nVisited++;
00330 
00331                                              streamFunction[curBPoints[pointI]]
00332                                                  =
00333                                                  currentBStream
00334                                                + phi.boundaryField()
00335                                                  [patchNo][faceNo]
00336                                                  *sign(nHat.x());
00337                                          }
00338                                          else if (edgeHat.y() < -VSMALL)
00339                                          {
00340                                              visitedPoint[curBPoints[pointI]] =
00341                                                  1;
00342                                              nVisited++;
00343 
00344                                              streamFunction[curBPoints[pointI]]
00345                                                  =
00346                                                  currentBStream
00347                                                - phi.boundaryField()
00348                                                  [patchNo][faceNo]
00349                                                  *sign(nHat.x());
00350                                          }
00351                                          else
00352                                          {
00353                                              if (edgeHat.x() > VSMALL)
00354                                              {
00355                                                  visitedPoint
00356                                                      [curBPoints[pointI]] = 1;
00357                                                  nVisited++;
00358 
00359                                                  streamFunction
00360                                                      [curBPoints[pointI]] =
00361                                                      currentBStream
00362                                                    + phi.boundaryField()
00363                                                      [patchNo][faceNo]
00364                                                      *sign(nHat.y());
00365                                              }
00366                                              else if (edgeHat.x() < -VSMALL)
00367                                              {
00368                                                  visitedPoint
00369                                                      [curBPoints[pointI]] = 1;
00370                                                  nVisited++;
00371 
00372                                                  streamFunction
00373                                                      [curBPoints[pointI]] =
00374                                                      currentBStream
00375                                                    - phi.boundaryField()
00376                                                      [patchNo][faceNo]
00377                                                      *sign(nHat.y());
00378                                              }
00379                                          }
00380                                      }
00381                                  }
00382                              }
00383                          }
00384                          else
00385                          {
00386                              finished = false;
00387                          }
00388                      }
00389 
00390                      for (label faceI=0; faceI<nInternalFaces; faceI++)
00391                      {
00392                          // Get the list of point labels for the face
00393                          const labelList& curPoints = faces[faceI];
00394 
00395                          bool pointFound = false;
00396 
00397                          scalar currentStream = 0.0;
00398                          point currentStreamPoint(0, 0, 0);
00399 
00400                          forAll (curPoints, pointI)
00401                          {
00402                              // Check if the point has been visited
00403                              if (visitedPoint[curPoints[pointI]] == 1)
00404                              {
00405                                  // The point has been visited
00406                                  currentStream =
00407                                      streamFunction[curPoints[pointI]];
00408                                  currentStreamPoint =
00409                                      points[curPoints[pointI]];
00410                                  pointFound = true;
00411 
00412                                  break;
00413                              }
00414                          }
00415 
00416                          if (pointFound)
00417                          {
00418                              // Sort out other points on the face
00419                              forAll (curPoints, pointI)
00420                              {
00421                                  // Check if the point has been visited
00422                                  if (visitedPoint[curPoints[pointI]] == 0)
00423                                  {
00424                                      vector edgeHat =
00425                                          points[curPoints[pointI]]
00426                                        - currentStreamPoint;
00427 
00428                                      edgeHat.replace(vector::Z, 0);
00429                                      edgeHat /= mag(edgeHat);
00430 
00431                                      vector nHat = unitAreas[faceI];
00432 
00433                                      if (edgeHat.y() > VSMALL)
00434                                      {
00435                                          visitedPoint[curPoints[pointI]] = 1;
00436                                          nVisited++;
00437 
00438                                          streamFunction[curPoints[pointI]] =
00439                                              currentStream
00440                                            + phi[faceI]*sign(nHat.x());
00441                                      }
00442                                      else if (edgeHat.y() < -VSMALL)
00443                                      {
00444                                          visitedPoint[curPoints[pointI]] = 1;
00445                                          nVisited++;
00446 
00447                                          streamFunction[curPoints[pointI]] =
00448                                              currentStream
00449                                            - phi[faceI]*sign(nHat.x());
00450                                      }
00451                                  }
00452                              }
00453                          }
00454                          else
00455                          {
00456                              finished = false;
00457                          }
00458                      }
00459 
00460                      Info<< ".";
00461 //                     Info<< "One pass, n visited = " << nVisited << endl;
00462 
00463                      if (nVisited == nVisitedOld)
00464                      {
00465                          // Find new seed.  This must be a
00466                          // multiply connected domain
00467                          Info<< nl << "Exhausted a seed. Looking for new seed "
00468                              << "(this is correct for multiply connected "
00469                              << "domains).";
00470 
00471                          break;
00472                      }
00473                      else
00474                      {
00475                          nVisitedOld = nVisited;
00476                      }
00477                 } while (!finished);
00478 
00479                 Info << endl;
00480             } while (!finished);
00481 
00482             streamFunction.boundaryField() = 0.0;
00483             streamFunction.write();
00484         }
00485         else
00486         {
00487             WarningIn(args.executable())
00488                 << "Flux field does not exist."
00489                 << " Stream function not calculated" << endl;
00490         }
00491     }
00492 
00493     Info<< "\nEnd\n" << endl;
00494 
00495     return 0;
00496 }
00497 
00498 
00499 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines