intersectionPatchToPatch.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021-2023 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::patchToPatches::intersection
26 
27 Description
28  Class to generate patchToPatch coupling geometry. A full geometric
29  intersection is done between a face and those opposite, and coupling
30  geometry is calculated accordingly.
31 
32 SourceFiles
33  intersection.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef intersectionPatchToPatch_H
38 #define intersectionPatchToPatch_H
39 
40 #include "patchToPatch.H"
41 #include "polygonTriangulate.H"
42 #include "triFaceList.H"
43 #include "triIntersectLocation.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 namespace patchToPatches
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class intersection Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class intersection
57 :
58  public patchToPatch
59 {
60 public:
61 
62  // Public Structures
63 
64  //- Structure to store the geometry associated with part of a patch
65  // face
66  struct part
67  {
68  //- The area of this part
69  vector area;
70 
71  //- The centre of this part
72  point centre;
73 
74  //- Default construct
75  part()
76  :
77  area(Zero),
79  {}
80 
81  //- Default construct
82  part(const zero&)
83  :
84  part()
85  {}
86 
87  //- Construct from an area and a centre
88  part(const vector& a, const point& c)
89  :
90  area(a),
91  centre(c)
92  {}
93 
94  //- Construct from a polygon
95  template<class PointField>
96  part(const PointField& ps)
97  {
98  const Tuple2<vector, point> areaAndCentre =
100 
101  area = areaAndCentre.first();
102  centre = areaAndCentre.second();
103  }
104 
105  //- Negate this part
106  part operator-() const
107  {
108  return part(- area, centre);
109  }
110 
111  //- Add another part to this one
112  void operator+=(const part& p)
113  {
114  const scalar magArea = mag(area);
115  const scalar magPArea = mag(p.area);
116 
117  area = area + p.area;
118 
119  centre =
120  magArea == 0 ? p.centre
121  : magPArea == 0 ? centre
122  : (magArea*centre + magPArea*p.centre)/(magArea + magPArea);
123  }
124 
125  //- Subtract another part from this one
126  void operator-=(const part& p)
127  {
128  this->operator+=(-p);
129  }
130 
131  //- Equality comparison
132  friend bool operator==(const part& a, const part& b)
133  {
134  return a.area == b.area && a.centre == b.centre;
135  }
136 
137  //- Inequality comparison
138  friend bool operator!=(const part& a, const part& b)
139  {
140  return !(a == b);
141  }
142 
143  //- Output stream operator
144  friend Ostream& operator<<(Ostream& os, const part& p)
145  {
146  return os << p.area << p.centre;
147  }
148 
149  //- Input stream operator
150  friend Istream& operator>>(Istream& is, part& p)
151  {
152  return is >> p.area >> p.centre;
153  }
154  };
155 
156  //- Structure to store the geometry associated with the coupling
157  // between parts of two patch faces
158  struct couple
159  :
160  public part
161  {
162  //- The neighbour part
163  part nbr;
164 
165  //- Default construct
166  couple()
167  :
168  part(),
169  nbr()
170  {}
171 
172  //- Construct from a coupled pair of parts
173  couple(const part& p, const part& nbrP)
174  :
175  part(p),
176  nbr(nbrP)
177  {}
178 
179  //- Equality comparison
180  friend bool operator==(const couple& a, const couple& b)
181  {
182  return
183  static_cast<const part&>(a) == static_cast<const part&>(b)
184  && a.nbr == b.nbr;
185  }
186 
187  //- Inequality comparison
188  friend bool operator!=(const couple& a, const couple& b)
189  {
190  return !(a == b);
191  }
192 
193  //- Output stream operator
194  friend Ostream& operator<<(Ostream& os, const couple& c)
195  {
196  return os << static_cast<const part&>(c) << c.nbr;
197  }
198 
199  //- Input stream operator
200  friend Istream& operator>>(Istream& is, couple& c)
201  {
202  return is >> static_cast<part&>(c) >> c.nbr;
203  }
204  };
205 
206 
207 private:
208 
209  // Private Member Data
210 
211  // Geometry
212 
213  //- The coupling geometry for each source face
214  List<DynamicList<couple>> srcCouples_;
215 
216  //- The proportion of the source faces that are coupled
217  List<scalar> srcCoverage_;
218 
219  //- The non-coupled geometry associated with each source edge
220  List<part> srcEdgeParts_;
221 
222  //- The non-coupled geometry associated with mismatch in each
223  // source face's couplings
224  List<part> srcErrorParts_;
225 
226  //- The coupling geometry for each target face
227  List<DynamicList<couple>> tgtCouples_;
228 
229  //- The proportion of the target faces that are coupled
230  List<scalar> tgtCoverage_;
231 
232 
233  //- Triangulation engine
234  mutable polygonTriangulate triEngine_;
235 
236 
237  // Workspace
238 
239  //- Source face triangulation points
240  mutable List<triFaceList> srcTriPoints_;
241 
242  //- Source face triangulation edges
243  mutable List<List<FixedList<label, 3>>> srcTriFaceEdges_;
244 
245  //- Target face triangulation points
246  mutable List<triFaceList> tgtTriPoints_;
247 
248  //- Target face triangulation edges
249  mutable List<List<FixedList<label, 3>>> tgtTriFaceEdges_;
250 
251  //- Source intersection points
252  mutable DynamicList<point> ictSrcPoints_;
253 
254  //- Source intersection point normals
255  mutable DynamicList<vector> ictSrcPointNormals_;
256 
257  //- Target intersection points
258  mutable DynamicList<point> ictTgtPoints_;
259 
260  //- Intersection locations
261  mutable DynamicList<triIntersect::location> ictPointLocations_;
262 
263  //- Source face edge parts
264  DynamicList<part> srcFaceEdgePart_;
265 
266  //- Target face edge parts
267  DynamicList<part> tgtFaceEdgePart_;
268 
269  //- Source face edge parts
270  List<List<part>> srcFaceEdgeParts_;
271 
272 
273  // Debugging
274 
275  //- Intersection points
276  mutable DynamicList<point> debugPoints_;
277 
278  //- Intersection faces
279  mutable DynamicList<labelList> debugFaces_;
280 
281  //- The source face associated with each intersection face
282  mutable DynamicList<label> debugFaceSrcFaces_;
283 
284  //- The target face associated with each intersection face
285  mutable DynamicList<label> debugFaceTgtFaces_;
286 
287  //- The side of the intersection each face is on; 1 is the source
288  // side, -1 is the target side, and 0 is a face that spans the
289  // intersection volume; e.g., an edge or error face.
290  mutable DynamicList<label> debugFaceSides_;
291 
292 
293  // Private Static Member Functions
294 
295  //- Return the values at the points of a tri face
296  template<class Type>
297  static FixedList<Type, 3> triPointValues
298  (
299  const triFace& t,
300  const UList<Type>& values
301  );
302 
303 
304  // Private Member Functions
305 
306  //- Get the bound box for a source face
307  virtual treeBoundBox srcBox
308  (
309  const face& srcFace,
310  const pointField& srcPoints,
311  const vectorField& srcPointNormals
312  ) const;
313 
314  //- Intersect two faces
315  virtual bool intersectFaces
316  (
317  const primitiveOldTimePatch& srcPatch,
318  const vectorField& srcPointNormals,
319  const vectorField& srcPointNormals0,
320  const primitiveOldTimePatch& tgtPatch,
321  const label srcFacei,
322  const label tgtFacei
323  );
324 
325  //- Initialise the workspace
326  virtual void initialise
327  (
328  const primitiveOldTimePatch& srcPatch,
329  const vectorField& srcPointNormals,
330  const vectorField& srcPointNormals0,
331  const primitiveOldTimePatch& tgtPatch
332  );
333 
334  //- Finalise the intersection locally. Trims the local target patch
335  // so that only parts that actually intersect the source remain.
336  // Returns new-to-old map for trimming target-associated data.
337  virtual labelList finaliseLocal
338  (
339  const primitiveOldTimePatch& srcPatch,
340  const vectorField& srcPointNormals,
341  const vectorField& srcPointNormals0,
342  const primitiveOldTimePatch& tgtPatch
343  );
344 
345  //- Send data that resulted from an intersection between the source
346  // patch and a distributed source-local-target patch back to the
347  // original target processes
348  virtual void rDistributeTgt(const primitiveOldTimePatch& tgtPatch);
349 
350  //- Finalise the intersection
351  virtual label finalise
352  (
353  const primitiveOldTimePatch& srcPatch,
354  const vectorField& srcPointNormals,
355  const vectorField& srcPointNormals0,
356  const primitiveOldTimePatch& tgtPatch,
357  const transformer& tgtToSrc
358  );
359 
360  //- For each source face, the coupled target weights
361  virtual tmpNrc<List<DynamicList<scalar>>> srcWeights() const;
362 
363  //- For each target face, the coupled source weights
364  virtual tmpNrc<List<DynamicList<scalar>>> tgtWeights() const;
365 
366 
367 public:
368 
369  //- Runtime type information
370  TypeName("intersection");
371 
372 
373  // Static Data Members
374 
375  //- Extra debugging for intersections between specific faces. Named
376  // "intersectionSrcFace" and "intersectionTgtFace" respectively.
377  static int debugSrcFacei, debugTgtFacei;
378 
379 
380  // Static Member Functions
381 
382  //- Get the bound box for a source face
384  (
385  const face& srcFace,
386  const pointField& srcPoints,
387  const vectorField& srcPointNormals
388  );
389 
390 
391  // Constructors
392 
393  //- Construct from components
394  intersection(const bool reverse);
395 
396 
397  //- Destructor
398  ~intersection();
399 
400 
401  // Member Functions
402 
403  //- For each source face, the source and target areas for each
404  // target coupling
405  inline const List<DynamicList<couple>>& srcCouples() const;
406 
407  //- For each source edge, the non-coupled geometry associated
408  // with its projection
409  inline const List<part>& srcEdgeParts() const;
410 
411  //- For each source face, the area associated with mismatch
412  // across the coupling
413  inline const List<part>& srcErrorParts() const;
414 
415  //- For each target face, the target and source areas for each
416  // source coupling
417  inline const List<DynamicList<couple>>& tgtCouples() const;
418 
419  //- Return the proportion of the source faces that are coupled
420  inline const List<scalar>& srcCoverage() const;
421 
422  //- Return the proportion of the target faces that are coupled
423  inline const List<scalar>& tgtCoverage() const;
424 };
425 
426 
427 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
428 
429 } // End namespace patchToPatches
430 } // End namespace Foam
431 
432 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
433 
435 
436 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
437 
438 #endif
439 
440 // ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type2 & second() const
Return second.
Definition: Tuple2.H:131
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
static Tuple2< vector, point > areaAndCentreStabilised(const PointField &)
Return vector area and centre point given face points. Stabilised.
Traits class for primitives.
Definition: pTraits.H:53
Class to generate coupling geometry between two primitive patches.
Definition: patchToPatch.H:56
tmp< Field< Type > > tgtToSrc(const Field< Type > &tgtFld) const
Interpolate a target patch field to the source with no left.
bool reverse() const
Flag to indicate that the two patches are co-directional and.
Definition: patchToPatchI.H:31
Class to generate patchToPatch coupling geometry. A full geometric intersection is done between a fac...
TypeName("intersection")
Runtime type information.
const List< DynamicList< couple > > & tgtCouples() const
For each target face, the target and source areas for each.
const List< part > & srcErrorParts() const
For each source face, the area associated with mismatch.
const List< DynamicList< couple > > & srcCouples() const
For each source face, the source and target areas for each.
const List< scalar > & srcCoverage() const
Return the proportion of the source faces that are coupled.
static treeBoundBox srcBoxStatic(const face &srcFace, const pointField &srcPoints, const vectorField &srcPointNormals)
Get the bound box for a source face.
static int debugSrcFacei
Extra debugging for intersections between specific faces. Named.
const List< part > & srcEdgeParts() const
For each source edge, the non-coupled geometry associated.
intersection(const bool reverse)
Construct from components.
const List< scalar > & tgtCoverage() const
Return the proportion of the target faces that are coupled.
Triangulation of three-dimensional polygons.
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:53
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:71
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
volScalarField & b
Definition: createFields.H:25
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
Structure to store the geometry associated with the coupling.
friend bool operator==(const couple &a, const couple &b)
Equality comparison.
friend bool operator!=(const couple &a, const couple &b)
Inequality comparison.
friend Istream & operator>>(Istream &is, couple &c)
Input stream operator.
friend Ostream & operator<<(Ostream &os, const couple &c)
Output stream operator.
Structure to store the geometry associated with part of a patch.
void operator-=(const part &p)
Subtract another part from this one.
friend Istream & operator>>(Istream &is, part &p)
Input stream operator.
friend Ostream & operator<<(Ostream &os, const part &p)
Output stream operator.
friend bool operator!=(const part &a, const part &b)
Inequality comparison.
friend bool operator==(const part &a, const part &b)
Equality comparison.
void operator+=(const part &p)
Add another part to this one.
Class to encapsulate the topology of a point within a triangle intersection.