RTK  2.7.0
Reconstruction Toolkit
rtkGgoFunctions.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 RTKGGOFUNCTIONS_H
20 // #define RTKGGOFUNCTIONS_H
21 
22 #ifndef rtkGgoFunctions_h
23 #define rtkGgoFunctions_h
24 
25 #include "rtkMacro.h"
26 #include "rtkConstantImageSource.h"
27 #include "rtkIOFactories.h"
28 #include "rtkProjectionsReader.h"
30 #include <itksys/RegularExpression.hxx>
31 
32 namespace rtk
33 {
34 
48 template <class TConstantImageSourceType, class TArgsInfo>
49 void
50 SetConstantImageSourceFromGgo(TConstantImageSourceType * source, const TArgsInfo & args_info)
51 {
52  using ImageType = typename TConstantImageSourceType::OutputImageType;
53 
54  const unsigned int Dimension = ImageType::GetImageDimension();
55 
56  typename ImageType::SizeType imageSize;
57  if (args_info.dimension_given)
58  {
59  itkGenericOutputMacro(
60  "Warning: --dimension is deprecated and will be removed in a future release. Use --size instead.");
61  imageSize.Fill(args_info.dimension_arg[0]);
62  for (unsigned int i = 0; i < std::min(args_info.dimension_given, Dimension); i++)
63  {
64  imageSize[i] = args_info.dimension_arg[i];
65  }
66  }
67  imageSize.Fill(args_info.size_arg[0]);
68  for (unsigned int i = 0; i < std::min(args_info.size_given, Dimension); i++)
69  imageSize[i] = args_info.size_arg[i];
70 
71  typename ImageType::SpacingType imageSpacing;
72  imageSpacing.Fill(args_info.spacing_arg[0]);
73  for (unsigned int i = 0; i < std::min(args_info.spacing_given, Dimension); i++)
74  imageSpacing[i] = args_info.spacing_arg[i];
75 
76  typename ImageType::PointType imageOrigin;
77  for (unsigned int i = 0; i < Dimension; i++)
78  imageOrigin[i] = imageSpacing[i] * (imageSize[i] - 1) * -0.5;
79  for (unsigned int i = 0; i < std::min(args_info.origin_given, Dimension); i++)
80  imageOrigin[i] = args_info.origin_arg[i];
81 
82  typename ImageType::DirectionType imageDirection;
83  if (args_info.direction_given)
84  for (unsigned int i = 0; i < Dimension; i++)
85  for (unsigned int j = 0; j < Dimension; j++)
86  imageDirection[i][j] = args_info.direction_arg[i * Dimension + j];
87  else
88  imageDirection.SetIdentity();
89 
90  source->SetOrigin(imageOrigin);
91  source->SetSpacing(imageSpacing);
92  source->SetDirection(imageDirection);
93  source->SetSize(imageSize);
94  source->SetConstant({});
95 
96  // Copy output image information from an existing file, if requested
97  // Overwrites parameters given in command line, if any
98  if (args_info.like_given)
99  {
100  using LikeReaderType = itk::ImageFileReader<ImageType>;
101  typename LikeReaderType::Pointer likeReader = LikeReaderType::New();
102  likeReader->SetFileName(args_info.like_arg);
103  TRY_AND_EXIT_ON_ITK_EXCEPTION(likeReader->UpdateOutputInformation());
104  source->SetInformationFromImage(likeReader->GetOutput());
105  }
106 
107  TRY_AND_EXIT_ON_ITK_EXCEPTION(source->UpdateOutputInformation());
108 }
109 
125 template <class TArgsInfo>
126 const std::vector<std::string>
127 GetProjectionsFileNamesFromGgo(const TArgsInfo & args_info)
128 {
129  // Generate file names
131  names->SetDirectory(args_info.path_arg);
132  names->SetNumericSort(args_info.nsort_flag);
133  names->SetRegularExpression(args_info.regexp_arg);
134  names->SetSubMatch(args_info.submatch_arg);
135 
136  if (args_info.verbose_flag)
137  std::cout << "Regular expression matches " << names->GetFileNames().size() << " file(s)..." << std::endl;
138 
139  // Check submatch in file names
140  if (args_info.submatch_given)
141  {
142  // Check that the submatch number returns something relevant
143  itksys::RegularExpression reg;
144  if (!reg.compile(args_info.regexp_arg))
145  {
146  itkGenericExceptionMacro(<< "Error compiling regular expression " << args_info.regexp_arg);
147  }
148 
149  // Store the full filename and the selected sub expression match
150  for (const std::string & name : names->GetFileNames())
151  {
152  reg.find(name);
153  if (reg.match(args_info.submatch_arg) == std::string(""))
154  {
155  itkGenericExceptionMacro(<< "Cannot find submatch " << args_info.submatch_arg << " in " << name
156  << " from regular expression " << args_info.regexp_arg);
157  }
158  }
159  }
160 
161  std::vector<std::string> fileNames = names->GetFileNames();
163  std::vector<size_t> idxtopop;
164  size_t i = 0;
165  for (const auto & fn : fileNames)
166  {
167  itk::ImageIOBase::Pointer imageio =
169 
170  if (imageio.IsNull())
171  {
172  std::cerr << "Ignoring file: " << fn << "\n";
173  idxtopop.push_back(i);
174  }
175  i++;
176  }
177  std::reverse(idxtopop.begin(), idxtopop.end());
178  for (const auto & id : idxtopop)
179  {
180  fileNames.erase(fileNames.begin() + id);
181  }
182 
183  return fileNames;
184 }
185 
186 template <class TProjectionsReaderType, class TArgsInfo>
187 void
188 SetProjectionsReaderFromGgo(TProjectionsReaderType * reader, const TArgsInfo & args_info)
189 {
190  const std::vector<std::string> fileNames = GetProjectionsFileNamesFromGgo(args_info);
191 
192  // Vector component extraction
193  if (args_info.component_given)
194  {
195  reader->SetVectorComponent(args_info.component_arg);
196  }
197 
198  // Change image information
199  const unsigned int Dimension = TProjectionsReaderType::OutputImageType::GetImageDimension();
200  typename TProjectionsReaderType::OutputImageDirectionType direction;
201  if (args_info.newdirection_given)
202  {
203  direction.Fill(args_info.newdirection_arg[0]);
204  for (unsigned int i = 0; i < args_info.newdirection_given; i++)
205  direction[i / Dimension][i % Dimension] = args_info.newdirection_arg[i];
206  reader->SetDirection(direction);
207  }
208  typename TProjectionsReaderType::OutputImageSpacingType spacing;
209  if (args_info.newspacing_given)
210  {
211  spacing.Fill(args_info.newspacing_arg[0]);
212  for (unsigned int i = 0; i < args_info.newspacing_given; i++)
213  spacing[i] = args_info.newspacing_arg[i];
214  reader->SetSpacing(spacing);
215  }
216  typename TProjectionsReaderType::OutputImagePointType origin;
217  if (args_info.neworigin_given)
218  {
219  origin.Fill(args_info.neworigin_arg[0]);
220  for (unsigned int i = 0; i < args_info.neworigin_given; i++)
221  origin[i] = args_info.neworigin_arg[i];
222  reader->SetOrigin(origin);
223  }
224 
225  // Crop boundaries
226  typename TProjectionsReaderType::OutputImageSizeType upperCrop, lowerCrop;
227  upperCrop.Fill(0);
228  lowerCrop.Fill(0);
229  for (unsigned int i = 0; i < args_info.lowercrop_given; i++)
230  lowerCrop[i] = args_info.lowercrop_arg[i];
231  reader->SetLowerBoundaryCropSize(lowerCrop);
232  for (unsigned int i = 0; i < args_info.uppercrop_given; i++)
233  upperCrop[i] = args_info.uppercrop_arg[i];
234  reader->SetUpperBoundaryCropSize(upperCrop);
235 
236  // Conditional median
237  typename TProjectionsReaderType::MedianRadiusType medianRadius;
238  medianRadius.Fill(0);
239  for (unsigned int i = 0; i < args_info.radius_given; i++)
240  medianRadius[i] = args_info.radius_arg[i];
241  reader->SetMedianRadius(medianRadius);
242  if (args_info.multiplier_given)
243  reader->SetConditionalMedianThresholdMultiplier(args_info.multiplier_arg);
244 
245  // Shrink / Binning
246  typename TProjectionsReaderType::ShrinkFactorsType binFactors;
247  binFactors.Fill(1);
248  for (unsigned int i = 0; i < args_info.binning_given; i++)
249  binFactors[i] = args_info.binning_arg[i];
250  reader->SetShrinkFactors(binFactors);
251 
252  // Boellaard scatter correction
253  if (args_info.spr_given)
254  reader->SetScatterToPrimaryRatio(args_info.spr_arg);
255  if (args_info.nonneg_given)
256  reader->SetNonNegativityConstraintThreshold(args_info.nonneg_arg);
257  if (args_info.airthres_given)
258  reader->SetAirThreshold(args_info.airthres_arg);
259 
260  // I0 and IDark
261  if (args_info.i0_given)
262  reader->SetI0(args_info.i0_arg);
263  reader->SetIDark(args_info.idark_arg);
264 
265  // Line integral flag
266  if (args_info.nolineint_flag)
267  reader->ComputeLineIntegralOff();
268 
269  // Water precorrection
270  if (args_info.wpc_given)
271  {
272  std::vector<double> coeffs;
273  coeffs.assign(args_info.wpc_arg, args_info.wpc_arg + args_info.wpc_given);
274  reader->SetWaterPrecorrectionCoefficients(coeffs);
275  }
276 
277  // Pass list to projections reader
278  reader->SetFileNames(fileNames);
279  TRY_AND_EXIT_ON_ITK_EXCEPTION(reader->UpdateOutputInformation());
280 }
281 
290 template <class TArgsInfo, class TIterativeReconstructionFilter>
291 void
292 SetBackProjectionFromGgo(const TArgsInfo & args_info, TIterativeReconstructionFilter * recon)
293 {
294  typename TIterativeReconstructionFilter::VolumeType::Pointer attenuationMap;
295  using VolumeType = typename TIterativeReconstructionFilter::VolumeType;
296  if (args_info.attenuationmap_given)
297  {
298  // Read an existing image to initialize the attenuation map
299  attenuationMap = itk::ReadImage<VolumeType>(args_info.attenuationmap_arg);
300  }
301 
302  switch (args_info.bp_arg)
303  {
304  case (-1): // bp__NULL, keep default
305  break;
306  case (0): // bp_arg_VoxelBasedBackProjection
307  recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_VOXELBASED);
308  break;
309  case (1): // bp_arg_Joseph
310  recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_JOSEPH);
311  break;
312  case (2): // bp_arg_CudaVoxelBased
313  recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_CUDAVOXELBASED);
314  break;
315  case (3): // bp_arg_CudaRayCast
316  recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_CUDARAYCAST);
317  if (args_info.step_given)
318  recon->SetStepSize(args_info.step_arg);
319  break;
320  case (4): // bp_arg_JosephAttenuated
321  recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_JOSEPHATTENUATED);
322  if (args_info.attenuationmap_given)
323  recon->SetAttenuationMap(attenuationMap);
324  break;
325  case (5): // bp_arg_RotationBased
326  recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_ZENG);
327  if (args_info.sigmazero_given)
328  recon->SetSigmaZero(args_info.sigmazero_arg);
329  if (args_info.alphapsf_given)
330  recon->SetAlphaPSF(args_info.alphapsf_arg);
331  if (args_info.attenuationmap_given)
332  recon->SetAttenuationMap(attenuationMap);
333  break;
334  }
335 }
336 
345 template <class TArgsInfo, class TIterativeReconstructionFilter>
346 void
347 SetForwardProjectionFromGgo(const TArgsInfo & args_info, TIterativeReconstructionFilter * recon)
348 {
349  typename TIterativeReconstructionFilter::VolumeType::Pointer attenuationMap;
350  using VolumeType = typename TIterativeReconstructionFilter::VolumeType;
351  if (args_info.attenuationmap_given)
352  {
353  // Read an existing image to initialize the attenuation map
354  attenuationMap = itk::ReadImage<VolumeType>(args_info.attenuationmap_arg);
355  }
357  typename ClipImageType::Pointer inferiorClipImage, superiorClipImage;
358  if (args_info.inferiorclipimage_given)
359  {
360  // Read an existing image to initialize the attenuation map
361  inferiorClipImage = itk::ReadImage<ClipImageType>(args_info.inferiorclipimage_arg);
362  }
363  if (args_info.superiorclipimage_given)
364  {
365  // Read an existing image to initialize the attenuation map
366  superiorClipImage = itk::ReadImage<ClipImageType>(args_info.superiorclipimage_arg);
367  }
368 
369  switch (args_info.fp_arg)
370  {
371  case (-1): // fp__NULL, keep default
372  break;
373  case (0): // fp_arg_Joseph
374  recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_JOSEPH);
375  if (args_info.superiorclipimage_given)
376  recon->SetSuperiorClipImage(superiorClipImage);
377  if (args_info.inferiorclipimage_given)
378  recon->SetInferiorClipImage(inferiorClipImage);
379  break;
380  case (1): // fp_arg_CudaRayCast
381  recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_CUDARAYCAST);
382  if (args_info.step_given)
383  recon->SetStepSize(args_info.step_arg);
384  break;
385  case (2): // fp_arg_JosephAttenuated
386  recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_JOSEPHATTENUATED);
387  if (args_info.superiorclipimage_given)
388  recon->SetSuperiorClipImage(superiorClipImage);
389  if (args_info.inferiorclipimage_given)
390  recon->SetInferiorClipImage(inferiorClipImage);
391  if (args_info.attenuationmap_given)
392  recon->SetAttenuationMap(attenuationMap);
393  break;
394  case (3): // fp_arg_RotationBased
395  recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_ZENG);
396  if (args_info.sigmazero_given)
397  recon->SetSigmaZero(args_info.sigmazero_arg);
398  if (args_info.alphapsf_given)
399  recon->SetAlphaPSF(args_info.alphapsf_arg);
400  if (args_info.attenuationmap_given)
401  recon->SetAttenuationMap(attenuationMap);
402  break;
403  }
404 }
405 
406 } // namespace rtk
407 
408 #endif // rtkGgoFunctions_h
static ImageIOBasePointer CreateImageIO(const char *path, IOFileModeEnum mode)
const std::vector< std::string > GetProjectionsFileNamesFromGgo(const TArgsInfo &args_info)
Read a stack of 2D projections from gengetopt specifications.
void SetBackProjectionFromGgo(const TArgsInfo &args_info, TIterativeReconstructionFilter *recon)
Set the correct RTK backprojection type from gengetopt Given a gengetopt bp_arg option from the rtkpr...
void SetForwardProjectionFromGgo(const TArgsInfo &args_info, TIterativeReconstructionFilter *recon)
Set the correct RTK forward projection type from gengetopt Given a gengetopt fp_arg option from the r...
void SetConstantImageSourceFromGgo(TConstantImageSourceType *source, const TArgsInfo &args_info)
Create 3D image from gengetopt specifications.
void SetProjectionsReaderFromGgo(TProjectionsReaderType *reader, const TArgsInfo &args_info)
void RTK_EXPORT RegisterIOFactories()
#define TRY_AND_EXIT_ON_ITK_EXCEPTION(execFunc)
Update a filter and catching/displaying exceptions.
Definition: rtkMacro.h:106