22 #ifndef rtkGgoFunctions_h 23 #define rtkGgoFunctions_h 30 #include <itksys/RegularExpression.hxx> 48 template <
class TConstantImageSourceType,
class TArgsInfo>
52 using ImageType =
typename TConstantImageSourceType::OutputImageType;
54 const unsigned int Dimension = ImageType::GetImageDimension();
56 typename ImageType::SizeType imageSize;
57 if (args_info.dimension_given)
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++)
64 imageSize[i] = args_info.dimension_arg[i];
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];
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];
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];
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];
88 imageDirection.SetIdentity();
90 source->SetOrigin(imageOrigin);
91 source->SetSpacing(imageSpacing);
92 source->SetDirection(imageDirection);
93 source->SetSize(imageSize);
94 source->SetConstant({});
98 if (args_info.like_given)
101 typename LikeReaderType::Pointer likeReader = LikeReaderType::New();
102 likeReader->SetFileName(args_info.like_arg);
104 source->SetInformationFromImage(likeReader->GetOutput());
125 template <
class TArgsInfo>
126 const std::vector<std::string>
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);
136 if (args_info.verbose_flag)
137 std::cout <<
"Regular expression matches " << names->GetFileNames().size() <<
" file(s)..." << std::endl;
140 if (args_info.submatch_given)
143 itksys::RegularExpression reg;
144 if (!reg.compile(args_info.regexp_arg))
146 itkGenericExceptionMacro(<<
"Error compiling regular expression " << args_info.regexp_arg);
150 for (
const std::string & name : names->GetFileNames())
153 if (reg.match(args_info.submatch_arg) == std::string(
""))
155 itkGenericExceptionMacro(<<
"Cannot find submatch " << args_info.submatch_arg <<
" in " << name
156 <<
" from regular expression " << args_info.regexp_arg);
161 std::vector<std::string> fileNames = names->GetFileNames();
163 std::vector<size_t> idxtopop;
165 for (
const auto & fn : fileNames)
170 if (imageio.IsNull())
172 std::cerr <<
"Ignoring file: " << fn <<
"\n";
173 idxtopop.push_back(i);
177 std::reverse(idxtopop.begin(), idxtopop.end());
178 for (
const auto &
id : idxtopop)
180 fileNames.erase(fileNames.begin() + id);
186 template <
class TProjectionsReaderType,
class TArgsInfo>
193 if (args_info.component_given)
195 reader->SetVectorComponent(args_info.component_arg);
199 const unsigned int Dimension = TProjectionsReaderType::OutputImageType::GetImageDimension();
200 typename TProjectionsReaderType::OutputImageDirectionType direction;
201 if (args_info.newdirection_given)
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);
208 typename TProjectionsReaderType::OutputImageSpacingType spacing;
209 if (args_info.newspacing_given)
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);
216 typename TProjectionsReaderType::OutputImagePointType origin;
217 if (args_info.neworigin_given)
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);
226 typename TProjectionsReaderType::OutputImageSizeType upperCrop, lowerCrop;
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);
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);
246 typename TProjectionsReaderType::ShrinkFactorsType binFactors;
248 for (
unsigned int i = 0; i < args_info.binning_given; i++)
249 binFactors[i] = args_info.binning_arg[i];
250 reader->SetShrinkFactors(binFactors);
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);
261 if (args_info.i0_given)
262 reader->SetI0(args_info.i0_arg);
263 reader->SetIDark(args_info.idark_arg);
266 if (args_info.nolineint_flag)
267 reader->ComputeLineIntegralOff();
270 if (args_info.wpc_given)
272 std::vector<double> coeffs;
273 coeffs.assign(args_info.wpc_arg, args_info.wpc_arg + args_info.wpc_given);
274 reader->SetWaterPrecorrectionCoefficients(coeffs);
278 reader->SetFileNames(fileNames);
290 template <
class TArgsInfo,
class TIterativeReconstructionFilter>
294 typename TIterativeReconstructionFilter::VolumeType::Pointer attenuationMap;
295 using VolumeType =
typename TIterativeReconstructionFilter::VolumeType;
296 if (args_info.attenuationmap_given)
299 attenuationMap = itk::ReadImage<VolumeType>(args_info.attenuationmap_arg);
302 switch (args_info.bp_arg)
307 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_VOXELBASED);
310 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_JOSEPH);
313 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_CUDAVOXELBASED);
316 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_CUDARAYCAST);
317 if (args_info.step_given)
318 recon->SetStepSize(args_info.step_arg);
321 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_JOSEPHATTENUATED);
322 if (args_info.attenuationmap_given)
323 recon->SetAttenuationMap(attenuationMap);
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);
345 template <
class TArgsInfo,
class TIterativeReconstructionFilter>
349 typename TIterativeReconstructionFilter::VolumeType::Pointer attenuationMap;
350 using VolumeType =
typename TIterativeReconstructionFilter::VolumeType;
351 if (args_info.attenuationmap_given)
354 attenuationMap = itk::ReadImage<VolumeType>(args_info.attenuationmap_arg);
357 typename ClipImageType::Pointer inferiorClipImage, superiorClipImage;
358 if (args_info.inferiorclipimage_given)
361 inferiorClipImage = itk::ReadImage<ClipImageType>(args_info.inferiorclipimage_arg);
363 if (args_info.superiorclipimage_given)
366 superiorClipImage = itk::ReadImage<ClipImageType>(args_info.superiorclipimage_arg);
369 switch (args_info.fp_arg)
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);
381 recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_CUDARAYCAST);
382 if (args_info.step_given)
383 recon->SetStepSize(args_info.step_arg);
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);
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);
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.