initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / autoMesh / autoHexMesh / meshRefinement / meshRefinementMerge.C
blobb1461cea77c00c17589eeb11153ebf2d415499bd
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 "meshRefinement.H"
28 #include "combineFaces.H"
29 #include "polyTopoChange.H"
30 #include "removePoints.H"
32 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
34 // Merge faces that are in-line.
35 Foam::label Foam::meshRefinement::mergePatchFaces
37     const scalar minCos,
38     const scalar concaveCos,
39     const labelList& patchIDs
42     // Patch face merging engine
43     combineFaces faceCombiner(mesh_);
45     const polyBoundaryMesh& patches = mesh_.boundaryMesh();
47     // Pick up all candidate cells on boundary
48     labelHashSet boundaryCells(mesh_.nFaces()-mesh_.nInternalFaces());
50     forAll(patchIDs, i)
51     {
52         label patchI = patchIDs[i];
54         const polyPatch& patch = patches[patchI];
56         if (!patch.coupled())
57         {
58             forAll(patch, i)
59             {
60                 boundaryCells.insert(mesh_.faceOwner()[patch.start()+i]);
61             }
62         }
63     }
65     // Get all sets of faces that can be merged
66     labelListList mergeSets
67     (
68         faceCombiner.getMergeSets
69         (
70             minCos,
71             concaveCos,
72             boundaryCells
73         )
74     );
76     label nFaceSets = returnReduce(mergeSets.size(), sumOp<label>());
78     Info<< "mergePatchFaces : Merging " << nFaceSets
79         << " sets of faces." << endl;
81     if (nFaceSets > 0)
82     {
83         // Topology changes container
84         polyTopoChange meshMod(mesh_);
86         // Merge all faces of a set into the first face of the set. Remove
87         // unused points.
88         faceCombiner.setRefinement(mergeSets, meshMod);
90         // Change the mesh (no inflation)
91         autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh_, false, true);
93         // Update fields
94         mesh_.updateMesh(map);
96         // Move mesh (since morphing does not do this)
97         if (map().hasMotionPoints())
98         {
99             mesh_.movePoints(map().preMotionPoints());
100         }
101         else
102         {
103             // Delete mesh volumes. No other way to do this?
104             mesh_.clearOut();
105         }
107         if (overwrite())
108         {
109             mesh_.setInstance(oldInstance());
110         }
112         faceCombiner.updateMesh(map);
114         // Get the kept faces that need to be recalculated.
115         // Merging two boundary faces might shift the cell centre
116         // (unless the faces are absolutely planar)
117         labelHashSet retestFaces(6*mergeSets.size());
119         forAll(mergeSets, setI)
120         {
121             label oldMasterI = mergeSets[setI][0];
123             label faceI = map().reverseFaceMap()[oldMasterI];
125             // faceI is always uncoupled boundary face
126             const cell& cFaces = mesh_.cells()[mesh_.faceOwner()[faceI]];
128             forAll(cFaces, i)
129             {
130                 retestFaces.insert(cFaces[i]);
131             }
132         }
133         updateMesh(map, retestFaces.toc());
134     }
137     return nFaceSets;
141 // Remove points not used by any face or points used by only two faces where
142 // the edges are in line
143 Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::mergeEdges
145     const scalar minCos
148     // Point removal analysis engine
149     removePoints pointRemover(mesh_);
151     // Count usage of points
152     boolList pointCanBeDeleted;
153     label nRemove = pointRemover.countPointUsage(minCos, pointCanBeDeleted);
155     Info<< "Removing " << nRemove
156         << " straight edge points." << endl;
158     autoPtr<mapPolyMesh> map;
160     if (nRemove > 0)
161     {
162         // Save my local faces that will change. These changed faces might
163         // cause a shift in the cell centre which needs to be retested.
164         // Have to do this before changing mesh since point will be removed.
165         labelHashSet retestOldFaces(nRemove / Pstream::nProcs());
167         {
168             const faceList& faces = mesh_.faces();
170             forAll(faces, faceI)
171             {
172                 const face& f = faces[faceI];
174                 forAll(f, fp)
175                 {
176                     if (pointCanBeDeleted[f[fp]])
177                     {
178                         retestOldFaces.insert(faceI);
179                         break;
180                     }
181                 }
182             }
183         }
185         // Topology changes container
186         polyTopoChange meshMod(mesh_);
188         pointRemover.setRefinement(pointCanBeDeleted, meshMod);
190         // Change the mesh (no inflation)
191         map = meshMod.changeMesh(mesh_, false, true);
193         // Update fields
194         mesh_.updateMesh(map);
196         // Move mesh (since morphing does not do this)
197         if (map().hasMotionPoints())
198         {
199             mesh_.movePoints(map().preMotionPoints());
200         }
201         else
202         {
203             // Delete mesh volumes. No other way to do this?
204             mesh_.clearOut();
205         }
207         if (overwrite())
208         {
209             mesh_.setInstance(oldInstance());
210         }
212         pointRemover.updateMesh(map);
214         // Get the kept faces that need to be recalculated.
215         labelHashSet retestFaces(6*retestOldFaces.size());
217         const cellList& cells = mesh_.cells();
219         forAllConstIter(labelHashSet, retestOldFaces, iter)
220         {
221             label faceI = map().reverseFaceMap()[iter.key()];
223             const cell& ownFaces = cells[mesh_.faceOwner()[faceI]];
225             forAll(ownFaces, i)
226             {
227                 retestFaces.insert(ownFaces[i]);
228             }
230             if (mesh_.isInternalFace(faceI))
231             {
232                 const cell& neiFaces = cells[mesh_.faceNeighbour()[faceI]];
234                 forAll(neiFaces, i)
235                 {
236                     retestFaces.insert(neiFaces[i]);
237                 }
238             }
239         }
240         updateMesh(map, retestFaces.toc());
241     }
243     return map;
247 // ************************************************************************* //