withGaps_searchableSurface.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) 2011-2025 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 "Time.H"
28 #include "ListOps.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace searchableSurfaces
36  {
38 
40  (
42  withGaps,
44  );
45 
47  (
49  withGaps,
50  dictionary,
51  searchableSurfaceWithGaps,
52  "searchableSurfaceWithGaps"
53  );
54  }
55 }
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
60 Foam::Pair<Foam::vector> Foam::searchableSurfaces::withGaps::offsetVecs
61 (
62  const point& start,
63  const point& end
64 ) const
65 {
66  Pair<vector> offsets(Zero, Zero);
67 
68  vector n(end-start);
69 
70  scalar magN = mag(n);
71 
72  if (magN > small)
73  {
74  n /= magN;
75 
76  // Do first offset vector. Is the coordinate axes with the smallest
77  // component along the vector n.
78  scalar minMag = great;
79  direction minCmpt = 0;
80 
81  for (direction cmpt = 0; cmpt < vector::nComponents; cmpt++)
82  {
83  if (mag(n[cmpt]) < minMag)
84  {
85  minMag = mag(n[cmpt]);
86  minCmpt = cmpt;
87  }
88  }
89 
90  offsets[0][minCmpt] = 1.0;
91  // Orthonormalise
92  offsets[0] -= n[minCmpt]*n;
93  offsets[0] /= mag(offsets[0]);
94  // Do second offset vector perp to original edge and first offset vector
95  offsets[1] = n ^ offsets[0];
96 
97  // Scale
98  offsets[0] *= gap_;
99  offsets[1] *= gap_;
100  }
101 
102  return offsets;
103 }
104 
105 
106 void Foam::searchableSurfaces::withGaps::offsetVecs
107 (
108  const pointField& start,
109  const pointField& end,
110  pointField& offset0,
111  pointField& offset1
112 ) const
113 {
114  offset0.setSize(start.size());
115  offset1.setSize(start.size());
116 
117  forAll(start, i)
118  {
119  const Pair<vector> offsets(offsetVecs(start[i], end[i]));
120  offset0[i] = offsets[0];
121  offset1[i] = offsets[1];
122  }
123 }
124 
125 
126 Foam::label Foam::searchableSurfaces::withGaps::countMisses
127 (
128  const List<pointIndexHit>& info,
129  labelList& missMap
130 )
131 {
132  label nMiss = 0;
133  forAll(info, i)
134  {
135  if (!info[i].hit())
136  {
137  nMiss++;
138  }
139  }
140 
141  missMap.setSize(nMiss);
142  nMiss = 0;
143 
144  forAll(info, i)
145  {
146  if (!info[i].hit())
147  {
148  missMap[nMiss++] = i;
149  }
150  }
151 
152  return nMiss;
153 }
154 
155 
156 // Anything not a hit in both counts as a hit
157 Foam::label Foam::searchableSurfaces::withGaps::countMisses
158 (
159  const List<pointIndexHit>& plusInfo,
160  const List<pointIndexHit>& minInfo,
161  labelList& missMap
162 )
163 {
164  label nMiss = 0;
165  forAll(plusInfo, i)
166  {
167  if (!plusInfo[i].hit() || !minInfo[i].hit())
168  {
169  nMiss++;
170  }
171  }
172 
173  missMap.setSize(nMiss);
174  nMiss = 0;
175 
176  forAll(plusInfo, i)
177  {
178  if (!plusInfo[i].hit() || !minInfo[i].hit())
179  {
180  missMap[nMiss++] = i;
181  }
182  }
183 
184  return nMiss;
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
189 
191 (
192  const IOobject& io,
193  const dictionary& dict
194 )
195 :
196  searchableSurface(io),
197  gap_(dict.lookup<scalar>("gap", dimLength)),
198  subGeom_(1)
199 {
200  const word subGeomName(dict.lookup("surface"));
201 
202  const searchableSurface& s =
203  io.db().lookupObject<searchableSurface>(subGeomName);
204 
205  subGeom_.set(0, &const_cast<searchableSurface&>(s));
206 
207  bounds() = subGeom_[0].bounds();
208 }
209 
210 
211 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
212 
214 {}
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
220 (
221  const pointField& start,
222  const pointField& end,
223  List<pointIndexHit>& info
224 ) const
225 {
226 
227  // Test with unperturbed vectors
228  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
229 
230  surface().findLine(start, end, info);
231 
232  // Count number of misses. Determine map
233  labelList compactMap;
234  label nMiss = countMisses(info, compactMap);
235 
236  if (returnReduce(nMiss, sumOp<label>()) > 0)
237  {
238  // Pout<< "** retesting with offset0 " << nMiss << " misses out of "
239  // << start.size() << endl;
240 
241  // extract segments according to map
242  pointField compactStart(start, compactMap);
243  pointField compactEnd(end, compactMap);
244 
245  // Calculate offset vector
246  pointField offset0, offset1;
247  offsetVecs
248  (
249  compactStart,
250  compactEnd,
251  offset0,
252  offset1
253  );
254 
255  // Test with offset0 perturbed vectors
256  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
257 
258  // test in pairs: only if both perturbations hit something
259  // do we accept the hit.
260 
261  const vectorField smallVec(1e-6*(compactEnd-compactStart));
262 
263  List<pointIndexHit> plusInfo;
264  surface().findLine
265  (
266  compactStart+offset0-smallVec,
267  compactEnd+offset0+smallVec,
268  plusInfo
269  );
270  List<pointIndexHit> minInfo;
271  surface().findLine
272  (
273  compactStart-offset0-smallVec,
274  compactEnd-offset0+smallVec,
275  minInfo
276  );
277 
278  // Extract any hits
279  forAll(plusInfo, i)
280  {
281  if (plusInfo[i].hit() && minInfo[i].hit())
282  {
283  info[compactMap[i]] = plusInfo[i];
284  info[compactMap[i]].rawPoint() -= offset0[i];
285  }
286  }
287 
288  labelList plusMissMap;
289  nMiss = countMisses(plusInfo, minInfo, plusMissMap);
290 
291  if (returnReduce(nMiss, sumOp<label>()) > 0)
292  {
293  // Pout<< "** retesting with offset1 " << nMiss << " misses out of "
294  // << start.size() << endl;
295 
296  // Test with offset1 perturbed vectors
297  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
298 
299  // Extract (inplace possible because of order)
300  forAll(plusMissMap, i)
301  {
302  label mapI = plusMissMap[i];
303  compactStart[i] = compactStart[mapI];
304  compactEnd[i] = compactEnd[mapI];
305  compactMap[i] = compactMap[mapI];
306  offset0[i] = offset0[mapI];
307  offset1[i] = offset1[mapI];
308  }
309  compactStart.setSize(plusMissMap.size());
310  compactEnd.setSize(plusMissMap.size());
311  compactMap.setSize(plusMissMap.size());
312  offset0.setSize(plusMissMap.size());
313  offset1.setSize(plusMissMap.size());
314 
315  const vectorField smallVec(1e-6*(compactEnd-compactStart));
316 
317  surface().findLine
318  (
319  compactStart+offset1-smallVec,
320  compactEnd+offset1+smallVec,
321  plusInfo
322  );
323  surface().findLine
324  (
325  compactStart-offset1-smallVec,
326  compactEnd-offset1+smallVec,
327  minInfo
328  );
329 
330  // Extract any hits
331  forAll(plusInfo, i)
332  {
333  if (plusInfo[i].hit() && minInfo[i].hit())
334  {
335  info[compactMap[i]] = plusInfo[i];
336  info[compactMap[i]].rawPoint() -= offset1[i];
337  }
338  }
339  }
340  }
341 }
342 
343 
345 (
346  const pointField& start,
347  const pointField& end,
348  List<pointIndexHit>& info
349 ) const
350 {
351  // To be done ...
352  findLine(start, end, info);
353 }
354 
355 
357 (
358  const pointField& start,
359  const pointField& end,
361 ) const
362 {
363  // To be done. Assume for now only one intersection.
364  List<pointIndexHit> nearestInfo;
365  findLine(start, end, nearestInfo);
366 
367  info.setSize(start.size());
368  forAll(info, pointi)
369  {
370  if (nearestInfo[pointi].hit())
371  {
372  info[pointi].setSize(1);
373  info[pointi][0] = nearestInfo[pointi];
374  }
375  else
376  {
377  info[pointi].clear();
378  }
379  }
380 }
381 
382 
383 // ************************************************************************* //
Various functions to operate on Lists.
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:309
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
void setSize(const label)
Reset size of List.
Definition: List.C:281
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:105
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
const boundBox & bounds() const
Return const reference to boundBox.
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
withGaps(const IOobject &io, const dictionary &dict)
Construct from dictionary (used by searchableSurface)
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
A class for handling words, derived from string.
Definition: word.H:62
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
addToRunTimeSelectionTable(searchableSurface, box, dictionary)
addBackwardCompatibleToRunTimeSelectionTable(searchableSurface, box, dictionary, searchableBox, "searchableBox")
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
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
const dimensionSet dimLength
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
uint8_t direction
Definition: direction.H:45
dictionary dict