patchToPatch.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) 2011-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::patchToPatch
26 
27 Description
28  Class to generate coupling geometry between two primitive patches
29 
30 SourceFiles
31  patchToPatch.C
32  patchToPatchParallelOps.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef patchToPatch_H
37 #define patchToPatch_H
38 
39 #include "distributionMap.H"
40 #include "primitivePatch.H"
41 #include "primitiveOldTimePatch.H"
42 #include "runTimeSelectionTables.H"
43 #include "treeBoundBox.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class patchToPatch Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 class patchToPatch
55 {
56 public:
57 
58  // Public Structures
59 
60  //- Structure to conveniently store processor and face indices
61  struct procFace
62  {
63  //- The processor index
64  label proci;
65 
66  //- The face index
67  label facei;
68 
69  //- Equality comparison
70  friend bool operator==(const procFace& a, const procFace& b)
71  {
72  return a.proci == b.proci && a.facei == b.facei;
73  }
74 
75  //- Inequality comparison
76  friend bool operator!=(const procFace& a, const procFace& b)
77  {
78  return !(a == b);
79  }
80 
81  //- Output stream operator
82  friend Ostream& operator<<(Ostream& os, const procFace& p)
83  {
84  return os << p.proci << token::SPACE << p.facei;
85  }
86 
87  //- Input stream operator
88  friend Istream& operator>>(Istream& is, procFace& p)
89  {
90  return is >> p.proci >> p.facei;
91  }
92  };
93 
94 
95 protected:
96 
97  // Private Data
98 
99  //- Flag to indicate that the two patches are co-directional and
100  // that the orientation of one should therefore be reversed
101  const bool reverse_;
102 
103  //- Index of the processor holding all faces of the patchToPatch, or -1
104  // if spread across multiple processors
106 
107  //- When running in parallel, a map from local source face index to
108  // source processor and face index
110 
111  //- When running in parallel, a map from local target face index to
112  // target processor and face index
114 
115  //- For each source face, the coupled local target faces
117 
118  //- For each target face, the coupled local source faces
120 
121 
122  // Private Member Functions
123 
124  // Indexing
125 
126  //- Transfer list-list b into list-list a
127  template<class SubListA, class SubListB>
128  static inline void transferListList
129  (
130  List<SubListA>& a,
132  );
133 
134  //- Reverse distribute a list-list given the map
135  template<class Type>
136  static inline void rDistributeListList
137  (
138  const label size,
139  const distributionMap& map,
141  );
142 
143  //- Reverse distribute a dynamic list-list given the map
144  template<class Type>
145  static inline void rDistributeListList
146  (
147  const label size,
148  const distributionMap& map,
150  );
151 
152  //- Map local faces to proc faces
154  (
155  const List<DynamicList<label>>& localFaces,
157  );
158 
159  //- Map proc faces to local faces
161  (
162  const List<List<procFace>>& procFaces,
163  const HashTable<label, procFace, Hash<procFace>>& map
164  );
165 
166 
167  // Searching
168 
169  //- Get the bound box for a source face
170  virtual treeBoundBox srcBox
171  (
172  const face& srcFace,
173  const pointField& srcPoints,
174  const vectorField& srcPointNormals
175  ) const = 0;
176 
177  //- Get the bound box for a source face
179  (
180  const primitiveOldTimePatch& srcPatch,
181  const vectorField& srcPointNormals,
182  const vectorField& srcPointNormals0,
183  const label srcFacei
184  ) const;
185 
186  //- Get the bound box for the source patch
188  (
189  const primitiveOldTimePatch& srcPatch,
190  const vectorField& srcPointNormals,
191  const vectorField& srcPointNormals0
192  ) const;
193 
194  //- Get the bound box for a target face
196  (
197  const primitiveOldTimePatch& tgtPatch,
198  const label tgtFacei
199  ) const;
200 
201  //- Get the bound box for the target patch
203  (
204  const primitiveOldTimePatch& tgtPatch
205  ) const;
206 
207 
208  // Intersection
209 
210  //- Intersect two faces
211  virtual bool intersectFaces
212  (
213  const primitiveOldTimePatch& srcPatch,
214  const vectorField& srcPointNormals,
215  const vectorField& srcPointNormals0,
216  const primitiveOldTimePatch& tgtPatch,
217  const label srcFacei,
218  const label tgtFacei
219  ) = 0;
220 
221  //- Intersect two faces
223  (
224  const primitiveOldTimePatch& srcPatch,
225  const vectorField& srcPointNormals,
226  const vectorField& srcPointNormals0,
227  const primitiveOldTimePatch& tgtPatch,
228  const label srcFacei,
229  const label tgtFacei
230  );
231 
232  //- Intersect a queue of source-target face pairs. Update completion
233  // lists and form a new queue of target-source face pairs.
235  (
236  const primitiveOldTimePatch& srcPatch,
237  const vectorField& srcPointNormals,
238  const vectorField& srcPointNormals0,
239  const primitiveOldTimePatch& tgtPatch,
240  const bool isSrc,
241  const DynamicList<labelPair>& queue,
242  labelList& faceComplete,
243  DynamicList<labelPair>& otherQueue,
244  const labelList& otherFaceComplete,
245  boolList& otherFaceQueued,
246  boolList& otherFaceVisited
247  );
248 
249  //- Intersect the patches
250  void intersectPatches
251  (
252  const primitiveOldTimePatch& srcPatch,
253  const vectorField& srcPointNormals,
254  const vectorField& srcPointNormals0,
255  const primitiveOldTimePatch& tgtPatch
256  );
257 
258 
259  // Parallel functionality
260 
261  //- Determine whether or not the intersection of the given patches
262  // is on a single process. Set singleProcess_ to that process if
263  // so, and to -1 if not.
264  void calcSingleProcess
265  (
266  const primitiveOldTimePatch& srcPatch,
267  const primitiveOldTimePatch& tgtPatch
268  );
269 
270  //- Determine which target faces need to be sent to the source.
271  // This is done before intersection. Bound boxes are used to
272  // estimate what faces will intersect.
274  (
275  const primitiveOldTimePatch& srcPatch,
276  const vectorField& srcPointNormals,
277  const vectorField& srcPointNormals0,
278  const primitiveOldTimePatch& tgtPatch
279  ) const;
280 
281  //- Determine which source faces need to be sent to the target.
282  // This is done after intersection, using the addressing generated
283  // by the intersection.
285 
286  //- Create a distribution map from a list-list of faces to be sent
287  // (i.e., the result of calling one of the above methods).
289  (
290  labelListList&& sendFaces
291  ) const;
292 
293  //- Distribute a patch given its distribution map. Just generate
294  // local-proc-face addressing; not the distributed patch itself.
295  void distributePatch
296  (
297  const distributionMap& map,
298  List<procFace>& localProcFaces
299  ) const;
300 
301  //- Distribute a patch given its distribution map
303  (
304  const distributionMap& map,
305  const primitiveOldTimePatch& patch,
306  List<procFace>& localProcFaces
307  ) const;
308 
309 
310  // Hooks
311 
312  //- Initialisation
313  virtual void initialise
314  (
315  const primitiveOldTimePatch& srcPatch,
316  const vectorField& srcPointNormals,
317  const vectorField& srcPointNormals0,
318  const primitiveOldTimePatch& tgtPatch
319  );
320 
321  //- Distribute the target patch so that enough is locally available
322  // for its intersection with the source patch can be computed.
323  // Happens before intersection. Bound boxes are used to determine
324  // what is needed where. Return the resulting local patch and the
325  // map from which it was calculated.
326  virtual
329  (
330  const primitiveOldTimePatch& srcPatch,
331  const vectorField& srcPointNormals,
332  const vectorField& srcPointNormals0,
333  const primitiveOldTimePatch& tgtPatch,
334  distributionMap& tgtMap
335  );
336 
337  //- Distribute the source patch so that everything the target
338  // intersects is locally available. Happens after intersection.
339  // The intersection addressing is used to determine what is needed
340  // where. Return the resulting local patch and the map from which
341  // it was calculated.
342  virtual
345  (
346  const primitiveOldTimePatch& srcPatch,
347  distributionMap& srcMap
348  );
349 
350  //- Send data that resulted from an intersection between the source
351  // patch and a distributed source-local-target patch back to the
352  // original target processes.
353  virtual void rDistributeTgt
354  (
355  const primitiveOldTimePatch& tgtPatch,
356  const distributionMap& tgtMap
357  );
358 
359  //- Finalising
360  virtual label finalise
361  (
362  const primitiveOldTimePatch& srcPatch,
363  const vectorField& srcPointNormals,
364  const vectorField& srcPointNormals0,
365  const primitiveOldTimePatch& tgtPatch,
366  const transformer& tgtToSrc
367  );
368 
369 
370 public:
371 
372  //- Runtime type information
373  TypeName("patchToPatch");
374 
375 
376  //- Declare runtime constructor selection table
378  (
379  autoPtr,
380  patchToPatch,
381  bool,
382  (const bool reverse),
383  (reverse)
384  );
385 
386 
387  // Constructors
388 
389  //- Construct from components
390  patchToPatch(const bool reverse);
391 
392  //- Disallow default bitwise copy construction
393  patchToPatch(const patchToPatch&) = delete;
394 
395 
396  //- Destructor
397  virtual ~patchToPatch();
398 
399 
400  // Selector
401 
402  //- Select from name
404  (
405  const word& patchToPatchType,
406  const bool reverse
407  );
408 
409 
410  // Member Functions
411 
412  // Access
413 
414  //- Flag to indicate that the two patches are co-directional and
415  // that the orientation of one should therefore be reversed
416  inline bool reverse() const;
417 
418  //- Index of the processor holding all faces of the patchToPatch,
419  // or -1 if spread across multiple processors
420  inline label singleProcess() const;
421 
422  //- Is this intersection on a single process?
423  inline bool isSingleProcess() const;
424 
425  //- For each source face, the coupled target procs and faces
426  inline List<List<procFace>> srcTgtProcFaces() const;
427 
428  //- For each target face, the coupled source procs and faces
429  inline List<List<procFace>> tgtSrcProcFaces() const;
430 
431  //- For each source face, the coupled target weights
433  (
434  const primitivePatch& srcPatch
435  ) const = 0;
436 
437  //- For each target face, the coupled source weights
439  (
440  const primitivePatch& tgtPatch
441  ) const = 0;
442 
443 
444  // Manipulation
445 
446  //- Update addressing and weights for the given patches
447  void update
448  (
449  const primitiveOldTimePatch& srcPatch,
450  const vectorField& srcPointNormals,
451  const vectorField& srcPointNormals0,
452  const primitiveOldTimePatch& tgtPatch,
453  const transformer& tgtToSrc = NullObjectRef<transformer>()
454  );
455 
456  //- Update addressing and weights for the given patches
457  void update
458  (
459  const primitivePatch& srcPatch,
460  const vectorField& srcPointNormals,
461  const primitivePatch& tgtPatch,
462  const transformer& tgtToSrc = NullObjectRef<transformer>()
463  );
464 
465 
466  // Member Operators
467 
468  //- Disallow default bitwise assignment
469  void operator=(const patchToPatch&) = delete;
470 };
471 
472 
473 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
474 
475 } // End namespace Foam
476 
477 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
478 
479 #include "patchToPatchI.H"
480 
481 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
482 
483 #endif
484 
485 // ************************************************************************* //
virtual treeBoundBox srcBox(const face &srcFace, const pointField &srcPoints, const vectorField &srcPointNormals) const =0
Get the bound box for a source face.
List< List< procFace > > srcTgtProcFaces() const
For each source face, the coupled target procs and faces.
const bool reverse_
Flag to indicate that the two patches are co-directional and.
Definition: patchToPatch.H:100
virtual void rDistributeTgt(const primitiveOldTimePatch &tgtPatch, const distributionMap &tgtMap)
Send data that resulted from an intersection between the source.
Definition: patchToPatch.C:775
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
patchToPatch(const bool reverse)
Construct from components.
Definition: patchToPatch.C:831
virtual tmpNrc< List< DynamicList< scalar > > > srcWeights(const primitivePatch &srcPatch) const =0
For each source face, the coupled target weights.
List< List< procFace > > tgtSrcProcFaces() const
For each target face, the coupled source procs and faces.
Structure to conveniently store processor and face indices.
Definition: patchToPatch.H:60
friend Istream & operator>>(Istream &is, procFace &p)
Input stream operator.
Definition: patchToPatch.H:87
label facei
The face index.
Definition: patchToPatch.H:66
label singleProcess_
Index of the processor holding all faces of the patchToPatch, or -1.
Definition: patchToPatch.H:104
friend bool operator!=(const procFace &a, const procFace &b)
Inequality comparison.
Definition: patchToPatch.H:75
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
label proci
The processor index.
Definition: patchToPatch.H:63
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
void intersectPatches(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
Intersect the patches.
Definition: patchToPatch.C:425
List< DynamicList< label > > srcLocalTgtFaces_
For each source face, the coupled local target faces.
Definition: patchToPatch.H:115
TypeName("patchToPatch")
Runtime type information.
static List< DynamicList< label > > procFacesToLocalFaces(const List< List< procFace >> &procFaces, const HashTable< label, procFace, Hash< procFace >> &map)
Map proc faces to local faces.
Definition: patchToPatch.C:123
virtual tmpNrc< PrimitiveOldTimePatch< faceList, pointField > > distributeSrc(const primitiveOldTimePatch &srcPatch, distributionMap &srcMap)
Distribute the source patch so that everything the target.
Definition: patchToPatch.C:751
static List< List< procFace > > localFacesToProcFaces(const List< DynamicList< label >> &localFaces, const List< procFace > &map=NullObjectRef< List< procFace >>())
Map local faces to proc faces.
Definition: patchToPatch.C:97
static void rDistributeListList(const label size, const distributionMap &map, List< List< Type >> &data)
Reverse distribute a list-list given the map.
Definition: patchToPatchI.H:47
label intersectPatchQueue(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const bool isSrc, const DynamicList< labelPair > &queue, labelList &faceComplete, DynamicList< labelPair > &otherQueue, const labelList &otherFaceComplete, boolList &otherFaceQueued, boolList &otherFaceVisited)
Intersect a queue of source-target face pairs. Update completion.
Definition: patchToPatch.C:283
A list of faces which address into the list of points.
static autoPtr< patchToPatch > New(const word &patchToPatchType, const bool reverse)
Select from name.
Definition: patchToPatch.C:851
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A class for handling words, derived from string.
Definition: word.H:59
virtual tmpNrc< PrimitiveOldTimePatch< faceList, pointField > > distributeTgt(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, distributionMap &tgtMap)
Distribute the target patch so that enough is locally available.
Definition: patchToPatch.C:713
declareRunTimeSelectionTable(autoPtr, patchToPatch, bool,(const bool reverse),(reverse))
Declare runtime constructor selection table.
virtual ~patchToPatch()
Destructor.
Definition: patchToPatch.C:844
An STL-conforming hash table.
Definition: HashTable.H:61
virtual bool intersectFaces(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const label srcFacei, const label tgtFacei)=0
Intersect two faces.
distributionMap patchDistributionMap(labelListList &&sendFaces) const
Create a distribution map from a list-list of faces to be sent.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
bool findOrIntersectFaces(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const label srcFacei, const label tgtFacei)
Intersect two faces.
Definition: patchToPatch.C:233
Database for solution and other reduced data.
Definition: data.H:51
void update(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const transformer &tgtToSrc=NullObjectRef< transformer >())
Update addressing and weights for the given patches.
Definition: patchToPatch.C:876
friend bool operator==(const procFace &a, const procFace &b)
Equality comparison.
Definition: patchToPatch.H:69
labelListList srcPatchSendFaces() const
Determine which source faces need to be sent to the target.
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:52
virtual tmpNrc< List< DynamicList< scalar > > > tgtWeights(const primitivePatch &tgtPatch) const =0
For each target face, the coupled source weights.
virtual label finalise(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const transformer &tgtToSrc)
Finalising.
Definition: patchToPatch.C:804
Class containing processor-to-processor mapping information.
autoPtr< List< procFace > > localSrcProcFacesPtr_
When running in parallel, a map from local source face index to.
Definition: patchToPatch.H:108
void calcSingleProcess(const primitiveOldTimePatch &srcPatch, const primitiveOldTimePatch &tgtPatch)
Determine whether or not the intersection of the given patches.
List< DynamicList< label > > tgtLocalSrcFaces_
For each target face, the coupled local source faces.
Definition: patchToPatch.H:118
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:52
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
treeBoundBox tgtBox(const primitiveOldTimePatch &tgtPatch, const label tgtFacei) const
Get the bound box for a target face.
Definition: patchToPatch.C:195
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
label singleProcess() const
Index of the processor holding all faces of the patchToPatch,.
Definition: patchToPatchI.H:93
Macros to ease declaration of run-time selection tables.
virtual void initialise(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
Initialisation.
Definition: patchToPatch.C:690
friend Ostream & operator<<(Ostream &os, const procFace &p)
Output stream operator.
Definition: patchToPatch.H:81
volScalarField & p
void operator=(const patchToPatch &)=delete
Disallow default bitwise assignment.
const T & NullObjectRef()
Return reference to the nullObject of type T.
Definition: nullObjectI.H:27
static void transferListList(List< SubListA > &a, List< SubListB > &b)
Transfer list-list b into list-list a.
Definition: patchToPatchI.H:32
Class to generate coupling geometry between two primitive patches.
Definition: patchToPatch.H:53
autoPtr< List< procFace > > localTgtProcFacesPtr_
When running in parallel, a map from local target face index to.
Definition: patchToPatch.H:112
Namespace for OpenFOAM.
void distributePatch(const distributionMap &map, List< procFace > &localProcFaces) const
Distribute a patch given its distribution map. Just generate.
bool isSingleProcess() const
Is this intersection on a single process?
Definition: patchToPatchI.H:99
labelListList tgtPatchSendFaces(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch) const
Determine which target faces need to be sent to the source.