Loading...
Searching...
No Matches
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 nInternalEdges_(-1)
43{}
44
45
46template<class FaceList, class PointField>
48(
49 FaceList&& faces,
50 const PointField& points
51)
52:
53 FaceList(std::move(faces)),
54 points_(points),
55 nInternalEdges_(-1)
56{}
57
58
59template<class FaceList, class PointField>
61(
62 FaceList& faces,
64 const bool reuse
65)
66:
67 FaceList(faces, reuse),
68 points_(points, reuse),
69 nInternalEdges_(-1)
70{}
71
72
73template<class FaceList, class PointField>
75(
77)
78:
79 FaceList(pp),
80 points_(pp.points_),
81 nInternalEdges_(-1)
82{}
83
84
85// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
86
87template<class FaceList, class PointField>
89{
90 clearOut();
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
96template<class FaceList, class PointField>
97void
99(
100 const Field<point_type>&
101)
102{
103 DebugInFunction << "Recalculating geometry following mesh motion" << endl;
105 clearGeom();
106}
107
108
109template<class FaceList, class PointField>
110const Foam::edgeList&
112{
113 if (!edgesPtr_)
114 {
115 calcAddressing();
116 }
118 return *edgesPtr_;
119}
120
121
122template<class FaceList, class PointField>
125{
126 const edgeList& allEdges = this->edges(); // Force demand-driven
127 return edgeList::subList(allEdges, nInternalEdges());
128}
129
130
131template<class FaceList, class PointField>
134{
135 const edgeList& allEdges = this->edges(); // Force demand-driven
136 return edgeList::subList(allEdges, nBoundaryEdges(), nInternalEdges());
137}
138
139
140template<class FaceList, class PointField>
141Foam::label
143{
144 if (!edgesPtr_)
145 {
146 calcAddressing();
147 }
149 return nInternalEdges_;
150}
151
152
153template<class FaceList, class PointField>
154Foam::label
156{
157 const edgeList& allEdges = this->edges(); // Force demand-driven
158 return (allEdges.size() - this->nInternalEdges());
159}
160
161
162template<class FaceList, class PointField>
163const Foam::labelList&
165{
166 if (!boundaryPointsPtr_)
167 {
168 calcBdryPoints();
169 }
171 return *boundaryPointsPtr_;
172}
173
174
175template<class FaceList, class PointField>
178{
179 if (!faceFacesPtr_)
180 {
181 calcAddressing();
182 }
184 return *faceFacesPtr_;
185}
186
187
188template<class FaceList, class PointField>
191{
192 if (!edgeFacesPtr_)
193 {
194 calcAddressing();
195 }
197 return *edgeFacesPtr_;
198}
199
200
201template<class FaceList, class PointField>
204{
205 if (!faceEdgesPtr_)
206 {
207 calcAddressing();
208 }
210 return *faceEdgesPtr_;
211}
212
213
214template<class FaceList, class PointField>
217{
218 if (!pointEdgesPtr_)
219 {
220 calcPointEdges();
221 }
223 return *pointEdgesPtr_;
224}
225
226
227template<class FaceList, class PointField>
230{
231 if (!pointFacesPtr_)
232 {
233 calcPointFaces();
234 }
235
236 return *pointFacesPtr_;
237}
239
240template<class FaceList, class PointField>
241const Foam::List
242<
244>&
246{
247 if (!localFacesPtr_)
248 {
249 calcMeshData();
250 }
252 return *localFacesPtr_;
253}
254
255
256template<class FaceList, class PointField>
257const Foam::labelList&
259{
260 if (!meshPointsPtr_)
261 {
262 calcMeshData();
263 }
265 return *meshPointsPtr_;
266}
267
268
269template<class FaceList, class PointField>
272{
273 if (!meshPointMapPtr_)
274 {
275 calcMeshPointMap();
276 }
277
278 return *meshPointMapPtr_;
279}
281
282template<class FaceList, class PointField>
283const Foam::Field
284<
286>&
288{
289 if (!localPointsPtr_)
290 {
291 calcLocalPoints();
292 }
294 return *localPointsPtr_;
295}
296
297
298template<class FaceList, class PointField>
299const Foam::labelList&
301{
302 if (!localPointOrderPtr_)
303 {
304 calcLocalPointOrder();
305 }
307 return *localPointOrderPtr_;
308}
309
310
311template<class FaceList, class PointField>
312Foam::label
314(
315 const label gp
316) const
317{
318 // The point found, or -1 if not found
319 return meshPointMap().lookup(gp, -1);
320}
322
323template<class FaceList, class PointField>
324const Foam::Field
325<
327>&
329{
330 if (!faceCentresPtr_)
331 {
332 calcFaceCentres();
333 }
334
335 return *faceCentresPtr_;
336}
338
339template<class FaceList, class PointField>
340const Foam::Field
341<
345{
346 if (!faceAreasPtr_)
347 {
348 calcFaceAreas();
349 }
351 return *faceAreasPtr_;
352}
353
354
355template<class FaceList, class PointField>
358{
359 if (!magFaceAreasPtr_)
360 {
361 calcMagFaceAreas();
362 }
363
364 return *magFaceAreasPtr_;
365}
367
368template<class FaceList, class PointField>
369const Foam::Field
370<
372>&
374{
375 if (!faceNormalsPtr_)
376 {
377 calcFaceNormals();
378 }
379
380 return *faceNormalsPtr_;
381}
383
384template<class FaceList, class PointField>
385const Foam::Field
386<
388>&
390{
391 if (!pointNormalsPtr_)
392 {
393 calcPointNormals();
394 }
395
396 return *pointNormalsPtr_;
398
399
400// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
401
402template<class FaceList, class PointField>
403void
405(
407)
408{
409 if (&rhs == this)
410 {
411 return;
412 }
413
414 clearOut();
415
416 FaceList::shallowCopy(rhs);
418 // Cannot copy assign points (could be const reference)
419}
420
421
422template<class FaceList, class PointField>
423void
425(
427)
428{
429 if (&rhs == this)
430 {
431 return;
432 }
433
434 clearOut();
435
436 FaceList::operator=(std::move(rhs));
437
438 // Cannot move assign points (could be const reference)
439}
440
441
442// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446#include "PrimitivePatchClear.C"
454#include "PrimitivePatchCheck.C"
455
456// ************************************************************************* //
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...
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
SubList< edge > subList
Definition List.H:129
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
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.
PrimitivePatch(const FaceList &faces, const PointField &points)
Construct from components.
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 > & points() const noexcept
Return reference to global points.
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.
virtual void movePoints(const Field< point_type > &)
Correct patch after moving points.
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.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const pointField & points
#define DebugInFunction
Report an information message using Foam::Info.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
GeometricField< Type, pointPatchField, pointMesh > PointField
A point field for a given type.