PrimitivePatch.C
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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2021 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
27\*---------------------------------------------------------------------------*/
28
29#include "Map.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class FaceList, class PointField>
35(
36 const FaceList& faces,
37 const PointField& points
38)
39:
40 FaceList(faces),
41 points_(points),
42 edgesPtr_(nullptr),
43 nInternalEdges_(-1),
44 boundaryPointsPtr_(nullptr),
45 faceFacesPtr_(nullptr),
46 edgeFacesPtr_(nullptr),
47 faceEdgesPtr_(nullptr),
48 pointEdgesPtr_(nullptr),
49 pointFacesPtr_(nullptr),
50 localFacesPtr_(nullptr),
51 meshPointsPtr_(nullptr),
52 meshPointMapPtr_(nullptr),
53 edgeLoopsPtr_(nullptr),
54 localPointsPtr_(nullptr),
55 localPointOrderPtr_(nullptr),
56 faceCentresPtr_(nullptr),
57 faceAreasPtr_(nullptr),
58 magFaceAreasPtr_(nullptr),
59 faceNormalsPtr_(nullptr),
60 pointNormalsPtr_(nullptr)
61{}
62
63
64template<class FaceList, class PointField>
66(
67 FaceList&& faces,
68 const PointField& points
69)
70:
71 FaceList(std::move(faces)),
72 points_(points),
73 edgesPtr_(nullptr),
74 nInternalEdges_(-1),
75 boundaryPointsPtr_(nullptr),
76 faceFacesPtr_(nullptr),
77 edgeFacesPtr_(nullptr),
78 faceEdgesPtr_(nullptr),
79 pointEdgesPtr_(nullptr),
80 pointFacesPtr_(nullptr),
81 localFacesPtr_(nullptr),
82 meshPointsPtr_(nullptr),
83 meshPointMapPtr_(nullptr),
84 edgeLoopsPtr_(nullptr),
85 localPointsPtr_(nullptr),
86 localPointOrderPtr_(nullptr),
87 faceCentresPtr_(nullptr),
88 faceAreasPtr_(nullptr),
89 magFaceAreasPtr_(nullptr),
90 faceNormalsPtr_(nullptr),
91 pointNormalsPtr_(nullptr)
92{}
93
94
95template<class FaceList, class PointField>
97(
98 FaceList& faces,
100 const bool reuse
101)
102:
103 FaceList(faces, reuse),
104 points_(points, reuse),
105 edgesPtr_(nullptr),
106 nInternalEdges_(-1),
107 boundaryPointsPtr_(nullptr),
108 faceFacesPtr_(nullptr),
109 edgeFacesPtr_(nullptr),
110 faceEdgesPtr_(nullptr),
111 pointEdgesPtr_(nullptr),
112 pointFacesPtr_(nullptr),
113 localFacesPtr_(nullptr),
114 meshPointsPtr_(nullptr),
115 meshPointMapPtr_(nullptr),
116 edgeLoopsPtr_(nullptr),
117 localPointsPtr_(nullptr),
118 localPointOrderPtr_(nullptr),
119 faceCentresPtr_(nullptr),
120 faceAreasPtr_(nullptr),
121 magFaceAreasPtr_(nullptr),
122 faceNormalsPtr_(nullptr),
123 pointNormalsPtr_(nullptr)
124{}
125
126
127template<class FaceList, class PointField>
129(
131)
132:
133 FaceList(pp),
134 points_(pp.points_),
135 edgesPtr_(nullptr),
136 nInternalEdges_(-1),
137 boundaryPointsPtr_(nullptr),
138 faceFacesPtr_(nullptr),
139 edgeFacesPtr_(nullptr),
140 faceEdgesPtr_(nullptr),
141 pointEdgesPtr_(nullptr),
142 pointFacesPtr_(nullptr),
143 localFacesPtr_(nullptr),
144 meshPointsPtr_(nullptr),
145 meshPointMapPtr_(nullptr),
146 edgeLoopsPtr_(nullptr),
147 localPointsPtr_(nullptr),
148 localPointOrderPtr_(nullptr),
149 faceCentresPtr_(nullptr),
150 faceAreasPtr_(nullptr),
151 magFaceAreasPtr_(nullptr),
152 faceNormalsPtr_(nullptr),
153 pointNormalsPtr_(nullptr)
154{}
155
156
157// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
158
159template<class FaceList, class PointField>
161{
162 clearOut();
163}
164
165
166// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167
168template<class FaceList, class PointField>
169void
171(
172 const Field<point_type>&
173)
174{
175 DebugInFunction << "Recalculating geometry following mesh motion" << endl;
176
177 clearGeom();
178}
179
180
181template<class FaceList, class PointField>
182const Foam::edgeList&
184{
185 if (!edgesPtr_)
186 {
187 calcAddressing();
188 }
189
190 return *edgesPtr_;
191}
192
193
194template<class FaceList, class PointField>
197{
198 const edgeList& allEdges = this->edges(); // Force demand-driven
199 return edgeList::subList(allEdges, nInternalEdges());
200}
201
202
203template<class FaceList, class PointField>
206{
207 const edgeList& allEdges = this->edges(); // Force demand-driven
208 return edgeList::subList(allEdges, nBoundaryEdges(), nInternalEdges());
209}
210
211
212template<class FaceList, class PointField>
213Foam::label
215{
216 if (!edgesPtr_)
217 {
218 calcAddressing();
219 }
220
221 return nInternalEdges_;
222}
223
224
225template<class FaceList, class PointField>
226Foam::label
228{
229 const edgeList& allEdges = this->edges(); // Force demand-driven
230 return (allEdges.size() - this->nInternalEdges());
231}
232
233
234template<class FaceList, class PointField>
235const Foam::labelList&
237{
238 if (!boundaryPointsPtr_)
239 {
240 calcBdryPoints();
241 }
242
243 return *boundaryPointsPtr_;
244}
245
246
247template<class FaceList, class PointField>
250{
251 if (!faceFacesPtr_)
252 {
253 calcAddressing();
254 }
255
256 return *faceFacesPtr_;
257}
258
259
260template<class FaceList, class PointField>
263{
264 if (!edgeFacesPtr_)
265 {
266 calcAddressing();
268
269 return *edgeFacesPtr_;
270}
271
272
273template<class FaceList, class PointField>
276{
277 if (!faceEdgesPtr_)
278 {
279 calcAddressing();
280 }
281
282 return *faceEdgesPtr_;
283}
284
285
286template<class FaceList, class PointField>
289{
290 if (!pointEdgesPtr_)
291 {
292 calcPointEdges();
293 }
294
295 return *pointEdgesPtr_;
296}
297
298
299template<class FaceList, class PointField>
302{
303 if (!pointFacesPtr_)
304 {
305 calcPointFaces();
306 }
307
308 return *pointFacesPtr_;
309}
310
311
312template<class FaceList, class PointField>
313const Foam::List
314<
316>&
318{
319 if (!localFacesPtr_)
320 {
321 calcMeshData();
322 }
323
324 return *localFacesPtr_;
325}
326
327
328template<class FaceList, class PointField>
329const Foam::labelList&
332 if (!meshPointsPtr_)
333 {
334 calcMeshData();
335 }
336
337 return *meshPointsPtr_;
338}
339
341template<class FaceList, class PointField>
344{
345 if (!meshPointMapPtr_)
346 {
347 calcMeshPointMap();
348 }
349
350 return *meshPointMapPtr_;
352
353
354template<class FaceList, class PointField>
355const Foam::Field
356<
358>&
361 if (!localPointsPtr_)
362 {
363 calcLocalPoints();
364 }
365
366 return *localPointsPtr_;
367}
368
370template<class FaceList, class PointField>
371const Foam::labelList&
373{
374 if (!localPointOrderPtr_)
375 {
376 calcLocalPointOrder();
377 }
378
379 return *localPointOrderPtr_;
380}
381
382
383template<class FaceList, class PointField>
384Foam::label
386(
387 const label gp
388) const
390 // The point found, or -1 if not found
391 return meshPointMap().lookup(gp, -1);
393
394
395template<class FaceList, class PointField>
396const Foam::Field
397<
401{
402 if (!faceCentresPtr_)
403 {
404 calcFaceCentres();
405 }
406
407 return *faceCentresPtr_;
408}
409
410
411template<class FaceList, class PointField>
412const Foam::Field
413<
415>&
417{
418 if (!faceAreasPtr_)
419 {
420 calcFaceAreas();
421 }
422
423 return *faceAreasPtr_;
424}
425
426
427template<class FaceList, class PointField>
430{
431 if (!magFaceAreasPtr_)
432 {
433 calcMagFaceAreas();
434 }
435
436 return *magFaceAreasPtr_;
437}
438
439
440template<class FaceList, class PointField>
441const Foam::Field
442<
444>&
446{
447 if (!faceNormalsPtr_)
448 {
449 calcFaceNormals();
450 }
451
452 return *faceNormalsPtr_;
453}
454
456template<class FaceList, class PointField>
457const Foam::Field
460>&
462{
463 if (!pointNormalsPtr_)
464 {
465 calcPointNormals();
466 }
467
468 return *pointNormalsPtr_;
469}
470
471
472// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
473
474template<class FaceList, class PointField>
475void
477(
479)
480{
481 if (&rhs == this)
482 {
483 return;
484 }
485
486 clearOut();
487
488 FaceList::shallowCopy(rhs);
489
490 // Cannot copy assign points (could be const reference)
491}
492
493
494template<class FaceList, class PointField>
495void
497(
499)
500{
501 if (&rhs == this)
502 {
503 return;
504 }
505
506 clearOut();
507
508 FaceList::operator=(std::move(rhs));
509
510 // Cannot move assign points (could be const reference)
511}
512
513
514// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
515
518#include "PrimitivePatchClear.C"
526#include "PrimitivePatchCheck.C"
527
528// ************************************************************************* //
This function calculates the list of patch edges, defined on the list of points supporting the patch....
Checks topology of the patch.
Create the list of loops of outside vertices. Goes wrong on multiply connected edges (loops will be u...
Orders the local points on the patch for most efficient search.
Point addressing on the patch: pointEdges and pointFaces.
For every point on the patch find the closest face on the target side. Return a target face label for...
Generic templated field type.
Definition: Field.H:82
SubList< edge > subList
Declare type of subList.
Definition: List.H:111
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
A list of faces which address into the list of points.
label nBoundaryEdges() const
Number of boundary edges == (nEdges() - nInternalEdges())
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const Map< label > & meshPointMap() const
Mesh point map.
const edgeList::subList internalEdges() const
Return sub-list of internal edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const labelList & localPointOrder() const
Return orders the local points for most efficient search.
const Field< point_type > & pointNormals() const
Return point normals for patch.
const edgeList::subList boundaryEdges() const
Return sub-list of boundary edges, address into LOCAL point list.
std::remove_reference< PointField >::type::value_type point_type
The point type.
const Field< point_type > & faceAreas() const
Return face area vectors for patch.
label whichPoint(const label gp) const
Given a global point index, return the local point index.
const labelList & boundaryPoints() const
Return list of boundary points, address into LOCAL point list.
const Field< point_type > & faceCentres() const
Return face centres for patch.
const labelListList & faceFaces() const
Return face-face addressing.
const labelListList & pointFaces() const
Return point-face addressing.
std::remove_reference< FaceList >::type::value_type face_type
The face type.
const labelListList & edgeFaces() const
Return edge-face addressing.
const labelListList & faceEdges() const
Return face-edge addressing.
const Field< scalar > & magFaceAreas() const
Return face area magnitudes for patch.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A List obtained as a section of another List.
Definition: SubList.H:70
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
void movePoints()
Update for new mesh geometry.
const pointField & points
#define DebugInFunction
Report an information message using Foam::Info.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372