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-2022 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),
78  centre(point::uniform(NaN))
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  area(face::area(ps)),
99  centre(face::centre(ps))
100  {}
101 
102  //- Negate this part
103  part operator-() const
104  {
105  return part(- area, centre);
106  }
107 
108  //- Add another part to this one
109  void operator+=(const part& p)
110  {
111  const scalar magArea = mag(area);
112  const scalar magPArea = mag(p.area);
113 
114  area = area + p.area;
115 
116  centre =
117  magArea == 0 ? p.centre
118  : magPArea == 0 ? centre
119  : (magArea*centre + magPArea*p.centre)/(magArea + magPArea);
120  }
121 
122  //- Subtract another part from this one
123  void operator-=(const part& p)
124  {
125  this->operator+=(-p);
126  }
127 
128  //- Equality comparison
129  friend bool operator==(const part& a, const part& b)
130  {
131  return a.area == b.area && a.centre == b.centre;
132  }
133 
134  //- Inequality comparison
135  friend bool operator!=(const part& a, const part& b)
136  {
137  return !(a == b);
138  }
139 
140  //- Output stream operator
141  friend Ostream& operator<<(Ostream& os, const part& p)
142  {
143  return os << p.area << p.centre;
144  }
145 
146  //- Input stream operator
147  friend Istream& operator>>(Istream& is, part& p)
148  {
149  return is >> p.area >> p.centre;
150  }
151  };
152 
153  //- Structure to store the geometry associated with the coupling
154  // between parts of two patch faces
155  struct couple
156  :
157  public part
158  {
159  //- The neighbour part
160  part nbr;
161 
162  //- Default construct
163  couple()
164  :
165  part(),
166  nbr()
167  {}
168 
169  //- Construct from a coupled pair of parts
170  couple(const part& p, const part& nbrP)
171  :
172  part(p),
173  nbr(nbrP)
174  {}
175 
176  //- Equality comparison
177  friend bool operator==(const couple& a, const couple& b)
178  {
179  return
180  static_cast<const part&>(a) == static_cast<const part&>(b)
181  && a.nbr == b.nbr;
182  }
183 
184  //- Inequality comparison
185  friend bool operator!=(const couple& a, const couple& b)
186  {
187  return !(a == b);
188  }
189 
190  //- Output stream operator
191  friend Ostream& operator<<(Ostream& os, const couple& c)
192  {
193  return os << static_cast<const part&>(c) << c.nbr;
194  }
195 
196  //- Input stream operator
197  friend Istream& operator>>(Istream& is, couple& c)
198  {
199  return is >> static_cast<part&>(c) >> c.nbr;
200  }
201  };
202 
203 
204 private:
205 
206  // Private Member Data
207 
208  // Geometry
209 
210  //- The coupling geometry for for each source face
211  List<DynamicList<couple>> srcCouples_;
212 
213  //- The non-coupled geometry associated with each source edge
214  List<part> srcEdgeParts_;
215 
216  //- The non-coupled geometry associated with mismatch in each
217  // source face's couplings
218  List<part> srcErrorParts_;
219 
220  //- The coupling geometry for for each target face
221  List<DynamicList<couple>> tgtCouples_;
222 
223 
224  //- Triangulation engine
225  mutable polygonTriangulate triEngine_;
226 
227 
228  // Workspace
229 
230  //- Source face triangulation points
231  mutable List<triFaceList> srcTriPoints_;
232 
233  //- Source face triangulation edges
234  mutable List<List<FixedList<label, 3>>> srcTriFaceEdges_;
235 
236  //- Target face triangulation points
237  mutable List<triFaceList> tgtTriPoints_;
238 
239  //- Target face triangulation edges
240  mutable List<List<FixedList<label, 3>>> tgtTriFaceEdges_;
241 
242  //- Source intersection points
243  mutable DynamicList<point> ictSrcPoints_;
244 
245  //- Source intersection point normals
246  mutable DynamicList<vector> ictSrcPointNormals_;
247 
248  //- Target intersection points
249  mutable DynamicList<point> ictTgtPoints_;
250 
251  //- Intersection locations
252  mutable DynamicList<triIntersect::location> ictPointLocations_;
253 
254  //- Source face edge parts
255  DynamicList<part> srcFaceEdgePart_;
256 
257  //- Target face edge parts
258  DynamicList<part> tgtFaceEdgePart_;
259 
260  //- Source face edge parts
261  List<List<part>> srcFaceEdgeParts_;
262 
263 
264  // Debugging
265 
266  //- Intersection points
267  mutable DynamicList<point> debugPoints_;
268 
269  //- Intersection faces
270  mutable DynamicList<labelList> debugFaces_;
271 
272  //- The source face associated with each intersection face
273  mutable DynamicList<label> debugFaceSrcFaces_;
274 
275  //- The target face associated with each intersection face
276  mutable DynamicList<label> debugFaceTgtFaces_;
277 
278  //- The side of the intersection each face is on; 1 is the source
279  // side, -1 is the target side, and 0 is a face that spans the
280  // intersection volume; e.g., an edge or error face.
281  mutable DynamicList<label> debugFaceSides_;
282 
283 
284  // Private Static Member Functions
285 
286  //- Return the values at the points of a tri face
287  template<class Type>
288  static FixedList<Type, 3> triPointValues
289  (
290  const triFace& t,
291  const UList<Type>& values
292  );
293 
294 
295  // Private Member Functions
296 
297  //- Get the bound box for a source face
298  virtual treeBoundBox srcBox
299  (
300  const face& srcFace,
301  const pointField& srcPoints,
302  const vectorField& srcPointNormals
303  ) const;
304 
305  //- Intersect two faces
306  virtual bool intersectFaces
307  (
308  const primitiveOldTimePatch& srcPatch,
309  const vectorField& srcPointNormals,
310  const vectorField& srcPointNormals0,
311  const primitiveOldTimePatch& tgtPatch,
312  const label srcFacei,
313  const label tgtFacei
314  );
315 
316  //- Initialise the workspace
317  virtual void initialise
318  (
319  const primitiveOldTimePatch& srcPatch,
320  const vectorField& srcPointNormals,
321  const vectorField& srcPointNormals0,
322  const primitiveOldTimePatch& tgtPatch
323  );
324 
325  //- Send data that resulted from an intersection between the source
326  // patch and a distributed source-local-target patch back to the
327  // original target processes.
328  virtual void rDistributeTgt
329  (
330  const primitiveOldTimePatch& tgtPatch,
331  const distributionMap& tgtMap
332  );
333 
334  //- Finalise the intersection
335  virtual label finalise
336  (
337  const primitiveOldTimePatch& srcPatch,
338  const vectorField& srcPointNormals,
339  const vectorField& srcPointNormals0,
340  const primitiveOldTimePatch& tgtPatch,
341  const transformer& tgtToSrc
342  );
343 
344 
345 public:
346 
347  //- Runtime type information
348  TypeName("intersection");
349 
350 
351  // Static Data Members
352 
353  //- Extra debugging for intersections between specific faces. Named
354  // "intersectionSrcFace" and "intersectionTgtFace" respectively.
355  static int debugSrcFacei, debugTgtFacei;
356 
357 
358  // Static Member Functions
359 
360  //- Get the bound box for a source face
362  (
363  const face& srcFace,
364  const pointField& srcPoints,
365  const vectorField& srcPointNormals
366  );
367 
368 
369  // Constructors
370 
371  //- Construct from components
372  intersection(const bool reverse);
373 
374 
375  //- Destructor
376  ~intersection();
377 
378 
379  // Member Functions
380 
381  //- For each source face, the source and target areas for each
382  // target coupling
383  inline const List<DynamicList<couple>>& srcCouples() const;
384 
385  //- For each source edge, the non-coupled geometry associated
386  // with its projection
387  inline const List<part>& srcEdgeParts() const;
388 
389  //- For each source face, the area associated with mismatch
390  // across the coupling
391  inline const List<part>& srcErrorParts() const;
392 
393  //- For each target face, the target and source areas for each
394  // source coupling
395  inline const List<DynamicList<couple>>& tgtCouples() const;
396 
397  //- For each source face, the coupled target weights
399  (
400  const primitivePatch& srcPatch
401  ) const;
402 
403  //- For each target face, the coupled source weights
405  (
406  const primitivePatch& tgtPatch
407  ) const;
408 };
409 
410 
411 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
412 
413 } // End namespace patchToPatches
414 } // End namespace Foam
415 
416 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
417 
419 
420 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
421 
422 #endif
423 
424 // ************************************************************************* //
const List< DynamicList< couple > > & srcCouples() const
For each source face, the source and target areas for each.
virtual tmpNrc< List< DynamicList< scalar > > > tgtWeights(const primitivePatch &tgtPatch) const
For each target face, the coupled source weights.
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
const List< part > & srcErrorParts() const
For each source face, the area associated with mismatch.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Class to encapsulate the topology of a point within a triangle intersection.
void operator-=(const part &p)
Subtract another part from this one.
const List< DynamicList< couple > > & tgtCouples() const
For each target face, the target and source areas for each.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
bool reverse() const
Flag to indicate that the two patches are co-directional and.
Definition: patchToPatchI.H:87
TypeName("intersection")
Runtime type information.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
friend bool operator!=(const part &a, const part &b)
Inequality comparison.
Structure to store the geometry associated with part of a patch.
Class to generate patchToPatch coupling geometry. A full geometric intersection is done between a fac...
const dimensionedScalar c
Speed of light in a vacuum.
friend bool operator==(const part &a, const part &b)
Equality comparison.
static int debugSrcFacei
Extra debugging for intersections between specific faces. Named.
A list of faces which address into the list of points.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
static const zero Zero
Definition: zero.H:97
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void operator+=(const part &p)
Add another part to this one.
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:52
virtual tmpNrc< List< DynamicList< scalar > > > srcWeights(const primitivePatch &srcPatch) const
For each source face, the coupled target weights.
Structure to store the geometry associated with the coupling.
Class containing processor-to-processor mapping information.
const List< part > & srcEdgeParts() const
For each source edge, the non-coupled geometry associated.
Triangulation of three-dimensional polygons.
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:49
dimensioned< scalar > mag(const dimensioned< Type > &)
static treeBoundBox srcBoxStatic(const face &srcFace, const pointField &srcPoints, const vectorField &srcPointNormals)
Get the bound box for a source face.
intersection(const bool reverse)
Construct from components.
volScalarField & p
Class to generate coupling geometry between two primitive patches.
Definition: patchToPatch.H:53
friend Istream & operator>>(Istream &is, part &p)
Input stream operator.
Namespace for OpenFOAM.
friend Ostream & operator<<(Ostream &os, const part &p)
Output stream operator.