RTK  2.4.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  * 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  typedef typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
97  itk::CudaImage<VectorForDVF, VolumeSeriesType::ImageDimension>>::type
99  typedef
100  typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
101  itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension - 1>,
102  itk::CudaImage<VectorForDVF, VolumeSeriesType::ImageDimension - 1>>::type DVFImageType;
103 #else
105  using DVFImageType = itk::Image<VectorForDVF, VolumeSeriesType::ImageDimension - 1>;
106 #endif
107 
111  itkNewMacro(Self);
112 
116 
118  void
120  {
121  itkExceptionMacro(<< "ForwardProjection cannot be changed");
122  }
123  void
124  SetBackProjectionFilter(BackProjectionType itkNotUsed(_arg)) override
125  {
126  itkExceptionMacro(<< "BackProjection cannot be changed");
127  }
129 
131  void
132  SetDisplacementField(const DVFSequenceImageType * DisplacementField);
133  void
134  SetInverseDisplacementField(const DVFSequenceImageType * InverseDisplacementField);
135  typename DVFSequenceImageType::ConstPointer
136  GetDisplacementField();
137  typename DVFSequenceImageType::ConstPointer
138  GetInverseDisplacementField();
140 
142  void
143  SetSignal(const std::vector<double> signal) override;
144 
145  // Sub filters type alias
147  using MCCGOperatorType =
149 
151  itkSetMacro(UseCudaCyclicDeformation, bool);
152  itkGetMacro(UseCudaCyclicDeformation, bool);
154 
155 protected:
158 
159  void
160  GenerateOutputInformation() override;
161  void
162  GenerateInputRequestedRegion() override;
163 
165 
166 }; // end of class
167 
168 } // end namespace rtk
169 
170 #ifndef ITK_MANUAL_INSTANTIATION
171 # include "rtkMotionCompensatedFourDConjugateGradientConeBeamReconstructionFilter.hxx"
172 #endif
173 
174 #endif
typename itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension > CPUVolumeSeriesType
#define itkSetMacro(name, type)
TOutputImage OutputImageType
Back projection part for motion compensated iterative 4D reconstruction.
Implements part of the 4D reconstruction by conjugate gradient.