RTK  2.0.1
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:
75 
81 
83  using InputImageType = VolumeSeriesType;
84  using OutputImageType = VolumeSeriesType;
85  using VolumeType = ProjectionStackType;
86  using VectorForDVF = itk::CovariantVector< typename VolumeSeriesType::ValueType, VolumeSeriesType::ImageDimension - 1>;
87 
90 
92  using CPUVolumeSeriesType = typename itk::Image< typename VolumeSeriesType::PixelType,
93  VolumeSeriesType::ImageDimension>;
94 #ifdef RTK_USE_CUDA
95  typedef typename std::conditional< std::is_same< VolumeSeriesType, CPUVolumeSeriesType >::value,
97  itk::CudaImage<VectorForDVF, VolumeSeriesType::ImageDimension> >::type
99  typedef typename std::conditional< std::is_same< VolumeSeriesType, CPUVolumeSeriesType >::value,
100  itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension - 1>,
101  itk::CudaImage<VectorForDVF, VolumeSeriesType::ImageDimension - 1> >::type
102  DVFImageType;
103 #else
104  using DVFSequenceImageType = itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension>;
105  using DVFImageType = itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension-1>;
106 #endif
107 
111  itkNewMacro(Self)
112 
113 
115 
117  void SetForwardProjectionFilter (ForwardProjectionType itkNotUsed(_arg)) override {itkExceptionMacro(<< "ForwardProjection cannot be changed");}
118  void SetBackProjectionFilter (BackProjectionType itkNotUsed(_arg)) override {itkExceptionMacro(<< "BackProjection cannot be changed");}
120 
122  void SetDisplacementField(const DVFSequenceImageType* DVFs);
123  void SetInverseDisplacementField(const DVFSequenceImageType* DVFs);
124  typename DVFSequenceImageType::ConstPointer GetDisplacementField();
125  typename DVFSequenceImageType::ConstPointer GetInverseDisplacementField();
127 
129  void SetSignal(const std::vector<double> signal) override;
130 
131  // Sub filters type alias
134 
136  itkSetMacro(UseCudaCyclicDeformation, bool)
137  itkGetMacro(UseCudaCyclicDeformation, bool)
138 
139 protected:
141  ~MotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter() override = default;
142 
143  void GenerateOutputInformation() override;
144  void GenerateInputRequestedRegion() override;
145 
146  bool m_UseCudaCyclicDeformation;
147 
148 }; // end of class
149 
150 } // end namespace rtk
151 
152 #ifndef ITK_MANUAL_INSTANTIATION
153 #include "rtkMotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter.hxx"
154 #endif
155 
156 #endif
typename 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.
#define itkSetMacro(name, type)