nearestPatchToPatch.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 "nearestPatchToPatch.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace patchToPatches
34 {
35  defineTypeNameAndDebug(nearest, 0);
36  addToRunTimeSelectionTable(patchToPatch, nearest, bool);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 Foam::treeBoundBox Foam::patchToPatches::nearest::srcBox
44 (
45  const face& srcFace,
46  const pointField& srcPoints,
47  const vectorField& srcPointNormals
48 ) const
49 {
50  const treeBoundBox bb(srcPoints, srcFace);
51 
52  const point c = bb.midpoint();
53  const scalar l = bb.maxDim();
54 
55  return treeBoundBox(c - l*vector::one, c + l*vector::one);
56 }
57 
58 
59 bool Foam::patchToPatches::nearest::intersectFaces
60 (
61  const primitivePatch& patch,
62  const primitivePatch& otherPatch,
63  const label facei,
64  const label otherFacei,
65  DynamicList<label>& faceOtherFaces,
66  scalar& faceDistance
67 ) const
68 {
69  auto closest = [&patch,&otherPatch]
70  (
71  const label facei,
72  const label otherFacei
73  )
74  {
75  const point& c = patch.faceCentres()[facei];
76  const point& otherC = otherPatch.faceCentres()[otherFacei];
77  const scalar distSqr = magSqr(c - otherC);
78 
79  forAll(otherPatch.faceEdges()[otherFacei], otherFaceEdgei)
80  {
81  const label otherEdgei =
82  otherPatch.faceEdges()[otherFacei][otherFaceEdgei];
83 
84  point otherNbrC;
85 
86  if (otherPatch.edgeFaces()[otherEdgei].size() == 2)
87  {
88  const label facej =
89  otherPatch.edgeFaces()[otherEdgei]
90  [otherPatch.edgeFaces()[otherEdgei][0] == otherFacei];
91 
92  otherNbrC = otherPatch.faceCentres()[facej];
93  }
94  else
95  {
96  const edge& e = otherPatch.edges()[otherEdgei];
97  const point& p = otherPatch.localPoints()[e[0]];
98  const vector dp = e.vec(otherPatch.localPoints());
99  const vector n = otherPatch.faceNormals()[otherFacei] ^ dp;
100 
101  otherNbrC = p + ((tensor::I - 2*sqr(n)) & (otherC - p));
102  }
103 
104  if (magSqr(c - otherNbrC) < distSqr)
105  {
106  return false;
107  }
108  }
109 
110  return true;
111  };
112 
113  if (closest(facei, otherFacei))
114  {
115  const point& c = patch.faceCentres()[facei];
116  const point& otherC = otherPatch.faceCentres()[otherFacei];
117  const scalar distSqr = magSqr(c - otherC);
118 
119  if (faceOtherFaces.empty() || faceDistance > distSqr)
120  {
121  faceOtherFaces.clear();
122  faceOtherFaces.append(otherFacei);
123  faceDistance = distSqr;
124  }
125 
126  return true;
127  }
128 
129  const labelList& otherFaceFaces = otherPatch.faceFaces()[otherFacei];
130  forAll(otherFaceFaces, otherFaceFacei)
131  {
132  if (closest(facei, otherFaceFaces[otherFaceFacei]))
133  {
134  return true;
135  }
136  }
137 
138  return false;
139 }
140 
141 
142 bool Foam::patchToPatches::nearest::intersectFaces
143 (
144  const primitiveOldTimePatch& srcPatch,
145  const vectorField& srcPointNormals,
146  const vectorField& srcPointNormals0,
147  const primitiveOldTimePatch& tgtPatch,
148  const label srcFacei,
149  const label tgtFacei
150 )
151 {
152  const bool srcCouples =
153  intersectFaces
154  (
155  srcPatch,
156  tgtPatch,
157  srcFacei,
158  tgtFacei,
159  srcLocalTgtFaces_[srcFacei],
160  srcDistances_[srcFacei]
161  );
162 
163  const bool tgtCouples =
164  intersectFaces
165  (
166  tgtPatch,
167  srcPatch,
168  tgtFacei,
169  srcFacei,
170  tgtLocalSrcFaces_[tgtFacei],
171  tgtDistances_[tgtFacei]
172  );
173 
174  return srcCouples || tgtCouples;
175 
176 }
177 
178 
179 void Foam::patchToPatches::nearest::initialise
180 (
181  const primitiveOldTimePatch& srcPatch,
182  const vectorField& srcPointNormals,
183  const vectorField& srcPointNormals0,
184  const primitiveOldTimePatch& tgtPatch
185 )
186 {
188  (
189  srcPatch,
190  srcPointNormals,
191  srcPointNormals0,
192  tgtPatch
193  );
194 
195  srcDistances_.resize(srcPatch.size());
196  srcDistances_ = vGreat;
197 
198  tgtDistances_.resize(tgtPatch.size());
199  tgtDistances_ = vGreat;
200 }
201 
202 
203 void Foam::patchToPatches::nearest::rDistributeTgt
204 (
205  const primitiveOldTimePatch& tgtPatch,
206  const distributionMap& tgtMap
207 )
208 {
209  // Create a list-list of distances to match the addressing
210  List<List<scalar>> tgtDistances(tgtLocalSrcFaces_.size());
211  forAll(tgtLocalSrcFaces_, tgtFacei)
212  {
213  if (!tgtLocalSrcFaces_[tgtFacei].empty())
214  {
215  tgtDistances[tgtFacei].resize(1, tgtDistances_[tgtFacei]);
216  }
217  }
218 
219  // Let the base class reverse distribute the addressing
220  patchToPatch::rDistributeTgt(tgtPatch, tgtMap);
221 
222  // Reverse distribute the distances
223  rDistributeListList(tgtPatch.size(), tgtMap, tgtDistances);
224 
225  // If there is more than one address, remove all but the closest
226  forAll(tgtLocalSrcFaces_, tgtFacei)
227  {
228  if (tgtLocalSrcFaces_[tgtFacei].size() > 1)
229  {
230  const label i = findMin(tgtDistances[tgtFacei]);
231 
232  const label srcFacei = tgtLocalSrcFaces_[tgtFacei][i];
233 
234  tgtLocalSrcFaces_[tgtFacei].resize(1);
235  tgtLocalSrcFaces_[tgtFacei][0] = srcFacei;
236 
237  tgtDistances_[tgtFacei] = tgtDistances[tgtFacei][i];
238  }
239  }
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
244 
246 :
247  patchToPatch(reverse)
248 {}
249 
250 
251 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
252 
254 {}
255 
256 
257 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
258 
261 {
263  (
265  );
266 
267  List<DynamicList<scalar>>& result = tResult.ref();
268 
269  forAll(srcLocalTgtFaces_, srcFacei)
270  {
271  if (!srcLocalTgtFaces_[srcFacei].empty())
272  {
273  result[srcFacei].resize(1, scalar(1));
274  }
275  }
276 
277  return tResult;
278 }
279 
280 
283 {
285  (
287  );
288 
289  List<DynamicList<scalar>>& result = tResult.ref();
290 
291  forAll(tgtLocalSrcFaces_, tgtFacei)
292  {
293  if (!tgtLocalSrcFaces_[tgtFacei].empty())
294  {
295  result[tgtFacei].resize(1, scalar(1));
296  }
297  }
298 
299  return tResult;
300 }
301 
302 
303 // ************************************************************************* //
addToRunTimeSelectionTable(patchToPatch, intersection, bool)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual void rDistributeTgt(const primitiveOldTimePatch &tgtPatch, const distributionMap &tgtMap)
Send data that resulted from an intersection between the source.
Definition: patchToPatch.C:775
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
label findMin(const ListType &, const label start=0)
Find index of min element (and less than given element).
List< DynamicList< label > > srcLocalTgtFaces_
For each source face, the coupled local target faces.
Definition: patchToPatch.H:115
T & ref()
Return non-const reference or generate a fatal error.
Definition: tmpNrcI.H:136
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
Macros for easy insertion into run-time selection tables.
virtual tmpNrc< List< DynamicList< scalar > > > srcWeights(const primitivePatch &srcPatch) const
For each source face, the coupled target weights.
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
defineTypeNameAndDebug(intersection, 0)
PrimitiveOldTimePatch< SubList< face >, const pointField & > primitiveOldTimePatch
Addressing for a faceList slice.
List< label > labelList
A List of labels.
Definition: labelList.H:56
nearest(const bool reverse)
Construct from components.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
static const Tensor I
Definition: Tensor.H:83
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
Addressing for a faceList slice.
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:52
vector point
Point is a vector.
Definition: point.H:41
List< DynamicList< label > > tgtLocalSrcFaces_
For each target face, the coupled local source faces.
Definition: patchToPatch.H:118
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual void initialise(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
Initialisation.
Definition: patchToPatch.C:690
Class to generate coupling geometry between two primitive patches.
Definition: patchToPatch.H:53
virtual tmpNrc< List< DynamicList< scalar > > > tgtWeights(const primitivePatch &tgtPatch) const
For each target face, the coupled source weights.
Namespace for OpenFOAM.