RTK  2.0.0
Reconstruction Toolkit
rtkMotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter.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 rtkMotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter_h
20 #define rtkMotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter_h
21 
25 #ifdef RTK_USE_CUDA
28 #endif
29 
30 namespace rtk
31 {
69 template< typename VolumeSeriesType, typename ProjectionStackType>
71  public rtk::FourDConjugateGradientConeBeamReconstructionFilter<VolumeSeriesType, ProjectionStackType>
72 {
73 public:
79 
81  typedef VolumeSeriesType InputImageType;
82  typedef VolumeSeriesType OutputImageType;
83  typedef ProjectionStackType VolumeType;
84  typedef itk::CovariantVector< typename VolumeSeriesType::ValueType, VolumeSeriesType::ImageDimension - 1> VectorForDVF;
85 
88 
90  typedef typename itk::Image< typename VolumeSeriesType::PixelType,
91  VolumeSeriesType::ImageDimension> CPUVolumeSeriesType;
92 #ifdef RTK_USE_CUDA
93  typedef typename std::conditional< std::is_same< VolumeSeriesType, CPUVolumeSeriesType >::value,
95  itk::CudaImage<VectorForDVF, VolumeSeriesType::ImageDimension> >::type
97  typedef typename std::conditional< std::is_same< VolumeSeriesType, CPUVolumeSeriesType >::value,
98  itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension - 1>,
99  itk::CudaImage<VectorForDVF, VolumeSeriesType::ImageDimension - 1> >::type
100  DVFImageType;
101 #else
102  typedef itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension> DVFSequenceImageType;
103  typedef itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension-1> DVFImageType;
104 #endif
105 
109  itkNewMacro(Self)
110 
111 
113 
115  void SetForwardProjectionFilter (ForwardProjectionType itkNotUsed(_arg)) ITK_OVERRIDE {itkExceptionMacro(<< "ForwardProjection cannot be changed");}
116  void SetBackProjectionFilter (BackProjectionType itkNotUsed(_arg)) ITK_OVERRIDE {itkExceptionMacro(<< "BackProjection cannot be changed");}
118 
120  void SetDisplacementField(const DVFSequenceImageType* DVFs);
121  void SetInverseDisplacementField(const DVFSequenceImageType* DVFs);
122  typename DVFSequenceImageType::ConstPointer GetDisplacementField();
123  typename DVFSequenceImageType::ConstPointer GetInverseDisplacementField();
125 
127  void SetSignal(const std::vector<double> signal) ITK_OVERRIDE;
128 
129  // Sub filters typedefs
132 
134  itkSetMacro(UseCudaCyclicDeformation, bool)
135  itkGetMacro(UseCudaCyclicDeformation, bool)
136 
137 protected:
139  virtual ~MotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter() ITK_OVERRIDE {}
140 
141  void GenerateOutputInformation() ITK_OVERRIDE;
142  void GenerateInputRequestedRegion() ITK_OVERRIDE;
143 
145 
146 private:
147  //purposely not implemented
149  void operator=(const Self&);
150 }; // end of class
151 
152 } // end namespace rtk
153 
154 #ifndef ITK_MANUAL_INSTANTIATION
155 #include "rtkMotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter.hxx"
156 #endif
157 
158 #endif
rtk::WarpProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType > MCProjStackToFourDType
rtk::MotionCompensatedFourDReconstructionConjugateGradientOperator< VolumeSeriesType, ProjectionStackType > MCCGOperatorType
itk::CovariantVector< typename VolumeSeriesType::ValueType, VolumeSeriesType::ImageDimension-1 > VectorForDVF
itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension > CPUVolumeSeriesType
Back projection part for motion compensated iterative 4D reconstruction.
Implements part of the 4D reconstruction by conjugate gradient.
FourDConjugateGradientConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType > Superclass
#define itkSetMacro(name, type)