distanceSurface.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 "distanceSurface.H"
30 #include "dictionary.H"
31 #include "volFields.H"
32 #include "volPointInterpolation.H"
34 #include "fvMesh.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(distanceSurface, 0);
41 }
42 
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 const Foam::Enum
47 <
48  Foam::distanceSurface::topologyFilterType
49 >
50 Foam::distanceSurface::topoFilterNames_
51 ({
52  { topologyFilterType::NONE, "none" },
53  { topologyFilterType::LARGEST_REGION, "largestRegion" },
54  { topologyFilterType::NEAREST_POINTS, "nearestPoints" },
55  { topologyFilterType::PROXIMITY, "proximity" },
56 });
57 
58 
59 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
60 
61 namespace Foam
62 {
63 
64 // Check that all point hits are valid
65 static inline void checkAllHits(const UList<pointIndexHit>& nearest)
66 {
67  label notHit = 0;
68  for (const pointIndexHit& pHit : nearest)
69  {
70  if (!pHit.hit())
71  {
72  ++notHit;
73  }
74  }
75 
76  if (notHit)
77  {
79  << "Had " << notHit << " from " << nearest.size()
80  << " without a point hit" << endl
81  << abort(FatalError);
82  }
83 }
84 
85 
86 // Normal distance from surface hit point to a point in the mesh
87 static inline scalar normalDistance_zero
88 (
89  const point& pt,
90  const pointIndexHit& pHit,
91  const vector& norm
92 )
93 {
94  const vector diff(pt - pHit.point());
95 
96  return (diff & norm);
97 }
98 
99 
100 // Signed distance from surface hit point to a point in the mesh,
101 // the sign is dictated by the normal
102 static inline scalar normalDistance_nonzero
103 (
104  const point& pt,
105  const pointIndexHit& pHit,
106  const vector& norm
107 )
108 {
109  const vector diff(pt - pHit.point());
110  const scalar normDist = (diff & norm);
111 
112  return Foam::sign(normDist) * Foam::mag(diff);
113 }
114 
115 
116 // Normal distance from surface hit point to a point in the mesh
117 static inline void calcNormalDistance_zero
118 (
120  const pointField& points,
121  const List<pointIndexHit>& nearest,
122  const vectorField& normals
123 )
124 {
125  forAll(nearest, i)
126  {
127  distance[i] =
128  normalDistance_zero(points[i], nearest[i], normals[i]);
129  }
130 }
131 
132 
133 // Signed distance from surface hit point -> point in the mesh,
134 // the sign is dictated by the normal
135 static inline void calcNormalDistance_nonzero
136 (
138  const pointField& points,
139  const List<pointIndexHit>& nearest,
140  const vectorField& normals
141 )
142 {
143  forAll(nearest, i)
144  {
145  distance[i] =
146  normalDistance_nonzero(points[i], nearest[i], normals[i]);
147  }
148 }
149 
150 
151 // Close to the surface: normal distance from surface hit point
152 // Far from surface: distance from surface hit point
153 //
154 // Note
155 // This switch may be helpful when working directly with
156 // distance/gradient fields. Has low overhead otherwise.
157 // May be replaced in the future (2020-11)
158 static inline void calcNormalDistance_filtered
159 (
161  const bitSet& ignoreLocation,
162  const pointField& points,
163  const List<pointIndexHit>& nearest,
164  const vectorField& normals
165 )
166 {
167  forAll(nearest, i)
168  {
169  if (ignoreLocation.test(i))
170  {
171  distance[i] =
172  normalDistance_nonzero(points[i], nearest[i], normals[i]);
173  }
174  else
175  {
176  distance[i] =
177  normalDistance_zero(points[i], nearest[i], normals[i]);
178  }
179  }
180 }
181 
182 
183 // Flat surfaces (eg, a plane) have an extreme change in
184 // the normal at the edge, which creates a zero-crossing
185 // extending to infinity.
186 //
187 // Ad hoc treatment: require that the surface hit
188 // point is within a somewhat generous bounding box
189 // for the cell
190 template<bool WantPointFilter = false>
191 static bitSet simpleGeometricFilter
192 (
193  bitSet& ignoreCells,
194  const List<pointIndexHit>& nearest,
195  const polyMesh& mesh
196 )
197 {
198  // A deny filter. Initially false (accept everything)
199  ignoreCells.resize(mesh.nCells());
200 
201  bitSet pointFilter;
202  if (WantPointFilter)
203  {
204  // Create as accept filter. Initially false (deny everything)
205  pointFilter.resize(mesh.nPoints());
206  }
207 
208  boundBox cellBb;
209 
210  forAll(nearest, celli)
211  {
212  const point& pt = nearest[celli].point();
213 
214  const labelList& cPoints = mesh.cellPoints(celli);
215 
216  cellBb.clear();
217  cellBb.add(mesh.points(), cPoints);
218 
219  // Expand slightly to catch corners
220  cellBb.inflate(0.1);
221 
222  if (!cellBb.contains(pt))
223  {
224  ignoreCells.set(celli);
225  }
226  else if (WantPointFilter)
227  {
228  // Good cell candidate, accept its points
229  pointFilter.set(cPoints);
230  }
231  }
232 
233  // Flip from accept to deny filter
234  pointFilter.flip();
235 
236  return pointFilter;
237 }
238 
239 
240 } // End namespace Foam
241 
242 
243 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
244 
246 (
247  const word& defaultSurfaceName,
248  const polyMesh& mesh,
249  const dictionary& dict
250 )
251 :
252  mesh_(mesh),
253  geometryPtr_
254  (
256  (
257  dict.get<word>("surfaceType"),
258  IOobject
259  (
260  dict.getOrDefault("surfaceName", defaultSurfaceName),
261  mesh.time().constant(), // directory
262  "triSurface", // instance
263  mesh.time(), // registry
264  IOobject::MUST_READ,
265  IOobject::NO_WRITE
266  ),
267  dict
268  )
269  ),
270  distance_(dict.getOrDefault<scalar>("distance", 0)),
271  withZeroDistance_(equal(distance_, 0)),
272  withSignDistance_
273  (
274  withZeroDistance_
275  || (distance_ < 0)
276  || dict.getOrDefault<bool>("signed", true)
277  ),
278  isoParams_
279  (
280  dict,
281  isoSurfaceParams::ALGO_TOPO,
282  isoSurfaceParams::filterType::DIAGCELL
283  ),
284  topoFilter_
285  (
286  topoFilterNames_.getOrDefault
287  (
288  "topology",
289  dict,
290  topologyFilterType::NONE
291  )
292  ),
293  nearestPoints_(),
294  maxDistanceSqr_(Foam::sqr(GREAT)),
295  absProximity_(dict.getOrDefault<scalar>("absProximity", 1e-5)),
296  cellDistancePtr_(nullptr),
297  pointDistance_(),
298  surface_(),
299  meshCells_(),
300  isoSurfacePtr_(nullptr)
301 {
302  if (topologyFilterType::NEAREST_POINTS == topoFilter_)
303  {
304  dict.readEntry("nearestPoints", nearestPoints_);
305  }
306 
307  if (dict.readIfPresent("maxDistance", maxDistanceSqr_))
308  {
309  maxDistanceSqr_ = Foam::sqr(maxDistanceSqr_);
310  }
311 }
312 
313 
315 (
316  const polyMesh& mesh,
317  const word& surfaceType,
318  const word& surfaceName,
319  const isoSurfaceParams& params,
320  const bool interpolate
321 )
322 :
324  (
325  mesh,
326  interpolate,
327  surfaceType,
328  surfaceName,
329  scalar(0),
330  true, // redundant - must be signed
331  params
332  )
333 {}
334 
335 
337 (
338  const polyMesh& mesh,
339  const bool interpolate,
340  const word& surfaceType,
341  const word& surfaceName,
342  const scalar distance,
343  const bool useSignedDistance,
344  const isoSurfaceParams& params
345 )
346 :
347  mesh_(mesh),
348  geometryPtr_
349  (
351  (
352  surfaceType,
353  IOobject
354  (
355  surfaceName, // name
356  mesh.time().constant(), // directory
357  "triSurface", // instance
358  mesh.time(), // registry
359  IOobject::MUST_READ,
360  IOobject::NO_WRITE
361  ),
362  dictionary()
363  )
364  ),
365  distance_(distance),
366  withZeroDistance_(equal(distance_, 0)),
367  withSignDistance_
368  (
369  withZeroDistance_
370  || (distance_ < 0)
371  || useSignedDistance
372  ),
373  isoParams_(params),
374  topoFilter_(topologyFilterType::NONE),
375  nearestPoints_(),
376  maxDistanceSqr_(Foam::sqr(GREAT)),
377  absProximity_(1e-5),
378  cellDistancePtr_(nullptr),
379  pointDistance_(),
380  surface_(),
381  meshCells_(),
382  isoSurfacePtr_(nullptr)
383 {}
384 
385 
386 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
387 
389 {
390  if (debug)
391  {
392  Pout<< "distanceSurface::createGeometry updating geometry." << endl;
393  }
394 
395  // Clear any previously stored topologies
396  isoSurfacePtr_.reset(nullptr);
397  surface_.clear();
398  meshCells_.clear();
399 
400  // Doing searches on this surface
401  const searchableSurface& geom = geometryPtr_();
402 
403  const fvMesh& fvmesh = static_cast<const fvMesh&>(mesh_);
404 
405  // Distance to cell centres
406  // ~~~~~~~~~~~~~~~~~~~~~~~~
407 
408  cellDistancePtr_.reset
409  (
410  new volScalarField
411  (
412  IOobject
413  (
414  "distanceSurface.cellDistance",
415  fvmesh.time().timeName(),
416  fvmesh.time(),
419  false
420  ),
421  fvmesh,
423  )
424  );
425  auto& cellDistance = *cellDistancePtr_;
426 
427 
428  // For distance = 0 we apply additional geometric filtering
429  // to limit the extent of open edges.
430  //
431  // Does not work with ALGO_POINT
432 
433  bitSet ignoreCells, ignoreCellPoints;
434 
435  const bool filterCells =
436  (
437  withZeroDistance_
438  && isoParams_.algorithm() != isoSurfaceParams::ALGO_POINT
439  );
440 
441 
442  // Internal field
443  {
444  const pointField& cc = fvmesh.C();
445  scalarField& fld = cellDistance.primitiveFieldRef();
446 
447  List<pointIndexHit> nearest;
448  geom.findNearest
449  (
450  cc,
451  // Use initialized field (GREAT) to limit search too
452  fld,
453  nearest
454  );
455  checkAllHits(nearest);
456 
457  // Geometric pre-filtering when distance == 0
458  if (filterCells)
459  {
460  ignoreCellPoints =
461  simpleGeometricFilter<false>(ignoreCells, nearest, fvmesh);
462  }
463 
464  if (withSignDistance_)
465  {
466  vectorField norms;
467  geom.getNormal(nearest, norms);
468 
469  if (filterCells)
470  {
471  // With inside/outside switching (see note above)
473  (
474  fld,
475  ignoreCells,
476  cc,
477  nearest,
478  norms
479  );
480  }
481  else if (withZeroDistance_)
482  {
483  calcNormalDistance_zero(fld, cc, nearest, norms);
484  }
485  else
486  {
487  calcNormalDistance_nonzero(fld, cc, nearest, norms);
488  }
489  }
490  else
491  {
492  calcAbsoluteDistance(fld, cc, nearest);
493  }
494  }
495 
496 
497  // Patch fields
498  {
499  forAll(fvmesh.C().boundaryField(), patchi)
500  {
501  const pointField& cc = fvmesh.C().boundaryField()[patchi];
502  scalarField& fld = cellDistance.boundaryFieldRef()[patchi];
503 
504  List<pointIndexHit> nearest;
505  geom.findNearest
506  (
507  cc,
508  scalarField(cc.size(), GREAT),
509  nearest
510  );
511  checkAllHits(nearest);
512 
513  if (withSignDistance_)
514  {
515  vectorField norms;
516  geom.getNormal(nearest, norms);
517 
518  if (withZeroDistance_)
519  {
520  // Slight inconsistency in boundary vs interior when
521  // cells are filtered, but the patch fields are only
522  // used by isoSurfacePoint, and filtering is disabled
523  // for that anyhow.
524 
525  calcNormalDistance_zero(fld, cc, nearest, norms);
526  }
527  else
528  {
529  calcNormalDistance_nonzero(fld, cc, nearest, norms);
530  }
531  }
532  else
533  {
534  calcAbsoluteDistance(fld, cc, nearest);
535  }
536  }
537  }
538 
539 
540  // On processor patches the mesh.C() will already be the cell centre
541  // on the opposite side so no need to swap cellDistance.
542 
543 
544  // Distance to points
545  pointDistance_.resize(fvmesh.nPoints());
546  pointDistance_ = GREAT;
547  {
548  const pointField& pts = fvmesh.points();
549  scalarField& fld = pointDistance_;
550 
551  List<pointIndexHit> nearest;
552  geom.findNearest
553  (
554  pts,
555  // Use initialized field (GREAT) to limit search too
556  pointDistance_,
557  nearest
558  );
559  checkAllHits(nearest);
560 
561  if (withSignDistance_)
562  {
563  vectorField norms;
564  geom.getNormal(nearest, norms);
565 
566  if (filterCells)
567  {
568  // With inside/outside switching (see note above)
570  (
571  fld,
572  ignoreCellPoints,
573  pts,
574  nearest,
575  norms
576  );
577  }
578  else if (withZeroDistance_)
579  {
580  calcNormalDistance_zero(fld, pts, nearest, norms);
581  }
582  else
583  {
584  calcNormalDistance_nonzero(fld, pts, nearest, norms);
585  }
586  }
587  else
588  {
589  calcAbsoluteDistance(fld, pts, nearest);
590  }
591  }
592 
593 
594  // Don't need ignoreCells if there is nothing to ignore.
595  if (ignoreCells.none())
596  {
597  ignoreCells.clearStorage();
598  }
599  else if (filterCells && topologyFilterType::NONE != topoFilter_)
600  {
601  // For refine blocked cells (eg, checking actual cells cut)
602  isoSurfaceBase isoCutter
603  (
604  mesh_,
605  cellDistance,
606  pointDistance_,
607  distance_
608  );
609 
610  if (topologyFilterType::LARGEST_REGION == topoFilter_)
611  {
612  refineBlockedCells(ignoreCells, isoCutter);
613  filterKeepLargestRegion(ignoreCells);
614  }
615  else if (topologyFilterType::NEAREST_POINTS == topoFilter_)
616  {
617  refineBlockedCells(ignoreCells, isoCutter);
618  filterKeepNearestRegions(ignoreCells);
619  }
620  }
621 
622  // Don't need point filter beyond this point
623  ignoreCellPoints.clearStorage();
624 
625 
626  if (debug)
627  {
628  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
629  cellDistance.write();
630  pointScalarField pDist
631  (
632  IOobject
633  (
634  "distanceSurface.pointDistance",
635  fvmesh.time().timeName(),
636  fvmesh.time(),
639  false
640  ),
641  pointMesh::New(fvmesh),
643  );
644  pDist.primitiveFieldRef() = pointDistance_;
645 
646  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
647  pDist.write();
648  }
649 
650  isoSurfacePtr_.reset
651  (
653  (
654  isoParams_,
655  cellDistance,
656  pointDistance_,
657  distance_,
658  ignoreCells
659  )
660  );
661 
662 
663  // ALGO_POINT still needs cell, point fields (for interpolate)
664  // The others can do straight transfer
665 
666  // But also flatten into a straight transfer for proximity filtering
667  if
668  (
670  || topologyFilterType::PROXIMITY == topoFilter_
671  )
672  {
673  surface_.transfer(static_cast<meshedSurface&>(*isoSurfacePtr_));
674  meshCells_.transfer(isoSurfacePtr_->meshCells());
675 
676  isoSurfacePtr_.reset(nullptr);
677  cellDistancePtr_.reset(nullptr);
678  pointDistance_.clear();
679  }
680 
681  if (topologyFilterType::PROXIMITY == topoFilter_)
682  {
684  }
685 
686  if (debug)
687  {
688  print(Pout);
689  Pout<< endl;
690  }
691 }
692 
693 
695 {
696  os << " surface:" << surfaceName()
697  << " distance:" << distance()
698  << " faces:" << surface().surfFaces().size()
699  << " points:" << surface().points().size();
700 }
701 
702 
703 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
volFields.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::isoSurfaceBase::New
static autoPtr< isoSurfaceBase > New(const isoSurfaceParams &params, const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const bitSet &ignoreCells=bitSet())
Create for specified algorithm type.
Definition: isoSurfaceBaseNew.C:37
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::isoSurfaceBase
Low-level components common to various iso-surface algorithms.
Definition: isoSurfaceBase.H:66
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::searchableSurface::getNormal
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const =0
From a set of points and indices get the normal.
Foam::primitiveMesh::cellPoints
const labelListList & cellPoints() const
Definition: primitiveMeshCellPoints.C:34
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::boundBox::inflate
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:175
Foam::isoSurfaceParams::algorithm
algorithmType algorithm() const noexcept
Get current algorithm.
Definition: isoSurfaceParams.H:190
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::distanceSurface::distanceSurface
distanceSurface(const word &defaultSurfaceName, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: distanceSurface.C:246
Foam::bitSet::set
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:574
Foam::calcNormalDistance_filtered
static void calcNormalDistance_filtered(scalarField &distance, const bitSet &ignoreLocation, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
Definition: distanceSurface.C:159
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::calcNormalDistance_zero
static void calcNormalDistance_zero(scalarField &distance, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
Definition: distanceSurface.C:118
Foam::simpleGeometricFilter
static bitSet simpleGeometricFilter(bitSet &ignoreCells, const List< pointIndexHit > &nearest, const polyMesh &mesh)
Definition: distanceSurface.C:192
Foam::MeshedSurface::clear
virtual void clear()
Clear all storage.
Definition: MeshedSurface.C:594
Foam::distanceSurface::filterKeepNearestRegions
void filterKeepNearestRegions(bitSet &ignoreCells) const
Keep region(s) closest to the nearest points.
Definition: distanceSurfaceFilter.C:186
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::bitSet::test
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:520
Foam::MatrixTools::equal
bool equal(const Matrix< Form1, Type > &A, const Matrix< Form2, Type > &B, const bool verbose=false, const label maxDiffs=10, const scalar relTol=1e-5, const scalar absTol=1e-8)
Compare matrix elements for absolute or relative equality.
Definition: MatrixTools.C:34
Foam::PackedList::resize
void resize(const label nElem, const unsigned int val=0u)
Reset addressable list size, does not shrink the allocated size.
Definition: PackedListI.H:399
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::distanceSurface::refineBlockedCells
bool refineBlockedCells(bitSet &ignoreCells, const isoSurfaceBase &isoContext) const
Re-filter the blocked cells based on the anticipated cuts.
Definition: distanceSurfaceFilter.C:36
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
Foam::distanceSurface
A surface defined by a distance from an input searchable surface. Uses an iso-surface algorithm (cell...
Definition: distanceSurface.H:198
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:52
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::PointIndexHit::hit
bool hit() const noexcept
Is there a hit?
Definition: PointIndexHit.H:130
Foam::normalDistance_zero
static scalar normalDistance_zero(const point &pt, const pointIndexHit &pHit, const vector &norm)
Definition: distanceSurface.C:88
Foam::MeshedSurface::transfer
void transfer(pointField &pointLst, List< Face > &faceLst)
Transfer the components.
Definition: MeshedSurface.C:1281
Foam::Field< scalar >
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:69
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:459
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
fld
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:341
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::bitSet::flip
void flip()
Invert all bits in the addressable region.
Definition: bitSetI.H:618
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::isoSurfaceParams
Preferences for controlling iso-surface algorithms.
Definition: isoSurfaceParams.H:91
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::isoSurfaceParams::ALGO_POINT
Definition: isoSurfaceParams.H:103
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::searchableSurface::findNearest
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const =0
distanceSurface.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::boundBox::clear
void clear()
Clear bounding box and make it an inverted box.
Definition: boundBoxI.H:184
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::boundBox::contains
bool contains(const point &pt) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:271
Foam::distanceSurface::print
void print(Ostream &os) const
Print information.
Definition: distanceSurface.C:694
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
volPointInterpolation.H
Foam::bitSet::none
bool none() const
True if no bits in this bitset are set.
Definition: bitSetI.H:487
Foam::checkAllHits
static void checkAllHits(const UList< pointIndexHit > &nearest)
Definition: distanceSurface.C:65
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
dictionary.H
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
Foam::PointIndexHit::point
const point_type & point() const noexcept
Return point, no checks.
Definition: PointIndexHit.H:142
Foam::PackedList::clearStorage
void clearStorage()
Clear the list and delete storage.
Definition: PackedListI.H:510
Foam::distanceSurface::filterKeepLargestRegion
void filterKeepLargestRegion(bitSet &ignoreCells) const
Keep region with the most cuts (after region split)
Definition: distanceSurfaceFilter.C:119
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::distanceSurface::createGeometry
void createGeometry()
Create/recreate the distance surface.
Definition: distanceSurface.C:388
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::distanceSurface::filterByProximity
void filterByProximity()
Adjust extracted iso-surface to remove far faces.
Definition: distanceSurfaceFilter.C:309
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::MeshedSurface< face >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::normalDistance_nonzero
static scalar normalDistance_nonzero(const point &pt, const pointIndexHit &pHit, const vector &norm)
Definition: distanceSurface.C:103
Foam::calcNormalDistance_nonzero
static void calcNormalDistance_nonzero(scalarField &distance, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
Definition: distanceSurface.C:136
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::boundBox::add
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191