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

mirrorFvMesh.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 \*---------------------------------------------------------------------------*/
00025 
00026 #include "mirrorFvMesh.H"
00027 #include <OpenFOAM/Time.H>
00028 #include <OpenFOAM/plane.H>
00029 
00030 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00031 
00032 Foam::mirrorFvMesh::mirrorFvMesh(const IOobject& io)
00033 :
00034     fvMesh(io),
00035     mirrorMeshDict_
00036     (
00037         IOobject
00038         (
00039             "mirrorMeshDict",
00040             time().system(),
00041             *this,
00042             IOobject::MUST_READ,
00043             IOobject::NO_WRITE
00044         )
00045     ),
00046     mirrorMeshPtr_(NULL)
00047 {
00048     plane mirrorPlane(mirrorMeshDict_);
00049 
00050     scalar planeTolerance
00051     (
00052         readScalar(mirrorMeshDict_.lookup("planeTolerance"))
00053     );
00054 
00055     const pointField& oldPoints = points();
00056     const faceList& oldFaces = faces();
00057     const cellList& oldCells = cells();
00058     const label nOldInternalFaces = nInternalFaces();
00059     const polyPatchList& oldPatches = boundaryMesh();
00060 
00061     // Mirror the points
00062     Info << "Mirroring points. Old points: " << oldPoints.size();
00063 
00064     pointField newPoints(2*oldPoints.size());
00065     label nNewPoints = 0;
00066 
00067     labelList mirrorPointLookup(oldPoints.size(), -1);
00068 
00069     // Grab the old points
00070     forAll (oldPoints, pointI)
00071     {
00072         newPoints[nNewPoints] = oldPoints[pointI];
00073         nNewPoints++;
00074     }
00075 
00076     forAll (oldPoints, pointI)
00077     {
00078         scalar alpha =
00079             mirrorPlane.normalIntersect
00080             (
00081                 oldPoints[pointI],
00082                 mirrorPlane.normal()
00083             );
00084 
00085         // Check plane on tolerance
00086         if (mag(alpha) > planeTolerance)
00087         {
00088             // The point gets mirrored
00089             newPoints[nNewPoints] =
00090                 oldPoints[pointI] + 2.0*alpha*mirrorPlane.normal();
00091 
00092             // remember the point correspondence
00093             mirrorPointLookup[pointI] = nNewPoints;
00094             nNewPoints++;
00095         }
00096         else
00097         {
00098             // The point is on the plane and does not get mirrored
00099             // Adjust plane location
00100             newPoints[nNewPoints] =
00101                 oldPoints[pointI] + alpha*mirrorPlane.normal();
00102 
00103             mirrorPointLookup[pointI] = pointI;
00104         }
00105     }
00106 
00107     // Reset the size of the point list
00108     Info << " New points: " << nNewPoints << endl;
00109     newPoints.setSize(nNewPoints);
00110 
00111     Info << "Mirroring faces. Old faces: " << oldFaces.size();
00112 
00113     // Algorithm:
00114     // During mirroring, the faces that were previously boundary faces
00115     // in the mirror plane may become ineternal faces. In order to
00116     // deal with the ordering of the faces, the algorithm is split
00117     // into two parts.  For original faces, the internal faces are
00118     // distributed to their owner cells.  Once all internal faces are
00119     // distributed, the boundary faces are visited and if they are in
00120     // the mirror plane they are added to the master cells (the future
00121     // boundary faces are not touched).  After the first phase, the
00122     // internal faces are collected in the cell order and numbering
00123     // information is added.  Then, the internal faces are mirrored
00124     // and the face numbering data is stored for the mirrored section.
00125     // Once all the internal faces are mirrored, the boundary faces
00126     // are added by mirroring the faces patch by patch.
00127 
00128     // Distribute internal faces
00129     labelListList newCellFaces(oldCells.size());
00130 
00131     const unallocLabelList& oldOwnerStart = lduAddr().ownerStartAddr();
00132 
00133     forAll (newCellFaces, cellI)
00134     {
00135         labelList& curFaces = newCellFaces[cellI];
00136 
00137         const label s = oldOwnerStart[cellI];
00138         const label e = oldOwnerStart[cellI + 1];
00139 
00140         curFaces.setSize(e - s);
00141 
00142         forAll (curFaces, i)
00143         {
00144             curFaces[i] = s + i;
00145         }
00146     }
00147 
00148     // Distribute boundary faces.  Remember the faces that have been inserted
00149     // as internal
00150     boolListList insertedBouFace(oldPatches.size());
00151 
00152     forAll (oldPatches, patchI)
00153     {
00154         const polyPatch& curPatch = oldPatches[patchI];
00155 
00156         if (curPatch.coupled())
00157         {
00158             WarningIn("mirrorFvMesh::mirrorFvMesh(const IOobject&)")
00159                 << "Found coupled patch " << curPatch.name() << endl
00160                 << "    Mirroring faces on coupled patches destroys"
00161                 << " the ordering. This might be fixed by running a dummy"
00162                 << " createPatch afterwards." << endl;
00163         }
00164 
00165         boolList& curInsBouFace = insertedBouFace[patchI];
00166 
00167         curInsBouFace.setSize(curPatch.size());
00168         curInsBouFace = false;
00169 
00170         // Get faceCells for face insertion
00171         const unallocLabelList& curFaceCells = curPatch.faceCells();
00172         const label curStart = curPatch.start();
00173 
00174         forAll (curPatch, faceI)
00175         {
00176             // Find out if the mirrored face is identical to the
00177             // original.  If so, the face needs to become internal and
00178             // added to its owner cell
00179             const face& origFace = curPatch[faceI];
00180 
00181             face mirrorFace(origFace.size());
00182             forAll (mirrorFace, pointI)
00183             {
00184                 mirrorFace[pointI] = mirrorPointLookup[origFace[pointI]];
00185             }
00186 
00187             if (origFace == mirrorFace)
00188             {
00189                 // The mirror is identical to current face.  This will
00190                 // become an internal face
00191                 const label oldSize = newCellFaces[curFaceCells[faceI]].size();
00192 
00193                 newCellFaces[curFaceCells[faceI]].setSize(oldSize + 1);
00194                 newCellFaces[curFaceCells[faceI]][oldSize] = curStart + faceI;
00195 
00196                 curInsBouFace[faceI] = true;
00197             }
00198         }
00199     }
00200 
00201     // Construct the new list of faces.  Boundary faces are added
00202     // last, cush that each patch is mirrored separately.  The
00203     // addressing is stored in two separate arrays: first for the
00204     // original cells (face order has changed) and then for the
00205     // mirrored cells.
00206     labelList masterFaceLookup(oldFaces.size(), -1);
00207     labelList mirrorFaceLookup(oldFaces.size(), -1);
00208 
00209     faceList newFaces(2*oldFaces.size());
00210     label nNewFaces = 0;
00211 
00212     // Insert original (internal) faces
00213     forAll (newCellFaces, cellI)
00214     {
00215         const labelList& curCellFaces = newCellFaces[cellI];
00216 
00217         forAll (curCellFaces, cfI)
00218         {
00219             newFaces[nNewFaces] = oldFaces[curCellFaces[cfI]];
00220             masterFaceLookup[curCellFaces[cfI]] = nNewFaces;
00221 
00222             nNewFaces++;
00223         }
00224     }
00225 
00226     // Mirror internal faces
00227     for (label faceI = 0; faceI < nOldInternalFaces; faceI++)
00228     {
00229         const face& oldFace = oldFaces[faceI];
00230         face& nf = newFaces[nNewFaces];
00231         nf.setSize(oldFace.size());
00232 
00233         nf[0] = mirrorPointLookup[oldFace[0]];
00234 
00235         for (label i = 1; i < oldFace.size(); i++)
00236         {
00237             nf[i] = mirrorPointLookup[oldFace[oldFace.size() - i]];
00238         }
00239 
00240         mirrorFaceLookup[faceI] = nNewFaces;
00241         nNewFaces++;
00242     }
00243 
00244     // Mirror boundary faces patch by patch
00245 
00246     wordList newPatchTypes(boundary().size());
00247     wordList newPatchNames(boundary().size());
00248     labelList newPatchSizes(boundary().size(), -1);
00249     labelList newPatchStarts(boundary().size(), -1);
00250     label nNewPatches = 0;
00251 
00252     forAll (boundaryMesh(), patchI)
00253     {
00254         const label curPatchSize = boundaryMesh()[patchI].size();
00255         const label curPatchStart = boundaryMesh()[patchI].start();
00256         const boolList& curInserted = insertedBouFace[patchI];
00257 
00258         newPatchStarts[nNewPatches] = nNewFaces;
00259 
00260         // Master side
00261         for (label faceI = 0; faceI < curPatchSize; faceI++)
00262         {
00263             // Check if the face has already been added.  If not, add it and
00264             // insert the numbering details.
00265             if (!curInserted[faceI])
00266             {
00267                 newFaces[nNewFaces] = oldFaces[curPatchStart + faceI];
00268 
00269                 masterFaceLookup[curPatchStart + faceI] = nNewFaces;
00270                 nNewFaces++;
00271             }
00272         }
00273 
00274         // Mirror side
00275         for (label faceI = 0; faceI < curPatchSize; faceI++)
00276         {
00277             // Check if the face has already been added.  If not, add it and
00278             // insert the numbering details.
00279             if (!curInserted[faceI])
00280             {
00281                 const face& oldFace = oldFaces[curPatchStart + faceI];
00282                 face& nf = newFaces[nNewFaces];
00283                 nf.setSize(oldFace.size());
00284 
00285                 nf[0] = mirrorPointLookup[oldFace[0]];
00286 
00287                 for (label i = 1; i < oldFace.size(); i++)
00288                 {
00289                     nf[i] = mirrorPointLookup[oldFace[oldFace.size() - i]];
00290                 }
00291 
00292                 mirrorFaceLookup[curPatchStart + faceI] = nNewFaces;
00293                 nNewFaces++;
00294             }
00295             else
00296             {
00297                 // Grab the index of the master face for the mirror side
00298                 mirrorFaceLookup[curPatchStart + faceI] =
00299                     masterFaceLookup[curPatchStart + faceI];
00300             }
00301         }
00302 
00303         // If patch exists, grab the name and type of the original patch
00304         if (nNewFaces > newPatchStarts[nNewPatches])
00305         {
00306             newPatchTypes[nNewPatches] = boundaryMesh()[patchI].type();
00307             newPatchNames[nNewPatches] = boundaryMesh()[patchI].name();
00308             newPatchSizes[nNewPatches] =
00309                 nNewFaces - newPatchStarts[nNewPatches];
00310 
00311             nNewPatches++;
00312         }
00313     }
00314 
00315     // Tidy up the lists
00316     newFaces.setSize(nNewFaces);
00317     Info << " New faces: " << nNewFaces << endl;
00318 
00319     newPatchTypes.setSize(nNewPatches);
00320     newPatchNames.setSize(nNewPatches);
00321     newPatchSizes.setSize(nNewPatches);
00322     newPatchStarts.setSize(nNewPatches);
00323 
00324     Info << "Mirroring patches. Old patches: " << boundary().size()
00325         << " New patches: " << nNewPatches << endl;
00326 
00327     Info<< "Mirroring cells.  Old cells: " << oldCells.size()
00328         << " New cells: " << 2*oldCells.size() << endl;
00329 
00330     cellList newCells(2*oldCells.size());
00331     label nNewCells = 0;
00332 
00333     // Grab the original cells.  Take care of face renumbering.
00334     forAll (oldCells, cellI)
00335     {
00336         const cell& oc = oldCells[cellI];
00337 
00338         cell& nc = newCells[nNewCells];
00339         nc.setSize(oc.size());
00340 
00341         forAll (oc, i)
00342         {
00343             nc[i] = masterFaceLookup[oc[i]];
00344         }
00345 
00346         nNewCells++;
00347     }
00348 
00349     // Mirror the cells
00350     forAll (oldCells, cellI)
00351     {
00352         const cell& oc = oldCells[cellI];
00353 
00354         cell& nc = newCells[nNewCells];
00355         nc.setSize(oc.size());
00356 
00357         forAll (oc, i)
00358         {
00359             nc[i] = mirrorFaceLookup[oc[i]];
00360         }
00361 
00362         nNewCells++;
00363     }
00364 
00365     // Mirror the cell shapes
00366     Info << "Mirroring cell shapes." << endl;
00367 
00368     Info << nl << "Creating new mesh" << endl;
00369     mirrorMeshPtr_ = new fvMesh
00370     (
00371         io,
00372         xferMove(newPoints),
00373         xferMove(newFaces),
00374         xferMove(newCells)
00375     );
00376 
00377     fvMesh& pMesh = *mirrorMeshPtr_;
00378 
00379     // Add the boundary patches
00380     List<polyPatch*> p(newPatchTypes.size());
00381 
00382     forAll (p, patchI)
00383     {
00384         p[patchI] = polyPatch::New
00385         (
00386             newPatchTypes[patchI],
00387             newPatchNames[patchI],
00388             newPatchSizes[patchI],
00389             newPatchStarts[patchI],
00390             patchI,
00391             pMesh.boundaryMesh()
00392         ).ptr();
00393     }
00394 
00395     pMesh.addPatches(p);
00396 }
00397 
00398 
00399 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00400 
00401 Foam::mirrorFvMesh::~mirrorFvMesh()
00402 {}
00403 
00404 
00405 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines