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

attachInterface.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 "attachDetach.H"
00027 #include <OpenFOAM/polyMesh.H>
00028 #include <OpenFOAM/primitiveMesh.H>
00029 #include <dynamicMesh/polyTopoChange.H>
00030 #include <dynamicMesh/polyTopoChanger.H>
00031 #include <dynamicMesh/polyRemovePoint.H>
00032 #include <dynamicMesh/polyRemoveFace.H>
00033 #include <dynamicMesh/polyModifyFace.H>
00034 
00035 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00036 
00037 const Foam::scalar Foam::attachDetach::positionDifference_ = 1e-8;
00038 
00039 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00040 
00041 void Foam::attachDetach::attachInterface
00042 (
00043     polyTopoChange& ref
00044 ) const
00045 {
00046     // Algorithm:
00047     // 1. Create the reverse patch out of the slave faces.
00048     // 2. Go through all the mesh points from the master and slave patch.
00049     //    If the point labels are different, insert them into the point
00050     //    renumbering list and remove them from the mesh.
00051     // 3. Remove all faces from the slave patch
00052     // 4. Modify all the faces from the master patch by making them internal
00053     //    between the faceCell cells for the two patches.  If the master owner
00054     //    is higher than the slave owner, turn the face around
00055     // 5. Get all the faces attached to the slave patch points.
00056     //    If they have not been removed, renumber them using the
00057     //    point renumbering list.
00058 
00059     if (debug)
00060     {
00061         Pout<< "void attachDetach::attachInterface("
00062             << "polyTopoChange& ref) const "
00063             << " for object " << name() << " : "
00064             << "Attaching interface" << endl;
00065     }
00066 
00067     const polyMesh& mesh = topoChanger().mesh();
00068     const faceList& faces = mesh.faces();
00069     const labelList& own = mesh.faceOwner();
00070     const labelList& nei = mesh.faceNeighbour();
00071 
00072     const polyPatch& masterPatch = mesh.boundaryMesh()[masterPatchID_.index()];
00073     const polyPatch& slavePatch = mesh.boundaryMesh()[slavePatchID_.index()];
00074 
00075     const label masterPatchStart = masterPatch.start();
00076     const label slavePatchStart = slavePatch.start();
00077 
00078     const labelList& slaveMeshPoints = slavePatch.meshPoints();
00079 
00080     const Map<label>& removedPointMap = pointMatchMap();
00081 
00082     const labelList removedPoints = removedPointMap.toc();
00083 
00084     forAll (removedPoints, pointI)
00085     {
00086         ref.setAction(polyRemovePoint(removedPoints[pointI]));
00087     }
00088 
00089 // Pout << "Points to be mapped: " << removedPoints << endl;
00090     // Remove all faces from the slave patch
00091     for (label i = 0; i < slavePatch.size(); i++)
00092     {
00093         ref.setAction(polyRemoveFace(i + slavePatchStart));
00094 // Pout << "Removing face " << i + slavePatchStart << endl;
00095     }
00096 
00097     // Modify the faces from the master patch
00098     const labelList& masterFaceCells = masterPatch.faceCells();
00099     const labelList& slaveFaceCells = slavePatch.faceCells();
00100 
00101     const boolList& mfFlip = mesh.faceZones()[faceZoneID_.index()].flipMap();
00102 
00103     forAll (masterFaceCells, faceI)
00104     {
00105         // If slave neighbour is greater than master, face does not need
00106         // turning.  Modify it to become internal
00107         if (masterFaceCells[faceI] < slaveFaceCells[faceI])
00108         {
00109             ref.setAction
00110             (
00111                 polyModifyFace
00112                 (
00113                     faces[masterPatchStart + faceI], // modified face
00114                     masterPatchStart + faceI,    // label of face being modified
00115                     masterFaceCells[faceI],          // owner
00116                     slaveFaceCells[faceI],           // neighbour
00117                     false,                           // face flip
00118                     -1,                              // patch for face
00119                     false,                           // remove from zone
00120                     faceZoneID_.index(),             // zone for face
00121                     mfFlip[faceI]                    // face flip in zone
00122                 )
00123             );
00124         }
00125         else
00126         {
00127             // Flip required
00128             ref.setAction
00129             (
00130                 polyModifyFace
00131                 (
00132                     faces[masterPatchStart + faceI].reverseFace(), // mod face
00133                     masterPatchStart + faceI,    // label of face being modified
00134                     slaveFaceCells[faceI],        // owner
00135                     masterFaceCells[faceI],       // neighbour
00136                     true,                         // face flip
00137                     -1,                           // patch for face
00138                     false,                        // remove from zone
00139                     faceZoneID_.index(),          // zone for face
00140                     !mfFlip[faceI]                // face flip in zone
00141                 )
00142             );
00143         }
00144     }
00145 
00146     // Renumber faces affected by point removal
00147 // Pout << "slaveMeshPoints: " << slaveMeshPoints << endl;
00148     // Make a map of faces that need to be renumbered
00149     labelHashSet facesToModifyMap
00150     (
00151         slaveMeshPoints.size()*primitiveMesh::facesPerPoint_
00152     );
00153 
00154     const labelListList& pf = mesh.pointFaces();
00155 
00156     // Grab all the faces off the points in the slave patch.  If the face has
00157     //  not been removed, add it to the map of faces to renumber
00158     forAll (slaveMeshPoints, pointI)
00159     {
00160         const labelList& curFaces = pf[slaveMeshPoints[pointI]];
00161 
00162         forAll (curFaces, faceI)
00163         {
00164             if (!ref.faceRemoved(curFaces[faceI]))
00165             {
00166                 facesToModifyMap.insert(curFaces[faceI]);
00167             }
00168         }
00169     }
00170 
00171     // Grab the faces to be renumbered
00172     const labelList ftm = facesToModifyMap.toc();
00173 
00174     forAll (ftm, faceI)
00175     {
00176         // For every face to modify, copy the face and re-map the vertices.
00177         // It is known all the faces will be changed since they hang off
00178         // re-mapped vertices
00179         label curFaceID = ftm[faceI];
00180 
00181         face newFace(faces[curFaceID]);
00182 
00183         forAll (newFace, pointI)
00184         {
00185             Map<label>::const_iterator rpmIter =
00186                 removedPointMap.find(newFace[pointI]);
00187 
00188             if (rpmIter != removedPointMap.end())
00189             {
00190                 // Point mapped. Replace it
00191                 newFace[pointI] = rpmIter();
00192             }
00193         }
00194 
00195 // Pout<< "face label: " << curFaceID << " old face: " << faces[curFaceID] << " new face: " << newFace << endl;
00196 
00197         // Get face zone and its flip
00198         label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
00199         bool modifiedFaceZoneFlip = false;
00200 
00201         if (modifiedFaceZone >= 0)
00202         {
00203             modifiedFaceZoneFlip =
00204                 mesh.faceZones()[modifiedFaceZone].flipMap()
00205                 [
00206                     mesh.faceZones()[modifiedFaceZone].whichFace(curFaceID)
00207                 ];
00208         }
00209 
00210 
00211         label patchID = mesh.boundaryMesh().whichPatch(curFaceID);
00212         label neiCell;
00213         if (patchID == -1)
00214         {
00215             neiCell = nei[curFaceID];
00216         }
00217         else
00218         {
00219             neiCell = -1;
00220         }
00221 
00222 
00223         // Modify the face
00224         ref.setAction
00225         (
00226             polyModifyFace
00227             (
00228                 newFace,                // modified face
00229                 curFaceID,              // label of face being modified
00230                 own[curFaceID],         // owner
00231                 neiCell,                // neighbour
00232                 false,                  // face flip
00233                 patchID,                // patch for face
00234                 false,                  // remove from zone
00235                 modifiedFaceZone,       // zone for face
00236                 modifiedFaceZoneFlip    // face flip in zone
00237             )
00238         );
00239     }
00240 
00241     if (debug)
00242     {
00243         Pout<< "void attachDetach::attachInterface("
00244             << "polyTopoChange& ref) const "
00245             << " for object " << name() << " : "
00246             << "Finished attaching interface" << endl;
00247     }
00248 }
00249 
00250 
00251 void Foam::attachDetach::modifyMotionPoints
00252 (
00253     pointField& motionPoints
00254 ) const
00255 {
00256     const Map<label>& removedPointMap = pointMatchMap();
00257 
00258     const labelList removedPoints = removedPointMap.toc();
00259 
00260     if (debug)
00261     {
00262         Pout<< "void attachDetach::modifyMotionPoints(" 
00263             << "pointField& motionPoints) const "
00264             << " for object " << name() << " : "
00265             << "Adjusting motion points." << endl;
00266 
00267         // Calculate the difference in motion point positions
00268         scalar pointDiff = 0;
00269 
00270         forAll (removedPoints, pointI)
00271         {
00272             pointDiff +=
00273                 mag
00274                 (
00275                     motionPoints[removedPoints[pointI]]
00276                   - motionPoints[removedPointMap.find(removedPoints[pointI])()]
00277                 );
00278         }
00279 
00280         if (pointDiff > removedPoints.size()*positionDifference_)
00281         {
00282             Pout<< "Point motion difference = " << pointDiff << endl;
00283         }
00284     }
00285 
00286     // Put the slave point on top of the master point
00287     forAll (removedPoints, pointI)
00288     {
00289         motionPoints[removedPoints[pointI]] =
00290             motionPoints[removedPointMap.find(removedPoints[pointI])()];
00291     }
00292 
00293 }
00294 
00295 
00296 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines