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

cellSplitter.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 
00026 \*---------------------------------------------------------------------------*/
00027 
00028 #include "cellSplitter.H"
00029 #include <OpenFOAM/polyMesh.H>
00030 #include <dynamicMesh/polyTopoChange.H>
00031 #include <dynamicMesh/polyAddCell.H>
00032 #include <dynamicMesh/polyAddFace.H>
00033 #include <dynamicMesh/polyAddPoint.H>
00034 #include <dynamicMesh/polyModifyFace.H>
00035 #include <OpenFOAM/mapPolyMesh.H>
00036 #include <meshTools/meshTools.H>
00037 
00038 
00039 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00040 
00041 namespace Foam
00042 {
00043 
00044 defineTypeNameAndDebug(cellSplitter, 0);
00045 
00046 }
00047 
00048 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00049 
00050 void Foam::cellSplitter::getFaceInfo
00051 (
00052     const label faceI,
00053     label& patchID,
00054     label& zoneID,
00055     label& zoneFlip
00056 ) const
00057 {
00058     patchID = -1;
00059 
00060     if (!mesh_.isInternalFace(faceI))
00061     {
00062         patchID = mesh_.boundaryMesh().whichPatch(faceI);
00063     }
00064 
00065     zoneID = mesh_.faceZones().whichZone(faceI);
00066 
00067     zoneFlip = false;
00068 
00069     if (zoneID >= 0)
00070     {
00071         const faceZone& fZone = mesh_.faceZones()[zoneID];
00072 
00073         zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
00074     }
00075 }
00076 
00077 
00078 // Find the new owner of faceI (since the original cell has been split into
00079 // newCells
00080 Foam::label Foam::cellSplitter::newOwner
00081 (
00082     const label faceI,
00083     const Map<labelList>& cellToCells
00084 ) const
00085 {
00086     label oldOwn = mesh_.faceOwner()[faceI];
00087 
00088     Map<labelList>::const_iterator fnd = cellToCells.find(oldOwn);
00089 
00090     if (fnd == cellToCells.end())
00091     {
00092         // Unsplit cell
00093         return oldOwn;
00094     }
00095     else
00096     {
00097         // Look up index of face in the cells' faces.
00098 
00099         const labelList& newCells = fnd();
00100 
00101         const cell& cFaces = mesh_.cells()[oldOwn];
00102 
00103         return newCells[findIndex(cFaces, faceI)];
00104     }
00105 }
00106 
00107 
00108 Foam::label Foam::cellSplitter::newNeighbour
00109 (
00110     const label faceI,
00111     const Map<labelList>& cellToCells
00112 ) const
00113 {
00114     label oldNbr = mesh_.faceNeighbour()[faceI];
00115 
00116     Map<labelList>::const_iterator fnd = cellToCells.find(oldNbr);
00117 
00118     if (fnd == cellToCells.end())
00119     {
00120         // Unsplit cell
00121         return oldNbr;
00122     }
00123     else
00124     {
00125         // Look up index of face in the cells' faces.
00126 
00127         const labelList& newCells = fnd();
00128 
00129         const cell& cFaces = mesh_.cells()[oldNbr];
00130 
00131         return newCells[findIndex(cFaces, faceI)];
00132     }
00133 }
00134 
00135 
00136 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00137 
00138 // Construct from components
00139 Foam::cellSplitter::cellSplitter(const polyMesh& mesh)
00140 :
00141     mesh_(mesh),
00142     addedPoints_()
00143 {}
00144 
00145 
00146 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00147 
00148 Foam::cellSplitter::~cellSplitter()
00149 {}
00150 
00151 
00152 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00153 
00154 void Foam::cellSplitter::setRefinement
00155 (
00156     const Map<point>& cellToMidPoint,
00157     polyTopoChange& meshMod
00158 )
00159 {
00160     addedPoints_.clear();
00161     addedPoints_.resize(cellToMidPoint.size());
00162 
00163 
00164     //
00165     // Introduce cellToMidPoints.
00166     //
00167 
00168     forAllConstIter(Map<point>, cellToMidPoint, iter)
00169     {
00170         label cellI = iter.key();
00171 
00172         label anchorPoint = mesh_.cellPoints()[cellI][0];
00173 
00174         label addedPointI =
00175             meshMod.setAction
00176             (
00177                 polyAddPoint
00178                 (
00179                     iter(),         // point
00180                     anchorPoint,    // master point
00181                     -1,             // zone for point
00182                     true            // supports a cell
00183                 )
00184             );
00185         addedPoints_.insert(cellI, addedPointI);
00186 
00187         //Pout<< "Added point " << addedPointI
00188         //    << iter() << " in cell " << cellI << " with centre "
00189         //    << mesh_.cellCentres()[cellI] << endl;
00190     }
00191 
00192 
00193     //
00194     // Add cells (first one is modified original cell)
00195     //
00196 
00197     Map<labelList> cellToCells(cellToMidPoint.size());
00198 
00199     forAllConstIter(Map<point>, cellToMidPoint, iter)
00200     {
00201         label cellI = iter.key();
00202 
00203         const cell& cFaces = mesh_.cells()[cellI];
00204 
00205         // Cells created for this cell.
00206         labelList newCells(cFaces.size());
00207 
00208         // First pyramid is the original cell
00209         newCells[0] = cellI;
00210 
00211         // Add other pyramids
00212         for (label i = 1; i < cFaces.size(); i++)
00213         {
00214             label addedCellI =
00215                 meshMod.setAction
00216                 (
00217                     polyAddCell
00218                     (
00219                         -1,     // master point
00220                         -1,     // master edge
00221                         -1,     // master face
00222                         cellI,  // master cell
00223                         -1      // zone
00224                     )
00225                 );
00226 
00227             newCells[i] = addedCellI;
00228         }
00229 
00230         cellToCells.insert(cellI, newCells);
00231 
00232         //Pout<< "Split cell " << cellI
00233         //    << " with centre " << mesh_.cellCentres()[cellI] << nl
00234         //    << " faces:" << cFaces << nl
00235         //    << " into :" << newCells << endl;
00236     }
00237 
00238 
00239     //
00240     // Introduce internal faces. These go from edges of the cell to the mid
00241     // point.
00242     //
00243 
00244     forAllConstIter(Map<point>, cellToMidPoint, iter)
00245     {
00246         label cellI = iter.key();
00247 
00248         label midPointI = addedPoints_[cellI];
00249 
00250         const cell& cFaces = mesh_.cells()[cellI];
00251 
00252         const labelList& cEdges = mesh_.cellEdges()[cellI];
00253 
00254         forAll(cEdges, i)
00255         {
00256             label edgeI = cEdges[i];
00257             const edge& e = mesh_.edges()[edgeI];
00258 
00259             // Get the faces on the cell using the edge
00260             label face0, face1;
00261             meshTools::getEdgeFaces(mesh_, cellI, edgeI, face0, face1);
00262 
00263             // Get the cells on both sides of the face by indexing into cFaces.
00264             // (since newly created cells are stored in cFaces order)
00265             const labelList& newCells = cellToCells[cellI];
00266 
00267             label cell0 = newCells[findIndex(cFaces, face0)];
00268             label cell1 = newCells[findIndex(cFaces, face1)];
00269 
00270             if (cell0 < cell1)
00271             {
00272                 // Construct face to midpoint that is pointing away from
00273                 // (pyramid split off from) cellI
00274 
00275                 const face& f0 = mesh_.faces()[face0];
00276 
00277                 label index = findIndex(f0, e[0]);
00278 
00279                 bool edgeInFaceOrder = (f0[f0.fcIndex(index)] == e[1]);
00280 
00281                 // Check if cellI is the face owner
00282 
00283                 face newF(3);
00284                 if (edgeInFaceOrder == (mesh_.faceOwner()[face0] == cellI))
00285                 {
00286                     // edge used in face order.
00287                     newF[0] = e[1];
00288                     newF[1] = e[0];
00289                     newF[2] = midPointI;
00290                 }
00291                 else
00292                 {
00293                     newF[0] = e[0];
00294                     newF[1] = e[1];
00295                     newF[2] = midPointI;
00296                 }
00297 
00298                 // Now newF points away from cell0
00299                 meshMod.setAction
00300                 (
00301                     polyAddFace
00302                     (
00303                         newF,                       // face
00304                         cell0,                      // owner
00305                         cell1,                      // neighbour
00306                         -1,                         // master point
00307                         -1,                         // master edge
00308                         face0,                      // master face for addition
00309                         false,                      // flux flip
00310                         -1,                         // patch for face
00311                         -1,                         // zone for face
00312                         false                       // face zone flip
00313                     )
00314                 );
00315             }
00316             else
00317             {
00318                 // Construct face to midpoint that is pointing away from
00319                 // (pyramid split off from) cellI
00320 
00321                 const face& f1 = mesh_.faces()[face1];
00322 
00323                 label index = findIndex(f1, e[0]);
00324 
00325                 bool edgeInFaceOrder = (f1[f1.fcIndex(index)] == e[1]);
00326 
00327                 // Check if cellI is the face owner
00328 
00329                 face newF(3);
00330                 if (edgeInFaceOrder == (mesh_.faceOwner()[face1] == cellI))
00331                 {
00332                     // edge used in face order.
00333                     newF[0] = e[1];
00334                     newF[1] = e[0];
00335                     newF[2] = midPointI;
00336                 }
00337                 else
00338                 {
00339                     newF[0] = e[0];
00340                     newF[1] = e[1];
00341                     newF[2] = midPointI;
00342                 }
00343 
00344                 // Now newF points away from cell1
00345                 meshMod.setAction
00346                 (
00347                     polyAddFace
00348                     (
00349                         newF,                       // face
00350                         cell1,                      // owner
00351                         cell0,                      // neighbour
00352                         -1,                         // master point
00353                         -1,                         // master edge
00354                         face0,                      // master face for addition
00355                         false,                      // flux flip
00356                         -1,                         // patch for face
00357                         -1,                         // zone for face
00358                         false                       // face zone flip
00359                     )
00360                 );
00361             }
00362         }
00363     }
00364 
00365 
00366     //
00367     // Update all existing faces for split owner or neighbour.
00368     //
00369 
00370 
00371     // Mark off affected face.
00372     boolList faceUpToDate(mesh_.nFaces(), true);
00373 
00374     forAllConstIter(Map<point>, cellToMidPoint, iter)
00375     {
00376         label cellI = iter.key();
00377 
00378         const cell& cFaces = mesh_.cells()[cellI];
00379 
00380         forAll(cFaces, i)
00381         {
00382             label faceI = cFaces[i];
00383 
00384             faceUpToDate[faceI] = false;
00385         }
00386     }
00387 
00388     forAll(faceUpToDate, faceI)
00389     {
00390         if (!faceUpToDate[faceI])
00391         {
00392             const face& f = mesh_.faces()[faceI];
00393 
00394             if (mesh_.isInternalFace(faceI))
00395             {
00396                 label newOwn = newOwner(faceI, cellToCells);
00397                 label newNbr = newNeighbour(faceI, cellToCells);
00398 
00399                 if (newOwn < newNbr)
00400                 {
00401                     meshMod.setAction
00402                     (
00403                         polyModifyFace
00404                         (
00405                             f,
00406                             faceI,
00407                             newOwn,         // owner
00408                             newNbr,         // neighbour
00409                             false,          // flux flip
00410                             -1,             // patch for face
00411                             false,          // remove from zone
00412                             -1,             // zone for face
00413                             false           // face zone flip
00414                         )
00415                     );
00416                 }
00417                 else
00418                 {
00419                     meshMod.setAction
00420                     (
00421                         polyModifyFace
00422                         (
00423                             f.reverseFace(),
00424                             faceI,
00425                             newNbr,         // owner
00426                             newOwn,         // neighbour
00427                             false,          // flux flip
00428                             -1,             // patch for face
00429                             false,          // remove from zone
00430                             -1,             // zone for face
00431                             false           // face zone flip
00432                         )
00433                     );
00434                 }
00435 
00436             }
00437             else
00438             {
00439                 label newOwn = newOwner(faceI, cellToCells);
00440 
00441                 label patchID, zoneID, zoneFlip;
00442                 getFaceInfo(faceI, patchID, zoneID, zoneFlip);
00443 
00444                 meshMod.setAction
00445                 (
00446                     polyModifyFace
00447                     (
00448                         mesh_.faces()[faceI],
00449                         faceI,
00450                         newOwn,         // owner
00451                         -1,             // neighbour
00452                         false,          // flux flip
00453                         patchID,        // patch for face
00454                         false,          // remove from zone
00455                         zoneID,         // zone for face
00456                         zoneFlip        // face zone flip
00457                     )
00458                 );
00459             }
00460 
00461             faceUpToDate[faceI] = true;
00462         }
00463     }
00464 }
00465 
00466 
00467 void Foam::cellSplitter::updateMesh(const mapPolyMesh& morphMap)
00468 {
00469     // Create copy since we're deleting entries. Only if both cell and added
00470     // point get mapped do they get inserted.
00471     Map<label> newAddedPoints(addedPoints_.size());
00472 
00473     forAllConstIter(Map<label>, addedPoints_, iter)
00474     {
00475         label oldCellI = iter.key();
00476 
00477         label newCellI = morphMap.reverseCellMap()[oldCellI];
00478 
00479         label oldPointI = iter();
00480 
00481         label newPointI = morphMap.reversePointMap()[oldPointI];
00482 
00483         if (newCellI >= 0 && newPointI >= 0)
00484         {
00485             newAddedPoints.insert(newCellI, newPointI);
00486         }
00487     }
00488 
00489     // Copy
00490     addedPoints_.transfer(newAddedPoints);
00491 }
00492 
00493 
00494 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines