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