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

fvMeshSubset.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 Description
00025     Post-processing mesh subset tool.  Given the original mesh and the
00026     list of selected cells, it creates the mesh consisting only of the
00027     desired cells, with the mapping list for points, faces, and cells.
00028 
00029     MJ 23/03/05 on coupled faces change the patch of the face to the
00030     oldInternalFaces patch.
00031 
00032 \*---------------------------------------------------------------------------*/
00033 
00034 #include "fvMeshSubset.H"
00035 #include <OpenFOAM/boolList.H>
00036 #include <OpenFOAM/Pstream.H>
00037 #include <OpenFOAM/emptyPolyPatch.H>
00038 #include <OpenFOAM/demandDrivenData.H>
00039 #include <OpenFOAM/cyclicPolyPatch.H>
00040 
00041 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00042 
00043 namespace Foam
00044 {
00045 
00046 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00047 
00048 bool Foam::fvMeshSubset::checkCellSubset() const
00049 {
00050     if (fvMeshSubsetPtr_.empty())
00051     {
00052         FatalErrorIn("bool fvMeshSubset::checkCellSubset() const")
00053             << "Mesh subset not set.  Please set the cell map using "
00054             << "void setCellSubset(const labelHashSet& cellsToSubset)" << endl
00055             << "before attempting to access subset data"
00056             << abort(FatalError);
00057 
00058         return false;
00059     }
00060     else
00061     {
00062         return true;
00063     }
00064 }
00065 
00066 
00067 void Foam::fvMeshSubset::markPoints
00068 (
00069     const labelList& curPoints,
00070     Map<label>& pointMap
00071 )
00072 {
00073     forAll (curPoints, pointI)
00074     {
00075         // Note: insert will only insert if not yet there.
00076         pointMap.insert(curPoints[pointI], 0);
00077     }
00078 }
00079 
00080 
00081 void Foam::fvMeshSubset::markPoints
00082 (
00083     const labelList& curPoints,
00084     labelList& pointMap
00085 )
00086 {
00087     forAll (curPoints, pointI)
00088     {
00089         pointMap[curPoints[pointI]] = 0;
00090     }
00091 }
00092 
00093 
00094 // Synchronize nCellsUsingFace on both sides of coupled patches. Marks
00095 // faces that become 'uncoupled' with 3.
00096 void Foam::fvMeshSubset::doCoupledPatches
00097 (
00098     const bool syncPar,
00099     labelList& nCellsUsingFace
00100 ) const
00101 {
00102     const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
00103 
00104     label nUncoupled = 0;
00105 
00106     if (syncPar && Pstream::parRun())
00107     {
00108         // Send face usage across processor patches
00109         forAll (oldPatches, oldPatchI)
00110         {
00111             const polyPatch& pp = oldPatches[oldPatchI];
00112 
00113             if (isA<processorPolyPatch>(pp))
00114             {
00115                 const processorPolyPatch& procPatch =
00116                     refCast<const processorPolyPatch>(pp);
00117 
00118                 OPstream toNeighbour
00119                 (
00120                     Pstream::blocking,
00121                     procPatch.neighbProcNo()
00122                 );
00123 
00124                 toNeighbour
00125                     << SubList<label>(nCellsUsingFace, pp.size(), pp.start());
00126             }
00127         }
00128 
00129 
00130         // Receive face usage count and check for faces that become uncoupled.
00131         forAll (oldPatches, oldPatchI)
00132         {
00133             const polyPatch& pp = oldPatches[oldPatchI];
00134 
00135             if (isA<processorPolyPatch>(pp))
00136             {
00137                 const processorPolyPatch& procPatch =
00138                     refCast<const processorPolyPatch>(pp);
00139 
00140                 IPstream fromNeighbour
00141                 (
00142                     Pstream::blocking,
00143                     procPatch.neighbProcNo()
00144                 );
00145 
00146                 labelList nbrCellsUsingFace(fromNeighbour);
00147 
00148                 // Combine with this side.
00149 
00150                 forAll (pp, i)
00151                 {
00152                     if
00153                     (
00154                         nCellsUsingFace[pp.start()+i] == 1
00155                      && nbrCellsUsingFace[i] == 0
00156                     )
00157                     {
00158                         // Face's neighbour is no longer there. Mark face off
00159                         // as coupled
00160                         nCellsUsingFace[pp.start()+i] = 3;
00161                         nUncoupled++;
00162                     }
00163                 }
00164             }
00165         }
00166     }
00167 
00168     // Do same for cyclics.
00169     forAll (oldPatches, oldPatchI)
00170     {
00171         const polyPatch& pp = oldPatches[oldPatchI];
00172 
00173         if (isA<cyclicPolyPatch>(pp))
00174         {
00175             const cyclicPolyPatch& cycPatch =
00176                 refCast<const cyclicPolyPatch>(pp);
00177 
00178             forAll (cycPatch, i)
00179             {
00180                 label thisFaceI = cycPatch.start() + i;
00181                 label otherFaceI = cycPatch.transformGlobalFace(thisFaceI);
00182 
00183                 if
00184                 (
00185                     nCellsUsingFace[thisFaceI] == 1
00186                  && nCellsUsingFace[otherFaceI] == 0
00187                 )
00188                 {
00189                     nCellsUsingFace[thisFaceI] = 3;
00190                     nUncoupled++;
00191                 }
00192             }
00193         }
00194     }
00195 
00196     if (syncPar)
00197     {
00198         reduce(nUncoupled, sumOp<label>());
00199     }
00200 
00201     if (nUncoupled > 0)
00202     {
00203         Info<< "Uncoupled " << nUncoupled << " faces on coupled patches. "
00204             << "(processorPolyPatch, cyclicPolyPatch)" << nl
00205             << "You might need to run couplePatches to restore the patch face"
00206             << " ordering." << nl << endl;
00207     }
00208 }
00209 
00210 
00211 labelList Foam::fvMeshSubset::subset
00212 (
00213     const label nElems,
00214     const labelList& selectedElements,
00215     const labelList& subsetMap
00216 )
00217 {
00218     // Mark selected elements.
00219     boolList selected(nElems, false);
00220     forAll(selectedElements, i)
00221     {
00222         selected[selectedElements[i]] = true;
00223     }
00224 
00225     // Count subset of selected elements
00226     label n = 0;
00227     forAll(subsetMap, i)
00228     {
00229         if (selected[subsetMap[i]])
00230         {
00231             n++;
00232         }
00233     }
00234 
00235     // Collect selected elements
00236     labelList subsettedElements(n);
00237     n = 0;
00238 
00239     forAll(subsetMap, i)
00240     {
00241         if (selected[subsetMap[i]])
00242         {
00243             subsettedElements[n++] = i;
00244         }
00245     }
00246 
00247     return subsettedElements;
00248 }
00249 
00250 
00251 void Foam::fvMeshSubset::subsetZones()
00252 {
00253     // Keep all zones, even if zero size.
00254 
00255     const pointZoneMesh& pointZones = baseMesh().pointZones();
00256 
00257     // PointZones
00258     List<pointZone*> pZonePtrs(pointZones.size());
00259 
00260     forAll(pointZones, i)
00261     {
00262         const pointZone& pz = pointZones[i];
00263 
00264         pZonePtrs[i] = new pointZone
00265         (
00266             pz.name(),
00267             subset(baseMesh().nPoints(), pz, pointMap()),
00268             i,
00269             fvMeshSubsetPtr_().pointZones()
00270         );
00271     }
00272 
00273 
00274     // FaceZones
00275 
00276     const faceZoneMesh& faceZones = baseMesh().faceZones();
00277 
00278 
00279     // Do we need to remove zones where the side we're interested in
00280     // no longer exists? Guess not.
00281     List<faceZone*> fZonePtrs(faceZones.size());
00282 
00283     forAll(faceZones, i)
00284     {
00285         const faceZone& fz = faceZones[i];
00286 
00287         // Expand faceZone to full mesh
00288         // +1 : part of faceZone, flipped
00289         // -1 :    ,,           , unflipped
00290         //  0 : not part of faceZone
00291         labelList zone(baseMesh().nFaces(), 0);
00292         forAll(fz, j)
00293         {
00294             if (fz.flipMap()[j])
00295             {
00296                 zone[fz[j]] = 1;
00297             }
00298             else
00299             {
00300                 zone[fz[j]] = -1;
00301             }
00302         }
00303 
00304         // Select faces
00305         label nSub = 0;
00306         forAll(faceMap(), j)
00307         {
00308             if (zone[faceMap()[j]] != 0)
00309             {
00310                 nSub++;
00311             }
00312         }
00313         labelList subAddressing(nSub);
00314         boolList subFlipStatus(nSub);
00315         nSub = 0;
00316         forAll(faceMap(), subFaceI)
00317         {
00318             label meshFaceI = faceMap()[subFaceI];
00319             if (zone[meshFaceI] != 0)
00320             {
00321                 subAddressing[nSub] = subFaceI;
00322                 label subOwner = subMesh().faceOwner()[subFaceI];
00323                 label baseOwner = baseMesh().faceOwner()[meshFaceI];
00324                 // If subowner is the same cell as the base keep the flip status
00325                 bool sameOwner = (cellMap()[subOwner] == baseOwner);
00326                 bool flip = (zone[meshFaceI] == 1);
00327                 subFlipStatus[nSub] = (sameOwner == flip);
00328 
00329                 nSub++;
00330             }
00331         }
00332 
00333         fZonePtrs[i] = new faceZone
00334         (
00335             fz.name(),
00336             subAddressing,
00337             subFlipStatus,
00338             i,
00339             fvMeshSubsetPtr_().faceZones()
00340         );
00341     }
00342 
00343 
00344     const cellZoneMesh& cellZones = baseMesh().cellZones();
00345 
00346     List<cellZone*> cZonePtrs(cellZones.size());
00347 
00348     forAll(cellZones, i)
00349     {
00350         const cellZone& cz = cellZones[i];
00351 
00352         cZonePtrs[i] = new cellZone
00353         (
00354             cz.name(),
00355             subset(baseMesh().nCells(), cz, cellMap()),
00356             i,
00357             fvMeshSubsetPtr_().cellZones()
00358         );
00359     }
00360 
00361 
00362     // Add the zones
00363     fvMeshSubsetPtr_().addZones(pZonePtrs, fZonePtrs, cZonePtrs);
00364 }
00365 
00366 
00367 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00368 
00369 // Construct from components
00370 Foam::fvMeshSubset::fvMeshSubset(const fvMesh& baseMesh)
00371 :
00372     baseMesh_(baseMesh),
00373     fvMeshSubsetPtr_(NULL),
00374     pointMap_(0),
00375     faceMap_(0),
00376     cellMap_(0),
00377     patchMap_(0)
00378 {}
00379 
00380 
00381 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00382 
00383 void Foam::fvMeshSubset::setCellSubset
00384 (
00385     const labelHashSet& globalCellMap,
00386     const label patchID
00387 )
00388 {
00389     // Initial check on patches before doing anything time consuming.
00390     const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
00391     const cellList& oldCells = baseMesh().cells();
00392     const faceList& oldFaces = baseMesh().faces();
00393     const pointField& oldPoints = baseMesh().points();
00394     const labelList& oldOwner = baseMesh().faceOwner();
00395     const labelList& oldNeighbour = baseMesh().faceNeighbour();
00396 
00397     label wantedPatchID = patchID;
00398 
00399     if (wantedPatchID == -1)
00400     {
00401         // No explicit patch specified. Put in oldInternalFaces patch.
00402         // Check if patch with this name already exists.
00403         wantedPatchID = oldPatches.findPatchID("oldInternalFaces");
00404     }
00405     else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
00406     {
00407         FatalErrorIn
00408         (
00409             "fvMeshSubset::setCellSubset(const labelHashSet&"
00410             ", const label patchID)"
00411         )   << "Non-existing patch index " << wantedPatchID << endl
00412             << "Should be between 0 and " << oldPatches.size()-1
00413             << abort(FatalError);
00414     }
00415 
00416 
00417     cellMap_ = globalCellMap.toc();
00418 
00419     // Sort the cell map in the ascending order
00420     sort(cellMap_);
00421 
00422     // Approximate sizing parameters for face and point lists
00423     const label avgNFacesPerCell = 6;
00424     const label avgNPointsPerFace = 4;
00425 
00426 
00427     label nCellsInSet = cellMap_.size();
00428 
00429     // Mark all used faces
00430 
00431     Map<label> facesToSubset(avgNFacesPerCell*nCellsInSet);
00432 
00433     forAll (cellMap_, cellI)
00434     {
00435         // Mark all faces from the cell
00436         const labelList& curFaces = oldCells[cellMap_[cellI]];
00437 
00438         forAll (curFaces, faceI)
00439         {
00440             if (!facesToSubset.found(curFaces[faceI]))
00441             {
00442                 facesToSubset.insert(curFaces[faceI], 1);
00443             }
00444             else
00445             {
00446                 facesToSubset[curFaces[faceI]]++;
00447             }
00448         }
00449     }
00450 
00451     // Mark all used points and make a global-to-local face map
00452     Map<label> globalFaceMap(facesToSubset.size());
00453 
00454     // Make a global-to-local point map
00455     Map<label> globalPointMap(avgNPointsPerFace*facesToSubset.size());
00456 
00457     // This is done in two goes, so that the boundary faces are last
00458     // in the list.  Because of this, I need to create the face map
00459     // along the way rather than just grab the table of contents.
00460     labelList facesToc = facesToSubset.toc();
00461     sort(facesToc);
00462     faceMap_.setSize(facesToc.size());
00463 
00464     // 1. Get all faces that will be internal to the submesh.
00465     forAll (facesToc, faceI)
00466     {
00467         if (facesToSubset[facesToc[faceI]] == 2)
00468         {
00469             // Mark face and increment number of points in set
00470             faceMap_[globalFaceMap.size()] = facesToc[faceI];
00471             globalFaceMap.insert(facesToc[faceI], globalFaceMap.size());
00472 
00473             // Mark all points from the face
00474             markPoints(oldFaces[facesToc[faceI]], globalPointMap);
00475         }
00476     }
00477 
00478     // These are all the internal faces in the mesh.
00479     label nInternalFaces = globalFaceMap.size();
00480 
00481 
00482     // Where to insert old internal faces.
00483     label oldPatchStart = labelMax;
00484     if (wantedPatchID != -1)
00485     {
00486         oldPatchStart = oldPatches[wantedPatchID].start();
00487     }
00488 
00489 
00490     label faceI = 0;
00491 
00492     // 2. Boundary faces up to where we want to insert old internal faces
00493     for (; faceI< facesToc.size(); faceI++)
00494     {
00495         if (facesToc[faceI] >= oldPatchStart)
00496         {
00497             break;
00498         }
00499         if
00500         (
00501             !baseMesh().isInternalFace(facesToc[faceI])
00502          && facesToSubset[facesToc[faceI]] == 1
00503         )
00504         {
00505             // Mark face and increment number of points in set
00506             faceMap_[globalFaceMap.size()] = facesToc[faceI];
00507             globalFaceMap.insert(facesToc[faceI], globalFaceMap.size());
00508 
00509             // Mark all points from the face
00510             markPoints(oldFaces[facesToc[faceI]], globalPointMap);
00511         }
00512     }
00513 
00514     // 3. old internal faces
00515     forAll(facesToc, intFaceI)
00516     {
00517         if
00518         (
00519             baseMesh().isInternalFace(facesToc[intFaceI])
00520          && facesToSubset[facesToc[intFaceI]] == 1
00521         )
00522         {
00523             // Mark face and increment number of points in set
00524             faceMap_[globalFaceMap.size()] = facesToc[intFaceI];
00525             globalFaceMap.insert(facesToc[intFaceI], globalFaceMap.size());
00526 
00527             // Mark all points from the face
00528             markPoints(oldFaces[facesToc[intFaceI]], globalPointMap);
00529         }
00530     }
00531 
00532     // 4. Remaining boundary faces
00533     for (; faceI< facesToc.size(); faceI++)
00534     {
00535         if
00536         (
00537             !baseMesh().isInternalFace(facesToc[faceI])
00538          && facesToSubset[facesToc[faceI]] == 1
00539         )
00540         {
00541             // Mark face and increment number of points in set
00542             faceMap_[globalFaceMap.size()] = facesToc[faceI];
00543             globalFaceMap.insert(facesToc[faceI], globalFaceMap.size());
00544 
00545             // Mark all points from the face
00546             markPoints(oldFaces[facesToc[faceI]], globalPointMap);
00547         }
00548     }
00549 
00550 
00551 
00552     // Grab the points map
00553     pointMap_ = globalPointMap.toc();
00554     sort(pointMap_);
00555 
00556     forAll (pointMap_, pointI)
00557     {
00558         globalPointMap[pointMap_[pointI]] = pointI;
00559     }
00560 
00561     Pout << "Number of cells in new mesh: " << nCellsInSet << endl;
00562     Pout << "Number of faces in new mesh: " << globalFaceMap.size() << endl;
00563     Pout << "Number of points in new mesh: " << globalPointMap.size() << endl;
00564 
00565     // Make a new mesh
00566     pointField newPoints(globalPointMap.size());
00567 
00568     label nNewPoints = 0;
00569 
00570     forAll (pointMap_, pointI)
00571     {
00572         newPoints[nNewPoints] = oldPoints[pointMap_[pointI]];
00573         nNewPoints++;
00574     }
00575 
00576     faceList newFaces(globalFaceMap.size());
00577 
00578     label nNewFaces = 0;
00579 
00580     // Make internal faces
00581     for (label faceI = 0; faceI < nInternalFaces; faceI++)
00582     {
00583         const face& oldF = oldFaces[faceMap_[faceI]];
00584 
00585         face newF(oldF.size());
00586 
00587         forAll (newF, i)
00588         {
00589             newF[i] = globalPointMap[oldF[i]];
00590         }
00591 
00592         newFaces[nNewFaces] = newF;
00593         nNewFaces++;
00594     }
00595 
00596     // Make boundary faces
00597 
00598     label nbSize = oldPatches.size();
00599     label oldInternalPatchID  = -1;
00600 
00601     if (wantedPatchID == -1)
00602     {
00603         // Create 'oldInternalFaces' patch at the end
00604         // and put all exposed internal faces in there.
00605         oldInternalPatchID = nbSize;
00606         nbSize++;
00607 
00608     }
00609     else
00610     {
00611         oldInternalPatchID = wantedPatchID;
00612     }
00613 
00614 
00615     // Grad size and start of each patch on the fly.  Because of the
00616     // structure of the underlying mesh, the patches will appear in the
00617     // ascending order
00618     labelList boundaryPatchSizes(nbSize, 0);
00619 
00620     // Assign boundary faces. Visited in order of faceMap_.
00621     for (label faceI = nInternalFaces; faceI < faceMap_.size(); faceI++)
00622     {
00623         label oldFaceI = faceMap_[faceI];
00624 
00625         face oldF = oldFaces[oldFaceI];
00626 
00627         // Turn the faces as necessary to point outwards
00628         if (baseMesh().isInternalFace(oldFaceI))
00629         {
00630             // Internal face. Possibly turned the wrong way round
00631             if
00632             (
00633                 !globalCellMap.found(oldOwner[oldFaceI])
00634              && globalCellMap.found(oldNeighbour[oldFaceI])
00635             )
00636             {
00637                 oldF = oldFaces[oldFaceI].reverseFace();
00638             }
00639 
00640             // Update count for patch
00641             boundaryPatchSizes[oldInternalPatchID]++;
00642         }
00643         else
00644         {
00645             // Boundary face. Increment the appropriate patch
00646             label patchOfFace = oldPatches.whichPatch(oldFaceI);
00647 
00648             // Update count for patch
00649             boundaryPatchSizes[patchOfFace]++;
00650         }
00651 
00652         face newF(oldF.size());
00653 
00654         forAll (newF, i)
00655         {
00656             newF[i] = globalPointMap[oldF[i]];
00657         }
00658 
00659         newFaces[nNewFaces] = newF;
00660         nNewFaces++;
00661     }
00662 
00663 
00664 
00665     // Create cells
00666     cellList newCells(nCellsInSet);
00667 
00668     label nNewCells = 0;
00669 
00670     forAll (cellMap_, cellI)
00671     {
00672         const labelList& oldC = oldCells[cellMap_[cellI]];
00673 
00674         labelList newC(oldC.size());
00675 
00676         forAll (newC, i)
00677         {
00678             newC[i] = globalFaceMap[oldC[i]];
00679         }
00680 
00681         newCells[nNewCells] = cell(newC);
00682         nNewCells++;
00683     }
00684 
00685 
00686     // Delete any old one
00687     fvMeshSubsetPtr_.clear();
00688     // Make a new mesh
00689     fvMeshSubsetPtr_.reset
00690     (
00691         new fvMesh
00692         (
00693             IOobject
00694             (
00695                 baseMesh().name() + "SubSet",
00696                 baseMesh().time().timeName(),
00697                 baseMesh().time(),
00698                 IOobject::NO_READ,
00699                 IOobject::NO_WRITE
00700             ),
00701             xferMove(newPoints),
00702             xferMove(newFaces),
00703             xferMove(newCells)
00704         )
00705     );
00706 
00707 
00708     // Add old patches
00709     List<polyPatch*> newBoundary(nbSize);
00710     patchMap_.setSize(nbSize);
00711     label nNewPatches = 0;
00712     label patchStart = nInternalFaces;
00713 
00714     forAll(oldPatches, patchI)
00715     {
00716         if (boundaryPatchSizes[patchI] > 0)
00717         {
00718             // Patch still exists. Add it
00719             newBoundary[nNewPatches] = oldPatches[patchI].clone
00720             (
00721                 fvMeshSubsetPtr_().boundaryMesh(),
00722                 nNewPatches,
00723                 boundaryPatchSizes[patchI],
00724                 patchStart
00725             ).ptr();
00726 
00727             patchStart += boundaryPatchSizes[patchI];
00728             patchMap_[nNewPatches] = patchI;
00729             nNewPatches++;
00730         }
00731     }
00732 
00733     if (wantedPatchID == -1)
00734     {
00735         // Newly created patch so is at end. Check if any faces in it.
00736         if (boundaryPatchSizes[oldInternalPatchID] > 0)
00737         {
00738             newBoundary[nNewPatches] = new emptyPolyPatch
00739             (
00740                 "oldInternalFaces",
00741                 boundaryPatchSizes[oldInternalPatchID],
00742                 patchStart,
00743                 nNewPatches,
00744                 fvMeshSubsetPtr_().boundaryMesh()
00745             );
00746 
00747             // The index for the first patch is -1 as it originates from
00748             // the internal faces
00749             patchMap_[nNewPatches] = -1;
00750             nNewPatches++;
00751         }
00752     }
00753 
00754     // Reset the patch lists
00755     newBoundary.setSize(nNewPatches);
00756     patchMap_.setSize(nNewPatches);
00757 
00758     // Add the fvPatches
00759     fvMeshSubsetPtr_().addFvPatches(newBoundary);
00760 
00761     // Subset and add any zones
00762     subsetZones();
00763 }
00764 
00765 
00766 void Foam::fvMeshSubset::setLargeCellSubset
00767 (
00768     const labelList& region,
00769     const label currentRegion,
00770     const label patchID,
00771     const bool syncPar
00772 )
00773 {
00774     const cellList& oldCells = baseMesh().cells();
00775     const faceList& oldFaces = baseMesh().faces();
00776     const pointField& oldPoints = baseMesh().points();
00777     const labelList& oldOwner = baseMesh().faceOwner();
00778     const labelList& oldNeighbour = baseMesh().faceNeighbour();
00779     const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
00780     const label oldNInternalFaces = baseMesh().nInternalFaces();
00781 
00782     // Initial checks
00783 
00784     if (region.size() != oldCells.size())
00785     {
00786         FatalErrorIn
00787         (
00788             "fvMeshSubset::setCellSubset(const labelList&"
00789             ", const label, const label, const bool)"
00790         )   << "Size of region " << region.size()
00791             << " is not equal to number of cells in mesh " << oldCells.size()
00792             << abort(FatalError);
00793     }
00794 
00795 
00796     label wantedPatchID = patchID;
00797 
00798     if (wantedPatchID == -1)
00799     {
00800         // No explicit patch specified. Put in oldInternalFaces patch.
00801         // Check if patch with this name already exists.
00802         wantedPatchID = oldPatches.findPatchID("oldInternalFaces");
00803     }
00804     else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
00805     {
00806         FatalErrorIn
00807         (
00808             "fvMeshSubset::setCellSubset(const labelList&"
00809             ", const label, const label, const bool)"
00810         )   << "Non-existing patch index " << wantedPatchID << endl
00811             << "Should be between 0 and " << oldPatches.size()-1
00812             << abort(FatalError);
00813     }
00814 
00815 
00816     // Get the cells for the current region.
00817     cellMap_.setSize(oldCells.size());
00818     label nCellsInSet = 0;
00819 
00820     forAll (region, oldCellI)
00821     {
00822         if (region[oldCellI] == currentRegion)
00823         {
00824             cellMap_[nCellsInSet++] = oldCellI;
00825         }
00826     }
00827     cellMap_.setSize(nCellsInSet);
00828 
00829 
00830     // Mark all used faces. Count number of cells using them
00831     // 0: face not used anymore
00832     // 1: face used by one cell, face becomes/stays boundary face
00833     // 2: face still used and remains internal face
00834     // 3: face coupled and used by one cell only (so should become normal,
00835     //    non-coupled patch face)
00836     //
00837     // Note that this is not really nessecary - but means we can size things
00838     // correctly. Also makes handling coupled faces much easier.
00839 
00840     labelList nCellsUsingFace(oldFaces.size(), 0);
00841 
00842     label nFacesInSet = 0;
00843     forAll (oldFaces, oldFaceI)
00844     {
00845         bool faceUsed = false;
00846 
00847         if (region[oldOwner[oldFaceI]] == currentRegion)
00848         {
00849             nCellsUsingFace[oldFaceI]++;
00850             faceUsed = true;
00851         }
00852 
00853         if
00854         (
00855             baseMesh().isInternalFace(oldFaceI)
00856          && (region[oldNeighbour[oldFaceI]] == currentRegion)
00857         )
00858         {
00859             nCellsUsingFace[oldFaceI]++;
00860             faceUsed = true;
00861         }
00862 
00863         if (faceUsed)
00864         {
00865             nFacesInSet++;
00866         }
00867     }
00868     faceMap_.setSize(nFacesInSet);
00869 
00870     // Handle coupled faces. Modifies patch faces to be uncoupled to 3.
00871     doCoupledPatches(syncPar, nCellsUsingFace);
00872 
00873 
00874     // See which patch to use for exposed internal faces.
00875     label oldInternalPatchID = 0;
00876 
00877     // Insert faces before which patch
00878     label nextPatchID = oldPatches.size();
00879 
00880     // old to new patches
00881     labelList globalPatchMap(oldPatches.size());
00882 
00883     // New patch size
00884     label nbSize = oldPatches.size();
00885 
00886     if (wantedPatchID == -1)
00887     {
00888         // Create 'oldInternalFaces' patch at the end (or before
00889         // processorPatches)
00890         // and put all exposed internal faces in there.
00891 
00892         forAll(oldPatches, patchI)
00893         {
00894             if (isA<processorPolyPatch>(oldPatches[patchI]))
00895             {
00896                 nextPatchID = patchI;
00897                 break;
00898             }
00899             oldInternalPatchID++;
00900         }
00901 
00902         nbSize++;
00903 
00904         // adapt old to new patches for inserted patch
00905         for (label oldPatchI = 0; oldPatchI < nextPatchID; oldPatchI++)
00906         {
00907             globalPatchMap[oldPatchI] = oldPatchI;
00908         }
00909         for
00910         (
00911             label oldPatchI = nextPatchID;
00912             oldPatchI < oldPatches.size();
00913             oldPatchI++
00914         )
00915         {
00916             globalPatchMap[oldPatchI] = oldPatchI+1;
00917         }
00918     }
00919     else
00920     {
00921         oldInternalPatchID = wantedPatchID;
00922         nextPatchID = wantedPatchID+1;
00923 
00924         // old to new patches
00925         globalPatchMap = identity(oldPatches.size());
00926     }
00927 
00928     labelList boundaryPatchSizes(nbSize, 0);
00929 
00930 
00931     // Make a global-to-local point map
00932     labelList globalPointMap(oldPoints.size(), -1);
00933 
00934     labelList globalFaceMap(oldFaces.size(), -1);
00935     label faceI = 0;
00936 
00937     // 1. Pick up all preserved internal faces.
00938     for (label oldFaceI = 0; oldFaceI < oldNInternalFaces; oldFaceI++)
00939     {
00940         if (nCellsUsingFace[oldFaceI] == 2)
00941         {
00942             globalFaceMap[oldFaceI] = faceI;
00943             faceMap_[faceI++] = oldFaceI;
00944 
00945             // Mark all points from the face
00946             markPoints(oldFaces[oldFaceI], globalPointMap);
00947         }
00948     }
00949 
00950     // These are all the internal faces in the mesh.
00951     label nInternalFaces = faceI;
00952 
00953     // 2. Boundary faces up to where we want to insert old internal faces
00954     for
00955     (
00956         label oldPatchI = 0;
00957         oldPatchI < oldPatches.size()
00958      && oldPatchI < nextPatchID;
00959         oldPatchI++
00960     )
00961     {
00962         const polyPatch& oldPatch = oldPatches[oldPatchI];
00963 
00964         label oldFaceI = oldPatch.start();
00965 
00966         forAll (oldPatch, i)
00967         {
00968             if (nCellsUsingFace[oldFaceI] == 1)
00969             {
00970                 // Boundary face is kept.
00971 
00972                 // Mark face and increment number of points in set
00973                 globalFaceMap[oldFaceI] = faceI;
00974                 faceMap_[faceI++] = oldFaceI;
00975 
00976                 // Mark all points from the face
00977                 markPoints(oldFaces[oldFaceI], globalPointMap);
00978 
00979                 // Increment number of patch faces
00980                 boundaryPatchSizes[globalPatchMap[oldPatchI]]++;
00981             }
00982             oldFaceI++;
00983         }
00984     }
00985 
00986     // 3a. old internal faces that have become exposed.
00987     for (label oldFaceI = 0; oldFaceI < oldNInternalFaces; oldFaceI++)
00988     {
00989         if (nCellsUsingFace[oldFaceI] == 1)
00990         {
00991             globalFaceMap[oldFaceI] = faceI;
00992             faceMap_[faceI++] = oldFaceI;
00993 
00994             // Mark all points from the face
00995             markPoints(oldFaces[oldFaceI], globalPointMap);
00996 
00997             // Increment number of patch faces
00998             boundaryPatchSizes[oldInternalPatchID]++;
00999         }
01000     }
01001 
01002     // 3b. coupled patch faces that have become uncoupled.
01003     for
01004     (
01005         label oldFaceI = oldNInternalFaces;
01006         oldFaceI < oldFaces.size();
01007         oldFaceI++
01008     )
01009     {
01010         if (nCellsUsingFace[oldFaceI] == 3)
01011         {
01012             globalFaceMap[oldFaceI] = faceI;
01013             faceMap_[faceI++] = oldFaceI;
01014 
01015             // Mark all points from the face
01016             markPoints(oldFaces[oldFaceI], globalPointMap);
01017 
01018             // Increment number of patch faces
01019             boundaryPatchSizes[oldInternalPatchID]++;
01020         }
01021     }
01022 
01023     // 4. Remaining boundary faces
01024     for
01025     (
01026         label oldPatchI = nextPatchID;
01027         oldPatchI < oldPatches.size();
01028         oldPatchI++
01029     )
01030     {
01031         const polyPatch& oldPatch = oldPatches[oldPatchI];
01032 
01033         label oldFaceI = oldPatch.start();
01034 
01035         forAll (oldPatch, i)
01036         {
01037             if (nCellsUsingFace[oldFaceI] == 1)
01038             {
01039                 // Boundary face is kept.
01040 
01041                 // Mark face and increment number of points in set
01042                 globalFaceMap[oldFaceI] = faceI;
01043                 faceMap_[faceI++] = oldFaceI;
01044 
01045                 // Mark all points from the face
01046                 markPoints(oldFaces[oldFaceI], globalPointMap);
01047 
01048                 // Increment number of patch faces
01049                 boundaryPatchSizes[globalPatchMap[oldPatchI]]++;
01050             }
01051             oldFaceI++;
01052         }
01053     }
01054 
01055     if (faceI != nFacesInSet)
01056     {
01057         FatalErrorIn
01058         (
01059             "fvMeshSubset::setCellSubset(const labelList&"
01060             ", const label, const label, const bool)"
01061         )   << "Problem" << abort(FatalError);
01062     }
01063 
01064 
01065     // Grab the points map
01066     label nPointsInSet = 0;
01067 
01068     forAll(globalPointMap, pointI)
01069     {
01070         if (globalPointMap[pointI] != -1)
01071         {
01072             nPointsInSet++;
01073         }
01074     }
01075     pointMap_.setSize(nPointsInSet);
01076 
01077     nPointsInSet = 0;
01078 
01079     forAll(globalPointMap, pointI)
01080     {
01081         if (globalPointMap[pointI] != -1)
01082         {
01083             pointMap_[nPointsInSet] = pointI;
01084             globalPointMap[pointI] = nPointsInSet;
01085             nPointsInSet++;
01086         }
01087     }
01088 
01089     //Pout<< "Number of cells in new mesh : " << cellMap_.size() << endl;
01090     //Pout<< "Number of faces in new mesh : " << faceMap_.size() << endl;
01091     //Pout<< "Number of points in new mesh: " << pointMap_.size() << endl;
01092 
01093     // Make a new mesh
01094     pointField newPoints(pointMap_.size());
01095 
01096     label nNewPoints = 0;
01097 
01098     forAll (pointMap_, pointI)
01099     {
01100         newPoints[nNewPoints] = oldPoints[pointMap_[pointI]];
01101         nNewPoints++;
01102     }
01103 
01104     faceList newFaces(faceMap_.size());
01105 
01106     label nNewFaces = 0;
01107 
01108     // Make internal faces
01109     for (label faceI = 0; faceI < nInternalFaces; faceI++)
01110     {
01111         const face& oldF = oldFaces[faceMap_[faceI]];
01112 
01113         face newF(oldF.size());
01114 
01115         forAll (newF, i)
01116         {
01117             newF[i] = globalPointMap[oldF[i]];
01118         }
01119 
01120         newFaces[nNewFaces] = newF;
01121         nNewFaces++;
01122     }
01123 
01124 
01125     // Make boundary faces. (different from internal since might need to be
01126     // flipped)
01127     for (label faceI = nInternalFaces; faceI < faceMap_.size(); faceI++)
01128     {
01129         label oldFaceI = faceMap_[faceI];
01130 
01131         face oldF = oldFaces[oldFaceI];
01132 
01133         // Turn the faces as necessary to point outwards
01134         if (baseMesh().isInternalFace(oldFaceI))
01135         {
01136             // Was internal face. Possibly turned the wrong way round
01137             if
01138             (
01139                 region[oldOwner[oldFaceI]] != currentRegion
01140              && region[oldNeighbour[oldFaceI]] == currentRegion
01141             )
01142             {
01143                 oldF = oldFaces[oldFaceI].reverseFace();
01144             }
01145         }
01146 
01147         // Relabel vertices of the (possibly turned) face.
01148         face newF(oldF.size());
01149 
01150         forAll (newF, i)
01151         {
01152             newF[i] = globalPointMap[oldF[i]];
01153         }
01154 
01155         newFaces[nNewFaces] = newF;
01156         nNewFaces++;
01157     }
01158 
01159 
01160 
01161     // Create cells
01162     cellList newCells(nCellsInSet);
01163 
01164     label nNewCells = 0;
01165 
01166     forAll (cellMap_, cellI)
01167     {
01168         const labelList& oldC = oldCells[cellMap_[cellI]];
01169 
01170         labelList newC(oldC.size());
01171 
01172         forAll (newC, i)
01173         {
01174             newC[i] = globalFaceMap[oldC[i]];
01175         }
01176 
01177         newCells[nNewCells] = cell(newC);
01178         nNewCells++;
01179     }
01180 
01181 
01182     // Make a new mesh
01183     // Note that mesh gets registered with same name as original mesh. This is
01184     // not proper but cannot be avoided since otherwise surfaceInterpolation
01185     // cannot find its fvSchemes (it will try to read e.g.
01186     // system/region0SubSet/fvSchemes)
01187     fvMeshSubsetPtr_.reset
01188     (
01189         new fvMesh
01190         (
01191             IOobject
01192             (
01193                 baseMesh().name(),
01194                 baseMesh().time().timeName(),
01195                 baseMesh().time(),
01196                 IOobject::NO_READ,
01197                 IOobject::NO_WRITE
01198             ),
01199             xferMove(newPoints),
01200             xferMove(newFaces),
01201             xferMove(newCells),
01202             syncPar           // parallel synchronisation
01203         )
01204     );
01205 
01206     // Add old patches
01207     List<polyPatch*> newBoundary(nbSize);
01208     patchMap_.setSize(nbSize);
01209     label nNewPatches = 0;
01210     label patchStart = nInternalFaces;
01211 
01212 
01213     // For parallel: only remove patch if none of the processors has it.
01214     // This only gets done for patches before the one being inserted
01215     // (so patches < nextPatchID)
01216 
01217     // Get sum of patch sizes. Zero if patch can be deleted.
01218     labelList globalPatchSizes(boundaryPatchSizes);
01219     globalPatchSizes.setSize(nextPatchID);
01220 
01221     if (syncPar && Pstream::parRun())
01222     {
01223         // Get patch names (up to nextPatchID)
01224         List<wordList> patchNames(Pstream::nProcs());
01225         patchNames[Pstream::myProcNo()] = oldPatches.names();
01226         patchNames[Pstream::myProcNo()].setSize(nextPatchID);
01227         Pstream::gatherList(patchNames);
01228         Pstream::scatterList(patchNames);
01229 
01230         // Get patch sizes (up to nextPatchID).
01231         // Note that up to nextPatchID the globalPatchMap is an identity so
01232         // no need to index through that.
01233         Pstream::listCombineGather(globalPatchSizes, plusEqOp<label>());
01234         Pstream::listCombineScatter(globalPatchSizes);
01235 
01236         // Now all processors have all the patchnames.
01237         // Decide: if all processors have the same patch names and size is zero
01238         // everywhere remove the patch.
01239         bool samePatches = true;
01240 
01241         for (label procI = 1; procI < patchNames.size(); procI++)
01242         {
01243             if (patchNames[procI] != patchNames[0])
01244             {
01245                 samePatches = false;
01246                 break;
01247             }
01248         }
01249 
01250         if (!samePatches)
01251         {
01252             // Patchnames not sync on all processors so disable removal of
01253             // zero sized patches.
01254             globalPatchSizes = labelMax;
01255         }
01256     }
01257 
01258 
01259     // Old patches
01260 
01261     for
01262     (
01263         label oldPatchI = 0;
01264         oldPatchI < oldPatches.size()
01265      && oldPatchI < nextPatchID;
01266         oldPatchI++
01267     )
01268     {
01269         label newSize = boundaryPatchSizes[globalPatchMap[oldPatchI]];
01270 
01271         // Clone (even if 0 size)
01272         newBoundary[nNewPatches] = oldPatches[oldPatchI].clone
01273         (
01274             fvMeshSubsetPtr_().boundaryMesh(),
01275             nNewPatches,
01276             newSize,
01277             patchStart
01278         ).ptr();
01279 
01280         patchStart += newSize;
01281         patchMap_[nNewPatches] = oldPatchI;    // compact patchMap
01282         nNewPatches++;
01283     }
01284 
01285     // Inserted patch
01286 
01287     if (wantedPatchID == -1)
01288     {
01289         label oldInternalSize = boundaryPatchSizes[oldInternalPatchID];
01290 
01291         if (syncPar)
01292         {
01293             reduce(oldInternalSize, sumOp<label>());
01294         }
01295 
01296         // Newly created patch so is at end. Check if any faces in it.
01297         if (oldInternalSize > 0)
01298         {
01299             newBoundary[nNewPatches] = new emptyPolyPatch
01300             (
01301                 "oldInternalFaces",
01302                 boundaryPatchSizes[oldInternalPatchID],
01303                 patchStart,
01304                 nNewPatches,
01305                 fvMeshSubsetPtr_().boundaryMesh()
01306             );
01307 
01308             //Pout<< "    oldInternalFaces : "
01309             //    << boundaryPatchSizes[oldInternalPatchID] << endl;
01310 
01311             // The index for the first patch is -1 as it originates from
01312             // the internal faces
01313             patchStart += boundaryPatchSizes[oldInternalPatchID];
01314             patchMap_[nNewPatches] = -1;
01315             nNewPatches++;
01316         }
01317     }
01318 
01319     // Old patches
01320 
01321     for
01322     (
01323         label oldPatchI = nextPatchID;
01324         oldPatchI < oldPatches.size();
01325         oldPatchI++
01326     )
01327     {
01328         label newSize = boundaryPatchSizes[globalPatchMap[oldPatchI]];
01329 
01330         // Patch still exists. Add it
01331         newBoundary[nNewPatches] = oldPatches[oldPatchI].clone
01332         (
01333             fvMeshSubsetPtr_().boundaryMesh(),
01334             nNewPatches,
01335             newSize,
01336             patchStart
01337         ).ptr();
01338 
01339         //Pout<< "    " << oldPatches[oldPatchI].name() << " : "
01340         //    << newSize << endl;
01341 
01342         patchStart += newSize;
01343         patchMap_[nNewPatches] = oldPatchI;    // compact patchMap
01344         nNewPatches++;
01345     }
01346 
01347 
01348     // Reset the patch lists
01349     newBoundary.setSize(nNewPatches);
01350     patchMap_.setSize(nNewPatches);
01351 
01352 
01353     // Add the fvPatches
01354     fvMeshSubsetPtr_().addFvPatches(newBoundary, syncPar);
01355 
01356     // Subset and add any zones
01357     subsetZones();
01358 }
01359 
01360 
01361 void Foam::fvMeshSubset::setLargeCellSubset
01362 (
01363     const labelHashSet& globalCellMap,
01364     const label patchID,
01365     const bool syncPar
01366 )
01367 {
01368     labelList region(baseMesh().nCells(), 0);
01369 
01370     forAllConstIter (labelHashSet, globalCellMap, iter)
01371     {
01372         region[iter.key()] = 1;
01373     }
01374     setLargeCellSubset(region, 1, patchID, syncPar);
01375 }
01376 
01377 
01378 const fvMesh& Foam::fvMeshSubset::subMesh() const
01379 {
01380     checkCellSubset();
01381 
01382     return fvMeshSubsetPtr_();
01383 }
01384 
01385 
01386 fvMesh& Foam::fvMeshSubset::subMesh()
01387 {
01388     checkCellSubset();
01389 
01390     return fvMeshSubsetPtr_();
01391 }
01392 
01393 
01394 const labelList& Foam::fvMeshSubset::pointMap() const
01395 {
01396     checkCellSubset();
01397 
01398     return pointMap_;
01399 }
01400 
01401 
01402 const labelList& Foam::fvMeshSubset::faceMap() const
01403 {
01404     checkCellSubset();
01405 
01406     return faceMap_;
01407 }
01408 
01409 
01410 const labelList& Foam::fvMeshSubset::cellMap() const
01411 {
01412     checkCellSubset();
01413 
01414     return cellMap_;
01415 }
01416 
01417 
01418 const labelList& Foam::fvMeshSubset::patchMap() const
01419 {
01420     checkCellSubset();
01421 
01422     return patchMap_;
01423 }
01424 
01425 
01426 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
01427 
01428 } // End namespace Foam
01429 
01430 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines