conformationSurfaces.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) 2012-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 
26 #include "conformationSurfaces.H"
27 #include "conformalVoronoiMesh.H"
28 #include "triSurface.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 defineTypeNameAndDebug(conformationSurfaces, 0);
36 }
37 
38 
39 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
40 
41 void Foam::conformationSurfaces::hasBoundedVolume
42 (
43  List<volumeType>& referenceVolumeTypes
44 ) const
45 {
46  vector sum(Zero);
47  label totalTriangles = 0;
48 
49  forAll(surfaces_, s)
50  {
51  const searchableSurface& surface(allGeometry_[surfaces_[s]]);
52 
53  if
54  (
55  surface.hasVolumeType()
56  && (
57  normalVolumeTypes_[regionOffset_[s]]
59  )
60  )
61  {
62  pointField pts(1, locationInMesh_);
63 
64  List<volumeType> vTypes
65  (
66  pts.size(),
68  );
69 
70  surface.getVolumeType(pts, vTypes);
71 
72  referenceVolumeTypes[s] = vTypes[0];
73 
74  Info<< " is "
75  << volumeType::names[referenceVolumeTypes[s]]
76  << " surface " << surface.name()
77  << endl;
78  }
79 
80  if (isA<triSurface>(surface))
81  {
82  const triSurface& triSurf = refCast<const triSurface>(surface);
83 
84  const pointField& surfPts = triSurf.points();
85 
86  Info<< " Checking " << surface.name() << endl;
87 
88  label nBaffles = 0;
89 
90  Info<< " Index = " << surfaces_[s] << endl;
91  Info<< " Offset = " << regionOffset_[s] << endl;
92 
93  forAll(triSurf, sI)
94  {
95  const label patchID =
96  triSurf[sI].region()
97  + regionOffset_[s];
98 
99  // Don't include baffle surfaces in the calculation
100  if
101  (
102  normalVolumeTypes_[patchID]
104  )
105  {
106  sum += triSurf[sI].normal(surfPts);
107  }
108  else
109  {
110  nBaffles++;
111  }
112  }
113  Info<< " has " << nBaffles << " baffles out of "
114  << triSurf.size() << " triangles" << endl;
115 
116  totalTriangles += triSurf.size();
117  }
118  }
119 
120  Info<< " Sum of all the surface normals (if near zero, surface is"
121  << " probably closed):" << nl
122  << " Note: Does not include baffle surfaces in calculation" << nl
123  << " Sum = " << sum/(totalTriangles + SMALL) << nl
124  << " mag(Sum) = " << mag(sum)/(totalTriangles + SMALL)
125  << endl;
126 }
127 
128 
129 void Foam::conformationSurfaces::readFeatures
130 (
131  const label surfI,
132  const dictionary& featureDict,
133  const word& surfaceName,
134  label& featureIndex
135 )
136 {
137  word featureMethod =
138  featureDict.lookupOrDefault<word>("featureMethod", "none");
139 
140  if (featureMethod == "extendedFeatureEdgeMesh")
141  {
142  fileName feMeshName(featureDict.lookup("extendedFeatureEdgeMesh"));
143 
144  Info<< " features: " << feMeshName << endl;
145 
146  features_.set
147  (
148  featureIndex,
149  new extendedFeatureEdgeMesh
150  (
151  IOobject
152  (
153  feMeshName,
154  runTime_.time().constant(),
155  "extendedFeatureEdgeMesh",
156  runTime_.time(),
159  )
160  )
161  );
162 
163  featureIndex++;
164  }
165  else if (featureMethod == "extractFeatures")
166  {
167  const searchableSurface& surface = allGeometry_[surfaces_[surfI]];
168 
169  Info<< " features: " << surface.name()
170  << " of type " << surface.type()
171  << ", id: " << featureIndex << endl;
172 
173  autoPtr<searchableSurfaceFeatures> ssFeatures
174  (
175  searchableSurfaceFeatures::New(surface, featureDict)
176  );
177 
178  if (ssFeatures().hasFeatures())
179  {
180  features_.set
181  (
182  featureIndex,
183  ssFeatures().features()
184  );
185 
186  featureIndex++;
187  }
188  else
189  {
191  << surface.name() << " of type "
192  << surface.type() << " does not have features"
193  << endl;
194  }
195  }
196  else if (featureMethod == "none")
197  {
198  // Currently nothing to do
199  }
200  else
201  {
203  << "No valid featureMethod found for surface " << surfaceName
204  << nl << "Use \"extendedFeatureEdgeMesh\" "
205  << "or \"extractFeatures\"."
206  << exit(FatalError);
207  }
208 }
209 
210 void Foam::conformationSurfaces::readFeatures
211 (
212  const dictionary& featureDict,
213  const word& surfaceName,
214  label& featureIndex
215 )
216 {
217  word featureMethod =
218  featureDict.lookupOrDefault<word>("featureMethod", "none");
219 
220  if (featureMethod == "extendedFeatureEdgeMesh")
221  {
222  fileName feMeshName(featureDict.lookup("extendedFeatureEdgeMesh"));
223 
224  Info<< " features: " << feMeshName << ", id: " << featureIndex
225  << endl;
226 
227  features_.set
228  (
229  featureIndex,
230  new extendedFeatureEdgeMesh
231  (
232  IOobject
233  (
234  feMeshName,
235  runTime_.time().constant(),
236  "extendedFeatureEdgeMesh",
237  runTime_.time(),
240  )
241  )
242  );
243 
244  featureIndex++;
245  }
246  else if (featureMethod == "none")
247  {
248  // Currently nothing to do
249  }
250  else
251  {
253  << "No valid featureMethod found for surface " << surfaceName
254  << nl << "Use \"extendedFeatureEdgeMesh\" "
255  << "or \"extractFeatures\"."
256  << exit(FatalError);
257  }
258 }
259 
260 
261 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
262 
263 Foam::conformationSurfaces::conformationSurfaces
264 (
265  const Time& runTime,
266  Random& rndGen,
267  const searchableSurfaces& allGeometry,
268  const dictionary& surfaceConformationDict
269 )
270 :
271  runTime_(runTime),
272  rndGen_(rndGen),
273  allGeometry_(allGeometry),
274  features_(),
275  locationInMesh_(surfaceConformationDict.lookup("locationInMesh")),
276  surfaces_(),
277  allGeometryToSurfaces_(),
278  normalVolumeTypes_(),
279  patchNames_(),
280  surfZones_(),
281  regionOffset_(),
282  patchInfo_(),
283  globalBounds_(),
284  referenceVolumeTypes_(0)
285 {
286  const dictionary& surfacesDict
287  (
288  surfaceConformationDict.subDict("geometryToConformTo")
289  );
290 
291  const dictionary& additionalFeaturesDict
292  (
293  surfaceConformationDict.subDict("additionalFeatures")
294  );
295 
296 
297  // Wildcard specification : loop over all surface, all regions
298  // and try to find a match.
299 
300  // Count number of surfaces.
301  label surfI = 0;
302  forAll(allGeometry.names(), geomI)
303  {
304  const word& geomName = allGeometry_.names()[geomI];
305 
306  if (surfacesDict.found(geomName))
307  {
308  surfI++;
309  }
310  }
311 
312  const label nAddFeat = additionalFeaturesDict.size();
313 
314  Info<< nl << "Reading geometryToConformTo" << endl;
315 
316  allGeometryToSurfaces_.setSize(allGeometry_.size(), -1);
317 
318  normalVolumeTypes_.setSize(surfI);
319  surfaces_.setSize(surfI);
320  surfZones_.setSize(surfI);
321 
322  // Features may be attached to host surfaces or independent
323  features_.setSize(surfI + nAddFeat);
324 
325  label featureI = 0;
326 
327  regionOffset_.setSize(surfI, 0);
328 
329  PtrList<dictionary> globalPatchInfo(surfI);
330  List<Map<autoPtr<dictionary>>> regionPatchInfo(surfI);
331  List<sideVolumeType> globalVolumeTypes(surfI);
332  List<Map<sideVolumeType>> regionVolumeTypes(surfI);
333 
334  HashSet<word> unmatchedKeys(surfacesDict.toc());
335 
336  surfI = 0;
337  forAll(allGeometry_.names(), geomI)
338  {
339  const word& geomName = allGeometry_.names()[geomI];
340 
341  const entry* ePtr = surfacesDict.lookupEntryPtr(geomName, false, true);
342 
343  if (ePtr)
344  {
345  const dictionary& dict = ePtr->dict();
346  unmatchedKeys.erase(ePtr->keyword());
347 
348  surfaces_[surfI] = geomI;
349 
350  const searchableSurface& surface = allGeometry_[surfaces_[surfI]];
351 
352  // Surface zones
353  if (dict.found("faceZone"))
354  {
355  surfZones_.set(surfI, new surfaceZonesInfo(surface, dict));
356  }
357 
358  allGeometryToSurfaces_[surfaces_[surfI]] = surfI;
359 
360  Info<< nl << " " << geomName << endl;
361 
362  const wordList& regionNames =
363  allGeometry_.regionNames()[surfaces_[surfI]];
364 
365  patchNames_.append(regionNames);
366 
367  globalVolumeTypes[surfI] =
368  (
370  [
371  dict.lookupOrDefault<word>
372  (
373  "meshableSide",
374  "inside"
375  )
376  ]
377  );
378 
379  if (!globalVolumeTypes[surfI])
380  {
381  if (!surface.hasVolumeType())
382  {
384  << "Non-baffle surface "
385  << surface.name()
386  << " does not allow inside/outside queries."
387  << " This usually is an error." << endl;
388  }
389  }
390 
391  // Load patch info
392  if (dict.found("patchInfo"))
393  {
394  globalPatchInfo.set
395  (
396  surfI,
397  dict.subDict("patchInfo").clone()
398  );
399  }
400 
401  readFeatures
402  (
403  surfI,
404  dict,
405  geomName,
406  featureI
407  );
408 
409  const wordList& rNames = surface.regions();
410 
411  if (dict.found("regions"))
412  {
413  const dictionary& regionsDict = dict.subDict("regions");
414 
415  forAll(rNames, regionI)
416  {
417  const word& regionName = rNames[regionI];
418 
419  if (regionsDict.found(regionName))
420  {
421  Info<< " region " << regionName << endl;
422 
423  // Get the dictionary for region
424  const dictionary& regionDict = regionsDict.subDict
425  (
426  regionName
427  );
428 
429  if (regionDict.found("patchInfo"))
430  {
431  regionPatchInfo[surfI].insert
432  (
433  regionI,
434  regionDict.subDict("patchInfo").clone()
435  );
436  }
437 
438  regionVolumeTypes[surfI].insert
439  (
440  regionI,
442  [
443  regionDict.lookupOrDefault<word>
444  (
445  "meshableSide",
446  "inside"
447  )
448  ]
449  );
450 
451  readFeatures(regionDict, regionName, featureI);
452  }
453  }
454  }
455 
456  surfI++;
457  }
458  }
459 
460 
461  if (unmatchedKeys.size() > 0)
462  {
464  (
465  surfacesDict
466  ) << "Not all entries in conformationSurfaces dictionary were used."
467  << " The following entries were not used : "
468  << unmatchedKeys.sortedToc()
469  << endl;
470  }
471 
472 
473  // Calculate local to global region offset
474  label nRegions = 0;
475 
476  forAll(surfaces_, surfI)
477  {
478  regionOffset_[surfI] = nRegions;
479 
480  const searchableSurface& surface = allGeometry_[surfaces_[surfI]];
481  nRegions += surface.regions().size();
482  }
483 
484  // Rework surface specific information into information per global region
485  patchInfo_.setSize(nRegions);
486  normalVolumeTypes_.setSize(nRegions);
487 
488  forAll(surfaces_, surfI)
489  {
490  const searchableSurface& surface = allGeometry_[surfaces_[surfI]];
491 
492  label nRegions = surface.regions().size();
493 
494  // Initialise to global (i.e. per surface)
495  for (label i = 0; i < nRegions; i++)
496  {
497  label globalRegionI = regionOffset_[surfI] + i;
498  normalVolumeTypes_[globalRegionI] = globalVolumeTypes[surfI];
499  if (globalPatchInfo.set(surfI))
500  {
501  patchInfo_.set
502  (
503  globalRegionI,
504  globalPatchInfo[surfI].clone()
505  );
506  }
507  }
508 
509  forAllConstIter(Map<sideVolumeType>, regionVolumeTypes[surfI], iter)
510  {
511  label globalRegionI = regionOffset_[surfI] + iter.key();
512 
513  normalVolumeTypes_[globalRegionI] =
514  regionVolumeTypes[surfI][iter.key()];
515  }
516 
517  const Map<autoPtr<dictionary>>& localInfo = regionPatchInfo[surfI];
518  forAllConstIter(Map<autoPtr<dictionary>>, localInfo, iter)
519  {
520  label globalRegionI = regionOffset_[surfI] + iter.key();
521 
522  patchInfo_.set(globalRegionI, iter()().clone());
523  }
524  }
525 
526 
527 
528  if (!additionalFeaturesDict.empty())
529  {
530  Info<< nl << "Reading additionalFeatures" << endl;
531  }
532 
533  forAllConstIter(dictionary, additionalFeaturesDict, iter)
534  {
535  word featureName = iter().keyword();
536 
537  Info<< nl << " " << iter().keyword() << endl;
538 
539  const dictionary& featureSubDict
540  (
541  additionalFeaturesDict.subDict(featureName)
542  );
543 
544  readFeatures(featureSubDict, featureName, featureI);
545  }
546 
547  // Remove unnecessary space from the features list
548  features_.setSize(featureI);
549 
550  globalBounds_ = treeBoundBox
551  (
552  searchableSurfacesQueries::bounds(allGeometry_, surfaces_)
553  );
554 
555  // Extend the global bounds to stop the bound box sitting on the surfaces
556  // to be conformed to
557  //globalBounds_ = globalBounds_.extend(rndGen_, 1e-4);
558 
559  vector newSpan = 1e-4*globalBounds_.span();
560 
561  globalBounds_.min() -= newSpan;
562  globalBounds_.max() += newSpan;
563 
564  // Look at all surfaces at determine whether the locationInMesh point is
565  // inside or outside each, to establish a signature for the domain to be
566  // meshed.
567 
568  referenceVolumeTypes_.setSize
569  (
570  surfaces_.size(),
572  );
573 
574  Info<< endl
575  << "Testing for locationInMesh " << locationInMesh_ << endl;
576 
577  hasBoundedVolume(referenceVolumeTypes_);
578 
579  if (debug)
580  {
581  Info<< "Names = " << allGeometry_.names() << endl;
582  Info<< "Surfaces = " << surfaces_ << endl;
583  Info<< "AllGeom to Surfaces = " << allGeometryToSurfaces_ << endl;
584  Info<< "Volume types = " << normalVolumeTypes_ << endl;
585  Info<< "Patch names = " << patchNames_ << endl;
586  Info<< "Region Offset = " << regionOffset_ << endl;
587 
588  forAll(features_, fI)
589  {
590  Info<< features_[fI].name() << endl;
591  }
592  }
593 }
594 
595 
596 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
597 
599 {}
600 
601 
602 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
603 
604 bool Foam::conformationSurfaces::overlaps(const treeBoundBox& bb) const
605 {
606  forAll(surfaces_, s)
607  {
608  if (allGeometry_[surfaces_[s]].overlaps(bb))
609  {
610  return true;
611  }
612  }
613 
614  return false;
615 }
616 
617 
619 (
620  const pointField& samplePts
621 ) const
622 {
623  return wellInside(samplePts, scalarField(samplePts.size(), 0.0));
624 }
625 
626 
628 (
629  const point& samplePt
630 ) const
631 {
632  return wellInside(pointField(1, samplePt), scalarField(1, 0))[0];
633 }
634 
635 
637 (
638  const pointField& samplePts
639 ) const
640 {
641  return wellOutside(samplePts, scalarField(samplePts.size(), 0.0));
642 }
643 
644 
646 (
647  const point& samplePt
648 ) const
649 {
650  return wellOutside(pointField(1, samplePt), scalarField(1, 0))[0];
651  //return !inside(samplePt);
652 }
653 
654 
656 (
657  const pointField& samplePts,
658  const scalarField& testDistSqr,
659  const bool testForInside
660 ) const
661 {
662  List<List<volumeType>> surfaceVolumeTests
663  (
664  surfaces_.size(),
665  List<volumeType>
666  (
667  samplePts.size(),
669  )
670  );
671 
672  // Get lists for the volumeTypes for each sample wrt each surface
673  forAll(surfaces_, s)
674  {
675  const searchableSurface& surface(allGeometry_[surfaces_[s]]);
676 
677  const label regionI = regionOffset_[s];
678 
679  if (normalVolumeTypes_[regionI] != extendedFeatureEdgeMesh::BOTH)
680  {
681  surface.getVolumeType(samplePts, surfaceVolumeTests[s]);
682  }
683  }
684 
685  // Compare the volumeType result for each point wrt to each surface with the
686  // reference value and if the points are inside the surface by a given
687  // distanceSquared
688 
689  // Assume that the point is wellInside until demonstrated otherwise.
690  Field<bool> insideOutsidePoint(samplePts.size(), testForInside);
691 
692  //Check if the points are inside the surface by the given distance squared
693 
694  labelList hitSurfaces;
695  List<pointIndexHit> hitInfo;
697  (
698  allGeometry_,
699  surfaces_,
700  samplePts,
701  testDistSqr,
702  hitSurfaces,
703  hitInfo
704  );
705 
706  forAll(samplePts, i)
707  {
708  const pointIndexHit& pHit = hitInfo[i];
709 
710  if (pHit.hit())
711  {
712  // If the point is within range of the surface, then it can't be
713  // well (in|out)side
714  insideOutsidePoint[i] = false;
715 
716  continue;
717  }
718 
719  forAll(surfaces_, s)
720  {
721  const label regionI = regionOffset_[s];
722 
723  if (normalVolumeTypes_[regionI] == extendedFeatureEdgeMesh::BOTH)
724  {
725  continue;
726  }
727 
728  const searchableSurface& surface(allGeometry_[surfaces_[s]]);
729 
730  if
731  (
732  !surface.hasVolumeType()
733  //&& surfaceVolumeTests[s][i] == volumeType::UNKNOWN
734  )
735  {
736  pointField sample(1, samplePts[i]);
737  scalarField nearestDistSqr(1, GREAT);
738  List<pointIndexHit> info;
739 
740  surface.findNearest(sample, nearestDistSqr, info);
741 
742  vector hitDir = info[0].rawPoint() - samplePts[i];
743  hitDir /= mag(hitDir) + SMALL;
744 
745  pointIndexHit surfHit;
746  label hitSurface;
747 
748  findSurfaceNearestIntersection
749  (
750  samplePts[i],
751  info[0].rawPoint() - 1e-3*mag(hitDir)*hitDir,
752  surfHit,
753  hitSurface
754  );
755 
756  if (surfHit.hit() && hitSurface != surfaces_[s])
757  {
758  continue;
759  }
760  }
761 
762  if (surfaceVolumeTests[s][i] == volumeType::OUTSIDE)
763  {
764  if
765  (
766  normalVolumeTypes_[regionI]
768  )
769  {
770  insideOutsidePoint[i] = !testForInside;
771  break;
772  }
773  }
774  else if (surfaceVolumeTests[s][i] == volumeType::INSIDE)
775  {
776  if
777  (
778  normalVolumeTypes_[regionI]
780  )
781  {
782  insideOutsidePoint[i] = !testForInside;
783  break;
784  }
785  }
786  }
787  }
788 
789  return insideOutsidePoint;
790 }
791 
792 
794 (
795  const pointField& samplePts,
796  const scalarField& testDistSqr
797 ) const
798 {
799  return wellInOutSide(samplePts, testDistSqr, true);
800 }
801 
802 
804 (
805  const point& samplePt,
806  scalar testDistSqr
807 ) const
808 {
809  return wellInside(pointField(1, samplePt), scalarField(1, testDistSqr))[0];
810 }
811 
812 
814 (
815  const pointField& samplePts,
816  const scalarField& testDistSqr
817 ) const
818 {
819  return wellInOutSide(samplePts, testDistSqr, false);
820 }
821 
822 
824 (
825  const point& samplePt,
826  scalar testDistSqr
827 ) const
828 {
829  return wellOutside(pointField(1, samplePt), scalarField(1, testDistSqr))[0];
830 }
831 
832 
834 (
835  const point& start,
836  const point& end
837 ) const
838 {
839  labelList hitSurfaces;
840  List<pointIndexHit> hitInfo;
841 
843  (
844  allGeometry_,
845  surfaces_,
846  pointField(1, start),
847  pointField(1, end),
848  hitSurfaces,
849  hitInfo
850  );
851 
852  return hitInfo[0].hit();
853 }
854 
855 
857 (
858  const point& start,
859  const point& end,
860  pointIndexHit& surfHit,
861  label& hitSurface
862 ) const
863 {
864  labelList hitSurfaces;
865  List<pointIndexHit> hitInfo;
866 
868  (
869  allGeometry_,
870  surfaces_,
871  pointField(1, start),
872  pointField(1, end),
873  hitSurfaces,
874  hitInfo
875  );
876 
877  surfHit = hitInfo[0];
878 
879  if (surfHit.hit())
880  {
881  // hitSurfaces has returned the index of the entry in surfaces_ that was
882  // found, not the index of the surface in allGeometry_, translating this
883  // to allGeometry_
884 
885  hitSurface = surfaces_[hitSurfaces[0]];
886  }
887 }
888 
889 
891 (
892  const point& start,
893  const point& end,
894  List<pointIndexHit>& surfHit,
895  labelList& hitSurface
896 ) const
897 {
898  labelListList hitSurfaces;
899  List<List<pointIndexHit>> hitInfo;
900 
902  (
903  allGeometry_,
904  surfaces_,
905  pointField(1, start),
906  pointField(1, end),
907  hitSurfaces,
908  hitInfo
909  );
910 
911  surfHit = hitInfo[0];
912 
913  hitSurface.setSize(hitSurfaces[0].size());
914 
915  forAll(hitSurfaces[0], surfI)
916  {
917  // hitSurfaces has returned the index of the entry in surfaces_ that was
918  // found, not the index of the surface in allGeometry_, translating this
919  // to allGeometry_
920 
921  hitSurface[surfI] = surfaces_[hitSurfaces[0][surfI]];
922  }
923 }
924 
925 
927 (
928  const point& start,
929  const point& end,
930  pointIndexHit& surfHit,
931  label& hitSurface
932 ) const
933 {
934  labelList hitSurfacesStart;
935  List<pointIndexHit> hitInfoStart;
936  labelList hitSurfacesEnd;
937  List<pointIndexHit> hitInfoEnd;
938 
940  (
941  allGeometry_,
942  surfaces_,
943  pointField(1, start),
944  pointField(1, end),
945  hitSurfacesStart,
946  hitInfoStart,
947  hitSurfacesEnd,
948  hitInfoEnd
949  );
950 
951  surfHit = hitInfoStart[0];
952 
953  if (surfHit.hit())
954  {
955  // hitSurfaces has returned the index of the entry in surfaces_ that was
956  // found, not the index of the surface in allGeometry_, translating this
957  // to allGeometry_
958 
959  hitSurface = surfaces_[hitSurfacesStart[0]];
960  }
961 }
962 
963 
965 (
966  const point& sample,
967  scalar nearestDistSqr,
968  pointIndexHit& surfHit,
969  label& hitSurface
970 ) const
971 {
972  labelList hitSurfaces;
973  List<pointIndexHit> surfaceHits;
974 
976  (
977  allGeometry_,
978  surfaces_,
979  pointField(1, sample),
980  scalarField(1, nearestDistSqr),
981  hitSurfaces,
982  surfaceHits
983  );
984 
985  surfHit = surfaceHits[0];
986 
987  if (surfHit.hit())
988  {
989  // hitSurfaces has returned the index of the entry in surfaces_ that was
990  // found, not the index of the surface in allGeometry_, translating this
991  // to allGeometry_
992 
993  hitSurface = surfaces_[hitSurfaces[0]];
994  }
995 }
996 
997 
999 (
1000  const pointField& samples,
1001  const scalarField& nearestDistSqr,
1002  List<pointIndexHit>& surfaceHits,
1003  labelList& hitSurfaces
1004 ) const
1005 {
1007  (
1008  allGeometry_,
1009  surfaces_,
1010  samples,
1011  nearestDistSqr,
1012  hitSurfaces,
1013  surfaceHits
1014  );
1015 
1016  forAll(surfaceHits, i)
1017  {
1018  if (surfaceHits[i].hit())
1019  {
1020  // hitSurfaces has returned the index of the entry in surfaces_ that
1021  // was found, not the index of the surface in allGeometry_,
1022  // translating this to the surface in allGeometry_.
1023 
1024  hitSurfaces[i] = surfaces_[hitSurfaces[i]];
1025  }
1026  }
1027 }
1028 
1029 
1031 (
1032  const point& sample,
1033  scalar nearestDistSqr,
1034  pointIndexHit& fpHit,
1035  label& featureHit
1036 ) const
1037 {
1038  // Work arrays
1039  scalar minDistSqr = nearestDistSqr;
1040  pointIndexHit hitInfo;
1041 
1042  forAll(features_, testI)
1043  {
1044  features_[testI].nearestFeaturePoint
1045  (
1046  sample,
1047  minDistSqr,
1048  hitInfo
1049  );
1050 
1051  if (hitInfo.hit())
1052  {
1053  minDistSqr = magSqr(hitInfo.hitPoint()- sample);
1054  fpHit = hitInfo;
1055  featureHit = testI;
1056  }
1057  }
1058 }
1059 
1060 
1062 (
1063  const point& sample,
1064  scalar nearestDistSqr,
1065  pointIndexHit& edgeHit,
1066  label& featureHit
1067 ) const
1068 {
1069  pointField samples(1, sample);
1070  scalarField nearestDistsSqr(1, nearestDistSqr);
1071 
1072  List<pointIndexHit> edgeHits;
1073  labelList featuresHit;
1074 
1075  findEdgeNearest
1076  (
1077  samples,
1078  nearestDistsSqr,
1079  edgeHits,
1080  featuresHit
1081  );
1082 
1083  edgeHit = edgeHits[0];
1084  featureHit = featuresHit[0];
1085 }
1086 
1087 
1089 (
1090  const pointField& samples,
1091  const scalarField& nearestDistsSqr,
1092  List<pointIndexHit>& edgeHits,
1093  labelList& featuresHit
1094 ) const
1095 {
1096  // Initialise
1097  featuresHit.setSize(samples.size());
1098  featuresHit = -1;
1099  edgeHits.setSize(samples.size());
1100 
1101  // Work arrays
1102  scalarField minDistSqr(nearestDistsSqr);
1103  List<pointIndexHit> hitInfo(samples.size());
1104 
1105  forAll(features_, testI)
1106  {
1107  features_[testI].nearestFeatureEdge
1108  (
1109  samples,
1110  minDistSqr,
1111  hitInfo
1112  );
1113 
1114  // Update minDistSqr and arguments
1115  forAll(hitInfo, pointi)
1116  {
1117  if (hitInfo[pointi].hit())
1118  {
1119  minDistSqr[pointi] = magSqr
1120  (
1121  hitInfo[pointi].hitPoint()
1122  - samples[pointi]
1123  );
1124  edgeHits[pointi] = hitInfo[pointi];
1125  featuresHit[pointi] = testI;
1126  }
1127  }
1128  }
1129 }
1130 
1131 
1133 (
1134  const point& sample,
1135  scalar nearestDistSqr,
1136  List<pointIndexHit>& edgeHits,
1137  List<label>& featuresHit
1138 ) const
1139 {
1140  // Initialise
1141  featuresHit.setSize(extendedFeatureEdgeMesh::nEdgeTypes);
1142  featuresHit = -1;
1143  edgeHits.setSize(extendedFeatureEdgeMesh::nEdgeTypes);
1144 
1145  // Work arrays
1146  scalarField minDistSqr(extendedFeatureEdgeMesh::nEdgeTypes, nearestDistSqr);
1147  List<pointIndexHit> hitInfo(extendedFeatureEdgeMesh::nEdgeTypes);
1148 
1149  forAll(features_, testI)
1150  {
1151  features_[testI].nearestFeatureEdgeByType
1152  (
1153  sample,
1154  minDistSqr,
1155  hitInfo
1156  );
1157 
1158  // Update minDistSqr and arguments
1159  forAll(hitInfo, typeI)
1160  {
1161  if (hitInfo[typeI].hit())
1162  {
1163  minDistSqr[typeI] = magSqr(hitInfo[typeI].hitPoint() - sample);
1164  edgeHits[typeI] = hitInfo[typeI];
1165  featuresHit[typeI] = testI;
1166  }
1167  }
1168  }
1169 }
1170 
1171 
1173 (
1174  const point& sample,
1175  const scalar searchRadiusSqr,
1176  List<List<pointIndexHit>>& edgeHitsByFeature,
1177  List<label>& featuresHit
1178 ) const
1179 {
1180  // Initialise
1181  //featuresHit.setSize(features_.size());
1182  //featuresHit = -1;
1183  //edgeHitsByFeature.setSize(features_.size());
1184 
1185  // Work arrays
1186  List<pointIndexHit> hitInfo(extendedFeatureEdgeMesh::nEdgeTypes);
1187 
1188  forAll(features_, testI)
1189  {
1190  features_[testI].allNearestFeatureEdges
1191  (
1192  sample,
1193  searchRadiusSqr,
1194  hitInfo
1195  );
1196 
1197  bool anyHit = false;
1198  forAll(hitInfo, hitI)
1199  {
1200  if (hitInfo[hitI].hit())
1201  {
1202  anyHit = true;
1203  }
1204  }
1205 
1206  if (anyHit)
1207  {
1208  edgeHitsByFeature.append(hitInfo);
1209  featuresHit.append(testI);
1210  }
1211  }
1212 }
1213 
1214 
1215 void Foam::conformationSurfaces::writeFeatureObj(const fileName& prefix) const
1216 {
1217  OFstream ftStr(runTime_.time().path()/prefix + "_allFeatures.obj");
1218 
1219  Pout<< nl << "Writing all features to " << ftStr.name() << endl;
1220 
1221  label verti = 0;
1222 
1223  forAll(features_, i)
1224  {
1225  const extendedFeatureEdgeMesh& fEM(features_[i]);
1226  const pointField pts(fEM.points());
1227  const edgeList eds(fEM.edges());
1228 
1229  ftStr << "g " << fEM.name() << endl;
1230 
1231  forAll(eds, j)
1232  {
1233  const edge& e = eds[j];
1234 
1235  meshTools::writeOBJ(ftStr, pts[e[0]]); verti++;
1236  meshTools::writeOBJ(ftStr, pts[e[1]]); verti++;
1237  ftStr << "l " << verti-1 << ' ' << verti << endl;
1238  }
1239  }
1240 }
1241 
1242 
1244 (
1245  const point& ptA,
1246  const point& ptB
1247 ) const
1248 {
1249  pointIndexHit surfHit;
1250  label hitSurface;
1251 
1252  findSurfaceAnyIntersection(ptA, ptB, surfHit, hitSurface);
1253 
1254  return getPatchID(hitSurface, surfHit);
1255 }
1256 
1257 
1259 {
1260  pointIndexHit surfHit;
1261  label hitSurface;
1262 
1263  findSurfaceNearest(pt, sqr(GREAT), surfHit, hitSurface);
1264 
1265  return getPatchID(hitSurface, surfHit);
1266 }
1267 
1268 
1270 (
1271  const label hitSurface,
1272  const pointIndexHit& surfHit
1273 ) const
1274 {
1275  if (!surfHit.hit())
1276  {
1277  return -1;
1278  }
1279 
1280  labelList surfLocalRegion;
1281 
1282  allGeometry_[hitSurface].getRegion
1283  (
1284  List<pointIndexHit>(1, surfHit),
1285  surfLocalRegion
1286  );
1287 
1288  const label patchID =
1289  surfLocalRegion[0]
1290  + regionOffset_[allGeometryToSurfaces_[hitSurface]];
1291 
1292  return patchID;
1293 }
1294 
1295 
1298 (
1299  const label hitSurface,
1300  const pointIndexHit& surfHit
1301 ) const
1302 {
1303  const label patchID = getPatchID(hitSurface, surfHit);
1304 
1305  if (patchID == -1)
1306  {
1308  }
1309 
1310  return normalVolumeTypes_[patchID];
1311 }
1312 
1313 
1315 (
1316  const label hitSurface,
1317  const List<pointIndexHit>& surfHit,
1318  vectorField& normal
1319 ) const
1320 {
1321  allGeometry_[hitSurface].getNormal(surfHit, normal);
1322 
1323  const label patchID = regionOffset_[allGeometryToSurfaces_[hitSurface]];
1324 
1325  // Now flip sign of normal depending on mesh side
1326  if (normalVolumeTypes_[patchID] == extendedFeatureEdgeMesh::OUTSIDE)
1327  {
1328  normal *= -1;
1329  }
1330 }
1331 
1332 
1333 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
Field< bool > outside(const pointField &samplePts) const
Check if points are outside surfaces to conform to.
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
Field< bool > wellOutside(const pointField &samplePts, const scalarField &testDistSqr) const
Check if point is outside surfaces to conform to by at least.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void findAllNearestEdges(const point &sample, const scalar searchRadiusSqr, List< List< pointIndexHit >> &edgeHitsByFeature, List< label > &featuresHit) const
Find the nearest points on each feature edge that is within.
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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.
static autoPtr< searchableSurfaceFeatures > New(const searchableSurface &surface, const dictionary &dict)
Return a reference to the selected searchableSurfaceFeatures.
Field< bool > wellInOutSide(const pointField &samplePts, const scalarField &testDistSqr, bool testForInside) const
Check if point is closer to the surfaces to conform to than.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const Foam::NamedEnum< sideVolumeType, 4 > sideVolumeTypeNames_
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool overlaps(const treeBoundBox &bb) const
Check if the supplied bound box overlaps any part of any of.
void writeFeatureObj(const fileName &prefix) const
Write all components of all the extendedFeatureEdgeMeshes as.
Field< bool > inside(const pointField &samplePts) const
Check if points are inside surfaces to conform to.
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
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.
void getNormal(const label hitSurface, const List< pointIndexHit > &surfHit, vectorField &normal) const
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
label getPatchID(const label hitSurface, const pointIndexHit &surfHit) const
Get the region number of a hit surface.
bool findSurfaceAnyIntersection(const point &start, const point &end) const
scalarField samples(nIntervals, 0)
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
stressControl lookup("compactNormalStress") >> compactNormalStress
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)
List< edge > edgeList
Definition: edgeList.H:38
Pre-declare SubField and related Field type.
Definition: Field.H:57
static label nEdgeTypes
Number of possible feature edge types (i.e. number of slices)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void findSurfaceNearest(const point &sample, scalar nearestDistSqr, pointIndexHit &surfHit, label &hitSurface) const
Find the nearest point to the sample and return it to the.
void findEdgeNearest(const point &sample, scalar nearestDistSqr, pointIndexHit &edgeHit, label &featureHit) const
Find the nearest point on any feature edge.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:97
sideVolumeType
Normals point to the outside.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
void findSurfaceAllIntersections(const point &start, const point &end, List< pointIndexHit > &surfHit, labelList &hitSurface) const
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.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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)
~conformationSurfaces()
Destructor.
label findPatch(const point &ptA, const point &ptB) const
Find which patch is intersected by the line from one point to.
List< word > wordList
A List of words.
Definition: fileName.H:54
void findSurfaceNearestIntersection(const point &start, const point &end, pointIndexHit &surfHit, label &hitSurface) const
Finding the nearestIntersection of the surface to start.
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
void findFeaturePointNearest(const point &sample, scalar nearestDistSqr, pointIndexHit &fpHit, label &featureHit) const
Find the nearest point on any feature edge.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
static const NamedEnum< volumeType, 4 > names
Definition: volumeType.H:80
Field< bool > wellInside(const pointField &samplePts, const scalarField &testDistSqr) const
Check if point is inside surfaces to conform to by at least.
void findEdgeNearestByType(const point &sample, scalar nearestDistSqr, List< pointIndexHit > &edgeHit, List< label > &featuresHit) const
Find the nearest point on each type of feature edge.
extendedFeatureEdgeMesh::sideVolumeType meshableSide(const label hitSurface, const pointIndexHit &surfHit) const
Is the surface a baffle.
static boundBox bounds(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest)
Find the boundBox of the selected surfaces.
Namespace for OpenFOAM.