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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::polyMeshFilter
29
30Description
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
36SourceFiles
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"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace Foam
55{
56
57class polyMesh;
58class fvMesh;
59class faceSet;
60class bitSet;
61
62/*---------------------------------------------------------------------------*\
63 Class polyMeshFilter Declaration
64\*---------------------------------------------------------------------------*/
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 template<class T>
96 static void updateSets(const mapPolyMesh& map);
97
98 template<class T>
99 static void copySets(const polyMesh& oldMesh, const polyMesh& newMesh);
100
101 label filterFacesLoop(const label nOriginalBadFaces);
102
103 label filterFaces
104 (
105 polyMesh& newMesh,
106 scalarField& newMeshFaceFilterFactor,
107 labelList& origToCurrentPointMap
108 );
109
110 label filterEdges
111 (
112 polyMesh& newMesh,
113 scalarField& newMeshMinEdgeLen,
114 labelList& origToCurrentPointMap
115 );
116
117 //- Increment pointErrorCount for points attached to a bad face
118 void updatePointErrorCount
119 (
120 const bitSet& isErrorPoint,
121 const labelList& oldToNewMesh,
122 labelList& pointErrorCount
123 ) const;
124
125
126 //- Given the new points that are part of bad faces, and a map from the
127 // old mesh points to the new mesh points, relax minEdgeLen_
128 void checkMeshEdgesAndRelaxEdges
129 (
130 const polyMesh& newMesh,
131 const labelList& oldToNewMesh,
132 const bitSet& isErrorPoint,
133 const labelList& pointErrorCount
134 );
135
136 //- Given the new points that are part of bad faces, and a map from the
137 // old mesh points to the new mesh points, relax faceFilterFactor_
138 void checkMeshFacesAndRelaxEdges
139 (
140 const polyMesh& newMesh,
141 const labelList& oldToNewMesh,
142 const bitSet& isErrorPoint,
143 const labelList& pointErrorCount
144 );
145
146 // Mark boundary points
147 // boundaryPoint:
148 // + -1 : point not on boundary
149 // + 0 : point on a real boundary
150 // + >0 : point on a processor patch with that ID
151 // TODO: Need to mark boundaryEdges as well, as an edge may have two
152 // boundary points but not itself lie on a boundary
153 void updatePointPriorities
154 (
155 const polyMesh& newMesh,
156 const labelList& pointMap
157 );
158
159 //- Print min/mean/max data for a field
160 void printScalarFieldStats
161 (
162 const string desc,
163 const scalarField& fld
164 ) const;
165
166 //- Update minEdgeLen_ for the new mesh based upon the movement of the
167 // old points to the new points
168 void mapOldMeshEdgeFieldToNewMesh
169 (
170 const polyMesh& newMesh,
171 const labelList& pointMap,
172 scalarField& newMeshMinEdgeLen
173 ) const;
174
175 //- Update faceFilterFactor_ for the new mesh based upon the movement
176 // of the old faces to the new faces
177 void mapOldMeshFaceFieldToNewMesh
178 (
179 const polyMesh& newMesh,
180 const labelList& faceMap,
181 scalarField& newMeshFaceFilterFactor
182 ) const;
183
184 //- Maintain a map of the original mesh points to the latest version of
185 // the filtered mesh.
186 void updateOldToNewPointMap
187 (
188 const labelList& currToNew,
189 labelList& origToCurrentPointMap
190 ) const;
191
192 //- No copy construct
193 polyMeshFilter(const polyMeshFilter&) = delete;
194
195 //- No copy assignment
196 void operator=(const polyMeshFilter&) = delete;
197
198
199public:
200
201 //- Runtime type information
202 ClassName("polyMeshFilter");
203
204
205 // Constructors
206
207 //- Construct from fvMesh
208 explicit polyMeshFilter(const fvMesh& mesh);
209
210 //- Construct from fvMesh and a label list of point priorities
212
213 //- Construct from fvMesh and a label list of point priorities
215 (
216 const fvMesh& mesh,
218 const dictionary& dict
219 );
220
221
222 //- Destructor
223 ~polyMeshFilter() = default;
224
225
226 // Member Functions
227
228 // Access
229
230 //- Return reference to the filtered mesh. Does not check if the
231 // mesh has actually been filtered.
232 const autoPtr<fvMesh>& filteredMesh() const;
233
234 //- Return the new pointPriority list.
235 const autoPtr<labelList>& pointPriority() const;
236
237
238 // Edit
239
240 //- Return a copy of an fvMesh
241 static autoPtr<fvMesh> copyMesh(const fvMesh& mesh);
242
243 //- Copy loaded topoSets from the old mesh to the new mesh
244 static void copySets
245 (
246 const polyMesh& oldMesh,
247 const polyMesh& newMesh
248 );
249
250 //- Update the loaded topoSets
251 static void updateSets(const mapPolyMesh& map);
252
253 //- Filter edges and faces
254 label filter(const label nOriginalBadFaces);
255
256 //- Filter all faces that are in the face set
257 label filter(const faceSet& fSet);
258
259 //- Filter edges only.
260 label filterEdges(const label nOriginalBadFaces);
261};
262
263
264// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265
266} // End namespace Foam
267
268// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269
270#ifdef NoRepository
272#endif
273
274// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275
276#endif
277
278// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A list of face labels.
Definition: faceSet.H:54
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Class to store the settings for the polyMeshFilter class.
Remove the edges and faces of a polyMesh whilst satisfying the given mesh quality criteria.
~polyMeshFilter()=default
Destructor.
const autoPtr< labelList > & pointPriority() const
Return the new pointPriority list.
label filter(const label nOriginalBadFaces)
Filter edges and faces.
static autoPtr< fvMesh > copyMesh(const fvMesh &mesh)
Return a copy of an fvMesh.
ClassName("polyMeshFilter")
Runtime type information.
const autoPtr< fvMesh > & filteredMesh() const
Return reference to the filtered mesh. Does not check if the.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition: className.H:67
dynamicFvMesh & mesh
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
dictionary dict