RTK  2.5.0
Reconstruction Toolkit
rtkProjectionStackToFourDImageFilter.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 rtkProjectionStackToFourDImageFilter_h
19 #define rtkProjectionStackToFourDImageFilter_h
20 
21 #include <itkExtractImageFilter.h>
22 #include <itkArray2D.h>
23 
26 #include "rtkConstantImageSource.h"
28 
29 #ifdef RTK_USE_CUDA
30 # include "rtkCudaSplatImageFilter.h"
33 #endif
34 
35 namespace rtk
36 {
104 template <typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
105 class ITK_TEMPLATE_EXPORT ProjectionStackToFourDImageFilter
106  : public itk::ImageToImageFilter<VolumeSeriesType, VolumeSeriesType>
107 {
108 public:
109  ITK_DISALLOW_COPY_AND_MOVE(ProjectionStackToFourDImageFilter);
110 
115 
117  using VolumeType = ProjectionStackType;
118 
120  itkNewMacro(Self);
121 
123 #ifdef itkOverrideGetNameOfClassMacro
124  itkOverrideGetNameOfClassMacro(ProjectionStackToFourDImageFilter);
125 #else
127 #endif
128 
129 
131  void
132  SetInputVolumeSeries(const VolumeSeriesType * VolumeSeries);
133  typename VolumeSeriesType::ConstPointer
134  GetInputVolumeSeries();
136 
138  void
139  SetInputProjectionStack(const ProjectionStackType * Projection);
140  typename ProjectionStackType::ConstPointer
141  GetInputProjectionStack();
143 
149 
151 
153  using CPUVolumeSeriesType =
155 #ifdef RTK_USE_CUDA
156  typedef typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
158  CudaSplatImageFilter>::type CudaSplatImageFilterType;
159  typedef typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
161  CudaConstantVolumeSource>::type CudaConstantVolumeSourceType;
162  typedef typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
164  CudaConstantVolumeSeriesSource>::type CudaConstantVolumeSeriesSourceType;
165 #else
169 #endif
170 
172  void
173  SetBackProjectionFilter(const typename BackProjectionFilterType::Pointer _arg);
174 
176  itkSetConstObjectMacro(Geometry, GeometryType);
177 
179  itkSetMacro(UseCudaSplat, bool);
180  itkGetMacro(UseCudaSplat, bool);
181  itkSetMacro(UseCudaSources, bool);
182  itkGetMacro(UseCudaSources, bool);
184 
186  itkGetMacro(Weights, itk::Array2D<float>);
189 
191  virtual void
192  SetSignal(const std::vector<double> signal);
193 
194 protected:
196  ~ProjectionStackToFourDImageFilter() override = default;
197 
199  void
200  VerifyPreconditions() ITKv5_CONST override;
201 
203  void
204  GenerateData() override;
205 
206  void
207  GenerateOutputInformation() override;
208 
209  void
210  GenerateInputRequestedRegion() override;
211 
212  void
213  InitializeConstantSource();
214 
216  typename SplatFilterType::Pointer m_SplatFilter;
217  typename BackProjectionFilterType::Pointer m_BackProjectionFilter;
218  typename ExtractFilterType::Pointer m_ExtractFilter;
219  typename ConstantVolumeSourceType::Pointer m_ConstantVolumeSource;
220  typename ConstantVolumeSeriesSourceType::Pointer m_ConstantVolumeSeriesSource;
221 
223  itk::Array2D<float> m_Weights;
225  bool m_UseCudaSplat;
226  bool m_UseCudaSources;
227  std::vector<double> m_Signal;
228 };
229 } // namespace rtk
230 
231 
232 #ifndef ITK_MANUAL_INSTANTIATION
233 # include "rtkProjectionStackToFourDImageFilter.hxx"
234 #endif
235 
236 #endif
Generate an n-dimensional image with constant pixel values.
STL namespace.
Splats (linearly) a 3D volume into a 3D+t sequence of volumes.
typename itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension > CPUVolumeSeriesType
Projection geometry for a source and a 2-D flat panel.
#define itkSetMacro(name, type)
Implements part of the 4D reconstruction by conjugate gradient.