ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkPolylineMaskImageFilter.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 #ifndef itkPolylineMaskImageFilter_h
19 #define itkPolylineMaskImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
36 template< typename TInputImage, typename TPolyline, typename TVector,
37  typename TOutputImage >
38 class ITK_TEMPLATE_EXPORT PolylineMaskImageFilter:public ImageToImageFilter< TInputImage, TOutputImage >
39 {
40 public:
41  ITK_DISALLOW_COPY_AND_ASSIGN(PolylineMaskImageFilter);
42 
48 
50  itkNewMacro(Self);
51 
54 
56  static constexpr unsigned int NDimensions = TInputImage::ImageDimension;
57 
58  static constexpr unsigned int InputDimension = 3;
59 
61  using InputImageType = TInputImage;
62  using InputImagePointer = typename InputImageType::ConstPointer;
64  using InputImagePixelType = typename InputImageType::PixelType;
67 
70 
72  using VectorType = TVector;
73 
75  using PolylineType = TPolyline;
76 
78  using OutputImageType = TOutputImage;
79  using OutputImagePointer = typename OutputImageType::Pointer;
81  using OutputImagePixelType = typename OutputImageType::PixelType;
82 
84  void SetInput1(const InputImageType *image);
85 
87  void SetInput2(const PolylineType *polyline);
88 
90  itkSetMacro(ViewVector, VectorType);
91  itkGetConstMacro(ViewVector, VectorType);
93 
95  itkSetMacro(UpVector, VectorType);
96  itkGetConstMacro(UpVector, VectorType);
98 
100  itkSetMacro(CameraCenterPoint, PointType);
101  itkGetConstMacro(CameraCenterPoint, PointType);
103 
105  itkSetMacro(FocalDistance, double);
106  itkGetConstMacro(FocalDistance, double);
108 
110  itkSetMacro(FocalPoint, ProjPlanePointType);
111  itkGetConstMacro(FocalPoint, ProjPlanePointType);
113 
115  void GenerateRotationMatrix();
116 
118  ProjPlanePointType TransformProjectPoint(PointType inputPoint);
119 
120 #ifdef ITK_USE_CONCEPT_CHECKING
121  // Begin concept checking
122  itkConceptMacro( VectorHasNumericTraitsCheck,
124  // End concept checking
125 #endif
126 
127 protected:
129  ~PolylineMaskImageFilter() override = default;
130  void PrintSelf(std::ostream & os, Indent indent) const override;
131 
132  void GenerateData() override;
133 
134 private:
137 
139 
141 
143 
144  double m_FocalDistance{ 0.0 };
145 };
146 } // end namespace itk
147 
148 #ifndef ITK_MANUAL_INSTANTIATION
149 #include "itkPolylineMaskImageFilter.hxx"
150 #endif
151 
152 #endif
typename OutputImageType::Pointer OutputImagePointer
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
Implements image masking operation constrained by a polyline on a plane perpendicular to certain view...
typename OutputImageType::PixelType OutputImagePixelType
typename InputImageType::PixelType InputImagePixelType
typename InputImageType::Pointer InputImagePointer
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
typename InputImageType::RegionType InputImageRegionType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
#define itkConceptMacro(name, concept)