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 && !args_info.size_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];
69 imageSize.Fill(args_info.size_arg[0]);
70 for (
unsigned int i = 0; i < std::min(args_info.size_given, Dimension); i++)
71 imageSize[i] = args_info.size_arg[i];
74 typename ImageType::SpacingType imageSpacing;
75 imageSpacing.Fill(args_info.spacing_arg[0]);
76 for (
unsigned int i = 0; i < std::min(args_info.spacing_given, Dimension); i++)
77 imageSpacing[i] = args_info.spacing_arg[i];
79 typename ImageType::PointType imageOrigin;
80 for (
unsigned int i = 0; i < Dimension; i++)
81 imageOrigin[i] = imageSpacing[i] * (imageSize[i] - 1) * -0.5;
82 for (
unsigned int i = 0; i < std::min(args_info.origin_given, Dimension); i++)
83 imageOrigin[i] = args_info.origin_arg[i];
85 typename ImageType::DirectionType imageDirection;
86 if (args_info.direction_given)
87 for (
unsigned int i = 0; i < Dimension; i++)
88 for (
unsigned int j = 0; j < Dimension; j++)
89 imageDirection[i][j] = args_info.direction_arg[i * Dimension + j];
91 imageDirection.SetIdentity();
93 source->SetOrigin(imageOrigin);
94 source->SetSpacing(imageSpacing);
95 source->SetDirection(imageDirection);
96 source->SetSize(imageSize);
97 source->SetConstant({});
101 if (args_info.like_given)
104 auto likeReader = LikeReaderType::New();
105 likeReader->SetFileName(args_info.like_arg);
107 source->SetInformationFromImage(likeReader->GetOutput());
128 template <
class TArgsInfo>
129 const std::vector<std::string>
134 names->SetDirectory(args_info.path_arg);
135 names->SetNumericSort(args_info.nsort_flag);
136 names->SetRegularExpression(args_info.regexp_arg);
137 names->SetSubMatch(args_info.submatch_arg);
139 if (args_info.verbose_flag)
140 std::cout <<
"Regular expression matches " << names->GetFileNames().size() <<
" file(s)..." << std::endl;
143 if (args_info.submatch_given)
146 itksys::RegularExpression reg;
147 if (!reg.compile(args_info.regexp_arg))
149 itkGenericExceptionMacro(<<
"Error compiling regular expression " << args_info.regexp_arg);
153 for (
const std::string & name : names->GetFileNames())
156 if (reg.match(args_info.submatch_arg) == std::string(
""))
158 itkGenericExceptionMacro(<<
"Cannot find submatch " << args_info.submatch_arg <<
" in " << name
159 <<
" from regular expression " << args_info.regexp_arg);
164 std::vector<std::string> fileNames = names->GetFileNames();
166 std::vector<size_t> idxtopop;
168 for (
const auto & fn : fileNames)
173 if (imageio.IsNull())
175 std::cerr <<
"Ignoring file: " << fn <<
"\n";
176 idxtopop.push_back(i);
180 std::reverse(idxtopop.begin(), idxtopop.end());
181 for (
const auto &
id : idxtopop)
183 fileNames.erase(fileNames.begin() + id);
189 template <
class TProjectionsReaderType,
class TArgsInfo>
196 if (args_info.component_given)
198 reader->SetVectorComponent(args_info.component_arg);
202 const unsigned int Dimension = TProjectionsReaderType::OutputImageType::GetImageDimension();
203 typename TProjectionsReaderType::OutputImageDirectionType direction;
204 if (args_info.newdirection_given)
206 direction.Fill(args_info.newdirection_arg[0]);
207 for (
unsigned int i = 0; i < args_info.newdirection_given; i++)
208 direction[i / Dimension][i % Dimension] = args_info.newdirection_arg[i];
209 reader->SetDirection(direction);
211 typename TProjectionsReaderType::OutputImageSpacingType spacing;
212 if (args_info.newspacing_given)
214 spacing.Fill(args_info.newspacing_arg[0]);
215 for (
unsigned int i = 0; i < args_info.newspacing_given; i++)
216 spacing[i] = args_info.newspacing_arg[i];
217 reader->SetSpacing(spacing);
219 typename TProjectionsReaderType::OutputImagePointType origin;
220 if (args_info.neworigin_given)
222 origin.Fill(args_info.neworigin_arg[0]);
223 for (
unsigned int i = 0; i < args_info.neworigin_given; i++)
224 origin[i] = args_info.neworigin_arg[i];
225 reader->SetOrigin(origin);
229 auto lowerCrop = itk::MakeFilled<typename TProjectionsReaderType::OutputImageSizeType>(0);
230 for (
unsigned int i = 0; i < args_info.lowercrop_given; i++)
231 lowerCrop[i] = args_info.lowercrop_arg[i];
232 if (args_info.lowercrop_given)
233 reader->SetLowerBoundaryCropSize(lowerCrop);
235 auto upperCrop = itk::MakeFilled<typename TProjectionsReaderType::OutputImageSizeType>(0);
236 for (
unsigned int i = 0; i < args_info.uppercrop_given; i++)
237 upperCrop[i] = args_info.uppercrop_arg[i];
238 if (args_info.uppercrop_given)
239 reader->SetUpperBoundaryCropSize(upperCrop);
242 typename TProjectionsReaderType::MedianRadiusType medianRadius;
243 medianRadius.Fill(0);
244 for (
unsigned int i = 0; i < args_info.radius_given; i++)
245 medianRadius[i] = args_info.radius_arg[i];
246 reader->SetMedianRadius(medianRadius);
247 if (args_info.multiplier_given)
248 reader->SetConditionalMedianThresholdMultiplier(args_info.multiplier_arg);
251 typename TProjectionsReaderType::ShrinkFactorsType binFactors;
253 for (
unsigned int i = 0; i < args_info.binning_given; i++)
254 binFactors[i] = args_info.binning_arg[i];
255 reader->SetShrinkFactors(binFactors);
258 if (args_info.spr_given)
259 reader->SetScatterToPrimaryRatio(args_info.spr_arg);
260 if (args_info.nonneg_given)
261 reader->SetNonNegativityConstraintThreshold(args_info.nonneg_arg);
262 if (args_info.airthres_given)
263 reader->SetAirThreshold(args_info.airthres_arg);
266 if (args_info.i0_given)
267 reader->SetI0(args_info.i0_arg);
268 reader->SetIDark(args_info.idark_arg);
271 if (args_info.nolineint_flag)
272 reader->ComputeLineIntegralOff();
275 if (args_info.wpc_given)
277 std::vector<double> coeffs;
278 coeffs.assign(args_info.wpc_arg, args_info.wpc_arg + args_info.wpc_given);
279 reader->SetWaterPrecorrectionCoefficients(coeffs);
283 reader->SetFileNames(fileNames);
295 template <
class TArgsInfo,
class TIterativeReconstructionFilter>
299 typename TIterativeReconstructionFilter::VolumeType::Pointer attenuationMap;
300 using VolumeType =
typename TIterativeReconstructionFilter::VolumeType;
301 if (args_info.attenuationmap_given)
304 attenuationMap = itk::ReadImage<VolumeType>(args_info.attenuationmap_arg);
307 switch (args_info.bp_arg)
312 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_VOXELBASED);
315 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_JOSEPH);
318 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_CUDAVOXELBASED);
321 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_CUDARAYCAST);
322 if (args_info.step_given)
323 recon->SetStepSize(args_info.step_arg);
326 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_JOSEPHATTENUATED);
327 if (args_info.attenuationmap_given)
328 recon->SetAttenuationMap(attenuationMap);
331 recon->SetBackProjectionFilter(TIterativeReconstructionFilter::BP_ZENG);
332 if (args_info.sigmazero_given)
333 recon->SetSigmaZero(args_info.sigmazero_arg);
334 if (args_info.alphapsf_given)
335 recon->SetAlphaPSF(args_info.alphapsf_arg);
336 if (args_info.attenuationmap_given)
337 recon->SetAttenuationMap(attenuationMap);
350 template <
class TArgsInfo,
class TIterativeReconstructionFilter>
354 typename TIterativeReconstructionFilter::VolumeType::Pointer attenuationMap;
355 using VolumeType =
typename TIterativeReconstructionFilter::VolumeType;
356 if (args_info.attenuationmap_given)
359 attenuationMap = itk::ReadImage<VolumeType>(args_info.attenuationmap_arg);
362 typename ClipImageType::Pointer inferiorClipImage, superiorClipImage;
363 if (args_info.inferiorclipimage_given)
366 inferiorClipImage = itk::ReadImage<ClipImageType>(args_info.inferiorclipimage_arg);
368 if (args_info.superiorclipimage_given)
371 superiorClipImage = itk::ReadImage<ClipImageType>(args_info.superiorclipimage_arg);
374 switch (args_info.fp_arg)
379 recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_JOSEPH);
380 if (args_info.superiorclipimage_given)
381 recon->SetSuperiorClipImage(superiorClipImage);
382 if (args_info.inferiorclipimage_given)
383 recon->SetInferiorClipImage(inferiorClipImage);
386 recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_CUDARAYCAST);
387 if (args_info.step_given)
388 recon->SetStepSize(args_info.step_arg);
391 recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_JOSEPHATTENUATED);
392 if (args_info.superiorclipimage_given)
393 recon->SetSuperiorClipImage(superiorClipImage);
394 if (args_info.inferiorclipimage_given)
395 recon->SetInferiorClipImage(inferiorClipImage);
396 if (args_info.attenuationmap_given)
397 recon->SetAttenuationMap(attenuationMap);
400 recon->SetForwardProjectionFilter(TIterativeReconstructionFilter::FP_ZENG);
401 if (args_info.sigmazero_given)
402 recon->SetSigmaZero(args_info.sigmazero_arg);
403 if (args_info.alphapsf_given)
404 recon->SetAlphaPSF(args_info.alphapsf_arg);
405 if (args_info.attenuationmap_given)
406 recon->SetAttenuationMap(attenuationMap);
413 #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.