ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkLevelSetTovtkImageData.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 
19 #ifndef itkLevelSetTovtkImageData_h
20 #define itkLevelSetTovtkImageData_h
21 
23 
24 #include "itkLevelSetDenseImage.h"
28 
29 #include "itkImage.h"
32 
34 
35 namespace itk
36 {
37 template< typename TLevelSet >
39  {};
40 
44 template< typename TImage >
46  public LevelSetTovtkImageDataBase< LevelSetDenseImage< TImage > >
47 {
48 public:
49  typedef TImage ImageType;
51 
56 
58  itkNewMacro(Self);
59 
62 
63  typedef typename LevelSetType::Pointer LevelSetPointer;
64 
65  vtkImageData* GetOutput() const;
66 
67 protected:
69  virtual ~LevelSetTovtkImageData();
70 
71  void GenerateData();
72 
73 private:
74  LevelSetTovtkImageData( const Self& );
75  void operator = ( const Self& );
76 
79 
81 };
82 
83 // -----------------------------------------------------------------------------
84 template< typename TOutput, unsigned int VDimension >
85 class LevelSetTovtkImageData< WhitakerSparseLevelSetImage< TOutput, VDimension > > :
86  public LevelSetTovtkImageDataBase< WhitakerSparseLevelSetImage< TOutput, VDimension > >
87 {
88 public:
90 
95 
97  itkNewMacro(Self);
98 
101 
102  typedef typename LevelSetType::Pointer LevelSetPointer;
103 
104  vtkImageData* GetOutput() const;
105 
106 protected:
108  virtual ~LevelSetTovtkImageData();
109 
110  void GenerateData();
111 
112 private:
113  LevelSetTovtkImageData( const Self& );
114  void operator = ( const Self& );
115 
118 
121 
124 };
125 
126 
127 // -----------------------------------------------------------------------------
128 template< unsigned int VDimension >
130  public LevelSetTovtkImageDataBase< ShiSparseLevelSetImage< VDimension > >
131 {
132 public:
134 
139 
141  itkNewMacro(Self);
142 
145 
146  typedef typename LevelSetType::Pointer LevelSetPointer;
147 
148  vtkImageData* GetOutput() const;
149 
150 protected:
152  virtual ~LevelSetTovtkImageData();
153 
154  void GenerateData();
155 
156 private:
157  LevelSetTovtkImageData( const Self& );
158  void operator = ( const Self& );
159 
162 
165 
168 
170 
173 
175 };
176 
177 // -----------------------------------------------------------------------------
178 template< unsigned int VDimension >
180  public LevelSetTovtkImageDataBase< MalcolmSparseLevelSetImage< VDimension > >
181 {
182 public:
184 
189 
191  itkNewMacro(Self);
192 
195 
196  typedef typename LevelSetType::Pointer LevelSetPointer;
197 
198  vtkImageData* GetOutput() const;
199 
200 protected:
202  virtual ~LevelSetTovtkImageData();
203 
204  void GenerateData();
205 
206 private:
207  LevelSetTovtkImageData( const Self& );
208  void operator = ( const Self& );
209 
212 
215 
218 
220 
223 
225 };
226 }
227 
228 #ifndef ITK_MANUAL_INSTANTIATION
229 #include "itkLevelSetTovtkImageData.hxx"
230 #endif
231 #endif // itkLevelSetTovtkImageData_h
Light weight base class for most itk classes.
Derived class for the shi representation of level-set function.
Superclass::LabelMapType LabelMapType
Derived class for the Malcolm representation of level-set function.
LabelMapToLabelImageFilter< LabelMapType, ImageType > LabelMapToLabelImageFilterType
Derived class for the sparse-field representation of level-set function.
Converts an ITK image into a VTK image and plugs a itk data pipeline to a VTK datapipeline.
Base class for the &quot;dense&quot; representation of a level-set function on one image.
LabelMapToLabelImageFilter< LabelMapType, ImageType > LabelMapToLabelImageFilterType
Templated n-dimensional image class.
Definition: itkImage.h:75
Converts a LabelMap to a labeled image.