RTK  2.7.0
Reconstruction Toolkit
rtkFourDROOSTERConeBeamReconstructionFilter.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 rtkFourDROOSTERConeBeamReconstructionFilter_h
19 #define rtkFourDROOSTERConeBeamReconstructionFilter_h
20 
23 #ifdef RTK_USE_CUDA
26 #else
29 #endif
35 
37 #include <itkSubtractImageFilter.h>
38 #include <itkAddImageFilter.h>
39 
40 #include <itkResampleImageFilter.h>
42 #include <itkIdentityTransform.h>
43 
44 
45 namespace rtk
46 {
206 template <typename VolumeSeriesType, typename ProjectionStackType>
208  : public rtk::IterativeConeBeamReconstructionFilter<VolumeSeriesType, ProjectionStackType>
209 {
210 public:
211  ITK_DISALLOW_COPY_AND_MOVE(FourDROOSTERConeBeamReconstructionFilter);
212 
217  using VolumeType = ProjectionStackType;
219  itk::CovariantVector<typename VolumeSeriesType::ValueType, VolumeSeriesType::ImageDimension - 1>;
222  using CPUVolumeSeriesType =
224 
226  typename VolumeSeriesType::template RebindImageType<CovariantVectorForSpatialGradient,
227  VolumeSeriesType::ImageDimension>;
229  typename VolumeSeriesType::template RebindImageType<CovariantVectorForTemporalGradient,
230  VolumeSeriesType::ImageDimension>;
231 
232  using DVFSequenceImageType =
233  typename VolumeSeriesType::template RebindImageType<DVFVectorType, VolumeSeriesType::ImageDimension>;
234  using DVFImageType =
235  typename VolumeSeriesType::template RebindImageType<DVFVectorType, VolumeSeriesType::ImageDimension - 1>;
236 
237 #ifdef RTK_USE_CUDA
238  using AverageOutOfROIFilterType = std::conditional_t<std::is_same_v<VolumeSeriesType, CPUVolumeSeriesType>,
242  std::conditional_t<std::is_same_v<VolumeSeriesType, CPUVolumeSeriesType>,
245 #else
249 #endif
250 
252  itkNewMacro(Self);
253 
255  itkOverrideGetNameOfClassMacro(FourDROOSTERConeBeamReconstructionFilter);
256 
258  void
259  SetInputVolumeSeries(const VolumeSeriesType * VolumeSeries);
260  typename VolumeSeriesType::ConstPointer
261  GetInputVolumeSeries();
263 
265  void
266  SetInputProjectionStack(const ProjectionStackType * Projection);
267  typename ProjectionStackType::Pointer
268  GetInputProjectionStack();
270 
272  void
273  SetMotionMask(const VolumeType * mask);
274  typename VolumeType::Pointer
275  GetMotionMask();
277 
279  void
280  SetDisplacementField(const DVFSequenceImageType * DVFs);
281  void
282  SetInverseDisplacementField(const DVFSequenceImageType * DVFs);
283  typename DVFSequenceImageType::Pointer
284  GetDisplacementField();
285  typename DVFSequenceImageType::Pointer
286  GetInverseDisplacementField();
288 
289  using FourDCGFilterType =
295  using WarpSequenceFilterType =
302  using TNVDenoisingFilterType =
304 
305  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
306  using BackProjectionType = typename Superclass::BackProjectionType;
307 
309  virtual void
310  SetWeights(const itk::Array2D<float> _arg);
311 
313  itkSetMacro(DisableDisplacedDetectorFilter, bool);
314  itkGetMacro(DisableDisplacedDetectorFilter, bool);
316 
317  // Regularization steps to perform
318  itkSetMacro(PerformPositivity, bool);
319  itkGetMacro(PerformPositivity, bool);
320  itkSetMacro(PerformMotionMask, bool);
321  itkGetMacro(PerformMotionMask, bool);
322  itkSetMacro(PerformTVSpatialDenoising, bool);
323  itkGetMacro(PerformTVSpatialDenoising, bool);
324  itkSetMacro(PerformWaveletsSpatialDenoising, bool);
325  itkGetMacro(PerformWaveletsSpatialDenoising, bool);
326  itkSetMacro(PerformWarping, bool);
327  itkGetMacro(PerformWarping, bool);
328  itkSetMacro(PerformTVTemporalDenoising, bool);
329  itkGetMacro(PerformTVTemporalDenoising, bool);
330  itkSetMacro(PerformL0TemporalDenoising, bool);
331  itkGetMacro(PerformL0TemporalDenoising, bool);
332  itkSetMacro(PerformTNVDenoising, bool);
333  itkGetMacro(PerformTNVDenoising, bool);
334  itkSetMacro(ComputeInverseWarpingByConjugateGradient, bool);
335  itkGetMacro(ComputeInverseWarpingByConjugateGradient, bool);
336  itkSetMacro(UseNearestNeighborInterpolationInWarping, bool);
337  itkGetMacro(UseNearestNeighborInterpolationInWarping, bool);
338  itkGetMacro(CudaConjugateGradient, bool);
339  itkSetMacro(CudaConjugateGradient, bool);
340 
342  itkSetMacro(UseCudaCyclicDeformation, bool);
343  itkGetMacro(UseCudaCyclicDeformation, bool);
345 
346  // Regularization parameters
347  itkSetMacro(GammaTVSpace, float);
348  itkGetMacro(GammaTVSpace, float);
349  itkSetMacro(GammaTVTime, float);
350  itkGetMacro(GammaTVTime, float);
351  itkSetMacro(GammaTNV, float);
352  itkGetMacro(GammaTNV, float);
353  itkSetMacro(LambdaL0Time, float);
354  itkGetMacro(LambdaL0Time, float);
355  itkSetMacro(SoftThresholdWavelets, float);
356  itkGetMacro(SoftThresholdWavelets, float);
357  itkSetMacro(PhaseShift, float);
358  itkGetMacro(PhaseShift, float);
359 
361  itkGetMacro(NumberOfLevels, unsigned int);
362  itkSetMacro(NumberOfLevels, unsigned int);
364 
366  itkGetMacro(Order, unsigned int);
367  itkSetMacro(Order, unsigned int);
369 
370  // Iterations
371  itkSetMacro(MainLoop_iterations, int);
372  itkGetMacro(MainLoop_iterations, int);
373  itkSetMacro(CG_iterations, int);
374  itkGetMacro(CG_iterations, int);
375  itkSetMacro(TV_iterations, int);
376  itkGetMacro(TV_iterations, int);
377  itkSetMacro(L0_iterations, int);
378  itkGetMacro(L0_iterations, int);
379 
380  // Geometry
381  itkSetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
382  itkGetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
383 
385  virtual void
386  SetSignal(const std::vector<double> signal);
387 
388 protected:
390  ~FourDROOSTERConeBeamReconstructionFilter() override = default;
391 
393  void
394  VerifyPreconditions() const override;
395 
397  void
398  GenerateData() override;
399 
400  void
401  GenerateOutputInformation() override;
402 
403  void
404  GenerateInputRequestedRegion() override;
405 
406  // Inputs are not supposed to occupy the same physical space,
407  // so there is nothing to verify
408  void
409  VerifyInputInformation() const override
410  {}
411 
416  typename AverageOutOfROIFilterType::Pointer m_AverageOutOfROIFilter;
420  typename TemporalTVDenoisingFilterType::Pointer m_TVDenoisingTime;
427 
428  // Booleans :
429  // should warping be performed ?
430  // should conjugate gradient be performed on GPU ?
431  // should wavelets replace TV in spatial denoising ?
441  bool m_UseNearestNeighborInterpolationInWarping; // Default is false, linear interpolation is used instead
445 
446  // Regularization parameters
449  float m_GammaTNV;
453  bool m_DimensionsProcessedForTVSpace[VolumeSeriesType::ImageDimension];
454  bool m_DimensionsProcessedForTVTime[VolumeSeriesType::ImageDimension];
455 
457 
459  unsigned int m_Order;
460  unsigned int m_NumberOfLevels;
461 
462  // Iterations
467 
468  // Geometry
470 
471  // Signal
472  std::vector<double> m_Signal;
473 };
474 } // namespace rtk
475 
476 
477 #ifndef ITK_MANUAL_INSTANTIATION
478 # include "rtkFourDROOSTERConeBeamReconstructionFilter.hxx"
479 #endif
480 
481 #endif
typename VolumeSeriesType::template RebindImageType< CovariantVectorForTemporalGradient, VolumeSeriesType::ImageDimension > TemporalGradientImageType
rtk::ThreeDCircularProjectionGeometry::ConstPointer m_Geometry
std::conditional_t< std::is_same_v< VolumeSeriesType, CPUVolumeSeriesType >, TotalVariationDenoisingBPDQImageFilter< VolumeSeriesType, TemporalGradientImageType >, CudaLastDimensionTVDenoisingImageFilter > TemporalTVDenoisingFilterType
Denoises along the last dimension, reducing the L0 norm of the gradient.
typename itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension > CPUVolumeSeriesType
Applies an N-D + time Motion Vector Field to an N-D + time sequence of images.
Finds the image sequence that, once warped, equals the input image sequence.
Applies 3D total variation denoising to a 3D + time sequence of images.
Applies a total variation denoising, only alm_SingularValueThresholdFilterong the dimensions specifie...
Implements 4D RecOnstructiOn using Spatial and TEmporal Regularization (short 4D ROOSTER) ...
Projection geometry for a source and a 2-D flat panel.
typename VolumeSeriesType::template RebindImageType< DVFVectorType, VolumeSeriesType::ImageDimension - 1 > DVFImageType
#define itkSetMacro(name, type)
Implements the AverageOutOfROIImageFilter on GPU.
std::conditional_t< std::is_same_v< VolumeSeriesType, CPUVolumeSeriesType >, AverageOutOfROIImageFilter< VolumeSeriesType, VolumeType >, CudaAverageOutOfROIImageFilter > AverageOutOfROIFilterType
Applies 3D Daubechies wavelets denoising to a 3D + time sequence of images.
SpatialWaveletsDenoisingFilterType::Pointer m_WaveletsDenoisingSpace
typename VolumeSeriesType::template RebindImageType< CovariantVectorForSpatialGradient, VolumeSeriesType::ImageDimension > SpatialGradientImageType
Mother class for cone beam reconstruction filters which need runtime selection of their forward and b...
itk::ImageToImageFilter< VolumeSeriesType, VolumeSeriesType >::Pointer m_DownstreamFilter
typename VolumeSeriesType::template RebindImageType< DVFVectorType, VolumeSeriesType::ImageDimension > DVFSequenceImageType
Implements the TotalVariationDenoisingBPDQImageFilter on GPU for a specific case : denoising only alo...
Implements part of the 4D reconstruction by conjugate gradient.
Averages along the last dimension if the pixel is outside ROI.