Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / dynamicMesh / attachDetach / attachInterface.C
blob74fb19211a58f9a689dcac4787468be7cf3e9ccc
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by
13     the Free Software Foundation, either version 3 of the License, or
14     (at your option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "attachDetach.H"
27 #include "polyMesh.H"
28 #include "primitiveMesh.H"
29 #include "polyTopoChange.H"
30 #include "polyTopoChanger.H"
31 #include "polyRemovePoint.H"
32 #include "polyRemoveFace.H"
33 #include "polyModifyFace.H"
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 const Foam::scalar Foam::attachDetach::positionDifference_ = 1e-8;
39 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
41 void Foam::attachDetach::attachInterface
43     polyTopoChange& ref
44 ) const
46     // Algorithm:
47     // 1. Create the reverse patch out of the slave faces.
48     // 2. Go through all the mesh points from the master and slave patch.
49     //    If the point labels are different, insert them into the point
50     //    renumbering list and remove them from the mesh.
51     // 3. Remove all faces from the slave patch
52     // 4. Modify all the faces from the master patch by making them internal
53     //    between the faceCell cells for the two patches.  If the master owner
54     //    is higher than the slave owner, turn the face around
55     // 5. Get all the faces attached to the slave patch points.
56     //    If they have not been removed, renumber them using the
57     //    point renumbering list.
59     if (debug)
60     {
61         Pout<< "void attachDetach::attachInterface("
62             << "polyTopoChange& ref) const "
63             << " for object " << name() << " : "
64             << "Attaching interface" << endl;
65     }
67     const polyMesh& mesh = topoChanger().mesh();
68     const faceList& faces = mesh.faces();
69     const labelList& own = mesh.faceOwner();
70     const labelList& nei = mesh.faceNeighbour();
72     const polyPatch& masterPatch = mesh.boundaryMesh()[masterPatchID_.index()];
73     const polyPatch& slavePatch = mesh.boundaryMesh()[slavePatchID_.index()];
75     const label masterPatchStart = masterPatch.start();
76     const label slavePatchStart = slavePatch.start();
78     const labelList& slaveMeshPoints = slavePatch.meshPoints();
80     const Map<label>& removedPointMap = pointMatchMap();
82     const labelList removedPoints = removedPointMap.toc();
84     forAll(removedPoints, pointI)
85     {
86         ref.setAction(polyRemovePoint(removedPoints[pointI]));
87     }
89 // Pout<< "Points to be mapped: " << removedPoints << endl;
90     // Remove all faces from the slave patch
91     forAll(slavePatch, i)
92     {
93         ref.setAction(polyRemoveFace(i + slavePatchStart));
94 // Pout<< "Removing face " << i + slavePatchStart << endl;
95     }
97     // Modify the faces from the master patch
98     const labelList& masterFaceCells = masterPatch.faceCells();
99     const labelList& slaveFaceCells = slavePatch.faceCells();
101     const boolList& mfFlip = mesh.faceZones()[faceZoneID_.index()].flipMap();
103     forAll(masterFaceCells, faceI)
104     {
105         // If slave neighbour is greater than master, face does not need
106         // turning.  Modify it to become internal
107         if (masterFaceCells[faceI] < slaveFaceCells[faceI])
108         {
109             ref.setAction
110             (
111                 polyModifyFace
112                 (
113                     faces[masterPatchStart + faceI], // modified face
114                     masterPatchStart + faceI,    // label of face being modified
115                     masterFaceCells[faceI],          // owner
116                     slaveFaceCells[faceI],           // neighbour
117                     false,                           // face flip
118                     -1,                              // patch for face
119                     false,                           // remove from zone
120                     faceZoneID_.index(),             // zone for face
121                     mfFlip[faceI]                    // face flip in zone
122                 )
123             );
124         }
125         else
126         {
127             // Flip required
128             ref.setAction
129             (
130                 polyModifyFace
131                 (
132                     faces[masterPatchStart + faceI].reverseFace(), // mod face
133                     masterPatchStart + faceI,    // label of face being modified
134                     slaveFaceCells[faceI],        // owner
135                     masterFaceCells[faceI],       // neighbour
136                     true,                         // face flip
137                     -1,                           // patch for face
138                     false,                        // remove from zone
139                     faceZoneID_.index(),          // zone for face
140                     !mfFlip[faceI]                // face flip in zone
141                 )
142             );
143         }
144     }
146     // Renumber faces affected by point removal
147 // Pout<< "slaveMeshPoints: " << slaveMeshPoints << endl;
148     // Make a map of faces that need to be renumbered
149     labelHashSet facesToModifyMap
150     (
151         slaveMeshPoints.size()*primitiveMesh::facesPerPoint_
152     );
154     const labelListList& pf = mesh.pointFaces();
156     // Grab all the faces off the points in the slave patch.  If the face has
157     //  not been removed, add it to the map of faces to renumber
158     forAll(slaveMeshPoints, pointI)
159     {
160         const labelList& curFaces = pf[slaveMeshPoints[pointI]];
162         forAll(curFaces, faceI)
163         {
164             if (!ref.faceRemoved(curFaces[faceI]))
165             {
166                 facesToModifyMap.insert(curFaces[faceI]);
167             }
168         }
169     }
171     // Grab the faces to be renumbered
172     const labelList ftm = facesToModifyMap.toc();
174     forAll(ftm, faceI)
175     {
176         // For every face to modify, copy the face and re-map the vertices.
177         // It is known all the faces will be changed since they hang off
178         // re-mapped vertices
179         label curFaceID = ftm[faceI];
181         face newFace(faces[curFaceID]);
183         forAll(newFace, pointI)
184         {
185             Map<label>::const_iterator rpmIter =
186                 removedPointMap.find(newFace[pointI]);
188             if (rpmIter != removedPointMap.end())
189             {
190                 // Point mapped. Replace it
191                 newFace[pointI] = rpmIter();
192             }
193         }
195         // Pout<< "face label: " << curFaceID
196         //     << " old face: " << faces[curFaceID]
197         //     << " new face: " << newFace
198         //     << endl;
200         // Get face zone and its flip
201         label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
202         bool modifiedFaceZoneFlip = false;
204         if (modifiedFaceZone >= 0)
205         {
206             modifiedFaceZoneFlip =
207                 mesh.faceZones()[modifiedFaceZone].flipMap()
208                 [
209                     mesh.faceZones()[modifiedFaceZone].whichFace(curFaceID)
210                 ];
211         }
214         label patchID = mesh.boundaryMesh().whichPatch(curFaceID);
215         label neiCell;
216         if (patchID == -1)
217         {
218             neiCell = nei[curFaceID];
219         }
220         else
221         {
222             neiCell = -1;
223         }
226         // Modify the face
227         ref.setAction
228         (
229             polyModifyFace
230             (
231                 newFace,                // modified face
232                 curFaceID,              // label of face being modified
233                 own[curFaceID],         // owner
234                 neiCell,                // neighbour
235                 false,                  // face flip
236                 patchID,                // patch for face
237                 false,                  // remove from zone
238                 modifiedFaceZone,       // zone for face
239                 modifiedFaceZoneFlip    // face flip in zone
240             )
241         );
242     }
244     if (debug)
245     {
246         Pout<< "void attachDetach::attachInterface("
247             << "polyTopoChange& ref) const "
248             << " for object " << name() << " : "
249             << "Finished attaching interface" << endl;
250     }
254 void Foam::attachDetach::modifyMotionPoints
256     pointField& motionPoints
257 ) const
259     const Map<label>& removedPointMap = pointMatchMap();
261     const labelList removedPoints = removedPointMap.toc();
263     if (debug)
264     {
265         Pout<< "void attachDetach::modifyMotionPoints("
266             << "pointField& motionPoints) const "
267             << " for object " << name() << " : "
268             << "Adjusting motion points." << endl;
270         // Calculate the difference in motion point positions
271         scalar pointDiff = 0;
273         forAll(removedPoints, pointI)
274         {
275             pointDiff +=
276                 mag
277                 (
278                     motionPoints[removedPoints[pointI]]
279                   - motionPoints[removedPointMap.find(removedPoints[pointI])()]
280                 );
281         }
283         if (pointDiff > removedPoints.size()*positionDifference_)
284         {
285             Pout<< "Point motion difference = " << pointDiff << endl;
286         }
287     }
289     // Put the slave point on top of the master point
290     forAll(removedPoints, pointI)
291     {
292         motionPoints[removedPoints[pointI]] =
293             motionPoints[removedPointMap.find(removedPoints[pointI])()];
294     }
299 // ************************************************************************* //