sampledMeshedSurface.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "sampledMeshedSurface.H"
30 #include "meshSearch.H"
31 #include "Tuple2.H"
32 #include "globalIndex.H"
33 #include "treeDataCell.H"
34 #include "treeDataFace.H"
35 #include "meshTools.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 const Foam::Enum
41 <
43 >
44 Foam::sampledMeshedSurface::samplingSourceNames_
45 ({
46  { samplingSource::cells, "cells" },
47  { samplingSource::insideCells, "insideCells" },
48  { samplingSource::boundaryFaces, "boundaryFaces" },
49 });
50 
51 
52 namespace Foam
53 {
54  defineTypeNameAndDebug(sampledMeshedSurface, 0);
55  // Use shorter name only
57  (
58  sampledSurface,
59  sampledMeshedSurface,
60  word,
62  );
63  // Compatibility name (1912)
65  (
66  sampledSurface,
67  sampledMeshedSurface,
68  word,
69  sampledTriSurfaceMesh
70  );
71 
72 } // End namespace Foam
73 
74 
75 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
76 
77 namespace Foam
78 {
79 
80 // The IOobject for reading
81 inline static IOobject selectReadIO(const word& name, const Time& runTime)
82 {
83  return IOobject
84  (
85  name,
86  runTime.constant(), // instance
87  "triSurface", // local
88  runTime, // registry
91  false // no register
92  );
93 }
94 
95 } // End namespace Foam
96 
97 
98 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
99 
100 void Foam::sampledMeshedSurface::setZoneMap()
101 {
102  // Populate zoneIds_ based on surfZone information
103 
104  const meshedSurface& s = static_cast<const meshedSurface&>(*this);
105 
106  const auto& zones = s.surfZones();
107 
108  zoneIds_.resize(s.size());
109 
110  // Trivial case
111  if (zoneIds_.empty() || zones.size() <= 1)
112  {
113  zoneIds_ = 0;
114  return;
115  }
116 
117 
118  label beg = 0;
119 
120  forAll(zones, zonei)
121  {
122  const label len = min(zones[zonei].size(), zoneIds_.size() - beg);
123 
124  // Assign sub-zone Ids
125  SubList<label>(zoneIds_, len, beg) = zonei;
126 
127  beg += len;
128  }
129 
130  // Anything remaining? Should not happen
131  {
132  const label len = (zoneIds_.size() - beg);
133 
134  if (len > 0)
135  {
136  SubList<label>(zoneIds_, len, beg) = max(0, zones.size()-1);
137  }
138  }
139 }
140 
141 
142 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
143 
144 bool Foam::sampledMeshedSurface::update(const meshSearch& meshSearcher)
145 {
146  // Global numbering for cells/faces
147  // - only used to uniquely identify local elements
148  globalIndex globalCells(onBoundary() ? mesh().nFaces() : mesh().nCells());
149 
150  // Find the cells the triangles of the surface are in.
151  // Does approximation by looking at the face centres only
152  const pointField& fc = surface_.faceCentres();
153 
154  // sqr(distance), global index
155  typedef Tuple2<scalar, label> nearInfo;
156  List<nearInfo> nearest(fc.size(), nearInfo(Foam::sqr(GREAT), labelMax));
157 
158  if (sampleSource_ == samplingSource::cells)
159  {
160  // Search for nearest cell
161 
162  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
163 
164  forAll(fc, facei)
165  {
166  const point& pt = fc[facei];
167  auto& near = nearest[facei];
168 
169  pointIndexHit info = cellTree.findNearest(pt, Foam::sqr(GREAT));
170 
171  if (info.hit())
172  {
173  near.first() = magSqr(info.hitPoint()-pt);
174  near.second() = globalCells.toGlobal(info.index());
175  }
176  }
177  }
178  else if (sampleSource_ == samplingSource::insideCells)
179  {
180  // Search for cell containing point
181 
182  const auto& cellTree = meshSearcher.cellTree();
183 
184  forAll(fc, facei)
185  {
186  const point& pt = fc[facei];
187  auto& near = nearest[facei];
188 
189  if (cellTree.bb().contains(pt))
190  {
191  const label index = cellTree.findInside(pt);
192  if (index != -1)
193  {
194  near.first() = 0;
195  near.second() = globalCells.toGlobal(index);
196  }
197  }
198  }
199  }
200  else // samplingSource::boundaryFaces
201  {
202  // Search for nearest boundary face
203  // on all non-coupled boundary faces
204 
205  const auto& bndTree = meshSearcher.nonCoupledBoundaryTree();
206 
207  forAll(fc, facei)
208  {
209  const point& pt = fc[facei];
210  auto& near = nearest[facei];
211 
212  pointIndexHit info = bndTree.findNearest(pt, sqr(GREAT));
213 
214  if (info.hit())
215  {
216  near.first() = magSqr(info.hitPoint()-pt);
217  near.second() =
218  globalCells.toGlobal
219  (
220  bndTree.shapes().faceLabels()[info.index()]
221  );
222  }
223  }
224  }
225 
226 
227  // See which processor has the nearest. Mark and subset
228  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
229 
230  Pstream::listCombineGather(nearest, minFirstEqOp<scalar>{});
232 
233  labelList cellOrFaceLabels(fc.size(), -1);
234 
235  bitSet facesToSubset(fc.size());
236 
237  forAll(nearest, facei)
238  {
239  const auto& near = nearest[facei];
240 
241  const label index = near.second();
242 
243  if (index == labelMax)
244  {
245  // Not found on any processor, or simply too far.
246  // How to map?
247  }
248  else if (globalCells.isLocal(index))
249  {
250  facesToSubset.set(facei);
251 
252  // Mark as special if too far away
253  cellOrFaceLabels[facei] =
254  (
255  (near.first() < maxDistanceSqr_)
256  ? globalCells.toLocal(index)
257  : -1
258  );
259  }
260  }
261 
262 
263  if (debug)
264  {
265  Pout<< "Local out of faces:" << cellOrFaceLabels.size()
266  << " keeping:" << facesToSubset.count() << endl;
267  }
268 
269 
270  // Subset the surface
271  meshedSurface& s = static_cast<meshedSurface&>(*this);
272 
273  labelList pointMap;
275 
276  s = surface_.subsetMesh(facesToSubset, pointMap, faceMap);
277 
278 
279  // Ensure zoneIds_ are indeed correct
280  setZoneMap();
281 
282 
283  // Subset cellOrFaceLabels (for compact faces)
284  cellOrFaceLabels = labelUIndList(cellOrFaceLabels, faceMap)();
285 
286 
287  // Collect the samplePoints and sampleElements
288  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
289 
291  {
292  // With point interpolation
293 
294  samplePoints_ = points();
295  sampleElements_.resize(pointMap.size(), -1);
296 
297  // Store any face per point (without using pointFaces())
298  labelList pointToFace(std::move(pointMap));
299 
300  forAll(s, facei)
301  {
302  const face& f = s[facei];
303 
304  for (const label labi : f)
305  {
306  pointToFace[labi] = facei;
307  }
308  }
309 
310 
311  if (sampleSource_ == samplingSource::cells)
312  {
313  // sampleElements_ : per surface point, the cell
314  // samplePoints_ : per surface point, a location inside the cell
315 
316  forAll(samplePoints_, pointi)
317  {
318  // Use _copy_ of point
319  const point pt = samplePoints_[pointi];
320 
321  const label celli = cellOrFaceLabels[pointToFace[pointi]];
322 
323  sampleElements_[pointi] = celli;
324 
325  if
326  (
327  celli >= 0
328  && !mesh().pointInCell(pt, celli, meshSearcher.decompMode())
329  )
330  {
331  // Point not inside cell
332  // Find nearest point on faces of cell
333 
334  scalar minDistSqr = VGREAT;
335 
336  for (const label facei : mesh().cells()[celli])
337  {
338  const face& f = mesh().faces()[facei];
339 
340  pointHit info =
341  f.nearestPoint
342  (
343  pt,
344  mesh().points()
345  );
346 
347  if (info.distance() < minDistSqr)
348  {
349  minDistSqr = info.distance();
350  samplePoints_[pointi] = info.rawPoint();
351  }
352  }
353  }
354  }
355  }
356  else if (sampleSource_ == samplingSource::insideCells)
357  {
358  // sampleElements_ : per surface point the cell
359  // samplePoints_ : per surface point a location inside the cell
360 
361  forAll(samplePoints_, pointi)
362  {
363  const label celli = cellOrFaceLabels[pointToFace[pointi]];
364 
365  sampleElements_[pointi] = celli;
366  }
367  }
368  else // samplingSource::boundaryFaces
369  {
370  // sampleElements_ : per surface point, the boundary face containing
371  // the location
372  // samplePoints_ : per surface point, a location on the boundary
373 
374  forAll(samplePoints_, pointi)
375  {
376  const point& pt = samplePoints_[pointi];
377 
378  const label facei = cellOrFaceLabels[pointToFace[pointi]];
379 
380  sampleElements_[pointi] = facei;
381 
382  if (facei >= 0)
383  {
384  samplePoints_[pointi] =
385  mesh().faces()[facei].nearestPoint
386  (
387  pt,
388  mesh().points()
389  ).rawPoint();
390  }
391  }
392  }
393  }
394  else
395  {
396  // if sampleSource_ == cells:
397  // sampleElements_ : per surface face, the cell
398  // samplePoints_ : n/a
399  // if sampleSource_ == insideCells:
400  // sampleElements_ : -1 or per surface face, the cell
401  // samplePoints_ : n/a
402  // if sampleSource_ == boundaryFaces:
403  // sampleElements_ : per surface face, the boundary face
404  // samplePoints_ : n/a
405 
406  sampleElements_.transfer(cellOrFaceLabels);
407  samplePoints_.clear();
408  }
409 
410 
411  if (debug)
412  {
413  this->clearOut();
414  OFstream str(mesh().time().path()/"surfaceToMesh.obj");
415  Info<< "Dumping correspondence from local surface (points or faces)"
416  << " to mesh (cells or faces) to " << str.name() << endl;
417 
418  const vectorField& centres =
419  (
420  onBoundary()
421  ? mesh().faceCentres()
422  : mesh().cellCentres()
423  );
424 
426  {
427  label vertI = 0;
428  forAll(samplePoints_, pointi)
429  {
430  meshTools::writeOBJ(str, points()[pointi]);
431  ++vertI;
432 
433  meshTools::writeOBJ(str, samplePoints_[pointi]);
434  ++vertI;
435 
436  label elemi = sampleElements_[pointi];
437  meshTools::writeOBJ(str, centres[elemi]);
438  ++vertI;
439  str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI << nl;
440  }
441  }
442  else
443  {
444  label vertI = 0;
445  forAll(sampleElements_, triI)
446  {
447  meshTools::writeOBJ(str, faceCentres()[triI]);
448  ++vertI;
449 
450  label elemi = sampleElements_[triI];
451  meshTools::writeOBJ(str, centres[elemi]);
452  ++vertI;
453  str << "l " << vertI-1 << ' ' << vertI << nl;
454  }
455  }
456  }
457 
458  needsUpdate_ = false;
459  return true;
460 }
461 
462 
463 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
464 
466 (
467  const word& name,
468  const polyMesh& mesh,
469  const word& surfaceName,
470  const samplingSource sampleSource
471 )
472 :
474  meshedSurface(),
475  surfaceName_(surfaceName),
476  surface_
477  (
478  selectReadIO(surfaceName, mesh.time()),
480  ),
481  sampleSource_(sampleSource),
482  needsUpdate_(true),
483  keepIds_(true),
484  zoneIds_(),
485  sampleElements_(),
486  samplePoints_(),
487  maxDistanceSqr_(Foam::sqr(GREAT)),
488  defaultValues_()
489 {}
490 
491 
493 (
494  const word& name,
495  const polyMesh& mesh,
496  const dictionary& dict
497 )
498 :
500  meshedSurface(),
501  surfaceName_
502  (
503  meshedSurface::findFile
504  (
505  selectReadIO(dict.get<word>("surface"), mesh.time()),
506  dict
507  ).name()
508  ),
509  surface_
510  (
511  selectReadIO(dict.get<word>("surface"), mesh.time()),
512  dict
513  ),
514  sampleSource_(samplingSourceNames_.get("source", dict)),
515  needsUpdate_(true),
516  keepIds_(dict.getOrDefault("keepIds", true)),
517  zoneIds_(),
518  sampleElements_(),
519  samplePoints_(),
520  maxDistanceSqr_(Foam::sqr(GREAT)),
521  defaultValues_(dict.subOrEmptyDict("defaultValue"))
522 {
523  if (dict.readIfPresent("maxDistance", maxDistanceSqr_))
524  {
525  // Info<< "Limit to maxDistance " << maxDistanceSqr_ << nl;
526  // if (defaultValues_.empty())
527  // {
528  // Info<< "defaultValues = zero" << nl;
529  // }
530  // else
531  // {
532  // defaultValues_.writeEntry(Info);
533  // }
534 
535  maxDistanceSqr_ = Foam::sqr(maxDistanceSqr_);
536  }
537 
538  wordRes includePatches;
539  dict.readIfPresent("patches", includePatches);
540  includePatches.uniq();
541 
542  // Could also shift this to the reader itself,
543  // but not yet necessary.
544 
545  if (!includePatches.empty())
546  {
547  Info<< "Subsetting surface " << surfaceName_
548  << " to patches: " << flatOutput(includePatches) << nl;
549 
550  const surfZoneList& zones = surface_.surfZones();
551 
552  const labelList zoneIndices
553  (
555  (
556  zones,
557  includePatches,
558  wordRes(),
560  )
561  );
562 
563  // Faces to subset
564  bitSet includeMap(surface_.size());
565 
566  for (const label zonei : zoneIndices)
567  {
568  const surfZone& zn = zones[zonei];
569  includeMap.set(zn.range());
570  }
571 
572  if (includeMap.none())
573  {
575  << "Patch selection results in an empty surface"
576  << " - ignoring" << nl;
577  }
578  else if (!includeMap.all())
579  {
580  meshedSurface subSurf(surface_.subsetMesh(includeMap));
581 
582  if (subSurf.empty())
583  {
585  << "Bad surface subset (empty)"
586  << " - skip and hope for the best" << nl;
587  }
588  else
589  {
590  // Replace
591  surface_.transfer(subSurf);
592  }
593  }
594  }
595 }
596 
597 
598 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
599 
601 {
602  return needsUpdate_;
603 }
604 
605 
607 {
608  // already marked as expired
609  if (needsUpdate_)
610  {
611  return false;
612  }
613 
616  zoneIds_.clear();
617 
618  sampleElements_.clear();
619  samplePoints_.clear();
620 
621  needsUpdate_ = true;
622  return true;
623 }
624 
625 
627 {
628  if (!needsUpdate_)
629  {
630  return false;
631  }
632 
633  // Calculate surface and mesh overlap bounding box
634  treeBoundBox bb(surface_.points(), surface_.meshPoints());
635 
636  // Check for overlap with (global!) mesh bb
637  const bool intersect = bb.intersect(mesh().bounds());
638 
639  if (!intersect)
640  {
641  // Surface and mesh do not overlap at all. Guarantee a valid
642  // bounding box so we don't get any 'invalid bounding box' errors.
643 
645  << "Surface " << surfaceName_
646  << " does not overlap bounding box of mesh " << mesh().bounds()
647  << endl;
648 
649  bb = treeBoundBox(mesh().bounds());
650  const vector span(bb.span());
651 
652  bb.min() += (0.5-1e-6)*span;
653  bb.max() -= (0.5-1e-6)*span;
654  }
655  else
656  {
657  // Extend a bit
658  const vector span(bb.span());
659  bb.min() -= 0.5*span;
660  bb.max() += 0.5*span;
661 
662  bb.inflate(1e-6);
663  }
664 
665  // Mesh search engine, no triangulation of faces.
666  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
667 
668  return update(meshSearcher);
669 }
670 
671 
673 {
674  if (!needsUpdate_)
675  {
676  return false;
677  }
678 
679  // Mesh search engine on subset, no triangulation of faces.
680  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
681 
682  return update(meshSearcher);
683 }
684 
685 
687 (
688  const interpolation<scalar>& sampler
689 ) const
690 {
691  return sampleOnFaces(sampler);
692 }
693 
694 
696 (
697  const interpolation<vector>& sampler
698 ) const
699 {
700  return sampleOnFaces(sampler);
701 }
702 
703 
705 (
706  const interpolation<sphericalTensor>& sampler
707 ) const
708 {
709  return sampleOnFaces(sampler);
710 }
711 
712 
714 (
715  const interpolation<symmTensor>& sampler
716 ) const
717 {
718  return sampleOnFaces(sampler);
719 }
720 
721 
723 (
724  const interpolation<tensor>& sampler
725 ) const
726 {
727  return sampleOnFaces(sampler);
728 }
729 
730 
732 (
733  const interpolation<scalar>& interpolator
734 ) const
735 {
736  return sampleOnPoints(interpolator);
737 }
738 
739 
741 (
742  const interpolation<vector>& interpolator
743 ) const
744 {
745  return sampleOnPoints(interpolator);
746 }
747 
749 (
750  const interpolation<sphericalTensor>& interpolator
751 ) const
752 {
753  return sampleOnPoints(interpolator);
754 }
755 
756 
758 (
759  const interpolation<symmTensor>& interpolator
760 ) const
761 {
762  return sampleOnPoints(interpolator);
763 }
764 
765 
767 (
768  const interpolation<tensor>& interpolator
769 ) const
770 {
771  return sampleOnPoints(interpolator);
772 }
773 
774 
776 {
777  os << "meshedSurface: " << name() << " :"
778  << " surface:" << surfaceName_;
779 
780  if (level)
781  {
782  os << " faces:" << faces().size()
783  << " points:" << points().size()
784  << " zoneids:" << zoneIds().size();
785  }
786 }
787 
788 
789 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
meshTools.H
Foam::sampledMeshedSurface::print
virtual void print(Ostream &os, int level=0) const
Print information.
Definition: sampledMeshedSurface.C:775
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::sampledMeshedSurface::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledMeshedSurface.C:600
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::labelMax
constexpr label labelMax
Definition: label.H:61
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:60
Foam::sampledMeshedSurface::sampledMeshedSurface
sampledMeshedSurface(const word &name, const polyMesh &mesh, const word &surfaceName, const samplingSource sampleSource)
Construct from components.
Definition: sampledMeshedSurface.C:466
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::List::resize
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
update
mesh update()
Tuple2.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::selectReadIO
static IOobject selectReadIO(const word &name, const Time &runTime)
Definition: sampledMeshedSurface.C:81
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::surfZone::range
labelRange range() const
The start/size range of this zone.
Definition: surfZone.H:152
globalIndex.H
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::MeshedSurface< face >::clear
virtual void clear()
Clear all storage.
Definition: MeshedSurface.C:598
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::sampledMeshedSurface::update
virtual bool update()
Update the surface as required.
Definition: sampledMeshedSurface.C:626
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:392
Foam::boundBox::intersect
bool intersect(const boundBox &bb)
Intersection (union) with the second box.
Definition: boundBox.C:191
treeDataFace.H
Foam::VectorSpace::min
static const Form min
Definition: VectorSpace.H:118
Foam::sampledSurface::interpolate
bool interpolate() const noexcept
Same as isPointData()
Definition: sampledSurface.H:598
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
Foam::nameOp
Extract name (as a word) from an object, typically using its name() method.
Definition: word.H:237
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:121
Foam::stringListOps::findMatching
labelList findMatching(const StringListType &input, const wordRes &allow, const wordRes &deny=wordRes(), AccessOp aop=noOp())
Return ids for items with matching names.
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::interpolation< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::pointHit
PointHit< point > pointHit
A PointIndexHit for 3D points.
Definition: pointHit.H:44
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::sampledSurface::isPointData
bool isPointData() const noexcept
Using interpolation to surface points.
Definition: sampledSurface.H:340
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:216
treeDataCell.H
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:450
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
meshSearch.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::VectorSpace::max
static const Form max
Definition: VectorSpace.H:117
Foam::meshedSurface
MeshedSurface< face > meshedSurface
Definition: MeshedSurfacesFwd.H:41
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< surfZone >
Foam::surfZone
A surface zone on a MeshedSurface.
Definition: surfZone.H:56
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:77
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
Definition: sampledSurface.C:53
Foam::List::set
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:341
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::polyMesh::FACE_PLANES
Definition: polyMesh.H:102
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::sampledMeshedSurface::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledMeshedSurface.C:606
Foam::sampledMeshedSurface::sample
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
Definition: sampledMeshedSurface.C:687
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::point
vector point
Point is a vector.
Definition: point.H:43
sampledMeshedSurface.H
Foam::wordRes::uniq
static wordRes uniq(const UList< wordRe > &input)
Return a wordRes with duplicate entries filtered out.
Definition: wordRes.C:32
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::MeshedSurface< face >
Foam::MeshedSurface< face >::size
label size() const
The surface size is the number of faces.
Definition: MeshedSurface.H:407
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::sampledMeshedSurface::samplingSource
samplingSource
Types of sampling regions.
Definition: sampledMeshedSurface.H:196
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58