searchableSurfaceWithGaps.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "Time.H"
32 #include "ListOps.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(searchableSurfaceWithGaps, 0);
40  (
41  searchableSurface,
42  searchableSurfaceWithGaps,
43  dict
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 Foam::Pair<Foam::vector> Foam::searchableSurfaceWithGaps::offsetVecs
51 (
52  const point& start,
53  const point& end
54 ) const
55 {
56  Pair<vector> offsets(Zero, Zero);
57 
58  vector n(end-start);
59 
60  scalar magN = mag(n);
61 
62  if (magN > SMALL)
63  {
64  n /= magN;
65 
66  // Do first offset vector. Is the coordinate axes with the smallest
67  // component along the vector n.
68  scalar minMag = GREAT;
69  direction minCmpt = 0;
70 
71  for (direction cmpt = 0; cmpt < vector::nComponents; cmpt++)
72  {
73  if (mag(n[cmpt]) < minMag)
74  {
75  minMag = mag(n[cmpt]);
76  minCmpt = cmpt;
77  }
78  }
79 
80  offsets[0][minCmpt] = 1.0;
81  // Orthonormalise
82  offsets[0] -= n[minCmpt]*n;
83  offsets[0] /= mag(offsets[0]);
84  // Do second offset vector perp to original edge and first offset vector
85  offsets[1] = n ^ offsets[0];
86 
87  // Scale
88  offsets[0] *= gap_;
89  offsets[1] *= gap_;
90  }
91 
92  return offsets;
93 }
94 
95 
96 void Foam::searchableSurfaceWithGaps::offsetVecs
97 (
98  const pointField& start,
99  const pointField& end,
100  pointField& offset0,
101  pointField& offset1
102 ) const
103 {
104  offset0.setSize(start.size());
105  offset1.setSize(start.size());
106 
107  forAll(start, i)
108  {
109  const Pair<vector> offsets(offsetVecs(start[i], end[i]));
110  offset0[i] = offsets[0];
111  offset1[i] = offsets[1];
112  }
113 }
114 
115 
116 Foam::label Foam::searchableSurfaceWithGaps::countMisses
117 (
118  const List<pointIndexHit>& info,
119  labelList& missMap
120 )
121 {
122  label nMiss = 0;
123  forAll(info, i)
124  {
125  if (!info[i].hit())
126  {
127  nMiss++;
128  }
129  }
130 
131  missMap.setSize(nMiss);
132  nMiss = 0;
133 
134  forAll(info, i)
135  {
136  if (!info[i].hit())
137  {
138  missMap[nMiss++] = i;
139  }
140  }
141 
142  return nMiss;
143 }
144 
145 
146 // Anything not a hit in both counts as a hit
147 Foam::label Foam::searchableSurfaceWithGaps::countMisses
148 (
149  const List<pointIndexHit>& plusInfo,
150  const List<pointIndexHit>& minInfo,
151  labelList& missMap
152 )
153 {
154  label nMiss = 0;
155  forAll(plusInfo, i)
156  {
157  if (!plusInfo[i].hit() || !minInfo[i].hit())
158  {
159  nMiss++;
160  }
161  }
162 
163  missMap.setSize(nMiss);
164  nMiss = 0;
165 
166  forAll(plusInfo, i)
167  {
168  if (!plusInfo[i].hit() || !minInfo[i].hit())
169  {
170  missMap[nMiss++] = i;
171  }
172  }
173 
174  return nMiss;
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
179 
180 Foam::searchableSurfaceWithGaps::searchableSurfaceWithGaps
181 (
182  const IOobject& io,
183  const dictionary& dict
184 )
185 :
187  gap_(dict.get<scalar>("gap")),
188  subGeom_(1)
189 {
190  const word subGeomName(dict.get<word>("surface"));
191 
192  subGeom_.set
193  (
194  0,
195  io.db().getObjectPtr<searchableSurface>(subGeomName)
196  );
197 
198  bounds() = subGeom_[0].bounds();
199 }
200 
201 
202 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
203 
205 (
206  const pointField& start,
207  const pointField& end,
208  List<pointIndexHit>& info
209 ) const
210 {
211 
212  // Test with unperturbed vectors
213  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
214 
215  surface().findLine(start, end, info);
216 
217  // Count number of misses. Determine map
218  labelList compactMap;
219  label nMiss = countMisses(info, compactMap);
220 
221  if (returnReduceOr(nMiss))
222  {
223  //Pout<< "** retesting with offset0 " << nMiss << " misses out of "
224  // << start.size() << endl;
225 
226  // extract segments according to map
227  pointField compactStart(start, compactMap);
228  pointField compactEnd(end, compactMap);
229 
230  // Calculate offset vector
231  pointField offset0, offset1;
232  offsetVecs
233  (
234  compactStart,
235  compactEnd,
236  offset0,
237  offset1
238  );
239 
240  // Test with offset0 perturbed vectors
241  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
242 
243  // test in pairs: only if both perturbations hit something
244  // do we accept the hit.
245 
246  const vectorField smallVec(1e-6*(compactEnd-compactStart));
247 
248  List<pointIndexHit> plusInfo;
249  surface().findLine
250  (
251  compactStart+offset0-smallVec,
252  compactEnd+offset0+smallVec,
253  plusInfo
254  );
255  List<pointIndexHit> minInfo;
256  surface().findLine
257  (
258  compactStart-offset0-smallVec,
259  compactEnd-offset0+smallVec,
260  minInfo
261  );
262 
263  // Extract any hits
264  forAll(plusInfo, i)
265  {
266  if (plusInfo[i].hit() && minInfo[i].hit())
267  {
268  info[compactMap[i]] = plusInfo[i];
269  info[compactMap[i]].point() -= offset0[i];
270  }
271  }
272 
273  labelList plusMissMap;
274  nMiss = countMisses(plusInfo, minInfo, plusMissMap);
275 
276  if (returnReduceOr(nMiss))
277  {
278  //Pout<< "** retesting with offset1 " << nMiss << " misses out of "
279  // << start.size() << endl;
280 
281  // Test with offset1 perturbed vectors
282  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
283 
284  // Extract (inplace possible because of order)
285  forAll(plusMissMap, i)
286  {
287  label mapI = plusMissMap[i];
288  compactStart[i] = compactStart[mapI];
289  compactEnd[i] = compactEnd[mapI];
290  compactMap[i] = compactMap[mapI];
291  offset0[i] = offset0[mapI];
292  offset1[i] = offset1[mapI];
293  }
294  compactStart.setSize(plusMissMap.size());
295  compactEnd.setSize(plusMissMap.size());
296  compactMap.setSize(plusMissMap.size());
297  offset0.setSize(plusMissMap.size());
298  offset1.setSize(plusMissMap.size());
299 
300  const vectorField smallVec(1e-6*(compactEnd-compactStart));
301 
302  surface().findLine
303  (
304  compactStart+offset1-smallVec,
305  compactEnd+offset1+smallVec,
306  plusInfo
307  );
308  surface().findLine
309  (
310  compactStart-offset1-smallVec,
311  compactEnd-offset1+smallVec,
312  minInfo
313  );
314 
315  // Extract any hits
316  forAll(plusInfo, i)
317  {
318  if (plusInfo[i].hit() && minInfo[i].hit())
319  {
320  info[compactMap[i]] = plusInfo[i];
321  info[compactMap[i]].point() -= offset1[i];
322  }
323  }
324  }
325  }
326 }
327 
328 
330 (
331  const pointField& start,
332  const pointField& end,
333  List<pointIndexHit>& info
334 ) const
335 {
336  // To be done ...
337  findLine(start, end, info);
338 }
339 
340 
342 (
343  const pointField& start,
344  const pointField& end,
346 ) const
347 {
348  // To be done. Assume for now only one intersection.
349  List<pointIndexHit> nearestInfo;
350  findLine(start, end, nearestInfo);
351 
352  info.setSize(start.size());
353  forAll(info, pointi)
354  {
355  if (nearestInfo[pointi].hit())
356  {
357  info[pointi].setSize(1);
358  info[pointi][0] = nearestInfo[pointi];
359  }
360  else
361  {
362  info[pointi].clear();
363  }
364  }
365 }
366 
367 
368 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
virtual const boundBox & bounds() const
Return const reference to boundBox.
Macros for easy insertion into run-time selection tables.
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
A class for handling words, derived from Foam::string.
Definition: word.H:63
Vector< scalar > vector
Definition: vector.H:57
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
defineTypeNameAndDebug(combustionModel, 0)
vector point
Point is a vector.
Definition: point.H:37
label n
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127