VMTK
vtkvmtkNonManifoldFastMarching.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3 Program: VMTK
4 Module: $RCSfile: vtkvmtkNonManifoldFastMarching.h,v $
5 Language: C++
6 
7  Copyright (c) Luca Antiga, David Steinman. All rights reserved.
8  See LICENSE file for details.
9 
10  Portions of this code are covered under the VTK copyright.
11  See VTKCopyright.txt or http://www.kitware.com/VTKCopyright.htm
12  for details.
13 
14  This software is distributed WITHOUT ANY WARRANTY; without even
15  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
16  PURPOSE. See the above copyright notices for more information.
17 
18 =========================================================================*/
37 #ifndef __vtkvmtkNonManifoldFastMarching_h
38 #define __vtkvmtkNonManifoldFastMarching_h
39 
40 #include "vtkPolyDataAlgorithm.h"
41 #include "vtkIdList.h"
42 #include "vtkIntArray.h"
43 #include "vtkPolyData.h"
44 #include "vtkvmtkMinHeap.h"
45 #include "vtkvmtkConstants.h"
46 //#include "vtkvmtkComputationalGeometryWin32Header.h"
47 #include "vtkvmtkWin32Header.h"
48 
49 const char VTK_VMTK_ACCEPTED_STATUS = 0x01;
50 const char VTK_VMTK_CONSIDERED_STATUS = 0x02;
51 const char VTK_VMTK_FAR_STATUS = 0x04;
52 
53 class vtkDoubleArray;
54 class vtkCharArray;
55 
56 class VTK_VMTK_COMPUTATIONAL_GEOMETRY_EXPORT vtkvmtkNonManifoldFastMarching : public vtkPolyDataAlgorithm
57 {
58  public:
59  vtkTypeMacro(vtkvmtkNonManifoldFastMarching,vtkPolyDataAlgorithm);
60  void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
61 
62  static vtkvmtkNonManifoldFastMarching *New();
63 
65 
66  vtkSetMacro(StopTravelTime,double);
67  vtkGetMacro(StopTravelTime,double);
69 
71 
73  vtkSetMacro(StopNumberOfPoints,int);
74  vtkGetMacro(StopNumberOfPoints,int);
76 
78 
79  vtkSetMacro(Regularization,double);
80  vtkGetMacro(Regularization,double);
82 
84 
85  vtkSetMacro(SeedsBoundaryConditions,int);
86  vtkGetMacro(SeedsBoundaryConditions,int);
87  vtkBooleanMacro(SeedsBoundaryConditions,int);
89 
91 
92  vtkSetMacro(PolyDataBoundaryConditions,int);
93  vtkGetMacro(PolyDataBoundaryConditions,int);
94  vtkBooleanMacro(PolyDataBoundaryConditions,int);
96 
98 
99  vtkSetObjectMacro(Seeds,vtkIdList);
100  vtkGetObjectMacro(Seeds,vtkIdList);
102 
104 
105  vtkSetObjectMacro(BoundaryPolyData,vtkPolyData);
106  vtkGetObjectMacro(BoundaryPolyData,vtkPolyData);
108 
110 
113  vtkSetStringMacro(IntersectedEdgesArrayName);
114  vtkGetStringMacro(IntersectedEdgesArrayName);
116 
118 
121  vtkSetMacro(InitializeFromScalars,int);
122  vtkGetMacro(InitializeFromScalars,int);
123  vtkBooleanMacro(InitializeFromScalars,int);
125 
127 
128  vtkSetStringMacro(InitializationArrayName);
129  vtkGetStringMacro(InitializationArrayName);
131 
133 
136  vtkSetMacro(UnitSpeed,int);
137  vtkGetMacro(UnitSpeed,int);
138  vtkBooleanMacro(UnitSpeed,int);
140 
142 
143  vtkSetStringMacro(CostFunctionArrayName);
144  vtkGetStringMacro(CostFunctionArrayName);
146 
148 
149  vtkSetStringMacro(SolutionArrayName);
150  vtkGetStringMacro(SolutionArrayName);
152 
153  protected:
156 
157  virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) VTK_OVERRIDE;
158 
159  void InitPropagation(vtkPolyData* input);
160 
161  void SolveQuadratic(double a, double b, double c, char &nSol, double &x0, double &x1);
162 
163  void GetNeighbors(vtkPolyData* input, vtkIdType pointId, vtkIdList* neighborIds);
164  double ComputeUpdateFromCellNeighbor(vtkPolyData* input, vtkIdType neighborId, vtkIdType* trianglePts);
165  void UpdateNeighbor(vtkPolyData* input, vtkIdType neighborId);
166  void UpdateNeighborhood(vtkPolyData* input, vtkIdType pointId);
167  void Propagate(vtkPolyData* input);
168 
169  static double Max(double a, double b)
170  { return a-b > VTK_VMTK_DOUBLE_TOL ? a : b; }
171 
172  static double Min(double a, double b)
173  { return a-b < - VTK_VMTK_DOUBLE_TOL ? a : b; }
174 
175  vtkDoubleArray* TScalars;
176  vtkCharArray* StatusScalars;
178 
179  vtkIdList* Seeds;
180  vtkPolyData* BoundaryPolyData;
181 
191 
194 
196 
199 
200  private:
202  void operator=(const vtkvmtkNonManifoldFastMarching&); // Not implemented.
203 };
204 
205 #endif
206 
const char VTK_VMTK_FAR_STATUS
Implementation of the Fast Marching Method on polygonal non-manifolds.
#define VTK_VMTK_DOUBLE_TOL
const char VTK_VMTK_ACCEPTED_STATUS
const char VTK_VMTK_CONSIDERED_STATUS
static double Min(double a, double b)
static double Max(double a, double b)
Implementation of the min heap data structure.