RTK  2.6.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  itkOverrideGetNameOfClassMacro(ProjectionStackToFourDImageFilter);
124 
126  void
127  SetInputVolumeSeries(const VolumeSeriesType * VolumeSeries);
128  typename VolumeSeriesType::ConstPointer
129  GetInputVolumeSeries();
131 
133  void
134  SetInputProjectionStack(const ProjectionStackType * Projection);
135  typename ProjectionStackType::ConstPointer
136  GetInputProjectionStack();
138 
144 
146 
148  using CPUVolumeSeriesType =
150 #ifdef RTK_USE_CUDA
151  typedef typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
153  CudaSplatImageFilter>::type CudaSplatImageFilterType;
154  typedef typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
156  CudaConstantVolumeSource>::type CudaConstantVolumeSourceType;
157  typedef typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
159  CudaConstantVolumeSeriesSource>::type CudaConstantVolumeSeriesSourceType;
160 #else
164 #endif
165 
167  void
168  SetBackProjectionFilter(BackProjectionFilterType * _arg);
169 
171  itkSetConstObjectMacro(Geometry, GeometryType);
172 
174  itkSetMacro(UseCudaSplat, bool);
175  itkGetMacro(UseCudaSplat, bool);
176  itkSetMacro(UseCudaSources, bool);
177  itkGetMacro(UseCudaSources, bool);
179 
181  itkGetMacro(Weights, itk::Array2D<float>);
184 
186  virtual void
187  SetSignal(const std::vector<double> signal);
188 
189 protected:
191  ~ProjectionStackToFourDImageFilter() override = default;
192 
194  void
195  VerifyPreconditions() const override;
196 
198  void
199  GenerateData() override;
200 
201  void
202  GenerateOutputInformation() override;
203 
204  void
205  GenerateInputRequestedRegion() override;
206 
207  void
208  InitializeConstantSource();
209 
216 
222  std::vector<double> m_Signal;
223 };
224 } // namespace rtk
225 
226 
227 #ifndef ITK_MANUAL_INSTANTIATION
228 # include "rtkProjectionStackToFourDImageFilter.hxx"
229 #endif
230 
231 #endif
Generate an n-dimensional image with constant pixel values.
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.
ConstantVolumeSeriesSourceType::Pointer m_ConstantVolumeSeriesSource