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-2020 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  //- Private class for finding nearest
73  // Comprising:
74  // - global index
75  // - sqr(distance)
77 
78  class nearestEqOp
79  {
80  public:
81 
82  void operator()(nearInfo& x, const nearInfo& y) const
83  {
84  if (y.first() < x.first())
85  {
86  x = y;
87  }
88  }
89  };
90 
91 } // End namespace Foam
92 
93 
94 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
95 
96 namespace Foam
97 {
98 
99 // The IOobject for reading
100 inline static IOobject selectReadIO(const word& name, const Time& runTime)
101 {
102  return IOobject
103  (
104  name,
105  runTime.constant(), // instance
106  "triSurface", // local
107  runTime, // registry
110  false // no register
111  );
112 }
113 
114 } // End namespace Foam
115 
116 
117 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
118 
119 void Foam::sampledMeshedSurface::setZoneMap()
120 {
121  // Populate zoneIds_ based on surfZone information
122 
123  const meshedSurface& s = static_cast<const meshedSurface&>(*this);
124 
125  const auto& zones = s.surfZones();
126 
127  zoneIds_.resize(s.size());
128 
129  // Trivial case
130  if (zoneIds_.empty() || zones.size() <= 1)
131  {
132  zoneIds_ = 0;
133  return;
134  }
135 
136 
137  label beg = 0;
138 
139  forAll(zones, zonei)
140  {
141  const label len = min(zones[zonei].size(), zoneIds_.size() - beg);
142 
143  // Assign sub-zone Ids
144  SubList<label>(zoneIds_, len, beg) = zonei;
145 
146  beg += len;
147  }
148 
149  // Anything remaining? Should not happen
150  {
151  const label len = (zoneIds_.size() - beg);
152 
153  if (len > 0)
154  {
155  SubList<label>(zoneIds_, len, beg) = max(0, zones.size()-1);
156  }
157  }
158 }
159 
160 
161 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
162 
163 bool Foam::sampledMeshedSurface::update(const meshSearch& meshSearcher)
164 {
165  // Global numbering for cells/faces
166  // - only used to uniquely identify local elements
167  globalIndex globalCells(onBoundary() ? mesh().nFaces() : mesh().nCells());
168 
169  // Find the cells the triangles of the surface are in.
170  // Does approximation by looking at the face centres only
171  const pointField& fc = surface_.faceCentres();
172 
173  List<nearInfo> nearest(fc.size(), nearInfo(GREAT, labelMax));
174 
175  if (sampleSource_ == cells)
176  {
177  // Search for nearest cell
178 
179  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
180 
181  forAll(fc, facei)
182  {
183  const point& pt = fc[facei];
184 
185  pointIndexHit info = cellTree.findNearest(pt, sqr(GREAT));
186 
187  if (info.hit())
188  {
189  nearest[facei].first() = magSqr(info.hitPoint()-pt);
190  nearest[facei].second() = globalCells.toGlobal(info.index());
191  }
192  }
193  }
194  else if (sampleSource_ == insideCells)
195  {
196  // Search for cell containing point
197 
198  const auto& cellTree = meshSearcher.cellTree();
199 
200  forAll(fc, facei)
201  {
202  const point& pt = fc[facei];
203 
204  if (cellTree.bb().contains(pt))
205  {
206  const label index = cellTree.findInside(pt);
207  if (index != -1)
208  {
209  nearest[facei].first() = 0;
210  nearest[facei].second() = globalCells.toGlobal(index);
211  }
212  }
213  }
214  }
215  else
216  {
217  // Search for nearest boundaryFace
218 
219  //- Search on all non-coupled boundary faces
220  const auto& bndTree = meshSearcher.nonCoupledBoundaryTree();
221 
222  forAll(fc, facei)
223  {
224  const point& pt = fc[facei];
225 
226  pointIndexHit info = bndTree.findNearest(pt, sqr(GREAT));
227 
228  if (info.hit())
229  {
230  nearest[facei].first() = magSqr(info.hitPoint()-pt);
231  nearest[facei].second() =
232  globalCells.toGlobal
233  (
234  bndTree.shapes().faceLabels()[info.index()]
235  );
236  }
237  }
238  }
239 
240 
241  // See which processor has the nearest. Mark and subset
242  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
243 
244  Pstream::listCombineGather(nearest, nearestEqOp());
246 
247  labelList cellOrFaceLabels(fc.size(), -1);
248 
249  bitSet facesToSubset(fc.size());
250 
251  forAll(nearest, facei)
252  {
253  const label index = nearest[facei].second();
254 
255  if (index == labelMax)
256  {
257  // Not found on any processor. How to map?
258  }
259  else if (globalCells.isLocal(index))
260  {
261  cellOrFaceLabels[facei] = globalCells.toLocal(index);
262  facesToSubset.set(facei);
263  }
264  }
265 
266 
267  if (debug)
268  {
269  Pout<< "Local out of faces:" << cellOrFaceLabels.size()
270  << " keeping:" << facesToSubset.count() << endl;
271  }
272 
273 
274  // Subset the surface
275  meshedSurface& s = static_cast<meshedSurface&>(*this);
276 
277  labelList pointMap;
279 
280  s = surface_.subsetMesh(facesToSubset, pointMap, faceMap);
281 
282 
283  // Ensure zoneIds_ are indeed correct
284  setZoneMap();
285 
286 
287  // Subset cellOrFaceLabels (for compact faces)
288  cellOrFaceLabels = labelUIndList(cellOrFaceLabels, faceMap)();
289 
290 
291  // Collect the samplePoints and sampleElements
292  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
293 
295  {
296  // With point interpolation
297 
298  samplePoints_.resize(pointMap.size());
299  sampleElements_.resize(pointMap.size(), -1);
300 
301  // Store any face per point (without using pointFaces())
302  labelList pointToFace(std::move(pointMap));
303 
304  forAll(s, facei)
305  {
306  const face& f = s[facei];
307 
308  for (const label labi : f)
309  {
310  pointToFace[labi] = facei;
311  }
312  }
313 
314 
315  if (sampleSource_ == cells)
316  {
317  // samplePoints_ : per surface point a location inside the cell
318  // sampleElements_ : per surface point the cell
319 
320  forAll(points(), pointi)
321  {
322  const point& pt = points()[pointi];
323 
324  const label celli = cellOrFaceLabels[pointToFace[pointi]];
325 
326  sampleElements_[pointi] = celli;
327 
328  // Check if point inside cell
329  if
330  (
331  mesh().pointInCell
332  (
333  pt,
334  sampleElements_[pointi],
335  meshSearcher.decompMode()
336  )
337  )
338  {
339  samplePoints_[pointi] = pt;
340  }
341  else
342  {
343  // Find nearest point on faces of cell
344 
345  scalar minDistSqr = VGREAT;
346 
347  for (const label facei : mesh().cells()[celli])
348  {
349  const face& f = mesh().faces()[facei];
350 
351  pointHit info = f.nearestPoint(pt, mesh().points());
352 
353  if (info.distance() < minDistSqr)
354  {
355  minDistSqr = info.distance();
356  samplePoints_[pointi] = info.rawPoint();
357  }
358  }
359  }
360  }
361  }
362  else if (sampleSource_ == insideCells)
363  {
364  // samplePoints_ : per surface point a location inside the cell
365  // sampleElements_ : per surface point the cell
366 
367  forAll(points(), pointi)
368  {
369  const point& pt = points()[pointi];
370 
371  const label celli = cellOrFaceLabels[pointToFace[pointi]];
372 
373  sampleElements_[pointi] = celli;
374  samplePoints_[pointi] = pt;
375  }
376  }
377  else
378  {
379  // samplePoints_ : per surface point a location on the boundary
380  // sampleElements_ : per surface point the boundary face containing
381  // the location
382 
383  forAll(points(), pointi)
384  {
385  const point& pt = points()[pointi];
386 
387  const label facei = cellOrFaceLabels[pointToFace[pointi]];
388 
389  sampleElements_[pointi] = facei;
390  samplePoints_[pointi] = mesh().faces()[facei].nearestPoint
391  (
392  pt,
393  mesh().points()
394  ).rawPoint();
395  }
396  }
397  }
398  else
399  {
400  // if sampleSource_ == cells:
401  // sampleElements_ : per surface triangle the cell
402  // samplePoints_ : n/a
403  // if sampleSource_ == insideCells:
404  // sampleElements_ : -1 or per surface triangle the cell
405  // samplePoints_ : n/a
406  // else:
407  // sampleElements_ : per surface triangle the boundary face
408  // samplePoints_ : n/a
409  sampleElements_.transfer(cellOrFaceLabels);
410  samplePoints_.clear();
411  }
412 
413 
414  if (debug)
415  {
416  this->clearOut();
417  OFstream str(mesh().time().path()/"surfaceToMesh.obj");
418  Info<< "Dumping correspondence from local surface (points or faces)"
419  << " to mesh (cells or faces) to " << str.name() << endl;
420 
421  const vectorField& centres =
422  (
423  onBoundary()
424  ? mesh().faceCentres()
425  : mesh().cellCentres()
426  );
427 
429  {
430  label vertI = 0;
431  forAll(samplePoints_, pointi)
432  {
433  meshTools::writeOBJ(str, points()[pointi]);
434  vertI++;
435 
436  meshTools::writeOBJ(str, samplePoints_[pointi]);
437  vertI++;
438 
439  label elemi = sampleElements_[pointi];
440  meshTools::writeOBJ(str, centres[elemi]);
441  vertI++;
442  str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI << nl;
443  }
444  }
445  else
446  {
447  label vertI = 0;
448  forAll(sampleElements_, triI)
449  {
450  meshTools::writeOBJ(str, faceCentres()[triI]);
451  vertI++;
452 
453  label elemi = sampleElements_[triI];
454  meshTools::writeOBJ(str, centres[elemi]);
455  vertI++;
456  str << "l " << vertI-1 << ' ' << vertI << nl;
457  }
458  }
459  }
460 
461  needsUpdate_ = false;
462  return true;
463 }
464 
465 
466 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
467 
469 (
470  const word& name,
471  const polyMesh& mesh,
472  const word& surfaceName,
473  const samplingSource sampleSource
474 )
475 :
477  MeshStorage(),
478  surfaceName_(surfaceName),
479  surface_
480  (
481  selectReadIO(surfaceName, mesh.time()),
483  ),
484  sampleSource_(sampleSource),
485  needsUpdate_(true),
486  keepIds_(true),
487  zoneIds_(),
488  sampleElements_(),
489  samplePoints_()
490 {}
491 
492 
494 (
495  const word& name,
496  const polyMesh& mesh,
497  const dictionary& dict
498 )
499 :
501  MeshStorage(),
502  surfaceName_
503  (
504  meshedSurface::findFile
505  (
506  selectReadIO(dict.get<word>("surface"), mesh.time()),
507  dict
508  ).name()
509  ),
510  surface_
511  (
512  selectReadIO(dict.get<word>("surface"), mesh.time()),
513  dict
514  ),
515  sampleSource_(samplingSourceNames_.get("source", dict)),
516  needsUpdate_(true),
517  keepIds_(dict.getOrDefault("keepIds", true)),
518  zoneIds_(),
519  sampleElements_(),
520  samplePoints_()
521 {
522  wordRes includePatches;
523  dict.readIfPresent("patches", includePatches);
524  includePatches.uniq();
525 
526  // Could also shift this to the reader itself,
527  // but not yet necessary.
528 
529  if (!includePatches.empty())
530  {
531  Info<< "Subsetting surface " << surfaceName_
532  << " to patches: " << flatOutput(includePatches) << nl;
533 
534  const surfZoneList& zones = surface_.surfZones();
535 
536  const labelList zoneIndices
537  (
539  (
540  zones,
541  includePatches,
542  wordRes(),
544  )
545  );
546 
547  // Faces to subset
548  bitSet includeMap(surface_.size());
549 
550  for (const label zonei : zoneIndices)
551  {
552  const surfZone& zn = zones[zonei];
553  includeMap.set(zn.range());
554  }
555 
556  if (includeMap.none())
557  {
559  << "Patch selection results in an empty surface"
560  << " - ignoring" << nl;
561  }
562  else if (!includeMap.all())
563  {
564  meshedSurface subSurf(surface_.subsetMesh(includeMap));
565 
566  if (subSurf.empty())
567  {
569  << "Bad surface subset (empty)"
570  << " - skip and hope for the best" << nl;
571  }
572  else
573  {
574  // Replace
575  surface_.transfer(subSurf);
576  }
577  }
578  }
579 }
580 
581 
582 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
583 
585 {
586  return needsUpdate_;
587 }
588 
589 
591 {
592  // already marked as expired
593  if (needsUpdate_)
594  {
595  return false;
596  }
597 
600  zoneIds_.clear();
601 
602  sampleElements_.clear();
603  samplePoints_.clear();
604 
605  needsUpdate_ = true;
606  return true;
607 }
608 
609 
611 {
612  if (!needsUpdate_)
613  {
614  return false;
615  }
616 
617  // Calculate surface and mesh overlap bounding box
618  treeBoundBox bb(surface_.points(), surface_.meshPoints());
619 
620  // Check for overlap with (global!) mesh bb
621  const bool intersect = bb.intersect(mesh().bounds());
622 
623  if (!intersect)
624  {
625  // Surface and mesh do not overlap at all. Guarantee a valid
626  // bounding box so we don't get any 'invalid bounding box' errors.
627 
629  << "Surface " << surfaceName_
630  << " does not overlap bounding box of mesh " << mesh().bounds()
631  << endl;
632 
633  bb = treeBoundBox(mesh().bounds());
634  const vector span(bb.span());
635 
636  bb.min() += (0.5-1e-6)*span;
637  bb.max() -= (0.5-1e-6)*span;
638  }
639  else
640  {
641  // Extend a bit
642  const vector span(bb.span());
643  bb.min() -= 0.5*span;
644  bb.max() += 0.5*span;
645 
646  bb.inflate(1e-6);
647  }
648 
649  // Mesh search engine, no triangulation of faces.
650  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
651 
652  return update(meshSearcher);
653 }
654 
655 
657 {
658  if (!needsUpdate_)
659  {
660  return false;
661  }
662 
663  // Mesh search engine on subset, no triangulation of faces.
664  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
665 
666  return update(meshSearcher);
667 }
668 
669 
671 (
672  const interpolation<scalar>& sampler
673 ) const
674 {
675  return sampleOnFaces(sampler);
676 }
677 
678 
680 (
681  const interpolation<vector>& sampler
682 ) const
683 {
684  return sampleOnFaces(sampler);
685 }
686 
687 
689 (
690  const interpolation<sphericalTensor>& sampler
691 ) const
692 {
693  return sampleOnFaces(sampler);
694 }
695 
696 
698 (
699  const interpolation<symmTensor>& sampler
700 ) const
701 {
702  return sampleOnFaces(sampler);
703 }
704 
705 
707 (
708  const interpolation<tensor>& sampler
709 ) const
710 {
711  return sampleOnFaces(sampler);
712 }
713 
714 
716 (
717  const interpolation<scalar>& interpolator
718 ) const
719 {
720  return sampleOnPoints(interpolator);
721 }
722 
723 
725 (
726  const interpolation<vector>& interpolator
727 ) const
728 {
729  return sampleOnPoints(interpolator);
730 }
731 
733 (
734  const interpolation<sphericalTensor>& interpolator
735 ) const
736 {
737  return sampleOnPoints(interpolator);
738 }
739 
740 
742 (
743  const interpolation<symmTensor>& interpolator
744 ) const
745 {
746  return sampleOnPoints(interpolator);
747 }
748 
749 
751 (
752  const interpolation<tensor>& interpolator
753 ) const
754 {
755  return sampleOnPoints(interpolator);
756 }
757 
758 
760 {
761  os << "meshedSurface: " << name() << " :"
762  << " surface:" << surfaceName_
763  << " faces:" << faces().size()
764  << " points:" << points().size()
765  << " zoneids:" << zoneIds().size();
766 }
767 
768 
769 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
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:104
meshTools.H
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:584
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:62
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:60
Foam::nearestEqOp
Definition: distributedTriSurfaceMesh.C:124
Foam::sampledMeshedSurface::samplingSource
samplingSource
Types of sampling regions.
Definition: sampledMeshedSurface.H:171
Foam::sampledMeshedSurface::sampledMeshedSurface
sampledMeshedSurface(const word &name, const polyMesh &mesh, const word &surfaceName, const samplingSource sampleSource)
Construct from components.
Definition: sampledMeshedSurface.C:469
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:64
update
mesh update()
Tuple2.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::selectReadIO
static IOobject selectReadIO(const word &name, const Time &runTime)
Definition: sampledMeshedSurface.C:100
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:87
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:350
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::sampledMeshedSurface::update
virtual bool update()
Update the surface as required.
Definition: sampledMeshedSurface.C:610
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::nearInfo
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
Definition: sampledMeshedSurface.C:76
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:385
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::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:459
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:238
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:120
Foam::stringListOps::findMatching
labelList findMatching(const StringListType &input, const wordRes &whitelist, const wordRes &blacklist=wordRes(), AccessOp aop=noOp())
Return ids for items with matching names.
Foam::sampledMeshedSurface::print
virtual void print(Ostream &os) const
Write.
Definition: sampledMeshedSurface.C:759
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::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
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:121
Foam::pointHit
PointHit< point > pointHit
Definition: pointHit.H:43
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
treeDataCell.H
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:441
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:175
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:45
meshSearch.H
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::flatOutput
FlatOutput< Container > flatOutput(const Container &obj, label len=0)
Global flatOutput function.
Definition: FlatOutput.H:85
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:144
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())
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
Definition: sampledSurface.C:55
Foam::nearestEqOp::operator()
void operator()(nearInfo &x, const nearInfo &y) const
Definition: sampledMeshedSurface.C:82
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:325
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:248
Foam::polyMesh::FACE_PLANES
Definition: polyMesh.H:103
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:590
Foam::sampledMeshedSurface::sample
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
Definition: sampledMeshedSurface.C:671
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: Tuple2.H:57
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:88
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::MeshedSurface< face >
Foam::MeshedSurface::size
label size() const
The surface size is the number of faces.
Definition: MeshedSurface.H:405
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation to nodes requested for surface.
Definition: sampledSurface.H:326
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:298
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58