RTK  1.4.0
Reconstruction Toolkit
rtkParkerShortScanImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright RTK 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 rtkParkerShortScanImageFilter_h
20 #define rtkParkerShortScanImageFilter_h
21 
22 #include <itkInPlaceImageFilter.h>
23 #include <mutex>
24 
26 #include "rtkConfiguration.h"
27 
28 namespace rtk
29 {
30 
47 template<class TInputImage, class TOutputImage=TInputImage>
48 class ITK_EXPORT ParkerShortScanImageFilter :
49  public itk::InPlaceImageFilter<TInputImage, TOutputImage>
50 {
51 public:
57 
59  typedef TInputImage InputImageType;
60  typedef TOutputImage OutputImageType;
61  typedef typename OutputImageType::RegionType OutputImageRegionType;
63 
66 
68  itkNewMacro(Self);
69 
72 
74  itkGetModifiableObjectMacro(Geometry, GeometryType);
75  itkSetObjectMacro(Geometry, GeometryType);
77 
78 protected:
79  ParkerShortScanImageFilter(){ this->SetInPlace(true); }
80  virtual ~ParkerShortScanImageFilter() ITK_OVERRIDE {}
81 
82 #if ITK_VERSION_MAJOR<5
83  void ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, ThreadIdType threadId) ITK_OVERRIDE;
84 #else
85  void DynamicThreadedGenerateData(const OutputImageRegionType& outputRegionForThread) ITK_OVERRIDE;
86 #endif
87 
88 private:
89  ParkerShortScanImageFilter(const Self&); //purposely not implemented
90  void operator=(const Self&); //purposely not implemented
91 
93  GeometryPointer m_Geometry;
94 
100 
101  std::mutex m_WarningMutex;
102 }; // end of class
103 
104 } // end namespace rtk
105 
106 #ifndef ITK_MANUAL_INSTANTIATION
107 #include "rtkParkerShortScanImageFilter.hxx"
108 #endif
109 
110 #endif
ThreeDCircularProjectionGeometry GeometryType
itk::Image< typename TOutputImage::PixelType, 1 > WeightImageType
Projection geometry for a source and a 2-D flat panel.
itk::SmartPointer< const Self > ConstPointer
unsigned int ThreadIdType
itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass
OutputImageType::RegionType OutputImageRegionType