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-2016 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 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 class polyMesh;
58 
59 /*---------------------------------------------------------------------------*\
60  Class isoSurfaceCell Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 class isoSurfaceCell
64 :
65  public triSurface
66 {
67  // Private data
68 
69  enum segmentCutType
70  {
71  NEARFIRST, // intersection close to e.first()
72  NEARSECOND, // ,, e.second()
73  NOTNEAR // no intersection
74  };
75 
76  enum cellCutType
77  {
78  NOTCUT, // not cut
79  SPHERE, // all edges to cell centre cut
80  CUT // normal cut
81  };
82 
83 
84  //- Reference to mesh
85  const polyMesh& mesh_;
86 
87  const scalarField& cVals_;
88 
89  const scalarField& pVals_;
90 
91  //- isoSurfaceCell value
92  const scalar iso_;
93 
94  //- When to merge points
95  const scalar mergeDistance_;
96 
97  //- Whether cell might be cut
98  List<cellCutType> cellCutType_;
99 
100  //- Estimated number of cut cells
101  label nCutCells_;
102 
103  //- For every triangle the original cell in mesh
104  labelList meshCells_;
105 
106  //- For every unmerged triangle point the point in the triSurface
107  labelList triPointMergeMap_;
108 
109 
110  // Private Member Functions
111 
112  //- Get location of iso value as fraction inbetween s0,s1
113  scalar isoFraction
114  (
115  const scalar s0,
116  const scalar s1
117  ) const;
118 
119  //- Does any edge of triangle cross iso value?
120  bool isTriCut
121  (
122  const triFace& tri,
123  const scalarField& pointValues
124  ) const;
125 
126  //- Determine whether cell is cut
127  cellCutType calcCutType
128  (
129  const PackedBoolList&,
130  const scalarField& cellValues,
131  const scalarField& pointValues,
132  const label
133  ) const;
134 
135  void calcCutTypes
136  (
137  const PackedBoolList&,
138  const scalarField& cellValues,
139  const scalarField& pointValues
140  );
141 
142  //- Return the two common points between two triangles
143  static labelPair findCommonPoints
144  (
145  const labelledTri&,
146  const labelledTri&
147  );
148 
149  //- Caculate centre of surface.
150  static point calcCentre(const triSurface&);
151 
152  //- Replace surface (localPoints, localTris) with single point.
153  // Returns point. Destroys arguments.
154  pointIndexHit collapseSurface
155  (
156  const label celli,
159  ) const;
160 
161  //- Determine per cc whether all near cuts can be snapped to single
162  // point.
163  void calcSnappedCc
164  (
165  const PackedBoolList& isTet,
166  const scalarField& cVals,
167  const scalarField& pVals,
168  DynamicList<point>& triPoints,
169  labelList& snappedCc
170  ) const;
171 
172  //- Generate triangles for face connected to pointi
173  void genPointTris
174  (
175  const scalarField& cellValues,
176  const scalarField& pointValues,
177  const label pointi,
178  const label facei,
179  const label celli,
180  DynamicList<point, 64>& localTriPoints
181  ) const;
182 
183  //- Generate triangles for tet connected to pointi
184  void genPointTris
185  (
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  //- Determine per point whether all near cuts can be snapped to single
194  // point.
195  void calcSnappedPoint
196  (
197  const PackedBoolList& isTet,
198  const scalarField& cVals,
199  const scalarField& pVals,
200  DynamicList<point>& triPoints,
201  labelList& snappedPoint
202  ) const;
203 
204  //- Generate single point by interpolation or snapping
205  template<class Type>
206  Type generatePoint
207  (
208  const DynamicList<Type>& snappedPoints,
209  const scalar s0,
210  const Type& p0,
211  const label p0Index,
212  const scalar s1,
213  const Type& p1,
214  const label p1Index
215  ) const;
216 
217  template<class Type>
218  void generateTriPoints
219  (
220  const DynamicList<Type>& snapped,
221  const scalar s0,
222  const Type& p0,
223  const label p0Index,
224  const scalar s1,
225  const Type& p1,
226  const label p1Index,
227  const scalar s2,
228  const Type& p2,
229  const label p2Index,
230  const scalar s3,
231  const Type& p3,
232  const label p3Index,
234  ) const;
235 
236  template<class Type>
237  void generateTriPoints
238  (
239  const scalarField& cVals,
240  const scalarField& pVals,
241 
242  const Field<Type>& cCoords,
243  const Field<Type>& pCoords,
244 
245  const DynamicList<Type>& snappedPoints,
246  const labelList& snappedCc,
247  const labelList& snappedPoint,
248 
249  DynamicList<Type>& triPoints,
250  DynamicList<label>& triMeshCells
251  ) const;
252 
253  triSurface stitchTriPoints
254  (
255  const bool checkDuplicates,
256  const List<point>& triPoints,
257  labelList& triPointReverseMap, // unmerged to merged point
258  labelList& triMap // merged to unmerged triangle
259  ) const;
260 
261  //- Check single triangle for (topological) validity
262  static bool validTri(const triSurface&, const label);
263 
264  //- Determine edge-face addressing
265  void calcAddressing
266  (
267  const triSurface& surf,
269  labelList& edgeFace0,
270  labelList& edgeFace1,
271  Map<labelList>& edgeFacesRest
272  ) const;
273 
274  //- Is triangle (given by 3 edges) not fully connected?
275  static bool danglingTriangle
276  (
277  const FixedList<label, 3>& fEdges,
278  const labelList& edgeFace1
279  );
280 
281  //- Mark all non-fully connected triangles
282  static label markDanglingTriangles
283  (
285  const labelList& edgeFace0,
286  const labelList& edgeFace1,
287  const Map<labelList>& edgeFacesRest,
288  boolList& keepTriangles
289  );
290 
291  static triSurface subsetMesh
292  (
293  const triSurface& s,
294  const labelList& newToOldFaces,
295  labelList& oldToNewPoints,
296  labelList& newToOldPoints
297  );
298 
299 
300 public:
301 
302  //- Runtime type information
303  TypeName("isoSurfaceCell");
304 
305 
306  // Constructors
307 
308  //- Construct from dictionary
310  (
311  const polyMesh& mesh,
312  const scalarField& cellValues,
313  const scalarField& pointValues,
314  const scalar iso,
315  const bool regularise,
316  const scalar mergeTol = 1e-6 // fraction of bounding box
317  );
318 
319 
320  // Member Functions
321 
322  //- For every face original cell in mesh
323  const labelList& meshCells() const
324  {
325  return meshCells_;
326  }
327 
328  //- Interpolates cCoords,pCoords.
329  template<class Type>
331  (
332  const Field<Type>& cCoords,
333  const Field<Type>& pCoords
334  ) const;
335 };
336 
337 
338 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
339 
340 } // End namespace Foam
341 
342 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
343 
344 #ifdef NoRepository
345  #include "isoSurfaceCellTemplates.C"
346 #endif
347 
348 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
349 
350 #endif
351 
352 // ************************************************************************* //
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
TypeName("isoSurfaceCell")
Runtime type information.
tmp< Field< Type > > interpolate(const Field< Type > &cCoords, const Field< Type > &pCoords) const
Interpolates cCoords,pCoords.
const double e
Elementary charge.
Definition: doubleFloat.H:78
const labelList & meshCells() const
For every face original cell in mesh.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
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))
dynamicFvMesh & mesh
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
const Field< point > & points() const
Return reference to global points.
Triangle with additional region number.
Definition: labelledTri.H:57
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.
const Field< point > & localPoints() const
Return pointField of points in patch.
isoSurfaceCell(const polyMesh &mesh, const scalarField &cellValues, const scalarField &pointValues, const scalar iso, const bool regularise, const scalar mergeTol=1e-6)
Construct from dictionary.
A bit-packed bool list.
const labelListList & faceEdges() const
Return face-edge addressing.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
Triangulated surface description with patch information.
Definition: triSurface.H:65
Namespace for OpenFOAM.
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49