searchableSurfaceCollection.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-2020 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 "SortableList.H"
29 #include "Time.H"
30 #include "ListOps.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(searchableSurfaceCollection, 0);
39 (
40  searchableSurface,
41  searchableSurfaceCollection,
42  dict
43 );
44 
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::searchableSurfaceCollection::findNearest
51 (
52  const pointField& samples,
53  scalarField& minDistSqr,
54  List<pointIndexHit>& nearestInfo,
55  labelList& nearestSurf
56 ) const
57 {
58  // Initialise
59  nearestInfo.setSize(samples.size());
60  nearestInfo = pointIndexHit();
61  nearestSurf.setSize(samples.size());
62  nearestSurf = -1;
63 
64  List<pointIndexHit> hitInfo(samples.size());
65 
66  const scalarField localMinDistSqr(samples.size(), great);
67 
68  forAll(subGeom_, surfI)
69  {
70  subGeom_[surfI].findNearest
71  (
72  cmptDivide // Transform then divide
73  (
74  transform_[surfI].localPosition(samples),
75  scale_[surfI]
76  ),
77  localMinDistSqr,
78  hitInfo
79  );
80 
81  forAll(hitInfo, pointi)
82  {
83  if (hitInfo[pointi].hit())
84  {
85  // Rework back into global coordinate sys. Multiply then
86  // transform
87  point globalPt = transform_[surfI].globalPosition
88  (
90  (
91  hitInfo[pointi].rawPoint(),
92  scale_[surfI]
93  )
94  );
95 
96  scalar distSqr = magSqr(globalPt - samples[pointi]);
97 
98  if (distSqr < minDistSqr[pointi])
99  {
100  minDistSqr[pointi] = distSqr;
101  nearestInfo[pointi].setPoint(globalPt);
102  nearestInfo[pointi].setHit();
103  nearestInfo[pointi].setIndex
104  (
105  hitInfo[pointi].index()
106  + indexOffset_[surfI]
107  );
108  nearestSurf[pointi] = surfI;
109  }
110  }
111  }
112  }
113 }
114 
115 
116 // Sort hits into per-surface bins. Misses are rejected. Maintains map back
117 // to position
118 void Foam::searchableSurfaceCollection::sortHits
119 (
120  const List<pointIndexHit>& info,
121  List<List<pointIndexHit>>& surfInfo,
122  labelListList& infoMap
123 ) const
124 {
125  // Count hits per surface.
126  labelList nHits(subGeom_.size(), 0);
127 
128  forAll(info, pointi)
129  {
130  if (info[pointi].hit())
131  {
132  label index = info[pointi].index();
133  label surfI = findLower(indexOffset_, index+1);
134  nHits[surfI]++;
135  }
136  }
137 
138  // Per surface the hit
139  surfInfo.setSize(subGeom_.size());
140  // Per surface the original position
141  infoMap.setSize(subGeom_.size());
142 
143  forAll(surfInfo, surfI)
144  {
145  surfInfo[surfI].setSize(nHits[surfI]);
146  infoMap[surfI].setSize(nHits[surfI]);
147  }
148  nHits = 0;
149 
150  forAll(info, pointi)
151  {
152  if (info[pointi].hit())
153  {
154  label index = info[pointi].index();
155  label surfI = findLower(indexOffset_, index+1);
156 
157  // Store for correct surface and adapt indices back to local
158  // ones
159  label localI = nHits[surfI]++;
160  surfInfo[surfI][localI] = pointIndexHit
161  (
162  info[pointi].hit(),
163  info[pointi].rawPoint(),
164  index-indexOffset_[surfI]
165  );
166  infoMap[surfI][localI] = pointi;
167  }
168  }
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
173 
175 (
176  const IOobject& io,
177  const dictionary& dict
178 )
179 :
180  searchableSurface(io),
181  instance_(dict.size()),
182  scale_(dict.size()),
183  transform_(dict.size()),
184  subGeom_(dict.size()),
185  mergeSubRegions_(dict.lookup("mergeSubRegions")),
186  indexOffset_(dict.size()+1)
187 {
188  Info<< "SearchableCollection : " << name() << endl;
189 
190  label surfI = 0;
191  label startIndex = 0;
192  forAllConstIter(dictionary, dict, iter)
193  {
194  if (dict.isDict(iter().keyword()))
195  {
196  instance_[surfI] = iter().keyword();
197 
198  const dictionary& subDict = dict.subDict(instance_[surfI]);
199 
200  scale_[surfI] = subDict.lookup("scale");
201  transform_.set
202  (
203  surfI,
205  (
206  io.db(),
207  subDict.subDict("transform")
208  )
209  );
210 
211  const word subGeomName(subDict.lookup("surface"));
212  // Pout<< "Trying to find " << subGeomName << endl;
213 
214  const searchableSurface& s =
215  io.db().lookupObject<searchableSurface>(subGeomName);
216 
217  // I don't know yet how to handle the globalSize combined with
218  // regionOffset. Would cause non-consecutive indices locally
219  // if all indices offset by globalSize() of the local region...
220  if (s.size() != s.globalSize())
221  {
223  << "Cannot use a distributed surface in a collection."
224  << exit(FatalError);
225  }
226 
227  subGeom_.set(surfI, &const_cast<searchableSurface&>(s));
228 
229  indexOffset_[surfI] = startIndex;
230  startIndex += subGeom_[surfI].size();
231 
232  Info<< " instance : " << instance_[surfI] << endl;
233  Info<< " surface : " << s.name() << endl;
234  Info<< " scale : " << scale_[surfI] << endl;
235  Info<< " coordsys : " << transform_[surfI] << endl;
236 
237  surfI++;
238  }
239  }
240  indexOffset_[surfI] = startIndex;
241 
242  instance_.setSize(surfI);
243  scale_.setSize(surfI);
244  transform_.setSize(surfI);
245  subGeom_.setSize(surfI);
246  indexOffset_.setSize(surfI+1);
247 
248  // Bounds is the overall bounds
249  bounds() = boundBox(point::max, point::min);
250 
251  forAll(subGeom_, surfI)
252  {
253  const boundBox& surfBb = subGeom_[surfI].bounds();
254 
255  // Transform back to global coordinate sys.
256  const point surfBbMin = transform_[surfI].globalPosition
257  (
259  (
260  surfBb.min(),
261  scale_[surfI]
262  )
263  );
264  const point surfBbMax = transform_[surfI].globalPosition
265  (
267  (
268  surfBb.max(),
269  scale_[surfI]
270  )
271  );
272 
273  bounds().min() = min(bounds().min(), surfBbMin);
274  bounds().max() = max(bounds().max(), surfBbMax);
275  }
276 }
277 
278 
279 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
280 
282 {}
283 
284 
285 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
286 
288 {
289  if (regions_.size() == 0)
290  {
291  regionOffset_.setSize(subGeom_.size());
292 
293  DynamicList<word> allRegions;
294  forAll(subGeom_, surfI)
295  {
296  regionOffset_[surfI] = allRegions.size();
297 
298  if (mergeSubRegions_)
299  {
300  // Single name regardless how many regions subsurface has
301  allRegions.append(instance_[surfI] + "_" + Foam::name(surfI));
302  }
303  else
304  {
305  const wordList& subRegions = subGeom_[surfI].regions();
306 
307  forAll(subRegions, i)
308  {
309  allRegions.append(instance_[surfI] + "_" + subRegions[i]);
310  }
311  }
312  }
313  regions_.transfer(allRegions.shrink());
314  }
315  return regions_;
316 }
317 
318 
320 {
321  return indexOffset_.last();
322 }
323 
324 
327 {
328  tmp<pointField> tCtrs = tmp<pointField>(new pointField(size()));
329  pointField& ctrs = tCtrs.ref();
330 
331  // Append individual coordinates
332  label coordI = 0;
333 
334  forAll(subGeom_, surfI)
335  {
336  const pointField subCoords(subGeom_[surfI].coordinates());
337 
338  forAll(subCoords, i)
339  {
340  ctrs[coordI++] = transform_[surfI].globalPosition
341  (
343  (
344  subCoords[i],
345  scale_[surfI]
346  )
347  );
348  }
349  }
350 
351  return tCtrs;
352 }
353 
354 
356 (
357  pointField& centres,
358  scalarField& radiusSqr
359 ) const
360 {
361  centres.setSize(size());
362  radiusSqr.setSize(centres.size());
363 
364  // Append individual coordinates
365  label coordI = 0;
366 
367  forAll(subGeom_, surfI)
368  {
369  scalar maxScale = cmptMax(scale_[surfI]);
370 
371  pointField subCentres;
372  scalarField subRadiusSqr;
373  subGeom_[surfI].boundingSpheres(subCentres, subRadiusSqr);
374 
375  forAll(subCentres, i)
376  {
377  centres[coordI] = transform_[surfI].globalPosition
378  (
380  (
381  subCentres[i],
382  scale_[surfI]
383  )
384  );
385  radiusSqr[coordI] = maxScale*subRadiusSqr[i];
386  coordI++;
387  }
388  }
389 }
390 
391 
394 {
395  // Get overall size
396  label nPoints = 0;
397 
398  forAll(subGeom_, surfI)
399  {
400  nPoints += subGeom_[surfI].points()().size();
401  }
402 
403  tmp<pointField> tPts(new pointField(nPoints));
404  pointField& pts = tPts.ref();
405 
406  // Append individual coordinates
407  nPoints = 0;
408 
409  forAll(subGeom_, surfI)
410  {
411  const pointField subCoords(subGeom_[surfI].points());
412 
413  forAll(subCoords, i)
414  {
415  pts[nPoints++] = transform_[surfI].globalPosition
416  (
418  (
419  subCoords[i],
420  scale_[surfI]
421  )
422  );
423  }
424  }
425 
426  return tPts;
427 }
428 
429 
430 void Foam::searchableSurfaceCollection::findNearest
431 (
432  const pointField& samples,
433  const scalarField& nearestDistSqr,
434  List<pointIndexHit>& nearestInfo
435 ) const
436 {
437  // How to scale distance?
438  scalarField minDistSqr(nearestDistSqr);
439 
440  labelList nearestSurf;
441  findNearest
442  (
443  samples,
444  minDistSqr,
445  nearestInfo,
446  nearestSurf
447  );
448 }
449 
450 
452 (
453  const pointField& start,
454  const pointField& end,
455  List<pointIndexHit>& info
456 ) const
457 {
458  info.setSize(start.size());
459  info = pointIndexHit();
460 
461  // Current nearest (to start) intersection
462  pointField nearest(end);
463 
464  List<pointIndexHit> hitInfo(start.size());
465 
466  forAll(subGeom_, surfI)
467  {
468  // Starting point
470  (
471  transform_[surfI].localPosition
472  (
473  start
474  ),
475  scale_[surfI]
476  );
477 
478  // Current best end point
480  (
481  transform_[surfI].localPosition
482  (
483  nearest
484  ),
485  scale_[surfI]
486  );
487 
488  subGeom_[surfI].findLine(e0, e1, hitInfo);
489 
490  forAll(hitInfo, pointi)
491  {
492  if (hitInfo[pointi].hit())
493  {
494  // Transform back to global coordinate sys.
495  nearest[pointi] = transform_[surfI].globalPosition
496  (
498  (
499  hitInfo[pointi].rawPoint(),
500  scale_[surfI]
501  )
502  );
503  info[pointi] = hitInfo[pointi];
504  info[pointi].rawPoint() = nearest[pointi];
505  info[pointi].setIndex
506  (
507  hitInfo[pointi].index()
508  + indexOffset_[surfI]
509  );
510  }
511  }
512  }
513 
514 
515  // Debug check
516  if (false)
517  {
518  forAll(info, pointi)
519  {
520  if (info[pointi].hit())
521  {
522  vector n(end[pointi] - start[pointi]);
523  scalar magN = mag(n);
524 
525  if (magN > small)
526  {
527  n /= mag(n);
528 
529  scalar s = ((info[pointi].rawPoint()-start[pointi])&n);
530 
531  if (s < 0 || s > 1)
532  {
534  << "point:" << info[pointi]
535  << " s:" << s
536  << " outside vector "
537  << " start:" << start[pointi]
538  << " end:" << end[pointi]
539  << abort(FatalError);
540  }
541  }
542  }
543  }
544  }
545 }
546 
547 
549 (
550  const pointField& start,
551  const pointField& end,
552  List<pointIndexHit>& info
553 ) const
554 {
555  // To be done ...
556  findLine(start, end, info);
557 }
558 
559 
561 (
562  const pointField& start,
563  const pointField& end,
565 ) const
566 {
567  // To be done. Assume for now only one intersection.
568  List<pointIndexHit> nearestInfo;
569  findLine(start, end, nearestInfo);
570 
571  info.setSize(start.size());
572  forAll(info, pointi)
573  {
574  if (nearestInfo[pointi].hit())
575  {
576  info[pointi].setSize(1);
577  info[pointi][0] = nearestInfo[pointi];
578  }
579  else
580  {
581  info[pointi].clear();
582  }
583  }
584 }
585 
586 
588 (
589  const List<pointIndexHit>& info,
590  labelList& region
591 ) const
592 {
593  if (subGeom_.size() == 0)
594  {}
595  else if (subGeom_.size() == 1)
596  {
597  if (mergeSubRegions_)
598  {
599  region.setSize(info.size());
600  region = regionOffset_[0];
601  }
602  else
603  {
604  subGeom_[0].getRegion(info, region);
605  }
606  }
607  else
608  {
609  // Multiple surfaces. Sort by surface.
610 
611  // Per surface the hit
612  List<List<pointIndexHit>> surfInfo;
613  // Per surface the original position
614  List<List<label>> infoMap;
615  sortHits(info, surfInfo, infoMap);
616 
617  region.setSize(info.size());
618  region = -1;
619 
620  // Do region tests
621 
622  if (mergeSubRegions_)
623  {
624  // Actually no need for surfInfo. Just take region for surface.
625  forAll(infoMap, surfI)
626  {
627  const labelList& map = infoMap[surfI];
628  forAll(map, i)
629  {
630  region[map[i]] = regionOffset_[surfI];
631  }
632  }
633  }
634  else
635  {
636  forAll(infoMap, surfI)
637  {
638  labelList surfRegion;
639  subGeom_[surfI].getRegion(surfInfo[surfI], surfRegion);
640 
641  const labelList& map = infoMap[surfI];
642  forAll(map, i)
643  {
644  region[map[i]] = regionOffset_[surfI] + surfRegion[i];
645  }
646  }
647  }
648  }
649 }
650 
651 
653 (
654  const List<pointIndexHit>& info,
655  vectorField& normal
656 ) const
657 {
658  if (subGeom_.size() == 0)
659  {}
660  else if (subGeom_.size() == 1)
661  {
662  subGeom_[0].getNormal(info, normal);
663  }
664  else
665  {
666  // Multiple surfaces. Sort by surface.
667 
668  // Per surface the hit
669  List<List<pointIndexHit>> surfInfo;
670  // Per surface the original position
671  List<List<label>> infoMap;
672  sortHits(info, surfInfo, infoMap);
673 
674  normal.setSize(info.size());
675 
676  // Do region tests
677  forAll(surfInfo, surfI)
678  {
679  vectorField surfNormal;
680  subGeom_[surfI].getNormal(surfInfo[surfI], surfNormal);
681 
682  // Transform back to global coordinate sys.
683  surfNormal = transform_[surfI].globalVector(surfNormal);
684 
685  const labelList& map = infoMap[surfI];
686  forAll(map, i)
687  {
688  normal[map[i]] = surfNormal[i];
689  }
690  }
691  }
692 }
693 
694 
696 (
697  const pointField& points,
698  List<volumeType>& volType
699 ) const
700 {
702  << "Volume type not supported for collection."
703  << exit(FatalError);
704 }
705 
706 
708 (
709  const List<treeBoundBox>& bbs,
710  const bool keepNonLocal,
711  autoPtr<mapDistribute>& faceMap,
712  autoPtr<mapDistribute>& pointMap
713 )
714 {
715  forAll(subGeom_, surfI)
716  {
717  // Note:Transform the bounding boxes? Something like
718  // pointField bbPoints =
719  // cmptDivide
720  // (
721  // transform_[surfI].localPosition
722  // (
723  // bbs[i].points()
724  // ),
725  // scale_[surfI]
726  // );
727  // treeBoundBox newBb(bbPoints);
728 
729  // Note: what to do with faceMap, pointMap from multiple surfaces?
730  subGeom_[surfI].distribute
731  (
732  bbs,
733  keepNonLocal,
734  faceMap,
735  pointMap
736  );
737  }
738 }
739 
740 
742 {
743  forAll(subGeom_, surfI)
744  {
745  subGeom_[surfI].setField
746  (
747  static_cast<const labelList&>
748  (
750  (
751  values,
752  subGeom_[surfI].size(),
753  indexOffset_[surfI]
754  )
755  )
756  );
757  }
758 }
759 
760 
762 (
763  const List<pointIndexHit>& info,
764  labelList& values
765 ) const
766 {
767  if (subGeom_.size() == 0)
768  {}
769  else if (subGeom_.size() == 1)
770  {
771  subGeom_[0].getField(info, values);
772  }
773  else
774  {
775  // Multiple surfaces. Sort by surface.
776 
777  // Per surface the hit
778  List<List<pointIndexHit>> surfInfo;
779  // Per surface the original position
780  List<List<label>> infoMap;
781  sortHits(info, surfInfo, infoMap);
782 
783  // Do surface tests
784  forAll(surfInfo, surfI)
785  {
786  labelList surfValues;
787  subGeom_[surfI].getField(surfInfo[surfI], surfValues);
788 
789  if (surfValues.size())
790  {
791  // Size values only when we have a surface that supports it.
792  values.setSize(info.size());
793 
794  const labelList& map = infoMap[surfI];
795  forAll(map, i)
796  {
797  values[map[i]] = surfValues[i];
798  }
799  }
800  }
801  }
802 }
803 
804 
805 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual void distribute(const List< treeBoundBox > &, const bool keepNonLocal, autoPtr< mapDistribute > &faceMap, autoPtr< mapDistribute > &pointMap)
Set bounds of surface. Bounds currently set as list of.
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
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
searchableSurfaceCollection(const IOobject &io, const dictionary &dict)
Construct from dictionary (used by searchableSurface)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
virtual void getVolumeType(const pointField &, List< volumeType > &) const
Determine type (inside/outside/mixed) for point. unknown if.
Macros for easy insertion into run-time selection tables.
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:936
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Various functions to operate on Lists.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
A List obtained as a section of another List.
Definition: SubList.H:53
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
const pointField & points
virtual tmp< pointField > points() const
Get the points that define the surface.
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
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
virtual label size() const
Range of local indices that can be returned.
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void setField(const labelList &values)
WIP. Store element-wise field.
void setSize(const label)
Reset size of List.
Definition: List.C:281
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
static autoPtr< coordinateSystem > New(const objectRegistry &obr, const dictionary &dict)
Select constructed from dictionary and objectRegistry.
vector point
Point is a vector.
Definition: point.H:41
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
messageStream Info
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 > &)
label n
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
A class for managing temporary objects.
Definition: PtrList.H:53
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
virtual const wordList & regions() const
Names of regions.