raysPatchToPatch.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "raysPatchToPatch.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace patchToPatches
35 {
36  defineTypeNameAndDebug(rays, 0);
37  addToRunTimeSelectionTable(patchToPatch, rays, bool);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::treeBoundBox Foam::patchToPatches::rays::srcBox
45 (
46  const face& srcFace,
47  const pointField& srcPoints,
48  const vectorField& srcPointNormals
49 ) const
50 {
51  return intersection::srcBoxStatic(srcFace, srcPoints, srcPointNormals);
52 }
53 
54 
55 bool Foam::patchToPatches::rays::intersectFaces
56 (
57  const primitiveOldTimePatch& srcPatch,
58  const vectorField& srcPointNormals,
59  const vectorField& srcPointNormals0,
60  const primitiveOldTimePatch& tgtPatch,
61  const label srcFacei,
62  const label tgtFacei
63 )
64 {
65  const treeBoundBox srcFaceBox =
67  (
68  srcPatch,
69  srcPointNormals,
70  srcPointNormals0,
71  srcFacei
72  );
73  const treeBoundBox tgtFaceBox =
74  patchToPatch::tgtBox(tgtPatch, tgtFacei);
75 
76  if (srcFaceBox.overlaps(tgtFaceBox))
77  {
78  srcLocalTgtFaces_[srcFacei].append(tgtFacei);
79  tgtLocalSrcFaces_[tgtFacei].append(srcFacei);
80 
81  return true;
82  }
83  else
84  {
85  return false;
86  }
87 }
88 
89 
91 Foam::patchToPatches::rays::distributeTgt
92 (
93  const primitiveOldTimePatch& srcPatch,
94  const vectorField& srcPointNormals,
95  const vectorField& srcPointNormals0,
96  const primitiveOldTimePatch& tgtPatch,
97  distributionMap& tgtMap
98 )
99 {
100  // Intercept generation of the local patch. Store it. Return a const
101  // reference tmp instead of a pointer.
102 
103  tmpNrc<PrimitiveOldTimePatch<faceList, pointField>> localTgtPatchPtr =
105  (
106  srcPatch,
107  srcPointNormals,
108  srcPointNormals0,
109  tgtPatch,
110  tgtMap
111  );
112 
113  localTgtPatchPtr_.reset(localTgtPatchPtr.ptr());
114 
115  return
116  tmpNrc<PrimitiveOldTimePatch<faceList, pointField>>
117  (
118  localTgtPatchPtr_()
119  );
120 }
121 
122 
124 Foam::patchToPatches::rays::distributeSrc
125 (
126  const primitiveOldTimePatch& srcPatch,
127  distributionMap& srcMap
128 )
129 {
130  // Intercept generation of the local patch. Store it. Return a const
131  // reference tmp instead of a pointer.
132 
133  tmpNrc<PrimitiveOldTimePatch<faceList, pointField>> localSrcPatchPtr =
134  patchToPatch::distributeSrc(srcPatch, srcMap);
135 
136  localSrcPatchPtr_.reset(localSrcPatchPtr.ptr());
137 
138  return
139  tmpNrc<PrimitiveOldTimePatch<faceList, pointField>>
140  (
141  localSrcPatchPtr_()
142  );
143 }
144 
145 
146 Foam::label Foam::patchToPatches::rays::finalise
147 (
148  const primitiveOldTimePatch& srcPatch,
149  const vectorField& srcPointNormals,
150  const vectorField& srcPointNormals0,
151  const primitiveOldTimePatch& tgtPatch,
152  const transformer& tgtToSrc
153 )
154 {
155  const label nCouples =
157  (
158  srcPatch,
159  srcPointNormals,
160  srcPointNormals0,
161  tgtPatch,
162  tgtToSrc
163  );
164 
165  // Transform the source-local target patch back to the target
166  if (!isSingleProcess() && !isNull(tgtToSrc))
167  {
168  autoPtr<PrimitiveOldTimePatch<faceList, pointField>>
169  localTgtPatchPtr(localTgtPatchPtr_.ptr());
170 
171  localTgtPatchPtr_.set
172  (
173  new PrimitiveOldTimePatch<faceList, pointField>
174  (
175  localTgtPatchPtr(),
176  tgtToSrc.invTransformPosition(localTgtPatchPtr().points()),
177  isNull(localTgtPatchPtr().points0())
178  ? NullObjectRef<pointField>()
179  : tgtToSrc.invTransformPosition(localTgtPatchPtr().points0())()
180  )
181  );
182  }
183 
184  return nCouples;
185 }
186 
187 
188 Foam::patchToPatch::procFace Foam::patchToPatches::rays::ray
189 (
190  const primitiveOldTimePatch& outPatch,
191  const autoPtr<PrimitiveOldTimePatch<faceList, pointField>>&
192  localOutPatchPtr,
193  const autoPtr<List<procFace>>& localOutProcFacesPtr,
194  const List<DynamicList<label>>& inLocalOutFaces,
195  const scalar fraction,
196  const label inFacei,
197  const point& inP,
198  const vector& inN,
199  point& outP
200 ) const
201 {
202  forAll(inLocalOutFaces[inFacei], i)
203  {
204  const label localOutFacei = inLocalOutFaces[inFacei][i];
205 
206  const face& outF =
207  isSingleProcess()
208  ? outPatch[localOutFacei]
209  : localOutPatchPtr()[localOutFacei];
210 
211  const pointField& outPoints =
212  isSingleProcess()
213  ? outPatch.points()
214  : localOutPatchPtr().points();
215  const pointField& outPoints0 =
216  isSingleProcess()
217  ? outPatch.points0()
218  : localOutPatchPtr().points0();
219 
220  const pointField outPoly
221  (
222  (1 - fraction)*outF.points(outPoints0)
223  + fraction*outF.points(outPoints)
224  );
225 
226  const pointHit ray =
227  face(identity(outPoly.size()))
228  .ray(inP, inN, outPoly, Foam::intersection::algorithm::visible);
229 
230  if (ray.hit())
231  {
232  outP = ray.rawPoint();
233 
234  return
235  isSingleProcess()
236  ? procFace({Pstream::myProcNo(), localOutFacei})
237  : localOutProcFacesPtr()[localOutFacei];
238  }
239  }
240 
241  return procFace({-1, -1});
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
246 
248 (
249  const bool reverse
250 )
251 :
252  patchToPatch(reverse)
253 {}
254 
255 
256 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
257 
259 {}
260 
261 
262 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
263 
266 {
268  return tmpNrc<List<DynamicList<scalar>>>(nullptr);
269 }
270 
271 
274 {
276  return tmpNrc<List<DynamicList<scalar>>>(nullptr);
277 }
278 
279 
281 (
282  const primitiveOldTimePatch& tgtPatch,
283  const scalar fraction,
284  const label srcFacei,
285  const vector& srcP,
286  const vector& srcN,
287  point& tgtP
288 ) const
289 {
290  return
291  ray
292  (
293  tgtPatch,
294  localTgtPatchPtr_,
295  localTgtProcFacesPtr_,
296  srcLocalTgtFaces_,
297  fraction,
298  srcFacei,
299  srcP,
300  srcN,
301  tgtP
302  );
303 }
304 
305 
307 (
308  const primitiveOldTimePatch& srcPatch,
309  const scalar fraction,
310  const label tgtFacei,
311  const vector& tgtP,
312  const vector& tgtN,
313  point& srcP
314 ) const
315 {
316  return
317  ray
318  (
319  srcPatch,
320  localSrcPatchPtr_,
321  localSrcProcFacesPtr_,
322  tgtLocalSrcFaces_,
323  fraction,
324  tgtFacei,
325  tgtP,
326  tgtN,
327  srcP
328  );
329 }
330 
331 
332 // ************************************************************************* //
addToRunTimeSelectionTable(patchToPatch, intersection, bool)
virtual treeBoundBox srcBox(const face &srcFace, const pointField &srcPoints, const vectorField &srcPointNormals) const =0
Get the bound box for a source face.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Structure to conveniently store processor and face indices.
Definition: patchToPatch.H:60
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
virtual tmpNrc< List< DynamicList< scalar > > > tgtWeights(const primitivePatch &tgtPatch) const
For each target face, the coupled source weights.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Macros for easy insertion into run-time selection tables.
virtual tmpNrc< PrimitiveOldTimePatch< faceList, pointField > > distributeSrc(const primitiveOldTimePatch &srcPatch, distributionMap &srcMap)
Distribute the source patch so that everything the target.
Definition: patchToPatch.C:751
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:46
defineTypeNameAndDebug(intersection, 0)
const pointField & points
PrimitiveOldTimePatch< SubList< face >, const pointField & > primitiveOldTimePatch
Addressing for a faceList slice.
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
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
procFace tgtToSrcRay(const primitiveOldTimePatch &srcPatch, const scalar fraction, const label tgtFacei, const vector &tgtP, const vector &tgtN, point &srcP) const
Compute a ray intersection from the target side to the source.
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:52
virtual label finalise(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const transformer &tgtToSrc)
Finalising.
Definition: patchToPatch.C:804
vector point
Point is a vector.
Definition: point.H:41
virtual tmpNrc< List< DynamicList< scalar > > > srcWeights(const primitivePatch &srcPatch) const
For each source face, the coupled target weights.
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
Field< vector > vectorField
Specialisation of Field<T> for vector.
procFace srcToTgtRay(const primitiveOldTimePatch &tgtPatch, const scalar fraction, const label srcFacei, const vector &srcP, const vector &srcN, point &tgtP) const
Compute a ray intersection from the source side to the target.
static treeBoundBox srcBoxStatic(const face &srcFace, const pointField &srcPoints, const vectorField &srcPointNormals)
Get the bound box for a source face.
PointHit< point > pointHit
Definition: pointHit.H:41
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
Class to generate coupling geometry between two primitive patches.
Definition: patchToPatch.H:53
rays(const bool reverse)
Construct from components.
Namespace for OpenFOAM.