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