RTK  2.7.0
Reconstruction Toolkit
rtkMotionCompensatedFourDReconstructionConjugateGradientOperator.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 #ifndef rtkMotionCompensatedFourDReconstructionConjugateGradientOperator_h
19 #define rtkMotionCompensatedFourDReconstructionConjugateGradientOperator_h
20 
24 
25 #ifdef RTK_USE_CUDA
29 #endif
30 
31 namespace rtk
32 {
99 template <typename VolumeSeriesType, typename ProjectionStackType>
101  : public FourDReconstructionConjugateGradientOperator<VolumeSeriesType, ProjectionStackType>
102 {
103 public:
105 
110 
112  using VolumeType = ProjectionStackType;
113  using VectorForDVF = itk::CovariantVector<typename VolumeSeriesType::ValueType, VolumeSeriesType::ImageDimension - 1>;
114 
116  using CPUVolumeSeriesType =
118  using DVFSequenceImageType =
119  typename VolumeSeriesType::template RebindImageType<VectorForDVF, VolumeSeriesType::ImageDimension>;
120  using DVFImageType =
121  typename VolumeSeriesType::template RebindImageType<VectorForDVF, VolumeSeriesType::ImageDimension - 1>;
122 #ifdef RTK_USE_CUDA
124  typename std::conditional_t<std::is_same_v<VolumeSeriesType, CPUVolumeSeriesType>,
128  typename std::conditional_t<std::is_same_v<VolumeSeriesType, CPUVolumeSeriesType>,
131 #else
135 #endif
137 #ifdef RTK_USE_CUDA
139  typename std::conditional_t<std::is_same_v<VolumeSeriesType, CPUVolumeSeriesType>,
142 #else
144 #endif
145 
147  itkNewMacro(Self);
148 
151 
153  void
154  SetForwardProjectionFilter(typename Superclass::ForwardProjectionFilterType * itkNotUsed(_arg))
155  {
156  itkExceptionMacro(<< "ForwardProjection cannot be changed");
157  }
158  void
159  SetBackProjectionFilter(typename Superclass::BackProjectionFilterType * itkNotUsed(_arg))
160  {
161  itkExceptionMacro(<< "BackProjection cannot be changed");
162  }
164 
166  void
167  SetDisplacementField(const DVFSequenceImageType * DisplacementField);
168  void
169  SetInverseDisplacementField(const DVFSequenceImageType * InverseDisplacementField);
170  typename DVFSequenceImageType::ConstPointer
171  GetInverseDisplacementField();
172  typename DVFSequenceImageType::ConstPointer
173  GetDisplacementField();
175 
177  void
178  SetSignal(const std::vector<double> signal) override;
179 
181  itkSetMacro(UseCudaCyclicDeformation, bool);
182  itkGetMacro(UseCudaCyclicDeformation, bool);
184 
185 protected:
188 
190  void
191  GenerateOutputInformation() override;
192 
194  void
195  VerifyInputInformation() const override
196  {}
197 
199  void
200  GenerateData() override;
201 
205  std::vector<double> m_Signal;
207 };
208 } // namespace rtk
209 
210 
211 #ifndef ITK_MANUAL_INSTANTIATION
212 # include "rtkMotionCompensatedFourDReconstructionConjugateGradientOperator.hxx"
213 #endif
214 
215 #endif
Return 3D deformation vector field according to input 4D vector field, phase signal and frame number...
typename VolumeSeriesType::template RebindImageType< VectorForDVF, VolumeSeriesType::ImageDimension - 1 > DVFImageType
typename itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension > CPUVolumeSeriesType
Voxel-based backprojection into warped volume implemented in CUDA.
typename VolumeSeriesType::template RebindImageType< VectorForDVF, VolumeSeriesType::ImageDimension > DVFSequenceImageType
#define itkSetMacro(name, type)
Implements part of the 4D reconstruction by conjugate gradient.
GPU version of the temporal DVF interpolator.
typename std::conditional_t< std::is_same_v< VolumeSeriesType, CPUVolumeSeriesType >, JosephForwardProjectionImageFilter< ProjectionStackType, ProjectionStackType >, CudaWarpForwardProjectionImageFilter > WarpForwardProjectionImageFilterType
typename std::conditional_t< std::is_same_v< VolumeSeriesType, CPUVolumeSeriesType >, BackProjectionImageFilter< VolumeType, VolumeType >, CudaWarpBackProjectionImageFilter > WarpBackProjectionImageFilterType
typename std::conditional_t< std::is_same_v< VolumeSeriesType, CPUVolumeSeriesType >, CPUDVFInterpolatorType, CudaCyclicDeformationImageFilter > CudaCyclicDeformationImageFilterType