RTK  2.7.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  * https://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 
88  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
89  using BackProjectionType = typename Superclass::BackProjectionType;
90 
92  using CPUVolumeSeriesType =
94 #ifdef RTK_USE_CUDA
95  using DVFSequenceImageType =
96  typename std::conditional_t<std::is_same_v<VolumeSeriesType, CPUVolumeSeriesType>,
98  itk::CudaImage<VectorForDVF, VolumeSeriesType::ImageDimension>>;
99  using DVFImageType = typename std::conditional_t<std::is_same_v<VolumeSeriesType, CPUVolumeSeriesType>,
100  itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension - 1>,
101  itk::CudaImage<VectorForDVF, VolumeSeriesType::ImageDimension - 1>>;
102 #else
104  using DVFImageType = itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension - 1>;
105 #endif
106 
110  itkNewMacro(Self);
111 
114 
116  void
118  {
119  itkExceptionMacro(<< "ForwardProjection cannot be changed");
120  }
121  void
122  SetBackProjectionFilter(BackProjectionType itkNotUsed(_arg)) override
123  {
124  itkExceptionMacro(<< "BackProjection cannot be changed");
125  }
127 
129  void
130  SetDisplacementField(const DVFSequenceImageType * DisplacementField);
131  void
132  SetInverseDisplacementField(const DVFSequenceImageType * InverseDisplacementField);
133  typename DVFSequenceImageType::ConstPointer
134  GetDisplacementField();
135  typename DVFSequenceImageType::ConstPointer
136  GetInverseDisplacementField();
138 
140  void
141  SetSignal(const std::vector<double> signal) override;
142 
143  // Sub filters type alias
145  using MCCGOperatorType =
147 
149  itkSetMacro(UseCudaCyclicDeformation, bool);
150  itkGetMacro(UseCudaCyclicDeformation, bool);
152 
153 protected:
156 
157  void
158  GenerateOutputInformation() override;
159  void
160  GenerateInputRequestedRegion() override;
161 
163 
164 }; // end of class
165 
166 } // end namespace rtk
167 
168 #ifndef ITK_MANUAL_INSTANTIATION
169 # include "rtkMotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter.hxx"
170 #endif
171 
172 #endif
typename itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension > CPUVolumeSeriesType
#define itkSetMacro(name, type)
typename std::conditional_t< std::is_same_v< VolumeSeriesType, CPUVolumeSeriesType >, itk::Image< VectorForDVF, VolumeSeriesType::ImageDimension - 1 >, itk::CudaImage< VectorForDVF, VolumeSeriesType::ImageDimension - 1 > > DVFImageType
typename std::conditional_t< std::is_same_v< VolumeSeriesType, CPUVolumeSeriesType >, itk::Image< VectorForDVF, VolumeSeriesType::ImageDimension >, itk::CudaImage< VectorForDVF, VolumeSeriesType::ImageDimension > > DVFSequenceImageType
Back projection part for motion compensated iterative 4D reconstruction.
Implements part of the 4D reconstruction by conjugate gradient.