searchableSurfacesQueries.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 | Copyright (C) 2011-2016 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 "ListOps.H"
28 #include "OFstream.H"
29 #include "meshTools.H"
30 #include "DynamicField.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(searchableSurfacesQueries, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 Foam::pointIndexHit Foam::searchableSurfacesQueries::tempFindNearest
43 (
44  const searchableSurface& surf,
45  const point& pt,
46  const scalar initDistSqr
47 )
48 {
49  pointField onePoint(1, pt);
50  scalarField oneDist(1, initDistSqr);
51  List<pointIndexHit> oneHit(1);
52  surf.findNearest(onePoint, oneDist, oneHit);
53  return oneHit[0];
54 }
55 
56 
57 // Calculate sum of distance to surfaces.
58 Foam::scalar Foam::searchableSurfacesQueries::sumDistSqr
59 (
60  const PtrList<searchableSurface>& allSurfaces,
61  const labelList& surfacesToTest,
62  const scalar initDistSqr,
63  const point& pt
64 )
65 {
66  scalar sum = 0;
67 
68  forAll(surfacesToTest, testI)
69  {
70  label surfI = surfacesToTest[testI];
71 
72  pointIndexHit hit
73  (
74  tempFindNearest(allSurfaces[surfI], pt, initDistSqr)
75  );
76 
77  // Note: make it fall over if not hit.
78  sum += magSqr(hit.hitPoint()-pt);
79  }
80  return sum;
81 }
82 
83 
84 // Reflects the point furthest away around the triangle centre by a factor fac.
85 // (triangle centre is the average of all points but the ihi. pSum is running
86 // sum of all points)
87 Foam::scalar Foam::searchableSurfacesQueries::tryMorphTet
88 (
89  const PtrList<searchableSurface>& allSurfaces,
90  const labelList& surfacesToTest,
91  const scalar initDistSqr,
92  List<vector>& p,
93  List<scalar>& y,
94  vector& pSum,
95  const label ihi,
96  const scalar fac
97 )
98 {
99  scalar fac1 = (1.0-fac)/vector::nComponents;
100  scalar fac2 = fac1-fac;
101 
102  vector ptry = pSum*fac1-p[ihi]*fac2;
103 
104  scalar ytry = sumDistSqr(allSurfaces, surfacesToTest, initDistSqr, ptry);
105 
106  if (ytry < y[ihi])
107  {
108  y[ihi] = ytry;
109  pSum += ptry - p[ihi];
110  p[ihi] = ptry;
111  }
112  return ytry;
113 }
114 
115 
116 bool Foam::searchableSurfacesQueries::morphTet
117 (
118  const PtrList<searchableSurface>& allSurfaces,
119  const labelList& surfacesToTest,
120  const scalar initDistSqr,
121  const scalar convergenceDistSqr,
122  const label maxIter,
123  List<vector>& p,
124  List<scalar>& y
125 )
126 {
127  vector pSum = sum(p);
128 
129  autoPtr<OFstream> str;
130  label vertI = 0;
131  if (debug)
132  {
133  wordList names(surfacesToTest.size());
134  forAll(surfacesToTest, i)
135  {
136  names[i] = allSurfaces[surfacesToTest[i]].name();
137  }
138  Pout<< "searchableSurfacesQueries::morphTet : intersection of "
139  << names << " starting from points:" << p << endl;
140  str.reset(new OFstream("track.obj"));
141  meshTools::writeOBJ(str(), p[0]);
142  vertI++;
143  }
144 
145  for (label iter = 0; iter < maxIter; iter++)
146  {
147  // Get the indices of lowest, highest and second-highest values.
148  label ilo, ihi, inhi;
149  {
150  labelList sortedIndices;
151  sortedOrder(y, sortedIndices);
152  ilo = sortedIndices[0];
153  ihi = sortedIndices[sortedIndices.size()-1];
154  inhi = sortedIndices[sortedIndices.size()-2];
155  }
156 
157  if (debug)
158  {
159  Pout<< "Iteration:" << iter
160  << " lowest:" << y[ilo] << " highest:" << y[ihi]
161  << " points:" << p << endl;
162 
163  meshTools::writeOBJ(str(), p[ilo]);
164  vertI++;
165  str()<< "l " << vertI-1 << ' ' << vertI << nl;
166  }
167 
168  if (y[ihi] < convergenceDistSqr)
169  {
170  // Get point on 0th surface.
171  Swap(p[0], p[ilo]);
172  Swap(y[0], y[ilo]);
173  return true;
174  }
175 
176  // Reflection: point furthest away gets reflected.
177  scalar ytry = tryMorphTet
178  (
179  allSurfaces,
180  surfacesToTest,
181  10*y[ihi], // search box.
182  p,
183  y,
184  pSum,
185  ihi,
186  -1.0
187  );
188 
189  if (ytry <= y[ilo])
190  {
191  // If in right direction (y lower) expand by two.
192  ytry = tryMorphTet
193  (
194  allSurfaces,
195  surfacesToTest,
196  10*y[ihi],
197  p,
198  y,
199  pSum,
200  ihi,
201  2.0
202  );
203  }
204  else if (ytry >= y[inhi])
205  {
206  // If inside tet try contraction.
207 
208  scalar ysave = y[ihi];
209 
210  ytry = tryMorphTet
211  (
212  allSurfaces,
213  surfacesToTest,
214  10*y[ihi],
215  p,
216  y,
217  pSum,
218  ihi,
219  0.5
220  );
221 
222  if (ytry >= ysave)
223  {
224  // Contract around lowest point.
225  forAll(p, i)
226  {
227  if (i != ilo)
228  {
229  p[i] = 0.5*(p[i] + p[ilo]);
230  y[i] = sumDistSqr
231  (
232  allSurfaces,
233  surfacesToTest,
234  y[ihi],
235  p[i]
236  );
237  }
238  }
239  pSum = sum(p);
240  }
241  }
242  }
243 
244  if (debug)
245  {
246  meshTools::writeOBJ(str(), p[0]);
247  vertI++;
248  str()<< "l " << vertI-1 << ' ' << vertI << nl;
249  }
250 
251  // Failure to converge. Return best guess so far.
252  label ilo = findMin(y);
253  Swap(p[0], p[ilo]);
254  Swap(y[0], y[ilo]);
255  return false;
256 }
257 
258 
260 //void Foam::searchableSurfacesQueries::findAllIntersections
261 //(
262 // const searchableSurface& s,
263 // const pointField& start,
264 // const pointField& end,
265 // const vectorField& smallVec,
266 // List<List<pointIndexHit>>& surfaceHitInfo
267 //)
268 //{
269 // surfaceHitInfo.setSize(start.size());
270 //
271 // // Current start point of vector
272 // pointField p0(start);
273 //
274 // List<pointIndexHit> intersectInfo(start.size());
275 //
276 // // For test whether finished doing vector.
277 // const vectorField dirVec(end-start);
278 // const scalarField magSqrDirVec(magSqr(dirVec));
279 //
280 // while (true)
281 // {
282 // // Find first intersection. Synced.
283 // s.findLine(p0, end, intersectInfo);
284 //
285 // label nHits = 0;
286 //
287 // forAll(intersectInfo, i)
288 // {
289 // if (intersectInfo[i].hit())
290 // {
291 // nHits++;
292 //
293 // label sz = surfaceHitInfo[i].size();
294 // surfaceHitInfo[i].setSize(sz+1);
295 // surfaceHitInfo[i][sz] = intersectInfo[i];
296 //
297 // p0[i] = intersectInfo[i].hitPoint() + smallVec[i];
298 //
299 // // If beyond endpoint set to endpoint so as not to pick up
300 // // any intersections. Could instead just filter out hits.
301 // if (((p0[i]-start[i])&dirVec[i]) > magSqrDirVec[i])
302 // {
303 // p0[i] = end[i];
304 // }
305 // }
306 // else
307 // {
308 // // Set to endpoint to stop intersection test. See above.
309 // p0[i] = end[i];
310 // }
311 // }
312 //
313 // // returnReduce(nHits) ?
314 // if (nHits == 0)
315 // {
316 // break;
317 // }
318 // }
319 //}
320 
321 
322 // Given current set of hits (allSurfaces, allInfo) merge in those coming from
323 // surface surfI.
324 void Foam::searchableSurfacesQueries::mergeHits
325 (
326  const point& start,
327 
328  const label testI, // index of surface
329  const List<pointIndexHit>& surfHits, // hits on surface
330 
331  labelList& allSurfaces,
332  List<pointIndexHit>& allInfo,
333  scalarList& allDistSqr
334 )
335 {
336  // Precalculate distances
337  scalarList surfDistSqr(surfHits.size());
338  forAll(surfHits, i)
339  {
340  surfDistSqr[i] = magSqr(surfHits[i].hitPoint() - start);
341  }
342 
343  forAll(surfDistSqr, i)
344  {
345  label index = findLower(allDistSqr, surfDistSqr[i]);
346 
347  // Check if equal to lower.
348  if (index >= 0)
349  {
350  // Same. Do not count.
351  //Pout<< "point:" << surfHits[i].hitPoint()
352  // << " considered same as:" << allInfo[index].hitPoint()
353  // << " within tol:" << mergeDist
354  // << endl;
355  }
356  else
357  {
358  // Check if equal to higher
359  label next = index + 1;
360 
361  if (next < allDistSqr.size())
362  {
363  //Pout<< "point:" << surfHits[i].hitPoint()
364  // << " considered same as:" << allInfo[next].hitPoint()
365  // << " within tol:" << mergeDist
366  // << endl;
367  }
368  else
369  {
370  // Insert after index
371  label sz = allSurfaces.size();
372  allSurfaces.setSize(sz+1);
373  allInfo.setSize(allSurfaces.size());
374  allDistSqr.setSize(allSurfaces.size());
375  // Make space.
376  for (label j = sz-1; j > index; --j)
377  {
378  allSurfaces[j+1] = allSurfaces[j];
379  allInfo[j+1] = allInfo[j];
380  allDistSqr[j+1] = allDistSqr[j];
381  }
382  // Insert new value
383  allSurfaces[index+1] = testI;
384  allInfo[index+1] = surfHits[i];
385  allDistSqr[index+1] = surfDistSqr[i];
386  }
387  }
388  }
389 }
390 
391 
392 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
393 
394 // Find any intersection
396 (
397  const PtrList<searchableSurface>& allSurfaces,
398  const labelList& surfacesToTest,
399  const pointField& start,
400  const pointField& end,
401  labelList& hitSurfaces,
402  List<pointIndexHit>& hitInfo
403 )
404 {
405  hitSurfaces.setSize(start.size());
406  hitSurfaces = -1;
407  hitInfo.setSize(start.size());
408 
409  // Work arrays
410  labelList hitMap(identity(start.size()));
411  pointField p0(start);
412  pointField p1(end);
413  List<pointIndexHit> intersectInfo(start.size());
414 
415  forAll(surfacesToTest, testI)
416  {
417  // Do synchronised call to all surfaces.
418  allSurfaces[surfacesToTest[testI]].findLineAny(p0, p1, intersectInfo);
419 
420  // Copy all hits into arguments, continue with misses
421  label newI = 0;
422  forAll(intersectInfo, i)
423  {
424  if (intersectInfo[i].hit())
425  {
426  hitInfo[hitMap[i]] = intersectInfo[i];
427  hitSurfaces[hitMap[i]] = testI;
428  }
429  else
430  {
431  if (i != newI)
432  {
433  hitMap[newI] = hitMap[i];
434  p0[newI] = p0[i];
435  p1[newI] = p1[i];
436  }
437  newI++;
438  }
439  }
440 
441  // All done? Note that this decision should be synchronised
442  if (newI == 0)
443  {
444  break;
445  }
446 
447  // Trim and continue
448  hitMap.setSize(newI);
449  p0.setSize(newI);
450  p1.setSize(newI);
451  intersectInfo.setSize(newI);
452  }
453 }
454 
455 
457 (
458  const PtrList<searchableSurface>& allSurfaces,
459  const labelList& surfacesToTest,
460  const pointField& start,
461  const pointField& end,
462  labelListList& hitSurfaces,
463  List<List<pointIndexHit>>& hitInfo
464 )
465 {
466  // Note: maybe move the single-surface all intersections test into
467  // searchable surface? Some of the tolerance issues might be
468  // lessened.
469 
470  // 2. Currently calling searchableSurface::findLine with start==end
471  // is expected to find no intersection. Problem if it does.
472 
473  hitSurfaces.setSize(start.size());
474  hitInfo.setSize(start.size());
475 
476  if (surfacesToTest.empty())
477  {
478  return;
479  }
480 
481  // Test first surface
482  allSurfaces[surfacesToTest[0]].findLineAll(start, end, hitInfo);
483 
484  // Set hitSurfaces and distance
485  List<scalarList> hitDistSqr(hitInfo.size());
486  forAll(hitInfo, pointi)
487  {
488  const List<pointIndexHit>& pHits = hitInfo[pointi];
489 
490  labelList& pSurfaces = hitSurfaces[pointi];
491  pSurfaces.setSize(pHits.size());
492  pSurfaces = 0;
493 
494  scalarList& pDistSqr = hitDistSqr[pointi];
495  pDistSqr.setSize(pHits.size());
496  forAll(pHits, i)
497  {
498  pDistSqr[i] = magSqr(pHits[i].hitPoint() - start[pointi]);
499  }
500  }
501 
502 
503  if (surfacesToTest.size() > 1)
504  {
505  // Test the other surfaces and merge (according to distance from start).
506  for (label testI = 1; testI < surfacesToTest.size(); testI++)
507  {
508  List<List<pointIndexHit>> surfHits;
509  allSurfaces[surfacesToTest[testI]].findLineAll
510  (
511  start,
512  end,
513  surfHits
514  );
515 
516  forAll(surfHits, pointi)
517  {
518  mergeHits
519  (
520  start[pointi], // Current segment
521 
522  testI, // Surface and its hits
523  surfHits[pointi],
524 
525  hitSurfaces[pointi], // Merge into overall hit info
526  hitInfo[pointi],
527  hitDistSqr[pointi]
528  );
529  }
530  }
531  }
532 }
533 
534 
535 //Find intersections of edge nearest to both endpoints.
537 (
538  const PtrList<searchableSurface>& allSurfaces,
539  const labelList& surfacesToTest,
540  const pointField& start,
541  const pointField& end,
542  labelList& surface1,
543  List<pointIndexHit>& hit1,
544  labelList& surface2,
545  List<pointIndexHit>& hit2
546 )
547 {
548  // 1. intersection from start to end
549  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
550 
551  // Initialize arguments
552  surface1.setSize(start.size());
553  surface1 = -1;
554  hit1.setSize(start.size());
555 
556  // Current end of segment to test.
557  pointField nearest(end);
558  // Work array
559  List<pointIndexHit> nearestInfo(start.size());
560 
561  forAll(surfacesToTest, testI)
562  {
563  // See if any intersection between start and current nearest
564  allSurfaces[surfacesToTest[testI]].findLine
565  (
566  start,
567  nearest,
568  nearestInfo
569  );
570 
571  forAll(nearestInfo, pointi)
572  {
573  if (nearestInfo[pointi].hit())
574  {
575  hit1[pointi] = nearestInfo[pointi];
576  surface1[pointi] = testI;
577  nearest[pointi] = hit1[pointi].hitPoint();
578  }
579  }
580  }
581 
582 
583  // 2. intersection from end to last intersection
584  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
585 
586  // Find the nearest intersection from end to start. Note that we
587  // initialize to the first intersection (if any).
588  surface2 = surface1;
589  hit2 = hit1;
590 
591  // Set current end of segment to test.
592  forAll(nearest, pointi)
593  {
594  if (hit1[pointi].hit())
595  {
596  nearest[pointi] = hit1[pointi].hitPoint();
597  }
598  else
599  {
600  // Disable testing by setting to end.
601  nearest[pointi] = end[pointi];
602  }
603  }
604 
605  forAll(surfacesToTest, testI)
606  {
607  // See if any intersection between end and current nearest
608  allSurfaces[surfacesToTest[testI]].findLine(end, nearest, nearestInfo);
609 
610  forAll(nearestInfo, pointi)
611  {
612  if (nearestInfo[pointi].hit())
613  {
614  hit2[pointi] = nearestInfo[pointi];
615  surface2[pointi] = testI;
616  nearest[pointi] = hit2[pointi].hitPoint();
617  }
618  }
619  }
620 }
621 
622 
623 // Find nearest. Return -1 or nearest point
625 (
626  const PtrList<searchableSurface>& allSurfaces,
627  const labelList& surfacesToTest,
628  const pointField& samples,
629  const scalarField& nearestDistSqr,
630  labelList& nearestSurfaces,
631  List<pointIndexHit>& nearestInfo
632 )
633 {
634  // Initialise
635  nearestSurfaces.setSize(samples.size());
636  nearestSurfaces = -1;
637  nearestInfo.setSize(samples.size());
638 
639  // Work arrays
640  scalarField minDistSqr(nearestDistSqr);
641  List<pointIndexHit> hitInfo(samples.size());
642 
643  forAll(surfacesToTest, testI)
644  {
645  allSurfaces[surfacesToTest[testI]].findNearest
646  (
647  samples,
648  minDistSqr,
649  hitInfo
650  );
651 
652  // Update minDistSqr and arguments
653  forAll(hitInfo, pointi)
654  {
655  if (hitInfo[pointi].hit())
656  {
657  minDistSqr[pointi] = magSqr
658  (
659  hitInfo[pointi].hitPoint()
660  - samples[pointi]
661  );
662  nearestInfo[pointi] = hitInfo[pointi];
663  nearestSurfaces[pointi] = testI;
664  }
665  }
666  }
667 }
668 
669 
670 // Find nearest. Return -1 or nearest point
672 (
673  const PtrList<searchableSurface>& allSurfaces,
674  const labelList& surfacesToTest,
675  const pointField& samples,
676  const scalarField& nearestDistSqr,
677  const labelList& regionIndices,
678  labelList& nearestSurfaces,
679  List<pointIndexHit>& nearestInfo
680 )
681 {
682  if (regionIndices.empty())
683  {
684  findNearest
685  (
686  allSurfaces,
687  surfacesToTest,
688  samples,
689  nearestDistSqr,
690  nearestSurfaces,
691  nearestInfo
692  );
693  }
694 
695  // Initialise
696  nearestSurfaces.setSize(samples.size());
697  nearestSurfaces = -1;
698  nearestInfo.setSize(samples.size());
699 
700  // Work arrays
701  scalarField minDistSqr(nearestDistSqr);
702  List<pointIndexHit> hitInfo(samples.size());
703 
704  forAll(surfacesToTest, testI)
705  {
706  allSurfaces[surfacesToTest[testI]].findNearest
707  (
708  samples,
709  minDistSqr,
710  regionIndices,
711  hitInfo
712  );
713 
714  // Update minDistSqr and arguments
715  forAll(hitInfo, pointi)
716  {
717  if (hitInfo[pointi].hit())
718  {
719  minDistSqr[pointi] = magSqr
720  (
721  hitInfo[pointi].hitPoint()
722  - samples[pointi]
723  );
724  nearestInfo[pointi] = hitInfo[pointi];
725  nearestSurfaces[pointi] = testI;
726  }
727  }
728  }
729 }
730 
731 
733 (
734  const PtrList<searchableSurface>& allSurfaces,
735  const labelList& surfacesToTest,
736  const pointField& samples,
737  const scalarField& nearestDistSqr,
738  const volumeType illegalHandling,
739  labelList& nearestSurfaces,
740  scalarField& distance
741 )
742 {
743  // Initialise
744  distance.setSize(samples.size());
745  distance = -GREAT;
746 
747  // Find nearest
748  List<pointIndexHit> nearestInfo;
749  findNearest
750  (
751  allSurfaces,
752  surfacesToTest,
753  samples,
754  nearestDistSqr,
755  nearestSurfaces,
756  nearestInfo
757  );
758 
759  // Determine sign of nearest. Sort by surface to do this.
760  DynamicField<point> surfPoints(samples.size());
761  DynamicList<label> surfIndices(samples.size());
762 
763  forAll(surfacesToTest, testI)
764  {
765  // Extract samples on this surface
766  surfPoints.clear();
767  surfIndices.clear();
768  forAll(nearestSurfaces, i)
769  {
770  if (nearestSurfaces[i] == testI)
771  {
772  surfPoints.append(samples[i]);
773  surfIndices.append(i);
774  }
775  }
776 
777  // Calculate sideness of these surface points
778  List<volumeType> volType;
779  allSurfaces[surfacesToTest[testI]].getVolumeType(surfPoints, volType);
780 
781  // Push back to original
782  forAll(volType, i)
783  {
784  label pointi = surfIndices[i];
785  scalar dist = mag(samples[pointi] - nearestInfo[pointi].hitPoint());
786 
787  volumeType vT = volType[i];
788 
789  if (vT == volumeType::OUTSIDE)
790  {
791  distance[pointi] = dist;
792  }
793  else if (vT == volumeType::INSIDE)
794  {
795  distance[i] = -dist;
796  }
797  else
798  {
799  switch (illegalHandling)
800  {
801  case volumeType::OUTSIDE:
802  {
803  distance[pointi] = dist;
804  break;
805  }
806  case volumeType::INSIDE:
807  {
808  distance[pointi] = -dist;
809  break;
810  }
811  default:
812  {
814  << "getVolumeType failure,"
815  << " neither INSIDE or OUTSIDE."
816  << " point:" << surfPoints[i]
817  << " surface:"
818  << allSurfaces[surfacesToTest[testI]].name()
819  << " volType:"
820  << volumeType::names[vT]
821  << exit(FatalError);
822  break;
823  }
824  }
825  }
826  }
827  }
828 }
829 
830 
832 (
833  const PtrList<searchableSurface>& allSurfaces,
834  const labelList& surfacesToTest
835 )
836 {
837  pointField bbPoints(2*surfacesToTest.size());
838 
839  forAll(surfacesToTest, testI)
840  {
841  const searchableSurface& surface(allSurfaces[surfacesToTest[testI]]);
842 
843  bbPoints[2*testI] = surface.bounds().min();
844 
845  bbPoints[2*testI + 1] = surface.bounds().max();
846  }
847 
848  return boundBox(bbPoints);
849 }
850 
851 
853 (
854  const PtrList<searchableSurface>& allSurfaces,
855  const labelList& surfacesToTest,
856  const scalar initDistSqr,
857  const scalar convergenceDistSqr,
858  const point& start
859 )
860 {
861  // Get four starting points. Take these as the projection of the
862  // starting point onto the surfaces and the mid point
863  List<point> nearest(surfacesToTest.size()+1);
864 
865  point sumNearest = Zero;
866 
867  forAll(surfacesToTest, i)
868  {
869  pointIndexHit hit
870  (
871  tempFindNearest(allSurfaces[surfacesToTest[i]], start, initDistSqr)
872  );
873 
874  if (hit.hit())
875  {
876  nearest[i] = hit.hitPoint();
877  sumNearest += nearest[i];
878  }
879  else
880  {
882  << "Did not find point within distance "
883  << initDistSqr << " of starting point " << start
884  << " on surface "
885  << allSurfaces[surfacesToTest[i]].IOobject::name()
886  << abort(FatalError);
887  }
888  }
889 
890  nearest.last() = sumNearest / surfacesToTest.size();
891 
892 
893  // Get the sum of distances (initial evaluation)
894  List<scalar> nearestDist(nearest.size());
895 
896  forAll(nearestDist, i)
897  {
898  nearestDist[i] = sumDistSqr
899  (
900  allSurfaces,
901  surfacesToTest,
902  initDistSqr,
903  nearest[i]
904  );
905  }
906 
907 
908  //- Downhill Simplex method
909 
910  bool converged = morphTet
911  (
912  allSurfaces,
913  surfacesToTest,
914  initDistSqr,
915  convergenceDistSqr,
916  2000,
917  nearest,
918  nearestDist
919  );
920 
921 
923 
924  if (converged)
925  {
926  // Project nearest onto 0th surface.
927  intersection = tempFindNearest
928  (
929  allSurfaces[surfacesToTest[0]],
930  nearest[0],
931  nearestDist[0]
932  );
933  }
934 
935  //if (!intersection.hit())
936  //{
937  // // Restart
938  // scalar smallDist = Foam::sqr(convergenceDistSqr);
939  // nearest[0] = intersection.hitPoint();
940  // nearest[1] = nearest[0];
941  // nearest[1].x() += smallDist;
942  // nearest[2] = nearest[0];
943  // nearest[2].y() += smallDist;
944  // nearest[3] = nearest[0];
945  // nearest[3].z() += smallDist;
946  //
947  // forAll(nearestDist, i)
948  // {
949  // nearestDist[i] = sumDistSqr
950  // (
951  // surfacesToTest,
952  // initDistSqr,
953  // nearest[i]
954  // );
955  // }
956  //
957  // intersection = morphTet
958  // (
959  // allSurfaces,
960  // surfacesToTest,
961  // initDistSqr,
962  // convergenceDistSqr,
963  // 1000,
964  // nearest,
965  // nearestDist
966  // );
967  //}
968 
969  return intersection;
970 }
971 
972 
973 
974 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
const boundBox & bounds() const
Return const reference to boundBox.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
static void findNearest(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &)
Find nearest. Return -1 (and a miss()) or surface and nearest.
bool hit() const
Is there a hit.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
label findMin(const ListType &, const label start=0)
Find index of min element (and less than given element).
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
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
static void findAllIntersections(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelListList &surfaces, List< List< pointIndexHit >> &surfaceHits)
Find all intersections in order from start to end. Returns for.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static void signedDistance(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest, const pointField &samples, const scalarField &nearestDistSqr, const volumeType illegalHandling, labelList &nearestSurfaces, scalarField &distance)
Find signed distance to nearest surface. Outside is positive.
Various functions to operate on Lists.
Foam::intersection.
Definition: intersection.H:49
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:96
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void Swap(T &a, T &b)
Definition: Swap.H:43
static void findNearestIntersection(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelList &surface1, List< pointIndexHit > &hit1, labelList &surface2, List< pointIndexHit > &hit2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Dynamically sized Field.
Definition: DynamicField.H:49
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
static void findAnyIntersection(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &)
Find any intersection. Return hit point information and.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< scalar > magSqr(const dimensioned< Type > &)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
static const char nl
Definition: Ostream.H:262
defineTypeNameAndDebug(combustionModel, 0)
static pointIndexHit facesIntersection(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest, const scalar initDistSqr, const scalar convergenceDistSqr, const point &start)
Calculate point which is on a set of surfaces. WIP.
List< word > wordList
A List of words.
Definition: fileName.H:54
void setSize(const label)
Reset size of List.
Definition: List.C:295
vector point
Point is a vector.
Definition: point.H:41
const Point & hitPoint() const
Return hit point.
label findLower(const ListType &, typename ListType::const_reference, const label stary, const BinaryOp &bop)
Find last element < given value in sorted list and return index,.
dimensioned< scalar > mag(const dimensioned< Type > &)
static const NamedEnum< volumeType, 4 > names
Definition: volumeType.H:80
static boundBox bounds(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest)
Find the boundBox of the selected surfaces.
Namespace for OpenFOAM.