All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
inverseDistancePatchToPatch.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 
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace patchToPatches
35 {
36  defineTypeNameAndDebug(inverseDistance, 0);
37  addToRunTimeSelectionTable(patchToPatch, inverseDistance, bool);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::treeBoundBox Foam::patchToPatches::inverseDistance::srcBox
45 (
46  const face& srcFace,
47  const pointField& srcPoints,
48  const vectorField& srcPointNormals
49 ) const
50 {
51  const treeBoundBox bb(srcPoints, srcFace);
52 
53  const point c = bb.midpoint();
54  const scalar l = bb.maxDim();
55 
56  return treeBoundBox(c - l*vector::one, c + l*vector::one);
57 }
58 
59 
60 bool Foam::patchToPatches::inverseDistance::inside
61 (
62  const face& f,
63  const pointField& ps,
64  const point& p,
65  const vector& r
66 ) const
67 {
68  using namespace constant::mathematical;
69 
70  const tensor T = tensor::I - sqr(r);
71 
72  scalar angle = 0;
73 
74  forAll(f, i)
75  {
76  const vector& a = T & (ps[f[i]] - p);
77  const vector& b = T & (ps[f[f.fcIndex(i)]] - p);
78  const scalar magAB = sqrt(magSqr(a)*magSqr(b));
79  angle -= sign(r & (a ^ b))*acos((a & b)/magAB);
80  }
81 
82  return pi < angle && angle < 3*pi;
83 }
84 
85 
86 bool Foam::patchToPatches::inverseDistance::intersectFaces
87 (
88  const primitivePatch& patch,
89  const primitivePatch& otherPatch,
90  const label facei,
91  const label otherFacei,
92  DynamicList<label>& faceOtherFaces,
93  DynamicList<scalar>& faceWeights
94 ) const
95 {
96  const face& f = otherPatch[otherFacei];
97  const pointField& ps = otherPatch.points();
98 
99  const point& p = patch.faceCentres()[facei];
100  const vector& r = patch.faceNormals()[facei];
101 
102  bool intersectsOther = inside(f, ps, p, r);
103 
104  if (!intersectsOther)
105  {
106  forAll(otherPatch.faceFaces()[otherFacei], otherFaceFacei)
107  {
108  const label otherFacej =
109  otherPatch.faceFaces()[otherFacei][otherFaceFacei];
110 
111  const face& g = otherPatch[otherFacej];
112 
113  if (inside(g, ps, p, r))
114  {
115  intersectsOther = true;
116  break;
117  }
118  }
119  }
120 
121  if (intersectsOther)
122  {
123  faceOtherFaces.append(otherFacei);
124  faceWeights.append
125  (
126  1/max(mag(p - otherPatch.faceCentres()[otherFacei]), vSmall)
127  );
128  }
129 
130  return intersectsOther;
131 }
132 
133 
134 bool Foam::patchToPatches::inverseDistance::intersectFaces
135 (
136  const primitiveOldTimePatch& srcPatch,
137  const vectorField& srcPointNormals,
138  const vectorField& srcPointNormals0,
139  const primitiveOldTimePatch& tgtPatch,
140  const label srcFacei,
141  const label tgtFacei
142 )
143 {
144  const bool srcCouples =
145  intersectFaces
146  (
147  srcPatch,
148  tgtPatch,
149  srcFacei,
150  tgtFacei,
151  srcLocalTgtFaces_[srcFacei],
152  srcWeights_[srcFacei]
153  );
154 
155  const bool tgtCouples =
156  intersectFaces
157  (
158  tgtPatch,
159  srcPatch,
160  tgtFacei,
161  srcFacei,
162  tgtLocalSrcFaces_[tgtFacei],
163  tgtWeights_[tgtFacei]
164  );
165 
166  return srcCouples || tgtCouples;
167 }
168 
169 
170 void Foam::patchToPatches::inverseDistance::initialise
171 (
172  const primitiveOldTimePatch& srcPatch,
173  const vectorField& srcPointNormals,
174  const vectorField& srcPointNormals0,
175  const primitiveOldTimePatch& tgtPatch
176 )
177 {
179  (
180  srcPatch,
181  srcPointNormals,
182  srcPointNormals0,
183  tgtPatch
184  );
185 
186  srcWeights_.resize(srcPatch.size());
187  tgtWeights_.resize(tgtPatch.size());
188 }
189 
190 
191 void Foam::patchToPatches::inverseDistance::rDistributeTgt
192 (
193  const primitiveOldTimePatch& tgtPatch,
194  const distributionMap& tgtMap
195 )
196 {
197  patchToPatch::rDistributeTgt(tgtPatch, tgtMap);
198 
199  rDistributeListList(tgtPatch.size(), tgtMap, tgtWeights_);
200 }
201 
202 
203 Foam::label Foam::patchToPatches::inverseDistance::finalise
204 (
205  const primitiveOldTimePatch& srcPatch,
206  const vectorField& srcPointNormals,
207  const vectorField& srcPointNormals0,
208  const primitiveOldTimePatch& tgtPatch,
209  const transformer& tgtToSrc
210 )
211 {
212  const label nCouples =
214  (
215  srcPatch,
216  srcPointNormals,
217  srcPointNormals0,
218  tgtPatch,
219  tgtToSrc
220  );
221 
222  forAll(srcWeights_, srcFacei)
223  {
224  const scalar w = sum(srcWeights_[srcFacei]);
225 
226  forAll(srcWeights_[srcFacei], i)
227  {
228  srcWeights_[srcFacei][i] /= max(w, vSmall);
229  }
230  }
231 
232  forAll(tgtWeights_, tgtFacei)
233  {
234  const scalar w = sum(tgtWeights_[tgtFacei]);
235 
236  forAll(tgtWeights_[tgtFacei], i)
237  {
238  tgtWeights_[tgtFacei][i] /= max(w, vSmall);
239  }
240  }
241 
242  return nCouples;
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
247 
249 :
250  patchToPatch(reverse),
251  srcWeights_(),
252  tgtWeights_()
253 {}
254 
255 
256 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
257 
259 {}
260 
261 
262 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
263 
266 (
267  const primitivePatch& srcPatch
268 ) const
269 {
270  return srcWeights_;
271 }
272 
273 
276 (
277  const primitivePatch& tgtPatch
278 ) const
279 {
280  return tgtWeights_;
281 }
282 
283 
284 // ************************************************************************* //
addToRunTimeSelectionTable(patchToPatch, intersection, bool)
dimensionedScalar sign(const dimensionedScalar &ds)
inverseDistance(const bool reverse)
Construct from components.
dimensionedScalar acos(const dimensionedScalar &ds)
#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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
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.
virtual tmpNrc< List< DynamicList< scalar > > > srcWeights(const primitivePatch &srcPatch) const
For each source face, the coupled target weights.
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
virtual tmpNrc< List< DynamicList< scalar > > > tgtWeights(const primitivePatch &tgtPatch) const
For each target face, the coupled source weights.
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
Addressing for a faceList slice.
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
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
dimensioned< scalar > mag(const dimensioned< Type > &)
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
volScalarField & p
Class to generate coupling geometry between two primitive patches.
Definition: patchToPatch.H:53
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Namespace for OpenFOAM.