Merge branch 'master' of ssh://opencfd@repo.or.cz/srv/git/OpenFOAM-1.5.x
[OpenFOAM-1.5.x.git] / src / dynamicMesh / attachDetach / attachInterface.C
blob714ee40567eacf607dc3882c330e84adc63775a9
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "attachDetach.H"
28 #include "polyMesh.H"
29 #include "primitiveMesh.H"
30 #include "polyTopoChange.H"
31 #include "polyTopoChanger.H"
32 #include "polyRemovePoint.H"
33 #include "polyRemoveFace.H"
34 #include "polyModifyFace.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 const Foam::scalar Foam::attachDetach::positionDifference_ = 1e-8;
40 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
42 void Foam::attachDetach::attachInterface
44     polyTopoChange& ref
45 ) const
47     // Algorithm:
48     // 1. Create the reverse patch out of the slave faces.
49     // 2. Go through all the mesh points from the master and slave patch.
50     //    If the point labels are different, insert them into the point
51     //    renumbering list and remove them from the mesh.
52     // 3. Remove all faces from the slave patch
53     // 4. Modify all the faces from the master patch by making them internal
54     //    between the faceCell cells for the two patches.  If the master owner
55     //    is higher than the slave owner, turn the face around
56     // 5. Get all the faces attached to the slave patch points.
57     //    If they have not been removed, renumber them using the
58     //    point renumbering list.
60     if (debug)
61     {
62         Pout<< "void attachDetach::attachInterface("
63             << "polyTopoChange& ref) const "
64             << " for object " << name() << " : "
65             << "Attaching interface" << endl;
66     }
68     const polyMesh& mesh = topoChanger().mesh();
69     const faceList& faces = mesh.faces();
70     const labelList& own = mesh.faceOwner();
71     const labelList& nei = mesh.faceNeighbour();
73     const polyPatch& masterPatch = mesh.boundaryMesh()[masterPatchID_.index()];
74     const polyPatch& slavePatch = mesh.boundaryMesh()[slavePatchID_.index()];
76     const label masterPatchStart = masterPatch.start();
77     const label slavePatchStart = slavePatch.start();
79     const labelList& slaveMeshPoints = slavePatch.meshPoints();
81     const Map<label>& removedPointMap = pointMatchMap();
83     const labelList removedPoints = removedPointMap.toc();
85     forAll (removedPoints, pointI)
86     {
87         ref.setAction(polyRemovePoint(removedPoints[pointI]));
88     }
90 // Pout << "Points to be mapped: " << removedPoints << endl;
91     // Remove all faces from the slave patch
92     for (label i = 0; i < slavePatch.size(); i++)
93     {
94         ref.setAction(polyRemoveFace(i + slavePatchStart));
95 // Pout << "Removing face " << i + slavePatchStart << endl;
96     }
98     // Modify the faces from the master patch
99     const labelList& masterFaceCells = masterPatch.faceCells();
100     const labelList& slaveFaceCells = slavePatch.faceCells();
102     const boolList& mfFlip = mesh.faceZones()[faceZoneID_.index()].flipMap();
104     forAll (masterFaceCells, faceI)
105     {
106         // If slave neighbour is greater than master, face does not need
107         // turning.  Modify it to become internal
108         if (masterFaceCells[faceI] < slaveFaceCells[faceI])
109         {
110             ref.setAction
111             (
112                 polyModifyFace
113                 (
114                     faces[masterPatchStart + faceI], // modified face
115                     masterPatchStart + faceI,    // label of face being modified
116                     masterFaceCells[faceI],          // owner
117                     slaveFaceCells[faceI],           // neighbour
118                     false,                           // face flip
119                     -1,                              // patch for face
120                     false,                           // remove from zone
121                     faceZoneID_.index(),             // zone for face
122                     mfFlip[faceI]                    // face flip in zone
123                 )
124             );
125         }
126         else
127         {
128             // Flip required
129             ref.setAction
130             (
131                 polyModifyFace
132                 (
133                     faces[masterPatchStart + faceI].reverseFace(), // mod face
134                     masterPatchStart + faceI,    // label of face being modified
135                     slaveFaceCells[faceI],        // owner
136                     masterFaceCells[faceI],       // neighbour
137                     true,                         // face flip
138                     -1,                           // patch for face
139                     false,                        // remove from zone
140                     faceZoneID_.index(),          // zone for face
141                     !mfFlip[faceI]                // face flip in zone
142                 )
143             );
144         }
145     }
147     // Renumber faces affected by point removal
148 // Pout << "slaveMeshPoints: " << slaveMeshPoints << endl;
149     // Make a map of faces that need to be renumbered
150     labelHashSet facesToModifyMap
151     (
152         slaveMeshPoints.size()*primitiveMesh::facesPerPoint_
153     );
155     const labelListList& pf = mesh.pointFaces();
157     // Grab all the faces off the points in the slave patch.  If the face has
158     //  not been removed, add it to the map of faces to renumber
159     forAll (slaveMeshPoints, pointI)
160     {
161         const labelList& curFaces = pf[slaveMeshPoints[pointI]];
163         forAll (curFaces, faceI)
164         {
165             if (!ref.faceRemoved(curFaces[faceI]))
166             {
167                 facesToModifyMap.insert(curFaces[faceI]);
168             }
169         }
170     }
172     // Grab the faces to be renumbered
173     const labelList ftm = facesToModifyMap.toc();
175     forAll (ftm, faceI)
176     {
177         // For every face to modify, copy the face and re-map the vertices.
178         // It is known all the faces will be changed since they hang off
179         // re-mapped vertices
180         label curFaceID = ftm[faceI];
182         face newFace(faces[curFaceID]);
184         forAll (newFace, pointI)
185         {
186             Map<label>::const_iterator rpmIter =
187                 removedPointMap.find(newFace[pointI]);
189             if (rpmIter != removedPointMap.end())
190             {
191                 // Point mapped. Replace it
192                 newFace[pointI] = rpmIter();
193             }
194         }
196 // Pout<< "face label: " << curFaceID << " old face: " << faces[curFaceID] << " new face: " << newFace << endl;
198         // Get face zone and its flip
199         label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
200         bool modifiedFaceZoneFlip = false;
202         if (modifiedFaceZone >= 0)
203         {
204             modifiedFaceZoneFlip =
205                 mesh.faceZones()[modifiedFaceZone].flipMap()
206                 [
207                     mesh.faceZones()[modifiedFaceZone].whichFace(curFaceID)
208                 ];
209         }
210             
211         // Modify the face
212         ref.setAction
213         (
214             polyModifyFace
215             (
216                 newFace,                // modified face
217                 curFaceID,              // label of face being modified
218                 own[curFaceID],         // owner
219                 nei[curFaceID],         // neighbour
220                 false,                  // face flip
221                 mesh.boundaryMesh().whichPatch(curFaceID),// patch for face
222                 false,                  // remove from zone
223                 modifiedFaceZone,       // zone for face
224                 modifiedFaceZoneFlip    // face flip in zone
225             )
226         );
227     }
229     if (debug)
230     {
231         Pout<< "void attachDetach::attachInterface("
232             << "polyTopoChange& ref) const "
233             << " for object " << name() << " : "
234             << "Finished attaching interface" << endl;
235     }
239 void Foam::attachDetach::modifyMotionPoints
241     pointField& motionPoints
242 ) const
244     const Map<label>& removedPointMap = pointMatchMap();
246     const labelList removedPoints = removedPointMap.toc();
248     if (debug)
249     {
250         Pout<< "void attachDetach::modifyMotionPoints(" 
251             << "pointField& motionPoints) const "
252             << " for object " << name() << " : "
253             << "Adjusting motion points." << endl;
255         // Calculate the difference in motion point positions
256         scalar pointDiff = 0;
258         forAll (removedPoints, pointI)
259         {
260             pointDiff +=
261                 mag
262                 (
263                     motionPoints[removedPoints[pointI]]
264                   - motionPoints[removedPointMap.find(removedPoints[pointI])()]
265                 );
266         }
268         if (pointDiff > removedPoints.size()*positionDifference_)
269         {
270             Pout<< "Point motion difference = " << pointDiff << endl;
271         }
272     }
274     // Put the slave point on top of the master point
275     forAll (removedPoints, pointI)
276     {
277         motionPoints[removedPoints[pointI]] =
278             motionPoints[removedPointMap.find(removedPoints[pointI])()];
279     }
284 // ************************************************************************* //