isoSurfaceCell.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::isoSurfaceCell
26 
27 Description
28  A surface formed by the iso value.
29  After "Polygonising A Scalar Field Using Tetrahedrons", Paul Bourke
30  (http://paulbourke.net/geometry/polygonise) and
31  "Regularised Marching Tetrahedra: improved iso-surface extraction",
32  G.M. Treece, R.W. Prager and A.H. Gee.
33 
34  See isoSurface. This is a variant which does tetrahedrisation from
35  triangulation of face to cell centre instead of edge of face to two
36  neighbouring cell centres. This gives much lower quality triangles
37  but they are local to a cell.
38 
39 SourceFiles
40  isoSurfaceCell.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef isoSurfaceCell_H
45 #define isoSurfaceCell_H
46 
47 #include "triSurface.H"
48 #include "labelPair.H"
49 #include "pointIndexHit.H"
50 #include "PackedBoolList.H"
51 #include "boundBox.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 class polyMesh;
59 class plane;
60 
61 /*---------------------------------------------------------------------------*\
62  Class isoSurfaceCell Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 class isoSurfaceCell
66 :
67  public triSurface
68 {
69  // Private data
70 
71  enum segmentCutType
72  {
73  NEARFIRST, // intersection close to e.first()
74  NEARSECOND, // ,, e.second()
75  NOTNEAR // no intersection
76  };
77 
78  enum cellCutType
79  {
80  NOTCUT, // not cut
81  SPHERE, // all edges to cell centre cut
82  CUT // normal cut
83  };
84 
85 
86  //- Reference to mesh
87  const polyMesh& mesh_;
88 
89  const scalarField& cVals_;
90 
91  const scalarField& pVals_;
92 
93  //- isoSurfaceCell value
94  const scalar iso_;
95 
96  //- Optional bounds
97  const boundBox bounds_;
98 
99  //- When to merge points
100  const scalar mergeDistance_;
101 
102  //- Whether cell might be cut
104 
105  //- Estimated number of cut cells
107 
108  //- For every triangle the original cell in mesh
110 
111  //- For every unmerged triangle point the point in the triSurface
113 
114  //- triSurface points that have weighted interpolation
116 
117  //- corresponding original, unmerged points
119 
120  //- corresponding weights
122 
123 
124  // Private Member Functions
125 
126  //- Get location of iso value as fraction inbetween s0,s1
127  scalar isoFraction
128  (
129  const scalar s0,
130  const scalar s1
131  ) const;
132 
133  //- Does any edge of triangle cross iso value?
134  bool isTriCut
135  (
136  const triFace& tri,
137  const scalarField& pointValues
138  ) const;
139 
140  //- Determine whether cell is cut
142  (
143  const PackedBoolList&,
144  const scalarField& cellValues,
145  const scalarField& pointValues,
146  const label
147  ) const;
148 
149  void calcCutTypes
150  (
151  const PackedBoolList&,
152  const scalarField& cellValues,
153  const scalarField& pointValues
154  );
155 
157  (
158  const labelledTri&,
159  const labelledTri&
160  );
161 
162  static point calcCentre(const triSurface&);
163 
165  (
166  const label cellI,
169  ) const;
170 
171  //- Determine per cc whether all near cuts can be snapped to single
172  // point.
173  void calcSnappedCc
174  (
175  const PackedBoolList& isTet,
176  const scalarField& cVals,
177  const scalarField& pVals,
179  labelList& snappedCc
180  ) const;
181 
182  //- Generate triangles for face connected to pointI
183  void genPointTris
184  (
185  const scalarField& cellValues,
186  const scalarField& pointValues,
187  const label pointI,
188  const label faceI,
189  const label cellI,
190  DynamicList<point, 64>& localTriPoints
191  ) const;
192 
193  //- Generate triangles for tet connected to pointI
194  void genPointTris
195  (
196  const scalarField& pointValues,
197  const label pointI,
198  const label faceI,
199  const label cellI,
200  DynamicList<point, 64>& localTriPoints
201  ) const;
202 
203  //- Determine per point whether all near cuts can be snapped to single
204  // point.
205  void calcSnappedPoint
206  (
207  const PackedBoolList& isTet,
208  const scalarField& cVals,
209  const scalarField& pVals,
211  labelList& snappedPoint
212  ) const;
213 
214  //- Generate single point by interpolation or snapping
215  template<class Type>
216  Type generatePoint
217  (
218  const DynamicList<Type>& snappedPoints,
219  const scalar s0,
220  const Type& p0,
221  const label p0Index,
222  const scalar s1,
223  const Type& p1,
224  const label p1Index
225  ) const;
226 
227  template<class Type>
228  void generateTriPoints
229  (
230  const DynamicList<Type>& snapped,
231  const scalar s0,
232  const Type& p0,
233  const label p0Index,
234  const scalar s1,
235  const Type& p1,
236  const label p1Index,
237  const scalar s2,
238  const Type& p2,
239  const label p2Index,
240  const scalar s3,
241  const Type& p3,
242  const label p3Index,
244  ) const;
245 
246  template<class Type>
247  void generateTriPoints
248  (
249  const scalarField& cVals,
250  const scalarField& pVals,
251 
252  const Field<Type>& cCoords,
253  const Field<Type>& pCoords,
254 
255  const DynamicList<Type>& snappedPoints,
256  const labelList& snappedCc,
257  const labelList& snappedPoint,
258 
260  DynamicList<label>& triMeshCells
261  ) const;
262 
264  (
265  const bool checkDuplicates,
266  const List<point>& triPoints,
267  labelList& triPointReverseMap, // unmerged to merged point
268  labelList& triMap // merged to unmerged triangle
269  ) const;
270 
271  //- Check single triangle for (topological) validity
272  static bool validTri(const triSurface&, const label);
273 
274  //- Determine edge-face addressing
275  void calcAddressing
276  (
277  const triSurface& surf,
279  labelList& edgeFace0,
280  labelList& edgeFace1,
281  Map<labelList>& edgeFacesRest
282  ) const;
283 
284  //- Is triangle (given by 3 edges) not fully connected?
285  static bool danglingTriangle
286  (
287  const FixedList<label, 3>& fEdges,
288  const labelList& edgeFace1
289  );
290 
291  //- Mark all non-fully connected triangles
293  (
295  const labelList& edgeFace0,
296  const labelList& edgeFace1,
297  const Map<labelList>& edgeFacesRest,
298  boolList& keepTriangles
299  );
300 
301  static triSurface subsetMesh
302  (
303  const triSurface& s,
304  const labelList& newToOldFaces,
305  labelList& oldToNewPoints,
306  labelList& newToOldPoints
307  );
308 
309 
310 public:
311 
312  //- Runtime type information
313  TypeName("isoSurfaceCell");
314 
315 
316  // Constructors
317 
318  //- Construct from dictionary
320  (
321  const polyMesh& mesh,
322  const scalarField& cellValues,
323  const scalarField& pointValues,
324  const scalar iso,
325  const bool regularise,
326  const boundBox& bounds = boundBox::greatBox,
327  const scalar mergeTol = 1e-6 // fraction of bounding box
328  );
329 
330 
331  // Member Functions
332 
333  //- For every face original cell in mesh
334  const labelList& meshCells() const
335  {
336  return meshCells_;
337  }
338 
339  //- Interpolates cCoords,pCoords.
340  template<class Type>
342  (
343  const Field<Type>& cCoords,
344  const Field<Type>& pCoords
345  ) const;
346 };
347 
348 
349 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
350 
351 } // End namespace Foam
352 
353 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354 
355 #ifdef NoRepository
356 # include "isoSurfaceCellTemplates.C"
357 #endif
358 
359 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
360 
361 #endif
362 
363 // ************************************************************************* //
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
pointIndexHit.H
Foam::isoSurfaceCell::genPointTris
void genPointTris(const scalarField &cellValues, const scalarField &pointValues, const label pointI, const label faceI, const label cellI, DynamicList< point, 64 > &localTriPoints) const
Generate triangles for face connected to pointI.
Definition: isoSurfaceCell.C:540
Foam::PrimitivePatch< labelledTri, List, pointField, point >::points
const Field< point > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::isoSurfaceCell::calcCutTypes
void calcCutTypes(const PackedBoolList &, const scalarField &cellValues, const scalarField &pointValues)
Definition: isoSurfaceCell.C:192
Foam::isoSurfaceCell::interpolatedPoints_
DynamicList< label > interpolatedPoints_
triSurface points that have weighted interpolation
Definition: isoSurfaceCell.H:114
isoSurfaceCellTemplates.C
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::DynamicList< label >
Foam::PrimitivePatch< labelledTri, List, pointField, point >::localPoints
const Field< point > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::isoSurfaceCell::markDanglingTriangles
static label markDanglingTriangles(const List< FixedList< label, 3 > > &faceEdges, const labelList &edgeFace0, const labelList &edgeFace1, const Map< labelList > &edgeFacesRest, boolList &keepTriangles)
Mark all non-fully connected triangles.
Definition: isoSurfaceCell.C:1257
Foam::isoSurfaceCell::NOTNEAR
@ NOTNEAR
Definition: isoSurfaceCell.H:74
Foam::isoSurfaceCell::interpolatedOldPoints_
DynamicList< FixedList< label, 3 > > interpolatedOldPoints_
corresponding original, unmerged points
Definition: isoSurfaceCell.H:117
Foam::isoSurfaceCell::isoFraction
scalar isoFraction(const scalar s0, const scalar s1) const
Get location of iso value as fraction inbetween s0,s1.
Definition: isoSurfaceCell.C:48
Foam::isoSurfaceCell::generateTriPoints
void generateTriPoints(const DynamicList< Type > &snapped, const scalar s0, const Type &p0, const label p0Index, const scalar s1, const Type &p1, const label p1Index, const scalar s2, const Type &p2, const label p2Index, const scalar s3, const Type &p3, const label p3Index, DynamicList< Type > &points) const
Definition: isoSurfaceCellTemplates.C:77
Foam::Map
A HashTable to objects of type <T> with a label key.
Definition: PrimitivePatchTemplate.H:68
Foam::isoSurfaceCell::cellCutType
cellCutType
Definition: isoSurfaceCell.H:77
Foam::isoSurfaceCell::NOTCUT
@ NOTCUT
Definition: isoSurfaceCell.H:79
Foam::PrimitivePatch< labelledTri, List, pointField, point >::faceEdges
const labelListList & faceEdges() const
Return face-edge addressing.
Definition: PrimitivePatchTemplate.C:312
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::isoSurfaceCell
A surface formed by the iso value. After "Polygonising A Scalar Field Using Tetrahedrons",...
Definition: isoSurfaceCell.H:64
Foam::isoSurfaceCell::TypeName
TypeName("isoSurfaceCell")
Runtime type information.
Foam::isoSurfaceCell::findCommonPoints
static labelPair findCommonPoints(const labelledTri &, const labelledTri &)
Definition: isoSurfaceCell.C:221
Foam::isoSurfaceCell::pVals_
const scalarField & pVals_
Definition: isoSurfaceCell.H:90
Foam::isoSurfaceCell::danglingTriangle
static bool danglingTriangle(const FixedList< label, 3 > &fEdges, const labelList &edgeFace1)
Is triangle (given by 3 edges) not fully connected?
Definition: isoSurfaceCell.C:1230
Foam::isoSurfaceCell::NEARFIRST
@ NEARFIRST
Definition: isoSurfaceCell.H:72
Foam::triPoints
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:49
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::isoSurfaceCell::generatePoint
Type generatePoint(const DynamicList< Type > &snappedPoints, const scalar s0, const Type &p0, const label p0Index, const scalar s1, const Type &p1, const label p1Index) const
Generate single point by interpolation or snapping.
Definition: isoSurfaceCellTemplates.C:35
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::isoSurfaceCell::calcSnappedPoint
void calcSnappedPoint(const PackedBoolList &isTet, const scalarField &cVals, const scalarField &pVals, DynamicList< point > &triPoints, labelList &snappedPoint) const
Determine per point whether all near cuts can be snapped to single.
Definition: isoSurfaceCell.C:694
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
Foam::isoSurfaceCell::meshCells
const labelList & meshCells() const
For every face original cell in mesh.
Definition: isoSurfaceCell.H:333
PackedBoolList.H
Foam::isoSurfaceCell::validTri
static bool validTri(const triSurface &, const label)
Check single triangle for (topological) validity.
Definition: isoSurfaceCell.C:1060
Foam::isoSurfaceCell::stitchTriPoints
triSurface stitchTriPoints(const bool checkDuplicates, const List< point > &triPoints, labelList &triPointReverseMap, labelList &triMap) const
Definition: isoSurfaceCell.C:912
Foam::isoSurfaceCell::nCutCells_
label nCutCells_
Estimated number of cut cells.
Definition: isoSurfaceCell.H:105
Foam::PrimitivePatch< labelledTri, List, pointField, point >::calcAddressing
void calcAddressing() const
Calculate addressing.
Definition: PrimitivePatchAddressing.C:52
Foam::isoSurfaceCell::cVals_
const scalarField & cVals_
Definition: isoSurfaceCell.H:88
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::isoSurfaceCell::interpolationWeights_
DynamicList< FixedList< scalar, 3 > > interpolationWeights_
corresponding weights
Definition: isoSurfaceCell.H:120
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::isoSurfaceCell::triPointMergeMap_
labelList triPointMergeMap_
For every unmerged triangle point the point in the triSurface.
Definition: isoSurfaceCell.H:111
Foam::isoSurfaceCell::isoSurfaceCell
isoSurfaceCell(const polyMesh &mesh, const scalarField &cellValues, const scalarField &pointValues, const scalar iso, const bool regularise, const boundBox &bounds=boundBox::greatBox, const scalar mergeTol=1e-6)
Construct from dictionary.
Definition: isoSurfaceCell.C:1379
Foam::isoSurfaceCell::calcCutType
cellCutType calcCutType(const PackedBoolList &, const scalarField &cellValues, const scalarField &pointValues, const label) const
Determine whether cell is cut.
Definition: isoSurfaceCell.C:81
Foam::isoSurfaceCell::calcSnappedCc
void calcSnappedCc(const PackedBoolList &isTet, const scalarField &cVals, const scalarField &pVals, DynamicList< point > &triPoints, labelList &snappedCc) const
Determine per cc whether all near cuts can be snapped to single.
Definition: isoSurfaceCell.C:365
boundBox.H
Foam::isoSurfaceCell::bounds_
const boundBox bounds_
Optional bounds.
Definition: isoSurfaceCell.H:96
Foam::isoSurfaceCell::iso_
const scalar iso_
isoSurfaceCell value
Definition: isoSurfaceCell.H:93
Foam::isoSurfaceCell::cellCutType_
List< cellCutType > cellCutType_
Whether cell might be cut.
Definition: isoSurfaceCell.H:102
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
Foam::boundBox::greatBox
static const boundBox greatBox
A very large boundBox: min/max == -/+ VGREAT.
Definition: boundBox.H:76
Foam::Vector< scalar >
Foam::List< cellCutType >
Foam::isoSurfaceCell::interpolate
tmp< Field< Type > > interpolate(const Field< Type > &cCoords, const Field< Type > &pCoords) const
Interpolates cCoords,pCoords.
Foam::isoSurfaceCell::isTriCut
bool isTriCut(const triFace &tri, const scalarField &pointValues) const
Does any edge of triangle cross iso value?
Definition: isoSurfaceCell.C:67
Foam::FixedList< label, 3 >
Foam::isoSurfaceCell::collapseSurface
pointIndexHit collapseSurface(const label cellI, pointField &localPoints, DynamicList< labelledTri, 64 > &localTris) const
Definition: isoSurfaceCell.C:272
Foam::isoSurfaceCell::NEARSECOND
@ NEARSECOND
Definition: isoSurfaceCell.H:73
Foam::isoSurfaceCell::mergeDistance_
const scalar mergeDistance_
When to merge points.
Definition: isoSurfaceCell.H:99
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::isoSurfaceCell::SPHERE
@ SPHERE
Definition: isoSurfaceCell.H:80
Foam::isoSurfaceCell::subsetMesh
static triSurface subsetMesh(const triSurface &s, const labelList &newToOldFaces, labelList &oldToNewPoints, labelList &newToOldPoints)
Definition: isoSurfaceCell.C:1305
Foam::isoSurfaceCell::meshCells_
labelList meshCells_
For every triangle the original cell in mesh.
Definition: isoSurfaceCell.H:108
Foam::isoSurfaceCell::CUT
@ CUT
Definition: isoSurfaceCell.H:81
Foam::isoSurfaceCell::mesh_
const polyMesh & mesh_
Reference to mesh.
Definition: isoSurfaceCell.H:86
Foam::isoSurfaceCell::segmentCutType
segmentCutType
Definition: isoSurfaceCell.H:70
Foam::isoSurfaceCell::calcCentre
static point calcCentre(const triSurface &)
Definition: isoSurfaceCell.C:257
labelPair.H