isoSurfacePoint.H
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-2022 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 Class
28  Foam::isoSurfacePoint
29 
30 Description
31  A surface formed by the iso value.
32  After "Regularised Marching Tetrahedra: improved iso-surface extraction",
33  G.M. Treece, R.W. Prager and A.H. Gee.
34 
35  Note:
36  - does tets without using cell centres/cell values. Not tested.
37  - regularisation can create duplicate triangles/non-manifold surfaces.
38  Current handling of those is bit ad-hoc for now and not perfect.
39  - regularisation does not do boundary points so as to preserve the
40  boundary perfectly.
41  - uses geometric merge with fraction of bounding box as distance.
42  - triangles can be between two cell centres so constant sampling
43  does not make sense.
44  - on empty patches behaves like zero gradient.
45  - does not do 2D correctly, creates non-flat iso surface.
46  - on processor boundaries might two overlapping (identical) triangles
47  (one from either side)
48 
49  The handling on coupled patches is a bit complex. All fields
50  (values and coordinates) get rewritten so
51  - empty patches get zerogradient (value) and facecentre (coordinates)
52  - separated processor patch faces get interpolate (value) and facecentre
53  (coordinates). (this is already the default for cyclics)
54  - non-separated processor patch faces keep opposite value and cell centre
55 
56  Now the triangle generation on non-separated processor patch faces
57  can use the neighbouring value. Any separated processor face or cyclic
58  face gets handled just like any boundary face.
59 
60 SourceFiles
61  isoSurfacePoint.C
62  isoSurfacePointTemplates.C
63 
64 \*---------------------------------------------------------------------------*/
65 
66 #ifndef Foam_isoSurfacePoint_H
67 #define Foam_isoSurfacePoint_H
68 
69 #include "bitSet.H"
70 #include "volFields.H"
71 #include "slicedVolFields.H"
72 #include "isoSurfaceBase.H"
73 
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
76 namespace Foam
77 {
78 
79 // Forward Declarations
80 class plane;
81 class treeBoundBox;
82 class triSurface;
83 
84 /*---------------------------------------------------------------------------*\
85  Class isoSurfacePoint Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 class isoSurfacePoint
89 :
90  public isoSurfaceBase
91 {
92  // Private Data
93 
94  //- Cell values.
95  //- Input volScalarField with separated coupled patches rewritten
97 
98  //- When to merge points
99  const scalar mergeDistance_;
100 
101  //- Whether face might be cut
102  List<cutType> faceCutType_;
103 
104  //- Whether cell might be cut
105  List<cutType> cellCutType_;
106 
107  //- Estimated number of cut cells
108  label nCutCells_;
109 
110  //- For every unmerged triangle point the point in the triSurface
111  labelList triPointMergeMap_;
112 
113  //- triSurface points that have weighted interpolation
114  DynamicList<label> interpolatedPoints_;
115 
116  //- corresponding original, unmerged points
117  DynamicList<FixedList<label, 3>> interpolatedOldPoints_;
118 
119  //- corresponding weights
120  DynamicList<FixedList<scalar, 3>> interpolationWeights_;
121 
122 
123  // Private Member Functions
124 
125  // Point synchronisation
126 
127  //- Per face whether is collocated
128  static bitSet collocatedFaces(const coupledPolyPatch& cpp);
129 
130  //- Synchronise points on all non-separated coupled patches
131  void syncUnseparatedPoints
132  (
133  pointField& collapsedPoint,
134  const point& nullValue
135  ) const;
136 
137 
138  //- Get location of iso value as fraction inbetween s0,s1
139  scalar isoFraction
140  (
141  const scalar s0,
142  const scalar s1
143  ) const;
144 
145  //- Check if any edge of a face is cut
146  bool isEdgeOfFaceCut
147  (
148  const scalarField& pVals,
149  const face& f,
150  const bool ownLower,
151  const bool neiLower
152  ) const;
153 
154  //- Get neighbour value and position.
155  void getNeighbour
156  (
157  const labelList& boundaryRegion,
158  const volVectorField& meshC,
159  const volScalarField& cVals,
160  const label celli,
161  const label facei,
162  scalar& nbrValue,
163  point& nbrPoint
164  ) const;
165 
166  //- Determine for every face/cell whether it (possibly) generates
167  // triangles.
168  void calcCutTypes
169  (
170  const labelList& boundaryRegion,
171  const volVectorField& meshC,
172  const volScalarField& cVals,
173  const scalarField& pVals
174  );
175 
176  static point calcCentre(const triSurface&);
177 
178  //- Determine per cc whether all near cuts can be snapped to single
179  // point.
180  void calcSnappedCc
181  (
182  const labelList& boundaryRegion,
183  const volVectorField& meshC,
184  const volScalarField& cVals,
185  const scalarField& pVals,
186  DynamicList<point>& snappedPoints,
187  labelList& snappedCc
188  ) const;
189 
190  //- Determine per point whether all near cuts can be snapped to single
191  // point.
192  void calcSnappedPoint
193  (
194  const bitSet& isBoundaryPoint,
195  const labelList& boundaryRegion,
196  const volVectorField& meshC,
197  const volScalarField& cVals,
198  const scalarField& pVals,
199  DynamicList<point>& snappedPoints,
200  labelList& snappedPoint
201  ) const;
202 
203 
204  //- Return input field with coupled (and empty) patch values rewritten
205  template<class Type>
207  adaptPatchFields(const VolumeField<Type>& fld) const;
208 
209  //- Generate single point by interpolation or snapping
210  template<class Type>
211  Type generatePoint
212  (
213  const scalar s0,
214  const Type& p0,
215  const bool hasSnap0,
216  const Type& snapP0,
217 
218  const scalar s1,
219  const Type& p1,
220  const bool hasSnap1,
221  const Type& snapP1
222  ) const;
223 
224 
225  //- Note: cannot use simpler isoSurfaceCell::generateTriPoints since
226  // the need here to sometimes pass in remote 'snappedPoints'
227  template<class Type>
228  void generateTriPoints
229  (
230  const scalar s0,
231  const Type& p0,
232  const bool hasSnap0,
233  const Type& snapP0,
234 
235  const scalar s1,
236  const Type& p1,
237  const bool hasSnap1,
238  const Type& snapP1,
239 
240  const scalar s2,
241  const Type& p2,
242  const bool hasSnap2,
243  const Type& snapP2,
244 
245  const scalar s3,
246  const Type& p3,
247  const bool hasSnap3,
248  const Type& snapP3,
249 
251  ) const;
252 
253  template<class Type>
254  label generateFaceTriPoints
255  (
256  const volScalarField& cVals,
257  const scalarField& pVals,
258 
259  const VolumeField<Type>& cCoords,
260  const Field<Type>& pCoords,
261 
262  const DynamicList<Type>& snappedPoints,
263  const labelList& snappedCc,
264  const labelList& snappedPoint,
265  const label facei,
266 
267  const scalar neiVal,
268  const Type& neiPt,
269  const bool hasNeiSnap,
270  const Type& neiSnapPt,
271 
273  DynamicList<label>& triMeshCells
274  ) const;
275 
276  template<class Type>
277  void generateTriPoints
278  (
279  const volScalarField& cVals,
280  const scalarField& pVals,
281 
282  const VolumeField<Type>& cCoords,
283  const Field<Type>& pCoords,
284 
285  const DynamicList<Type>& snappedPoints,
286  const labelList& snappedCc,
287  const labelList& snappedPoint,
288 
290  DynamicList<label>& triMeshCells
291  ) const;
292 
293  template<class Type>
294  static tmp<Field<Type>> interpolate
295  (
296  const label nPoints,
297  const labelList& triPointMergeMap,
298  const labelList& interpolatedPoints,
299  const List<FixedList<label, 3>>& interpolatedOldPoints,
301  const DynamicList<Type>& unmergedValues
302  );
303 
304  triSurface stitchTriPoints
305  (
306  const bool checkDuplicates,
307  const List<point>& triPoints,
308  labelList& triPointReverseMap, // unmerged to merged point
309  labelList& triMap // merged to unmerged triangle
310  ) const;
311 
312  //- Trim triangle to planes
313  static void trimToPlanes
314  (
315  const PtrList<plane>& planes,
316  const triPointRef& tri,
317  DynamicList<point>& newTriPoints
318  );
319 
320  //- Trim all triangles to box
321  static void trimToBox
322  (
323  const treeBoundBox& bb,
325  DynamicList<label>& triMeshCells
326  );
327 
328  //- Trim all triangles to box. Determine interpolation
329  // for existing and new points
330  static void trimToBox
331  (
332  const treeBoundBox& bb,
334  DynamicList<label>& triMap,
335  labelList& triPointMap,
336  labelList& interpolatedPoints,
337  List<FixedList<label, 3>>& interpolatedOldPoints,
339  );
340 
341  static triSurface subsetMesh
342  (
343  const triSurface&,
344  const labelList& newToOldFaces,
345  labelList& oldToNewPoints,
346  labelList& newToOldPoints
347  );
348 
349 
350  //- Interpolates cCoords, pCoords.
351  // Uses the references to the original fields used to create the
352  // iso surface.
353  template<class Type>
354  tmp<Field<Type>> interpolateTemplate
355  (
356  const VolumeField<Type>& cCoords,
357  const Field<Type>& pCoords
358  ) const;
359 
360 
361 public:
362 
363  //- Declare friendship to share some functionality
364  friend class isoSurfaceCell;
365  friend class isoSurfaceTopo;
366 
367 
368  //- Runtime type information
369  TypeName("isoSurfacePoint");
370 
371 
372  // Constructors
373 
374  //- Construct from cell values and point values.
375  // Uses boundaryField for boundary values.
376  // Holds reference to cellIsoVals and pointIsoVals.
377  //
378  // Control parameters include
379  // - bounds optional bounding box for trimming
380  // - mergeTol fraction of mesh bounding box for merging points
382  (
383  const volScalarField& cellValues,
384  const scalarField& pointValues,
385  const scalar iso,
386  const isoSurfaceParams& params = isoSurfaceParams(),
387  const bitSet& ignoreCells = bitSet()
388  );
389 
390 
391  //- Destructor
392  virtual ~isoSurfacePoint() = default;
393 
394 
395  // Member Functions
396 
397  // Sampling
398 
404 };
405 
406 
407 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
408 
409 } // End namespace Foam
410 
411 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
412 
413 #ifdef NoRepository
414  #include "isoSurfacePointTemplates.C"
415 #endif
417 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
418 
419 #endif
420 
421 // ************************************************************************* //
const scalarField & pointValues() const noexcept
The mesh point values used for creating the iso-surface.
label nPoints() const
Number of points supporting patch faces.
virtual ~isoSurfacePoint()=default
Destructor.
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
Marching tet iso surface algorithm with optional filtering to keep only points originating from mesh ...
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
Preferences for controlling iso-surface algorithms.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Abstract base class for interpolating in 1D.
TypeName("isoSurfacePoint")
Runtime type information.
isoSurfacePoint(const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const isoSurfaceParams &params=isoSurfaceParams(), const bitSet &ignoreCells=bitSet())
Construct from cell values and point values.
Low-level components common to various iso-surface algorithms.
isoSurfaceParams(const algorithmType algo=algorithmType::ALGO_DEFAULT, const filterType filter=filterType::DIAGCELL) noexcept
Default construct, or with specified algorithm.
#define declareIsoSurfaceInterpolateMethod(Type)
A surface formed by the iso value. After "Regularised Marching Tetrahedra: improved iso-surface extra...
Triangle point storage. Default constructable (triangle is not)
Definition: triangle.H:74
const scalarField & cellValues() const noexcept
The mesh cell values used for creating the iso-surface.
A surface formed by the iso value. After "Polygonising A Scalar Field Using Tetrahedrons", Paul Bourke (http://paulbourke.net/geometry/polygonise) and "Regularised Marching Tetrahedra: improved iso-surface extraction", G.M. Treece, R.W. Prager and A.H. Gee.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
labelList f(nPoints)
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))
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Triangulated surface description with patch information.
Definition: triSurface.H:71
The boundaryRegion persistent data saved as a Map<dictionary>.
Tensor of scalars, i.e. Tensor<scalar>.
const volScalarField & p0
Definition: EEqn.H:36
Namespace for OpenFOAM.
const pointField & pts