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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::polyMeshFilter
29 
30 Description
31  Remove the edges and faces of a polyMesh whilst satisfying the given mesh
32  quality criteria.
33 
34  Works on a copy of the mesh.
35 
36 SourceFiles
37  polyMeshFilter.C
38  polyMeshFilterTemplates.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef polyMeshFilter_H
43 #define polyMeshFilter_H
44 
45 #include "IOdictionary.H"
46 #include "Time.H"
47 #include "List.H"
48 #include "autoPtr.H"
49 #include "scalarField.H"
50 #include "polyMeshFilterSettings.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 class polyMesh;
58 class fvMesh;
59 class faceSet;
60 class bitSet;
61 
62 /*---------------------------------------------------------------------------*\
63  Class polyMeshFilter Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class polyMeshFilter
67 :
69 {
70  // Private data
71 
72  //- Reference to the original mesh
73  const fvMesh& mesh_;
74 
75  //- Copy of the original mesh to perform the filtering on
76  autoPtr<fvMesh> newMeshPtr_;
77 
78  //- Original point priorities. If a point has a higher priority than
79  // another point then the edge between them collapses towards the
80  // point with the higher priority. (e.g. 2 is a higher priority than 1)
81  labelList originalPointPriority_;
82 
83  //- Point priority associated with the new mesh
84  autoPtr<labelList> pointPriority_;
85 
86  //- The minimum edge length for each edge
87  scalarField minEdgeLen_;
88 
89  //- The face filter factor for each face
90  scalarField faceFilterFactor_;
91 
92 
93  // Private Member Functions
94 
95  //- Update all sets of the given type
96  template<class SetType>
97  static void updateSets(const mapPolyMesh& map);
98 
99  template<class SetType>
100  static void copySets(const polyMesh& oldMesh, const polyMesh& newMesh);
101 
102  label filterFacesLoop(const label nOriginalBadFaces);
103 
104  label filterFaces
105  (
106  polyMesh& newMesh,
107  scalarField& newMeshFaceFilterFactor,
108  labelList& origToCurrentPointMap
109  );
110 
111  label filterEdges
112  (
113  polyMesh& newMesh,
114  scalarField& newMeshMinEdgeLen,
115  labelList& origToCurrentPointMap
116  );
117 
118  //- Increment pointErrorCount for points attached to a bad face
119  void updatePointErrorCount
120  (
121  const bitSet& isErrorPoint,
122  const labelList& oldToNewMesh,
123  labelList& pointErrorCount
124  ) const;
125 
126 
127  //- Given the new points that are part of bad faces, and a map from the
128  // old mesh points to the new mesh points, relax minEdgeLen_
129  void checkMeshEdgesAndRelaxEdges
130  (
131  const polyMesh& newMesh,
132  const labelList& oldToNewMesh,
133  const bitSet& isErrorPoint,
134  const labelList& pointErrorCount
135  );
136 
137  //- Given the new points that are part of bad faces, and a map from the
138  // old mesh points to the new mesh points, relax faceFilterFactor_
139  void checkMeshFacesAndRelaxEdges
140  (
141  const polyMesh& newMesh,
142  const labelList& oldToNewMesh,
143  const bitSet& isErrorPoint,
144  const labelList& pointErrorCount
145  );
146 
147  // Mark boundary points
148  // boundaryPoint:
149  // + -1 : point not on boundary
150  // + 0 : point on a real boundary
151  // + >0 : point on a processor patch with that ID
152  // TODO: Need to mark boundaryEdges as well, as an edge may have two
153  // boundary points but not itself lie on a boundary
154  void updatePointPriorities
155  (
156  const polyMesh& newMesh,
157  const labelList& pointMap
158  );
159 
160  //- Print min/mean/max data for a field
161  void printScalarFieldStats
162  (
163  const string desc,
164  const scalarField& fld
165  ) const;
166 
167  //- Update minEdgeLen_ for the new mesh based upon the movement of the
168  // old points to the new points
169  void mapOldMeshEdgeFieldToNewMesh
170  (
171  const polyMesh& newMesh,
172  const labelList& pointMap,
173  scalarField& newMeshMinEdgeLen
174  ) const;
175 
176  //- Update faceFilterFactor_ for the new mesh based upon the movement
177  // of the old faces to the new faces
178  void mapOldMeshFaceFieldToNewMesh
179  (
180  const polyMesh& newMesh,
181  const labelList& faceMap,
182  scalarField& newMeshFaceFilterFactor
183  ) const;
184 
185  //- Maintain a map of the original mesh points to the latest version of
186  // the filtered mesh.
187  void updateOldToNewPointMap
188  (
189  const labelList& currToNew,
190  labelList& origToCurrentPointMap
191  ) const;
192 
193  //- No copy construct
194  polyMeshFilter(const polyMeshFilter&) = delete;
195 
196  //- No copy assignment
197  void operator=(const polyMeshFilter&) = delete;
198 
199 
200 public:
201 
202  //- Runtime type information
203  ClassName("polyMeshFilter");
204 
205 
206  // Constructors
207 
208  //- Construct from fvMesh
209  explicit polyMeshFilter(const fvMesh& mesh);
210 
211  //- Construct from fvMesh and a label list of point priorities
213 
214  //- Construct from fvMesh and a label list of point priorities
216  (
217  const fvMesh& mesh,
218  const labelList& pointPriority,
219  const dictionary& dict
220  );
221 
222 
223  //- Destructor
224  ~polyMeshFilter() = default;
225 
226 
227  // Member Functions
228 
229  // Access
230 
231  //- Return reference to the filtered mesh. Does not check if the
232  // mesh has actually been filtered.
233  const autoPtr<fvMesh>& filteredMesh() const;
234 
235  //- Return the new pointPriority list.
236  const autoPtr<labelList>& pointPriority() const;
237 
238 
239  // Edit
240 
241  //- Return a copy of an fvMesh
242  static autoPtr<fvMesh> copyMesh(const fvMesh& mesh);
243 
244  //- Copy loaded topoSets from the old mesh to the new mesh
245  static void copySets
246  (
247  const polyMesh& oldMesh,
248  const polyMesh& newMesh
249  );
250 
251  //- Update the loaded topoSets
252  static void updateSets(const mapPolyMesh& map);
253 
254  //- Filter edges and faces
255  label filter(const label nOriginalBadFaces);
256 
257  //- Filter all faces that are in the face set
258  label filter(const faceSet& fSet);
259 
260  //- Filter edges only.
261  label filterEdges(const label nOriginalBadFaces);
262 };
263 
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 } // End namespace Foam
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 #ifdef NoRepository
272  #include "polyMeshFilterTemplates.C"
273 #endif
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 #endif
278 
279 // ************************************************************************* //
dictionary dict
label filter(const label nOriginalBadFaces)
Filter edges and faces.
A list of face labels.
Definition: faceSet.H:47
const autoPtr< labelList > & pointPriority() const
Return the new pointPriority list.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static autoPtr< fvMesh > copyMesh(const fvMesh &mesh)
Return a copy of an fvMesh.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
dynamicFvMesh & mesh
Remove the edges and faces of a polyMesh whilst satisfying the given mesh quality criteria...
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
ClassName("polyMeshFilter")
Runtime type information.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
~polyMeshFilter()=default
Destructor.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const autoPtr< fvMesh > & filteredMesh() const
Return reference to the filtered mesh. Does not check if the.
Class to store the settings for the polyMeshFilter class.
Namespace for OpenFOAM.