RTK  2.5.0
Reconstruction Toolkit
rtkUnwarpSequenceImageFilter.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 rtkUnwarpSequenceImageFilter_h
20 #define rtkUnwarpSequenceImageFilter_h
21 
25 #include "rtkConstantImageSource.h"
26 
27 #ifdef RTK_USE_CUDA
30 #endif
31 
32 namespace rtk
33 {
73 template <typename TImageSequence,
74  typename TDVFImageSequence =
75  itk::Image<itk::CovariantVector<typename TImageSequence::ValueType, TImageSequence::ImageDimension - 1>,
76  TImageSequence::ImageDimension>,
77  typename TImage = itk::Image<typename TImageSequence::ValueType, TImageSequence::ImageDimension - 1>,
78  typename TDVFImage =
79  itk::Image<itk::CovariantVector<typename TImageSequence::ValueType, TImageSequence::ImageDimension - 1>,
80  TImageSequence::ImageDimension - 1>>
81 class ITK_TEMPLATE_EXPORT UnwarpSequenceImageFilter : public itk::ImageToImageFilter<TImageSequence, TImageSequence>
82 {
83 public:
84  ITK_DISALLOW_COPY_AND_MOVE(UnwarpSequenceImageFilter);
85 
90 
92  itkNewMacro(Self);
93 
95 #ifdef itkOverrideGetNameOfClassMacro
96  itkOverrideGetNameOfClassMacro(UnwarpSequenceImageFilter);
97 #else
99 #endif
100 
101 
102  using CGOperatorFilterType =
106 
109 #ifdef RTK_USE_CUDA
110  typedef typename std::conditional<std::is_same<TImageSequence, CPUImageSequence>::value,
112  CudaConstantVolumeSeriesSource>::type ConstantSourceType;
113  typedef typename std::conditional<std::is_same<TImageSequence, CPUImageSequence>::value,
115  CudaConjugateGradientImageFilter<TImageSequence>>::type CudaConjugateGradientType;
116 #else
119 #endif
120 
122  void
123  SetDisplacementField(const TDVFImageSequence * DVFs);
124 
126  typename TDVFImageSequence::Pointer
127  GetDisplacementField();
128 
130  itkSetMacro(NumberOfIterations, float);
131  itkGetMacro(NumberOfIterations, float);
133 
135  itkSetMacro(PhaseShift, float);
136  itkGetMacro(PhaseShift, float);
138 
139  itkSetMacro(UseNearestNeighborInterpolationInWarping, bool);
140  itkGetMacro(UseNearestNeighborInterpolationInWarping, bool);
141 
142  itkSetMacro(CudaConjugateGradient, bool);
143  itkGetMacro(CudaConjugateGradient, bool);
144 
146  itkSetMacro(UseCudaCyclicDeformation, bool);
147  itkGetMacro(UseCudaCyclicDeformation, bool);
149 
150 protected:
152  ~UnwarpSequenceImageFilter() override = default;
153 
155  void
156  GenerateData() override;
157 
163 
166 
170  void
171  VerifyInputInformation() const override
172  {}
173 
176  void
177  GenerateInputRequestedRegion() override;
178  void
179  GenerateOutputInformation() override;
181 
182  bool m_UseNearestNeighborInterpolationInWarping; // Default is false, linear interpolation is used instead
185 
186 private:
187  unsigned int m_NumberOfIterations;
188 };
189 } // namespace rtk
190 
191 
192 #ifndef ITK_MANUAL_INSTANTIATION
193 # include "rtkUnwarpSequenceImageFilter.hxx"
194 #endif
195 
196 #endif
Implements the operator A used in the conjugate gradient unwarp sequence filter.
CGOperatorFilterType::Pointer m_CGOperator
Generate an n-dimensional image with constant pixel values.
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.
#define itkSetMacro(name, type)
ConstantSourceType::Pointer m_ConstantSource
typename itk::Image< typename TImageSequence::PixelType, TImageSequence::ImageDimension > CPUImageSequence
ConjugateGradientFilterType::Pointer m_ConjugateGradientFilter
Solves AX = B by conjugate gradient.
WarpForwardFilterType::Pointer m_WarpForwardFilter