polyMeshFilter.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2012-2013 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::polyMeshFilter
26 
27 Description
28  Remove the edges and faces of a polyMesh whilst satisfying the given mesh
29  quality criteria.
30 
31  Works on a copy of the mesh.
32 
33 SourceFiles
34  polyMeshFilter.C
35  polyMeshFilterTemplates.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef polyMeshFilter_H
40 #define polyMeshFilter_H
41 
42 #include "IOdictionary.H"
43 #include "Time.H"
44 #include "List.H"
45 #include "autoPtr.H"
46 #include "scalarField.H"
47 #include "polyMeshFilterSettings.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 class polyMesh;
55 class fvMesh;
56 class PackedBoolList;
57 class faceSet;
58 
59 /*---------------------------------------------------------------------------*\
60  Class polyMeshFilter Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 class polyMeshFilter
64 :
66 {
67  // Private data
68 
69  //- Reference to the original mesh
70  const fvMesh& mesh_;
71 
72  //- Copy of the original mesh to perform the filtering on
74 
75  //- Original point priorities. If a point has a higher priority than
76  // another point then the edge between them collapses towards the
77  // point with the higher priority. (e.g. 2 is a higher priority than 1)
79 
80  //- Point priority associated with the new mesh
82 
83  //- The minimum edge length for each edge
85 
86  //- The face filter factor for each face
88 
89 
90  // Private Member Functions
91 
92  template<typename T>
93  static void updateSets(const mapPolyMesh& map);
94 
95  template<typename T>
96  static void copySets(const polyMesh& oldMesh, const polyMesh& newMesh);
97 
98  label filterFacesLoop(const label nOriginalBadFaces);
99 
101  (
102  polyMesh& newMesh,
103  scalarField& newMeshFaceFilterFactor,
104  labelList& origToCurrentPointMap
105  );
106 
108  (
109  polyMesh& newMesh,
110  scalarField& newMeshMinEdgeLen,
111  labelList& origToCurrentPointMap
112  );
113 
114  //- Increment pointErrorCount for points attached to a bad face
116  (
117  const PackedBoolList& isErrorPoint,
118  const labelList& oldToNewMesh,
119  labelList& pointErrorCount
120  ) const;
121 
122 
123  //- Given the new points that are part of bad faces, and a map from the
124  // old mesh points to the new mesh points, relax minEdgeLen_
126  (
127  const polyMesh& newMesh,
128  const labelList& oldToNewMesh,
129  const PackedBoolList& isErrorPoint,
130  const labelList& pointErrorCount
131  );
132 
133  //- Given the new points that are part of bad faces, and a map from the
134  // old mesh points to the new mesh points, relax faceFilterFactor_
136  (
137  const polyMesh& newMesh,
138  const labelList& oldToNewMesh,
139  const PackedBoolList& isErrorPoint,
140  const labelList& pointErrorCount
141  );
142 
143  // Mark boundary points
144  // boundaryPoint:
145  // + -1 : point not on boundary
146  // + 0 : point on a real boundary
147  // + >0 : point on a processor patch with that ID
148  // @todo Need to mark boundaryEdges as well, as an edge may have two
149  // boundary points but not itself lie on a boundary
151  (
152  const polyMesh& newMesh,
153  const labelList& pointMap
154  );
155 
156  //- Print min/mean/max data for a field
158  (
159  const string desc,
160  const scalarField& fld
161  ) const;
162 
163  //- Update minEdgeLen_ for the new mesh based upon the movement of the
164  // old points to the new points
166  (
167  const polyMesh& newMesh,
168  const labelList& pointMap,
169  scalarField& newMeshMinEdgeLen
170  ) const;
171 
172  //- Update faceFilterFactor_ for the new mesh based upon the movement
173  // of the old faces to the new faces
175  (
176  const polyMesh& newMesh,
177  const labelList& faceMap,
178  scalarField& newMeshFaceFilterFactor
179  ) const;
180 
181  //- Maintain a map of the original mesh points to the latest version of
182  // the filtered mesh.
184  (
185  const labelList& currToNew,
186  labelList& origToCurrentPointMap
187  ) const;
188 
189  //- Disallow default bitwise copy construct
191 
192  //- Disallow default bitwise assignment
193  void operator=(const polyMeshFilter&);
194 
195 
196 public:
197 
198  //- Runtime type information
199  ClassName("polyMeshFilter");
200 
201 
202  // Constructors
203 
204  //- Construct from fvMesh
205  explicit polyMeshFilter(const fvMesh& mesh);
206 
207  //- Construct from fvMesh and a label list of point priorities
209 
210 
211  //- Destructor
212  ~polyMeshFilter();
213 
214 
215  // Member Functions
216 
217  // Access
218 
219  //- Return reference to the filtered mesh. Does not check if the
220  // mesh has actually been filtered.
221  const autoPtr<fvMesh>& filteredMesh() const;
222 
223  //- Return the new pointPriority list.
224  const autoPtr<labelList>& pointPriority() const;
225 
226 
227  // Edit
228 
229  //- Return a copy of an fvMesh
230  static autoPtr<fvMesh> copyMesh(const fvMesh& mesh);
231 
232  //- Copy loaded topoSets from the old mesh to the new mesh
233  static void copySets
234  (
235  const polyMesh& oldMesh,
236  const polyMesh& newMesh
237  );
238 
239  //- Update the loaded topoSets
240  static void updateSets(const mapPolyMesh& map);
241 
242  //- Filter edges and faces
243  label filter(const label nOriginalBadFaces);
244 
245  //- Filter all faces that are in the face set
246  label filter(const faceSet& fSet);
247 
248  //- Filter edges only.
249  label filterEdges(const label nOriginalBadFaces);
250 };
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 } // End namespace Foam
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 #ifdef NoRepository
260 # include "polyMeshFilterTemplates.C"
261 #endif
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 #endif
266 
267 // ************************************************************************* //
Foam::polyMeshFilter::checkMeshFacesAndRelaxEdges
void checkMeshFacesAndRelaxEdges(const polyMesh &newMesh, const labelList &oldToNewMesh, const PackedBoolList &isErrorPoint, const labelList &pointErrorCount)
Given the new points that are part of bad faces, and a map from the.
Definition: polyMeshFilter.C:624
Foam::polyMeshFilter::mesh_
const fvMesh & mesh_
Reference to the original mesh.
Definition: polyMeshFilter.H:69
Foam::polyMeshFilter::newMeshPtr_
autoPtr< fvMesh > newMeshPtr_
Copy of the original mesh to perform the filtering on.
Definition: polyMeshFilter.H:72
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
List.H
Foam::polyMeshFilter::checkMeshEdgesAndRelaxEdges
void checkMeshEdgesAndRelaxEdges(const polyMesh &newMesh, const labelList &oldToNewMesh, const PackedBoolList &isErrorPoint, const labelList &pointErrorCount)
Given the new points that are part of bad faces, and a map from the.
Definition: polyMeshFilter.C:547
Foam::polyMeshFilter::ClassName
ClassName("polyMeshFilter")
Runtime type information.
scalarField.H
Foam::polyMeshFilter::updateOldToNewPointMap
void updateOldToNewPointMap(const labelList &currToNew, labelList &origToCurrentPointMap) const
Maintain a map of the original mesh points to the latest version of.
Definition: polyMeshFilter.C:860
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::polyMeshFilter::mapOldMeshFaceFieldToNewMesh
void mapOldMeshFaceFieldToNewMesh(const polyMesh &newMesh, const labelList &faceMap, scalarField &newMeshFaceFilterFactor) const
Update faceFilterFactor_ for the new mesh based upon the movement.
Definition: polyMeshFilter.C:833
Foam::polyMeshFilter::~polyMeshFilter
~polyMeshFilter()
Destructor.
Definition: polyMeshFilter.C:952
Foam::polyMeshFilter::operator=
void operator=(const polyMeshFilter &)
Disallow default bitwise assignment.
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::polyMeshFilter::faceFilterFactor_
scalarField faceFilterFactor_
The face filter factor for each face.
Definition: polyMeshFilter.H:86
polyMeshFilterSettings.H
Foam::polyMeshFilter::filterFacesLoop
label filterFacesLoop(const label nOriginalBadFaces)
Definition: polyMeshFilter.C:105
Foam::polyMeshFilter::filterEdges
label filterEdges(polyMesh &newMesh, scalarField &newMeshMinEdgeLen, labelList &origToCurrentPointMap)
Definition: polyMeshFilter.C:413
Foam::polyMeshFilter
Remove the edges and faces of a polyMesh whilst satisfying the given mesh quality criteria.
Definition: polyMeshFilter.H:62
Foam::polyMeshFilter::printScalarFieldStats
void printScalarFieldStats(const string desc, const scalarField &fld) const
Print min/mean/max data for a field.
Definition: polyMeshFilter.C:748
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::polyMeshFilter::polyMeshFilter
polyMeshFilter(const polyMeshFilter &)
Disallow default bitwise copy construct.
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::polyMeshFilter::updatePointErrorCount
void updatePointErrorCount(const PackedBoolList &isErrorPoint, const labelList &oldToNewMesh, labelList &pointErrorCount) const
Increment pointErrorCount for points attached to a bad face.
Definition: polyMeshFilter.C:530
Foam::polyMeshFilter::filteredMesh
const autoPtr< fvMesh > & filteredMesh() const
Return reference to the filtered mesh. Does not check if the.
Definition: polyMeshFilter.C:1107
Foam::polyMeshFilter::filter
label filter(const label nOriginalBadFaces)
Filter edges and faces.
Definition: polyMeshFilter.C:958
Foam::polyMeshFilter::updateSets
static void updateSets(const mapPolyMesh &map)
Definition: polyMeshFilterTemplates.C:34
Foam::polyMeshFilter::minEdgeLen_
scalarField minEdgeLen_
The minimum edge length for each edge.
Definition: polyMeshFilter.H:83
Foam::polyMeshFilter::originalPointPriority_
labelList originalPointPriority_
Original point priorities. If a point has a higher priority than.
Definition: polyMeshFilter.H:77
Foam::polyMeshFilterSettings
Class to store the settings for the polyMeshFilter class.
Definition: polyMeshFilterSettings.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::polyMeshFilter::pointPriority
const autoPtr< labelList > & pointPriority() const
Return the new pointPriority list.
Definition: polyMeshFilter.C:1114
Foam::polyMeshFilter::pointPriority_
autoPtr< labelList > pointPriority_
Point priority associated with the new mesh.
Definition: polyMeshFilter.H:80
Foam::polyMeshFilter::mapOldMeshEdgeFieldToNewMesh
void mapOldMeshEdgeFieldToNewMesh(const polyMesh &newMesh, const labelList &pointMap, scalarField &newMeshMinEdgeLen) const
Update minEdgeLen_ for the new mesh based upon the movement of the.
Definition: polyMeshFilter.C:797
IOdictionary.H
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
polyMeshFilterTemplates.C
Foam::polyMeshFilter::copySets
static void copySets(const polyMesh &oldMesh, const polyMesh &newMesh)
Definition: polyMeshFilterTemplates.C:71
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::polyMeshFilter::copyMesh
static autoPtr< fvMesh > copyMesh(const fvMesh &mesh)
Return a copy of an fvMesh.
Definition: polyMeshFilter.C:67
Foam::polyMeshFilter::updatePointPriorities
void updatePointPriorities(const polyMesh &newMesh, const labelList &pointMap)
Definition: polyMeshFilter.C:719
Foam::polyMeshFilter::filterFaces
label filterFaces(polyMesh &newMesh, scalarField &newMeshFaceFilterFactor, labelList &origToCurrentPointMap)
Definition: polyMeshFilter.C:297
autoPtr.H