porting changes
[OpenFOAM-1.5.x.git] / src / dynamicMesh / attachDetach / detachInterface.C
blob73d330a36d0d659d4b00e5bba97436bf7625d67c
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 "polyAddPoint.H"
33 #include "polyModifyFace.H"
34 #include "polyAddFace.H"
36 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
38 void Foam::attachDetach::detachInterface
40     polyTopoChange& ref
41 ) const
43     // Algorithm:
44     // 1. Create new points for points of the master face zone
45     // 2. Modify all faces of the master zone, by putting them into the master
46     //    patch (look for orientation) and their renumbered mirror images
47     //    into the slave patch
48     // 3. Create a point renumbering list, giving a new point index for original
49     //    points in the face patch
50     // 4. Grab all faces in cells on the master side and renumber them 
51     //    using the point renumbering list.  Exclude the ones that belong to
52     //    the master face zone
53     //
54     // Note on point creation:
55     // In order to take into account the issues related to partial
56     // blocking in an attach/detach mesh modifier, special treatment
57     // is required for the duplication of points on the edge of the
58     // face zone.  Points which are shared only by internal edges need
59     // not to be duplicated, as this would propagate the discontinuity
60     // in the mesh beyond the face zone.  Therefore, before creating
61     // the new points, check the external edge loop.  For each edge
62     // check if the edge is internal (i.e. does not belong to a
63     // patch); if so, exclude both of its points from duplication.
65     if (debug)
66     {
67         Pout<< "void attachDetach::detachInterface("
68             << "polyTopoChange& ref) const "
69             << " for object " << name() << " : "
70             << "Detaching interface" << endl;
71     }
73     const polyMesh& mesh = topoChanger().mesh();
74     const faceZoneMesh& zoneMesh = mesh.faceZones();
76     const primitiveFacePatch& masterFaceLayer = zoneMesh[faceZoneID_.index()]();
77     const pointField& points = mesh.points();
78     const labelListList& meshEdgeFaces = mesh.edgeFaces();
80     const labelList& mp = masterFaceLayer.meshPoints();
81     const edgeList& zoneLocalEdges = masterFaceLayer.edges();
83     const labelList& meshEdges = zoneMesh[faceZoneID_.index()].meshEdges();
85     // Create the points
87     labelList addedPoints(mp.size(), -1);
89     // Go through boundary edges of the master patch.  If all the faces from
90     // this patch are internal, mark the points in the addedPoints lookup
91     // with their original labels to stop duplication
92     label nIntEdges = masterFaceLayer.nInternalEdges();
94     for (label curEdgeID = nIntEdges; curEdgeID < meshEdges.size(); curEdgeID++)
95     {
96         const labelList& curFaces = meshEdgeFaces[meshEdges[curEdgeID]];
98         bool edgeIsInternal = true;
100         forAll (curFaces, faceI)
101         {
102             if (!mesh.isInternalFace(curFaces[faceI]))
103             {
104                 // The edge belongs to a boundary face
105                 edgeIsInternal = false;
106                 break;
107             }
108         }
110         if (edgeIsInternal)
111         {
112 // Pout<< "Internal edge found: (" << mp[zoneLocalEdges[curEdgeID].start()] << " " << mp[zoneLocalEdges[curEdgeID].end()] << ")" << endl;
114             // Reset the point creation
115             addedPoints[zoneLocalEdges[curEdgeID].start()] =
116                 mp[zoneLocalEdges[curEdgeID].start()];
118             addedPoints[zoneLocalEdges[curEdgeID].end()] =
119                 mp[zoneLocalEdges[curEdgeID].end()];
120         }
121     }
122 // Pout << "addedPoints before point creation: " << addedPoints << endl;
124     // Create new points for face zone
125     forAll (addedPoints, pointI)
126     {
127         if (addedPoints[pointI] < 0)
128         {
129             addedPoints[pointI] =
130                 ref.setAction
131                 (
132                     polyAddPoint
133                     (
134                         points[mp[pointI]],        // point
135                         mp[pointI],                // master point
136                         -1,                        // zone ID
137                         true                       // supports a cell
138                     )
139                 );
140 // Pout << "Adding point " << points[mp[pointI]] << " for original point " << mp[pointI] << endl;
141         }
142     }
144     // Modify faces in the master zone and duplicate for the slave zone
146     const labelList& mf = zoneMesh[faceZoneID_.index()];
147     const boolList& mfFlip = zoneMesh[faceZoneID_.index()].flipMap();
148     const faceList& zoneFaces = masterFaceLayer.localFaces();
150     const faceList& faces = mesh.faces();
151     const labelList& own = mesh.faceOwner();
152     const labelList& nei = mesh.faceNeighbour();
154     forAll (mf, faceI)
155     {
156         const label curFaceID = mf[faceI];
158         // Build the face for the slave patch by renumbering
159         const face oldFace = zoneFaces[faceI].reverseFace();
161         face newFace(oldFace.size());
163         forAll (oldFace, pointI)
164         {
165             newFace[pointI] = addedPoints[oldFace[pointI]];
166         }
168         if (mfFlip[faceI])
169         {
170             // Face needs to be flipped for the master patch
171             ref.setAction
172             (
173                 polyModifyFace
174                 (
175                     faces[curFaceID].reverseFace(), // modified face
176                     curFaceID,                   // label of face being modified
177                     nei[curFaceID],                 // owner
178                     -1,                             // neighbour
179                     true,                           // face flip
180                     masterPatchID_.index(),         // patch for face
181                     false,                          // remove from zone
182                     faceZoneID_.index(),            // zone for face
183                     !mfFlip[faceI]                  // face flip in zone
184                 )
185             );
187             // Add renumbered face into the slave patch
188             ref.setAction
189             (
190                 polyAddFace
191                 (
192                     newFace,                        // face
193                     own[curFaceID],                 // owner
194                     -1,                             // neighbour
195                     -1,                             // master point
196                     -1,                             // master edge
197                     curFaceID,                      // master face
198                     false,                          // flip flux
199                     slavePatchID_.index(),          // patch to add the face to
200                     -1,                             // zone for face
201                     false                           // zone flip
202                 )
203             );
204 // Pout << "Flip.  Modifying face: " << faces[curFaceID].reverseFace() << " next to cell: " << nei[curFaceID] << " and adding face: " << newFace << " next to cell: " << own[curFaceID] << endl;
205         }
206         else
207         {
208             // No flip
209             ref.setAction
210             (
211                 polyModifyFace
212                 (
213                     faces[curFaceID],         // modified face
214                     curFaceID,                // label of face being modified
215                     own[curFaceID],           // owner
216                     -1,                       // neighbour
217                     false,                    // face flip
218                     masterPatchID_.index(),   // patch for face
219                     false,                    // remove from zone
220                     faceZoneID_.index(),      // zone for face
221                     mfFlip[faceI]             // face flip in zone
222                 )
223             );
225             // Add renumbered face into the slave patch
226             ref.setAction
227             (
228                 polyAddFace
229                 (
230                     newFace,                        // face
231                     nei[curFaceID],                 // owner
232                     -1,                             // neighbour
233                     -1,                             // master point
234                     -1,                             // master edge
235                     curFaceID,                      // master face
236                     true,                           // flip flux
237                     slavePatchID_.index(),          // patch to add the face to
238                     -1,                             // zone for face
239                     false                           // face flip in zone
240                 )
241             );
242 // Pout << "No flip.  Modifying face: " << faces[curFaceID] << " next to cell: " << own[curFaceID] << " and adding face: " << newFace << " next to cell: " << nei[curFaceID] << endl;
243         }
244     }
246     // Modify the remaining faces of the master cells to reconnect to the new
247     // layer of faces.
249     // Algorithm: Go through all the cells of the master zone and make
250     // a map of faces to avoid duplicates.  Do not insert the faces in
251     // the master patch (as they have already been dealt with).  Make
252     // a master layer point renumbering map, which for every point in
253     // the master layer gives its new label. Loop through all faces in
254     // the map and attempt to renumber them using the master layer
255     // point renumbering map.  Once the face is renumbered, compare it
256     // with the original face; if they are the same, the face has not
257     // changed; if not, modify the face but keep all of its old
258     // attributes (apart from the vertex numbers).
260     // Create the map of faces in the master cell layer
261     const labelList& mc =
262         mesh.faceZones()[faceZoneID_.index()].masterCells();
264     labelHashSet masterCellFaceMap(6*mc.size());
266     const cellList& cells = mesh.cells();
268     forAll (mc, cellI)
269     {
270         const labelList& curFaces = cells[mc[cellI]];
272         forAll (curFaces, faceI)
273         {
274             // Check if the face belongs to the master patch; if not add it
275             if (zoneMesh.whichZone(curFaces[faceI]) != faceZoneID_.index())
276             {
277                 masterCellFaceMap.insert(curFaces[faceI]);
278             }
279         }
280     }
282     // Extend the map to include first neighbours of the master cells to
283     // deal with multiple corners.
284     { // Protection and memory management
285         // Make a map of master cells for quick reject
286         labelHashSet mcMap(2*mc.size());
288         forAll (mc, mcI)
289         {
290             mcMap.insert(mc[mcI]);
291         }
293         // Go through all the faces in the masterCellFaceMap.  If the
294         // cells around them are not already used, add all of their
295         // faces to the map
296         const labelList mcf = masterCellFaceMap.toc();
298         forAll (mcf, mcfI)
299         {
300             // Do the owner side
301             const label ownCell = own[mcf[mcfI]];
303             if (!mcMap.found(ownCell))
304             {
305                 // Cell not found. Add its faces to the map
306                 const cell& curFaces = cells[ownCell];
308                 forAll (curFaces, faceI)
309                 {
310                     masterCellFaceMap.insert(curFaces[faceI]);
311                 }
312             }
314             // Do the neighbour side if face is internal
315             if (mesh.isInternalFace(mcf[mcfI]))
316             {
317                 const label neiCell = nei[mcf[mcfI]];
319                 if (!mcMap.found(neiCell))
320                 {
321                     // Cell not found. Add its faces to the map
322                     const cell& curFaces = cells[neiCell];
324                     forAll (curFaces, faceI)
325                     {
326                         masterCellFaceMap.insert(curFaces[faceI]);
327                     }
328                 }
329             }
330         }
331     }
333     // Create the master layer point map
334     Map<label> masterLayerPointMap(2*mp.size());
336     forAll (mp, pointI)
337     {
338         masterLayerPointMap.insert
339         (
340             mp[pointI],
341             addedPoints[pointI]
342         );
343     }
345     // Grab the list of faces of the master layer
346     const labelList masterCellFaces = masterCellFaceMap.toc();
348     forAll (masterCellFaces, faceI)
349     {
350         // Attempt to renumber the face using the masterLayerPointMap.
351         // Missing point remain the same
353         const label curFaceID = masterCellFaces[faceI];
355         const face& oldFace = faces[curFaceID];
357         face newFace(oldFace.size());
359         bool changed = false;
361         forAll (oldFace, pointI)
362         {
363             if (masterLayerPointMap.found(oldFace[pointI]))
364             {
365                 changed = true;
367                 newFace[pointI] = masterLayerPointMap.find(oldFace[pointI])();
368             }
369             else
370             {
371                 newFace[pointI] = oldFace[pointI];
372             }
373         }
375         // If the face has changed, create a modification entry
376         if (changed)
377         {
378             if (mesh.isInternalFace(curFaceID))
379             {
380                 ref.setAction
381                 (
382                     polyModifyFace
383                     (
384                         newFace,                    // face
385                         curFaceID,                  // master face
386                         own[curFaceID],             // owner
387                         nei[curFaceID],             // neighbour
388                         false,                      // flip flux
389                         -1,                         // patch for face
390                         false,                      // remove from zone
391                         -1,                         // zone for face
392                         false                       // face zone flip
393                     )
394                 );
395 // Pout << "modifying stick-out face. Internal Old face: " << oldFace << " new face: " << newFace << " own: " << own[curFaceID] << " nei: " << nei[curFaceID] << endl;
396             }
397             else
398             {
399                 ref.setAction
400                 (
401                     polyModifyFace
402                     (
403                         newFace,                     // face
404                         curFaceID,                   // master face
405                         own[curFaceID],              // owner
406                         -1,                          // neighbour
407                         false,                       // flip flux
408                         mesh.boundaryMesh().whichPatch(curFaceID), // patch
409                         false,                        // remove from zone
410                         -1,                           // zone for face
411                         false                         // face zone flip
412                     )
413                 );   
414 // Pout << "modifying stick-out face. Boundary Old face: " << oldFace << " new face: " << newFace << " own: " << own[curFaceID] << " patch: " << mesh.boundaryMesh().whichPatch(curFaceID) << endl;
415             }                                                  
416         }
417     }
419     if (debug)
420     {
421         Pout<< "void attachDetach::detachInterface("
422             << "polyTopoChange& ref) const "
423             << " for object " << name() << " : "
424             << "Finished detaching interface" << endl;
425     }
429 // ************************************************************************* //