CASToR  1.1
Tomographic Reconstruction (PET/SPECT)
 All Classes Files Functions Variables Typedefs Macros Groups Pages
castor-proj.cc
Go to the documentation of this file.
1 
9 #include "gVariables.hh"
10 #include "gOptions.hh"
11 #include "oAnalyticProjection.hh"
14 #include "iDataFilePET.hh"
15 #include "iDataFileSPECT.hh"
16 #include "iDataFileTransmission.hh"
17 #include "sOutputManager.hh"
18 #include "sScannerManager.hh"
20 #include "iScannerPET.hh"
21 
22 
23 // =====================================================================
24 // ---------------------------------------------------------------------
25 // ---------------------------------------------------------------------
26 // =====================================================================
31 void ShowHelp(int a_returnCode, int a_mpiRank)
32 {
33  // Return when using MPI and mpi_rank is not 0
34  #ifdef CASTOR_MPI
35  int mpi_rank = 0;
36  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
37  if (mpi_rank!=0) return;
38  #endif
39 
40  // Show help
41  cout << endl;
42  cout << "Usage: castor-proj.exe -img path_to_file.hdr -sc name -(f/d)out output [settings]" << endl;
43  cout << endl;
44  cout << "[Main options]:" << endl;
45  cout << " -img path_to_img.hdr : Give the interfile header of the image to project (no default)." << endl;
46  cout << " -sc scanner_alias : Give the scanner model. It must correspond to the name of one of the file in the scanner repository (w/o file extension)" << endl;
47  cout << " -fout name : Give the root name for all output files (no default, alternative to -dout)" << endl;
48  cout << " -dout name : Give the name of the output directory where all output files will be written (no default, alternative to -fout)" << endl;
49  cout << " -fileType type : Give the type of events that will be generated : ('TYPE_PET' ou '0' = PET, 'TYPE_SPECT' ou '1' = SPECT) (default : TYPE_PET)" << endl;
50  cout << " -vb : Give the verbosity level, from 0 (no verbose) to above 5 (at the event level) (default: 1)." << endl;
51  cout << endl;
52  cout << "[Specific options]:" << endl;
53  cout << " -help-in : Print out specific help about input settings." << endl; // managed by main
54  cout << " -help-out : Print out specific help about output settings." << endl; // managed by main
55  cout << " -help-proj : Print out specific help about projector settings." << endl; // managed by main
56  cout << " -help-comp : Print out specific help about computing settings." << endl; // managed by main
57  cout << " -help-misc : Print out specific help about miscellaneous and verbose settings." << endl; // managed by main
58  cout << " -help-dynamic : Print out specific help about dynamic options." << endl; // managed by main
59  cout << " -help-pet : Print out specific help about PET settings for analytic projection." << endl; // managed by main
60  cout << " -help-spect : Print out specific help about SPECT settings for analytic projection." << endl; // managed by main
61  cout << endl;
62  cout << "[Implemented Modules]:" << endl;
63  cout << " -help-scan : Show the list of all scanners from the configuration directory." << endl; // managed by sScannerManager
64  cout << " -help-projm : Show the list and description of all implemented projectors." << endl; // managed by oProjectorManager
65  cout << endl;
66  cout << " --help,-h,-help : Print out this help page." << endl; // managed by main
67  cout << endl;
68  #ifdef CASTOR_MPI
69  cout << " Compiled with MPI" << endl;
70  #endif
71  #ifdef CASTOR_OMP
72  cout << " Compiled with OpenMP" << endl;
73  #endif
74  #ifdef CASTOR_GPU
75  cout << " Compiled for GPU" << endl;
76  #endif
77  #if defined(CASTOR_OMP) || defined(CASTOR_MPI) || defined(CASTOR_GPU)
78  cout << endl;
79  #endif
80  #ifdef BUILD_DATE
81  cout << " Build date: " << BUILD_DATE << endl;
82  cout << endl;
83  #endif
84  #ifdef CASTOR_VERSION
85  cout << " This program is part of the CASToR release version " << CASTOR_VERSION << "." << endl;
86  cout << endl;
87  #endif
88  Exit(a_returnCode);
89 }
90 
91 
92 
93 
94 // =====================================================================
95 // ---------------------------------------------------------------------
96 // ---------------------------------------------------------------------
97 // =====================================================================
103 {
104  // Return when using MPI and mpi_rank is not 0
105  #ifdef CASTOR_MPI
106  int mpi_rank = 0;
107  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
108  if (mpi_rank!=0) return;
109  #endif
110  // Show help
111  cout << endl;
112  cout << "[Input settings]:" << endl;
113  cout << endl;
114  cout << " -off x,y,z : Give the offset of the field-of-view in each dimension, in mm (default: 0.,0.,0.)." << endl;
115  cout << endl;
116  cout << " -atn path_to_img.hdr : Give the interfile header of an input attenuation image (unit has to be cm-1) (default: uniform value)." << endl;
117  cout << endl;
118  cout << " -help-in : Print out this help." << endl;
119  cout << endl;
120 }
121 
122 
123 
124 // =====================================================================
125 // ---------------------------------------------------------------------
126 // ---------------------------------------------------------------------
127 // =====================================================================
133 {
134  // Return when using MPI and mpi_rank is not 0
135  #ifdef CASTOR_MPI
136  int mpi_rank = 0;
137  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
138  if (mpi_rank!=0) return;
139  #endif
140  // Show help
141  cout << endl;
142  cout << "[Output settings]:" << endl;
143  cout << endl;
144  cout << " -fout name : Give the root name for all output files. All output files will be written as 'name_suffix.ext'." << endl;
145  cout << " So the provided name should not end with '.' or '/' character. (no default, alternative to -dout)" << endl;
146  cout << " -dout name : Give the name of the output directory where all output files will be written. All files will also" << endl;
147  cout << " be prefixed by the name of the directory. The provided name should not end with '.' or '/' character." << endl;
148  cout << " (no default, alternative to -fout)" << endl;
149  cout << endl;
150  cout << " -olut : If a scanner LUT (geometry information) is computed from a .geom file, it will be written on disk in the scanner repository" << endl;
151  cout << endl;
152  cout << " -help-out : Print out this help." << endl;
153  cout << endl;
154 }
155 
156 
157 
158 // =====================================================================
159 // ---------------------------------------------------------------------
160 // ---------------------------------------------------------------------
161 // =====================================================================
167 {
168  // Return when using MPI and mpi_rank is not 0
169  #ifdef CASTOR_MPI
170  int mpi_rank = 0;
171  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
172  if (mpi_rank!=0) return;
173  #endif
174  // Show help
175  cout << "[Projector settings]:" << endl;
176  cout << endl;
177  cout << " -proj param : Give the projector to be used for both forward and backward projections, along with a configuration file (proj:file.conf)" << endl;
178  cout << " or the list of parameters associated to the projector (proj,param1,param2,...). If the projector only is specified, the" << endl;
179  cout << " default configuration file is used. By default, the Siddon projector is used." << endl;
180  cout << endl;
181  cout << " -proj-comp : Give the strategy for projection line computation. Here are the three different strategies that can be used:" << endl;
182  cout << " 1 : Image-based system matrix elements storage: The voxels weights are added in a matrix representing the whole image, so" << endl;
183  cout << " the addition of a new line to the previous ones is straightforward only by adding the weights to the corresponding voxels." << endl;
184  cout << " As it is insanely long, it can possibly be used for example with extremely complex projectors that makes use of huge number" << endl;
185  cout << " of ray tracings for a single event, where the list of contributions can become longer than the number of voxels in the image." << endl;
186  cout << " 2 : Fixed-size list storage of system matrix elements: The voxels are added one by one in two separated lists, one containing voxel" << endl;
187  cout << " indices and the other voxel weights. When a voxel is added to the oProjectionLine, it is simply pilled-up to the list. The list" << endl;
188  cout << " has a fixed size which is provided by the EstimateMaxNumberOfVoxelsPerLine() function from the vProjector class. There are no" << endl;
189  cout << " ckecks at all for possible buffer overflows. This is the fastest strategy and default one." << endl;
190  cout << " 3 : Adaptative-size list storage of system matrix elements: This is the same as the fixed-size strategy except that the size can be" << endl;
191  cout << " upgraded if the current number of contributing voxels exceed the list's size. The first allocated size corresponds to the diagonal" << endl;
192  cout << " of the image. During the first iteration, this size will be upgraded until it will reach a size suitable for all events. Thus it" << endl;
193  cout << " is a bit slower than the fixed-list strategy during the first iteration, but is optimal with respect to RAM usage." << endl;
194  cout << endl;
195  cout << " -help-projm : Print out specific help about projector settings." << endl; // managed by oProjectorManager
196  cout << endl;
197  cout << " -help-proj : Print out this help." << endl; // managed by oProjectorManager
198  cout << endl;
199 }
200 
201 
202 
203 // =====================================================================
204 // ---------------------------------------------------------------------
205 // ---------------------------------------------------------------------
206 // =====================================================================
212 {
213  // Return when using MPI and mpi_rank is not 0
214  #ifdef CASTOR_MPI
215  int mpi_rank = 0;
216  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
217  if (mpi_rank!=0) return;
218  #endif
219  // Show help
220  cout << endl;
221  cout << "[Computation settings]:" << endl;
222  cout << endl;
223  #ifdef CASTOR_OMP
224  cout << " -th param : Set the number of threads for parallel computation (default: 1). If 0 is given, then the maximum number of available threads is automatically determined." << endl;
225  cout << " Can also give two parameters separated by a comma (e.g. 16,4), to distinguish between the number of threads for projection and image operations respectively." << endl;
226  cout << endl;
227  #endif
228  #ifdef CASTOR_GPU
229  cout << " -gpu : Flag to say that we want to use the GPU device (default: use the CPU only)." << endl;
230  cout << endl;
231  #endif
232  cout << " -conv param : Give an image convolver model to be used within the algorithm, along with a configuration file (conv:file.conf) or the" << endl;
233  cout << " list of parameters associated to the convolver (conv,param1,param2,...). If the convolver only is specified, its default" << endl;
234  cout << " configuration file is used. By default, no convolver is applied" << endl;
235  cout << endl;
236  cout << " -help-conv : Show the list and description of all implemented image convolvers." << endl; // managed by oImageConvolverManager
237  cout << endl;
238  cout << " -noise poisson : Give the noise model to be used. Only Poisson noise is implemented for the moment (default : no noise model)" << endl;
239  cout << endl;
240  cout << " -rng : Give a seed for the random number generator (should be >=0)" << endl;
241  cout << endl;
242  cout << " -help-comp : Print out this help." << endl;
243  cout << endl;
244  #ifdef CASTOR_MPI
245  cout << " Compiled with MPI" << endl;
246  #endif
247  #ifdef CASTOR_OMP
248  cout << " Compiled with OpenMP" << endl;
249  #endif
250  #ifdef CASTOR_GPU
251  cout << " Compiled for GPU" << endl;
252  #endif
253  #if defined(CASTOR_OMP) || defined(CASTOR_MPI) || defined(CASTOR_GPU)
254  cout << endl;
255  #endif
256 }
257 
258 
259 
260 // =====================================================================
261 // ---------------------------------------------------------------------
262 // ---------------------------------------------------------------------
263 // =====================================================================
269 {
270  // Return when using MPI and mpi_rank is not 0
271  #ifdef CASTOR_MPI
272  int mpi_rank = 0;
273  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
274  if (mpi_rank!=0) return;
275  #endif
276  // Show help
277  cout << endl;
278  cout << "[Dynamic settings]:" << endl;
279  cout << endl;
280  cout << " -frm list : Give the framing for the reconstruction where 'list' is a list of all frames durations in seconds separated" << endl;
281  cout << " by commas. To specify a dead frame, add ';0' concatenated to the frame duration. Milliseconds maximum precision." << endl;
282  cout << " (default: read from the 'image duration (sec)' key in the interfile image header" << endl;
283  cout << endl;
284  cout << " (about gated images) : An image contening several respiratory/cardiac 'gates' can be provided to generate one datafile containing histograms of these gated images." << endl;
285  cout << " The number of gates being directly read from the 'number of respiratory gates' and 'number of cardiac gates' keys in the interfile header" << endl;
286  cout << endl;
287  cout << " -help-dynamic : Print out this help." << endl;
288  cout << endl;
289 }
290 
291 
292 
293 // =====================================================================
294 // ---------------------------------------------------------------------
295 // ---------------------------------------------------------------------
296 // =====================================================================
302 {
303  // Return when using MPI and mpi_rank is not 0
304  #ifdef CASTOR_MPI
305  int mpi_rank = 0;
306  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
307  if (mpi_rank!=0) return;
308  #endif
309  // Show help
310  cout << "[SPECT settings]:" << endl;
311  cout << endl;
312  cout << " -SPECT_bins par : Give two integer parameters corresponding to the transaxial number of bins, separated by a comma (nbBinsX, nbBinsY). " << endl;
313  cout << " Default = transaxial number of pixels in the SPECT scanner file " << endl;
314  cout << endl;
315  cout << " -SPECT_ang par : Give the total number of projections, followed by the first angle and the angle step in float" << endl;
316  cout << " Format : nb_projs:first_angle,step_angle" << endl;
317  cout << endl;
318  cout << " -SPECT_c_ang par : SPECT custom angles : give the total number of projections, followed by the projections angles in float" << endl;
319  cout << " Format : nb_projs:angle1,angle2,... " << endl;
320  cout << endl;
321  cout << " -SPECT_radius par : Give a distance between center of rotation to detectors which will be used for all heads" << endl;
322  cout << endl;
323  cout << " -SPECT_rot par : Give the heads rotation direction. Accept two parameters: CW (clockwise, default), CCW (counter-clockwise)" << endl;
324  cout << endl;
325  cout << " -help-spect : Print out this specific help." << endl; // managed by oProjectorManager
326  cout << endl;
327 }
328 
329 
330 
331 // =====================================================================
332 // ---------------------------------------------------------------------
333 // ---------------------------------------------------------------------
334 // =====================================================================
340 {
341  // Return when using MPI and mpi_rank is not 0
342  #ifdef CASTOR_MPI
343  int mpi_rank = 0;
344  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
345  if (mpi_rank!=0) return;
346  #endif
347  // Show help
348  cout << "[PET settings]:" << endl;
349  cout << endl;
350  cout << " -PET_maxRingDiff param : Set the maximum number of differences between the ring ID of the two crystals (1 = 2D mode) (default = total number of rings in the scanner)" << endl;
351  cout << endl;
352 }
353 
354 
355 
356 // =====================================================================
357 // ---------------------------------------------------------------------
358 // ---------------------------------------------------------------------
359 // =====================================================================
365 {
366  // Return when using MPI and mpi_rank is not 0
367  #ifdef CASTOR_MPI
368  int mpi_rank = 0;
369  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
370  if (mpi_rank!=0) return;
371  #endif
372  // Show help
373  cout << endl;
374  cout << "[Miscellaneous settings]:" << endl;
375  cout << endl;
376  cout << " -vb : Give the general verbosity level, from 0 (no verbose) to 5 and above (at the event level) (default: 1)." << endl;
377  cout << " -vb-algo : Give the verbose level specific to the analytic projection algorithm (default: same as general verbose level)." << endl;
378  cout << " -vb-proj : Give the verbose level specific to the projector (default: same as general verbose level)." << endl;
379  cout << " -vb-conv : Give the verbose level specific to the image convolver (default: same as general verbose level)." << endl;
380  cout << " -vb-scan : Give the verbose level specific to the scanner (default: same as general verbose level)." << endl;
381  cout << " -vb-data : Give the verbose level specific to the data and image management (default: same as general verbose level)." << endl;
382  cout << endl;
383  cout << " -conf : Give the path to the CASToR config directory (default: located through the CASTOR_CONFIG environment variable)." << endl;
384  cout << endl;
385  cout << " -help-misc : Print out this help." << endl;
386  cout << endl;
387 }
388 
389 
390 
391 
392 
393 // =====================================================================
394 // ---------------------------------------------------------------------
395 // ---------------------------------------------------------------------
396 // =====================================================================
397 /*
398  Main program
399 */
400 int main(int argc, char** argv)
401 {
402  // ============================================================================================================
403  // MPI management
404  // ============================================================================================================
405  int mpi_rank = 0;
406  //int mpi_size = 1;
407  #ifdef CASTOR_MPI
408  int mpi_size = 1;
409  MPI_Init(&argc, &argv );
410  MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
411  MPI_Comm_size(MPI_COMM_WORLD, &mpi_size);
412  #endif
413 
414  // No argument, then show help
415  if (argc==1) ShowHelp(0,mpi_rank);
416  if (mpi_rank==0) cout << endl; // TODO : Most stuff is handled with processor 1 at first
417 
418  // ============================================================================================================
419  // Parameterized variables with their default values
420  // ============================================================================================================
421 
422  // Initialization of the voxel and field-of-view values in each spatial dimensions.
423  // Read from the interfile header of the image
424  int nb_VoxX=-1, nb_VoxY=-1, nb_VoxZ=-1;
425  FLTNB fov_SizeX=-1, fov_SizeY=-1, fov_SizeZ=-1;
426  FLTNB vox_SizeX=-1, vox_SizeY=-1, vox_SizeZ=-1;
427  // Initialization offset values of the field-of-view in each spatial dimensions
428  FLTNB offsetX = 0, offsetY = 0, offsetZ = 0;
429 
430  // Path to an image used as initialization. no default
431  string path_to_initial_img = "";
432  // Path to an image used for the attenuation. default : uniform image
433  string path_to_attenuation_img = "";
434  // Scanner name provided by the user
435  string scanner_name = "";
436 
437  // DataFile mode
438  int datafile_mode = MODE_HISTOGRAM;
439  int datafile_type = TYPE_PET;
440 
441  // Frame list descriptor
442  string frame_list = "";
443  // Number of resp gates
444  int nb_resp_gates = 1;
445  // Number of card gates
446  int nb_card_gates = 1;
447 
448  // Number of beds position in the data. default : 1
449  int nb_beds = 1;
450 
451  // Maximum ring difference between the crystals in a LOR
452  int max_ring_diff = -1;
453 
454  // Output directory name.
455  string path_dout = "";
456  // Or root name
457  string path_fout = "";
458 
459  // Write scanner LUT generated by a geom file on disk
460  bool save_LUT_flag = false;
461 
462  // String gathering the name of the used noise model (default : none)
463  string noise_model = "";
464 
465  // String gathering the name of the used projector for forward/backward projection (default : Siddon)
466  string options_projector = "incrementalSiddon";
467 
468  // Default projector computation strategy
469  int projector_computation_strategy = FIXED_LIST_COMPUTATION_STRATEGY;
470 
471  // String vector gathering the name of the image convolvers with specific options (default: no image convolver)
472  vector<string> options_image_convolver;
473 
474  // Using GPU (flag) ->NOTE : default : only CPU
475  bool gpu_flag = false;
476  // Size of the buffer used to read the data file (0 by default for projection as it is not used, but should be initialized to use the Datafile objects)
477  int data_file_percentage_load = 0;
478 
479  // Verbose level
480  int verbose_general = 1,
481  verbose_algo = -1,
482  verbose_proj = -1,
483  verbose_conv = -1,
484  verbose_scan = -1,
485  verbose_data = -1;
486 
487  // Path to config directory
488  string path_to_config_dir = "";
489 
490  // Initial seed for random number generator
491  int64_t seed_RNG = -1;
492  // Number of threads
493  string nb_threads = "1";
494 
495  // SPECT specific parameters
496 
497  // SPECT number of transaxial bins
498  uint16_t SPECT_nb_bins[2] = {0, 0};
499  // SPECT Number of projection angles
500  uint32_t SPECT_nb_projection_angles = 0;
501  // SPECT projections angles
502  FLTNB* SPECT_projection_angles = NULL;
503  // SPECT first angle and step angle for automatic initialization
504  FLTNB SPECT_first_angle= -1., SPECT_step_angle=-1.;
505  // SPECT distance from center of rotation to detectors
506  FLTNB SPECT_radius = -1.;
507  // SPECT head rotation orientation
508  string SPECT_rot = "CW";
509 
510  // ============================================================================================================
511  // Read command-line parameters
512  // ============================================================================================================
513 
514  for (int i=1; i<argc; i++)
515  {
516  string option = (string)argv[i];
517 
518  if (option=="-h" || option=="--help" || option=="-help") ShowHelp(0,mpi_rank);
519 
520  // Specific help for integrated scanners (managed by scanner manager)
521  else if (option=="-help-scan")
522  {
523  if(sScannerManager::GetInstance()->ShowScannersDescription())
524  Cerr("***** castor-proj() -> An error occured when trying to output the available scanners from the scanner repository !'" << endl;);
525  Exit(EXIT_SUCCESS);
526  }
527 
528 
529  // Specific help for projector settings (managed by projector children)
530  else if (option=="-help-projm")
531  {
532  // Call specific showHelp function from vProjector children
534  // Call the static showHelp function from vProjector
536  Exit(EXIT_SUCCESS);
537  }
538 
539  // Specific help for image convolver settings (managed by image convolver children)
540  else if (option=="-help-conv")
541  {
542  // Call specific showHelp function from vImageConvolver children
544  // Call the static showHelp function from oImageConvolverManager
545  // Call the static showHelp function from oImageConvolverManager
547  Exit(EXIT_SUCCESS);
548  }
549 
550  // Specific help for input settings (managed by main)
551  else if (option=="-help-in")
552  {
553  ShowHelpInput();
554  Exit(EXIT_SUCCESS);
555  }
556  // Specific help for output settings (managed by main)
557  else if (option=="-help-out")
558  {
559  ShowHelpOutput();
560  Exit(EXIT_SUCCESS);
561  }
562  // Specific help for projector settings (managed by main)
563  else if (option=="-help-proj")
564  {
565  ShowHelpProj();
566  Exit(EXIT_SUCCESS);
567  }
568  // Specific help for computation settings (managed by main)
569  else if (option=="-help-comp")
570  {
572  Exit(EXIT_SUCCESS);
573  }
574  // Specific help for dynamic settings (managed by main)
575  else if (option=="-help-dynamic")
576  {
577  ShowHelpDynamic();
578  Exit(EXIT_SUCCESS);
579  }
580  // Specific help for PET projection settings (managed by main)
581  else if (option=="-help-pet")
582  {
583  ShowHelpPET();
584  Exit(EXIT_SUCCESS);
585  }
586  // Specific help for SPECT projection settings (managed by main)
587  else if (option=="-help-spect")
588  {
589  ShowHelpSPECT();
590  Exit(EXIT_SUCCESS);
591  }
592  // Specific help for miscellaneous settings (managed by main)
593  else if (option=="-help-misc")
594  {
596  Exit(EXIT_SUCCESS);
597  }
598 
599  else if (option=="-off")
600  {
601  if (i>=argc-1)
602  {
603  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
604  Exit(EXIT_FAILURE);
605  }
606  else
607  {
608  FLTNB input[3];
609  if (ReadStringOption(argv[i+1], input, 3, ",", option))
610  {
611  Cerr("***** castor-proj() -> Invalid argument " << argv[i+1] << " for option " << option << " !" << endl);
612  Exit(EXIT_FAILURE);
613  }
614  offsetX = input[0];
615  offsetY = input[1];
616  offsetZ = input[2];
617  }
618  i++;
619  }
620 
621  else if (option=="-frm")
622  {
623  if (i>=argc-1)
624  {
625  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
626  Exit(EXIT_FAILURE);
627  }
628  frame_list = (string)argv[i+1];
629  i++;
630  }
631 
632  // Scanner name
633  else if (option=="-sc")
634  {
635  if (i>=argc-1)
636  {
637  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
638  Exit(EXIT_FAILURE);
639  }
640  scanner_name = argv[i+1];
641  i++;
642  }
643 
644  // Image for the initialisation of the algorithm : What should we do in case of multiple frames ?
645  else if (option=="-img")
646  {
647  if (i>=argc-1)
648  {
649  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
650  Exit(EXIT_FAILURE);
651  }
652  path_to_initial_img = argv[i+1];
653  i++;
654  }
655 
656 
657  else if (option=="-fileType")
658  {
659  if (i>=argc-1)
660  {
661  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
662  Exit(EXIT_FAILURE);
663  }
664 
665  string input_str = argv[i+1];
666 
667  if(!input_str.compare("TYPE_PET") )
668  datafile_type = TYPE_PET;
669  else if(!input_str.compare("TYPE_SPECT") )
670  datafile_type = TYPE_SPECT;
671  else if(!input_str.compare("TYPE_TRANSMISSION") )
672  datafile_type = TYPE_TRANSMISSION;
673  else // check for number
674  {
675  if(ConvertFromString(input_str , &datafile_type))
676  {
677  Cerr("***** castor-proj() -> Exception when trying to read provided entry '" << input_str << " for option: " << option << endl);
678  Exit(EXIT_FAILURE);
679  }
680  }
681 
682  i++;
683  }
684 
685  // Image for the attenuation
686  else if (option=="-atn")
687  {
688  if (i>=argc-1)
689  {
690  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
691  Exit(EXIT_FAILURE);
692  }
693  path_to_attenuation_img = argv[i+1];
694  i++;
695  }
696 
697  // Number of beds
698  else if (option=="-bed")
699  {
700  if (i>=argc-1)
701  {
702  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
703  Exit(EXIT_FAILURE);
704  }
705  nb_beds = atoi(argv[i+1]);
706  i++;
707  }
708 
709  // maximum number of differences between the ring ID of the two crystals
710  else if (option=="-PET_maxRingDiff")
711  {
712  if (i>=argc-1)
713  {
714  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
715  Exit(EXIT_FAILURE);
716  }
717  max_ring_diff = atoi(argv[i+1]);
718  i++;
719  }
720 
721  // Name of the output directory
722  else if (option=="-dout") // This is a mandatory option
723  {
724  if (i>=argc-1)
725  {
726  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
727  Exit(EXIT_FAILURE);
728  }
729  path_dout = argv[i+1];
730  i++;
731  }
732  // Base name of the output files
733  else if (option=="-fout") // This is a mandatory option
734  {
735  if (i>=argc-1)
736  {
737  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
738  Exit(EXIT_FAILURE);
739  }
740  path_fout = argv[i+1];
741  i++;
742  }
743  // Flag to say that we want to save time basis functions too
744  else if (option=="-olut")
745  {
746  save_LUT_flag = true;
747  }
748  else if (option=="-noise")
749  {
750  noise_model = argv[i+1];
751  i++;
752  }
753 
754  else if (option=="-proj")
755  {
756  if (i>=argc-1)
757  {
758  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
759  Exit(EXIT_FAILURE);
760  }
761  options_projector = (string)argv[i+1];
762  i++;
763  }
764 
765  // Projection line computation strategy
766  else if (option=="-proj-comp")
767  {
768  if (i>=argc-1)
769  {
770  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
771  Exit(EXIT_FAILURE);
772  }
773  projector_computation_strategy = atoi(argv[i+1]);
774  i++;
775  }
776 
777  // Image convolver settings
778  else if (option=="-conv")
779  {
780  if (i>=argc-1)
781  {
782  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
783  Exit(EXIT_FAILURE);
784  }
785  string convolver = (string)argv[i+1];
786  convolver.append("::forward"); // convolution only allowed before projection
787  options_image_convolver.push_back(convolver);
788  i++;
789  }
790 
791 
792  // --- Verbosities ---
793 
794  // General verbosity level
795  else if (option=="-vb")
796  {
797  if (i>=argc-1)
798  {
799  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
800  Exit(EXIT_FAILURE);
801  }
802  if (ConvertFromString(argv[i+1], &verbose_general))
803  {
804  Cerr("***** castor-proj() -> Exception when trying to read provided verbosity level '" << verbose_general << " for option: " << option << endl);
805  Exit(EXIT_FAILURE);
806  }
807  i++;
808  }
809 
810  // Algorithm verbosity level
811  else if (option=="-vb-algo")
812  {
813  if (i>=argc-1)
814  {
815  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
816  Exit(EXIT_FAILURE);
817  }
818  if (ConvertFromString(argv[i+1], &verbose_algo))
819  {
820  Cerr("***** castor-proj() -> Exception when trying to read provided verbosity level '" << verbose_algo << " for option: " << option << endl);
821  Exit(EXIT_FAILURE);
822  }
823  i++;
824  }
825 
826  // Projector verbosity level
827  else if (option=="-vb-proj")
828  {
829  if (i>=argc-1)
830  {
831  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
832  Exit(EXIT_FAILURE);
833  }
834  if (ConvertFromString(argv[i+1], &verbose_proj))
835  {
836  Cerr("***** castor-proj() -> Exception when trying to read provided verbosity level '" << verbose_proj << " for option: " << option << endl);
837  Exit(EXIT_FAILURE);
838  }
839  i++;
840  }
841 
842  // Image convolver verbosity level
843  else if (option=="-vb-conv")
844  {
845  if (i>=argc-1)
846  {
847  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
848  Exit(EXIT_FAILURE);
849  }
850  if (ConvertFromString(argv[i+1], &verbose_conv))
851  {
852  Cerr("***** castor-proj() -> Exception when trying to read provided verbosity level '" << verbose_conv << " for option: " << option << endl);
853  Exit(EXIT_FAILURE);
854  }
855  i++;
856  }
857 
858  // Scanner verbosity level
859  else if (option=="-vb-scan")
860  {
861  if (i>=argc-1)
862  {
863  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
864  Exit(EXIT_FAILURE);
865  }
866  if (ConvertFromString(argv[i+1], &verbose_scan))
867  {
868  Cerr("***** castor-proj() -> Exception when trying to read provided verbosity level '" << verbose_scan << " for option: " << option << endl);
869  Exit(EXIT_FAILURE);
870  }
871  i++;
872  }
873 
874  // Data and image verbosity level
875  else if (option=="-vb-data")
876  {
877  if (i>=argc-1)
878  {
879  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
880  Exit(EXIT_FAILURE);
881  }
882  if (ConvertFromString(argv[i+1], &verbose_data))
883  {
884  Cerr("***** castor-proj() -> Exception when trying to read provided verbosity level '" << verbose_data << " for option: " << option << endl);
885  Exit(EXIT_FAILURE);
886  }
887  i++;
888  }
889 
890 
891  // RNG seed
892  else if (option=="-rng")
893  {
894  if (i>=argc-1)
895  {
896  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
897  Exit(EXIT_FAILURE);
898  }
899  if(ConvertFromString(argv[i+1], &seed_RNG))
900  {
901  Cerr("***** castor-proj() -> Exception when trying to read provided number '" << seed_RNG << " for option: " << option << endl);
902  Exit(EXIT_FAILURE);
903  }
904  i++;
905  }
906 
907  // Path to config directory
908  else if (option=="-conf")
909  {
910  if (i>=argc-1)
911  {
912  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
913  Exit(EXIT_FAILURE);
914  }
915  path_to_config_dir = (string)argv[i+1];
916  i++;
917  }
918 
919  #ifdef CASTOR_GPU
920  else if (option=="-gpu")
921  {
922  gpu_flag = 1;
923  }
924  #endif
925  #ifdef CASTOR_OMP
926  else if (option=="-th")
927  {
928  if (i>=argc-1)
929  {
930  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
931  Exit(EXIT_FAILURE);
932  }
933  nb_threads = (string)argv[i+1];
934  i++;
935  }
936  #endif
937 
938  else if (option=="-load")
939  {
940  if (i>=argc-1)
941  {
942  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
943  Exit(EXIT_FAILURE);
944  }
945  data_file_percentage_load = atoi(argv[i+1]);
946  if (data_file_percentage_load>100 || data_file_percentage_load<0)
947  {
948  Cerr("***** castor-proj() -> Incorrect initialization of the size data buffer" << endl);
949  Cerr(" Number provided: " << argv[i+1] << " is not in the expected [0;100] interval" << endl);
950  Exit(EXIT_FAILURE);
951  }
952  i++;
953  }
954  else if (option=="-SPECT_bins")
955  {
956  if (i>=argc-1)
957  {
958  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
959  Exit(EXIT_FAILURE);
960  }
961 
962  if(ReadStringOption(argv[i+1], SPECT_nb_bins, 2, ",", option))
963  {
964  Cerr("***** castor-proj() -> Invalid argument " << argv[i+1] << " for option " << option << " !" << endl);
965  Cerr("***** castor-proj() -> Expected two uint16 parameters separated by a comma" << endl);
966  Exit(EXIT_FAILURE);
967  }
968 
969  i++;
970  }
971 
972  // Custom initialization of SPECT projection angles using a number of projections, and an equivalent number of pre-set projection angles
973  else if (option=="-SPECT_c_ang")
974  {
975  if (i>=argc-1)
976  {
977  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
978  Exit(EXIT_FAILURE);
979  }
980 
981  // parsing
982  string input_str = argv[i+1];
983  size_t colon = input_str.find_first_of(":");
984 
985  if (colon==string::npos)
986  {
987  Cerr("***** castor-proj() -> Parsing error : colon missing in option: " << option << endl);
988  Exit(EXIT_FAILURE);
989  }
990  else
991  {
992  if(ConvertFromString(input_str.substr(0,colon), &SPECT_nb_projection_angles))
993  {
994  Cerr("***** castor-proj() -> An error occured while trying to read the number of SPECT projection values " << input_str.substr(0,colon) << " for option " << option << " !" << endl);
995  Exit(EXIT_FAILURE);
996  }
997 
998  SPECT_projection_angles = new FLTNB[SPECT_nb_projection_angles];
999 
1000  if(ReadStringOption(input_str.substr(colon+1), SPECT_projection_angles, SPECT_nb_projection_angles, ",", option))
1001  {
1002  Cerr("***** castor-proj() -> Invalid argument " << argv[i+1] << " for option " << option << " !" << endl);
1003  Exit(EXIT_FAILURE);
1004  }
1005  }
1006  i++;
1007  }
1008 
1009  // Generic initialization of SPECT projection angles using a first angle, an angle step, and a number of projection
1010  else if (option=="-SPECT_ang")
1011  {
1012  if (i>=argc-1)
1013  {
1014  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
1015  Exit(EXIT_FAILURE);
1016  }
1017 
1018  // parsing
1019  string input_str = argv[i+1];
1020  size_t colon = input_str.find_first_of(":");
1021 
1022  if (colon==string::npos)
1023  {
1024  Cerr("***** castor-proj() -> Parsing error : colon missing in option: " << option << endl);
1025  Exit(EXIT_FAILURE);
1026  }
1027  else
1028  {
1029  if(ConvertFromString(input_str.substr(0,colon), &SPECT_nb_projection_angles))
1030  {
1031  Cerr("***** castor-proj() -> An error occured while trying to read the number of SPECT projection values " << input_str.substr(0,colon) << " for option " << option << " !" << endl);
1032  Exit(EXIT_FAILURE);
1033  }
1034 
1035  FLTNB input[2];
1036  if(ReadStringOption(input_str.substr(colon+1), input, 2, ",", option))
1037  {
1038  Cerr("***** castor-proj() -> Invalid argument " << argv[i+1] << " for option " << option << " !" << endl);
1039  Exit(EXIT_FAILURE);
1040  }
1041  else
1042  {
1043  SPECT_first_angle = input[0];
1044  SPECT_step_angle = input[1];
1045  }
1046  }
1047  i++;
1048  }
1049 
1050  else if (option=="-SPECT_radius")
1051  {
1052  if (i>=argc-1)
1053  {
1054  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
1055  Exit(EXIT_FAILURE);
1056  }
1057 
1058  else
1059  {
1060  if(ConvertFromString(argv[i+1], &SPECT_radius) )
1061  {
1062  Cerr("***** castor-proj() -> An error occured while trying to read the SPECT global radius " << argv[i+1] << " for option " << option << " !" << endl);
1063  Exit(EXIT_FAILURE);
1064  }
1065  }
1066  i++;
1067  }
1068 
1069  else if (option=="-SPECT_rot")
1070  {
1071  if (i>=argc-1)
1072  {
1073  Cerr("***** castor-proj() -> Argument missing for option: " << option << endl);
1074  Exit(EXIT_FAILURE);
1075  }
1076 
1077  else
1078  {
1079  SPECT_rot = argv[i+1];
1080 
1081  if(SPECT_rot != "CW" && SPECT_rot != "CCW")
1082  {
1083  Cerr("***** castor-proj() -> '" << SPECT_rot <<"' is unknown argument for option " << option << ".");
1084  Cerr(" SPECT rotation direction must be either 'CW' (clockwise) or 'CCW' (counter-clockwise) !" << endl);
1085  Exit(EXIT_FAILURE);
1086  }
1087  }
1088  i++;
1089  }
1090 
1091  else
1092  {
1093  if (mpi_rank==0) Cerr("***** castor-proj() -> Unknown option '" << option << "' !" << endl);
1094  Exit(EXIT_FAILURE);
1095  }
1096 
1097  }
1098 
1099  #ifdef CASTOR_MPI
1100  MPI_Barrier(MPI_COMM_WORLD);
1101  #endif
1102 
1103 
1104  // Affect specific verbose levels if not set
1105  if (verbose_algo==-1) verbose_algo = verbose_general;
1106  if (verbose_proj==-1) verbose_proj = verbose_general;
1107  if (verbose_conv==-1) verbose_conv = verbose_general;
1108  if (verbose_scan==-1) verbose_scan = verbose_general;
1109  if (verbose_data==-1) verbose_data = verbose_general;
1110 
1111 
1112 
1113  // ============================================================================================================
1114  // Mandatory checks: the idea is that the most likely problems are checked here, so that we do not apply a
1115  // too much secured approach in the rest of the classes in order to not overload the code
1116  // with checks everywhere, which could also possibly slows down its execution.
1117  // ============================================================================================================
1118 
1119  // Basic initialization checks (minimal initializations mandatory for the next steps)
1120 
1121  // scanner
1122  if (scanner_name.empty())
1123  {
1124  Cerr("***** castor-proj() -> Please provide a scanner name (format: Modality_Constuctor_Model) :" << endl);
1125  Cerr(" -sc scanner_alias: Give the scanner model." << endl);
1126  Cerr(" It must correspond to the name of one of the file in the scanner repository (without file extension)" << endl);
1127  Exit(EXIT_FAILURE);
1128  }
1129  else
1130  {
1131  if (verbose_general>=2) Cout(" selected scanner name: " << scanner_name << endl);
1132  }
1133 
1134  // Output files
1135  if (path_fout.empty() && path_dout.empty())
1136  {
1137  Cerr("***** castor-proj() -> Please provide an output option for output files (-fout or -dout) !" << endl);
1138  Exit(EXIT_FAILURE);
1139  }
1140  // Check that only one option has been provided
1141  if (!path_fout.empty() && !path_dout.empty())
1142  {
1143  Cerr("***** castor-proj() -> Please provide either output option -fout or -dout but not both !" << endl);
1144  Exit(EXIT_FAILURE);
1145  }
1146 
1147  // Image source
1148  if (path_to_initial_img.empty())
1149  {
1150  Cerr("***** castor-proj() -> -img path_to_img.hdr : Give an interfile header of the image to project" << endl);
1151  Exit(EXIT_FAILURE);
1152  }
1153  else
1154  {
1155  if (verbose_general>=2) Cout(" selected input image path: " << path_to_initial_img << endl);
1156  }
1157 
1158  Cout(endl);
1159 
1160 
1161 
1162 
1163  // ============================================================================================================
1164  // Singletons initialization: create here all the needed singletons.
1165  // ============================================================================================================
1166 
1167  // Get user endianness (for interfile)
1169 
1170  // OutputManager
1171  sOutputManager* p_outputManager;
1172  p_outputManager = sOutputManager::GetInstance();
1173  p_outputManager->SetVerbose(verbose_general);
1174 
1175  // Set path to the config directory
1176  if (p_outputManager->CheckConfigDir(path_to_config_dir))
1177  {
1178  Cerr("***** castor-proj() -> A problem occured while checking for the config directory path !" << endl);
1179  Exit(EXIT_FAILURE);
1180  }
1181 
1182  // Initialize output directory and base name
1183  if (p_outputManager->InitOutputDirectory(path_fout, path_dout))
1184  {
1185  Cerr("***** castor-proj() -> A problem occured while initializing output directory !" << endl);
1186  Exit(EXIT_FAILURE);
1187  }
1188  // Log command line
1189  if (p_outputManager->LogCommandLine(argc,argv))
1190  {
1191  Cerr("***** castor-proj() -> A problem occured while logging command line arguments !" << endl);
1192  Exit(EXIT_FAILURE);
1193  }
1194 
1195  // ScannerManager
1196  sScannerManager* p_scannermanager = sScannerManager::GetInstance();
1197  p_scannermanager->SetVerbose(verbose_scan);
1198  p_scannermanager->SetSaveLUTFlag(save_LUT_flag);
1199 
1200 
1201 
1202  // ============================================================================================================
1203  // Main part of the program: create all elements needed for the algorithm, create the algorithm, initialize it
1204  // and launch it.
1205  // ============================================================================================================
1206 
1207 
1208  // ============================================================================================================
1209  // Generate system geometry
1210  // ============================================================================================================
1211 
1212  scanner_name = (scanner_name.find(OS_SEP)) ?
1213  scanner_name.substr(scanner_name.find_last_of(OS_SEP)+1) :
1214  scanner_name;
1215 
1216  if(p_scannermanager->FindScannerSystem(scanner_name) )
1217  {
1218  Cerr("***** castor-proj() -> A problem occurred while searching for scanner system !" << endl);
1219  Exit(EXIT_FAILURE);
1220  }
1221 
1222  if(p_scannermanager->BuildScannerObject() ) //TODO verifier si correct d'envoyer uniquement le chemin du premier DataFile (pb si dataFile ne proviennent pas du même scan ou de la même acquisition)
1223  {
1224  Cerr("***** castor-proj() -> A problem occurred during scanner object construction !" << endl);
1225  Exit(EXIT_FAILURE);
1226  }
1227 
1228  if(p_scannermanager->InstantiateScanner() )
1229  {
1230  Cerr("***** castor-proj() -> A problem occurred while creating Scanner object !" << endl);
1231  Exit(EXIT_FAILURE);
1232  }
1233 
1234 
1235  //Set the specific mandatory geometry information of the modality using the SetXXXSpecificParameters functions of the scannerManager
1236  if (p_scannermanager->GetScannerType() == SCANNER_PET)
1237  {
1238  if(p_scannermanager->PROJ_SetPETSpecificParameters(max_ring_diff) )
1239  {
1240  Cerr("***** castor-proj() -> A problem occurred while setting PET geometric parameters !" << endl);
1241  Exit(EXIT_FAILURE);
1242  }
1243  }
1244  else if (p_scannermanager->GetScannerType() == SCANNER_SPECT_CONVERGENT || p_scannermanager->GetScannerType() == SCANNER_SPECT_PINHOLE )
1245  {
1246  if(p_scannermanager->PROJ_SetSPECTSpecificParameters(SPECT_nb_bins,
1247  SPECT_nb_projection_angles,
1248  SPECT_first_angle,
1249  SPECT_step_angle,
1250  SPECT_projection_angles,
1251  SPECT_radius,
1252  SPECT_rot) )
1253  {
1254  Cerr("***** castor-proj() -> A problem occurred while setting SPECT geometric parameters !" << endl);
1255  Exit(EXIT_FAILURE);
1256  }
1257  }
1258 
1259 
1260  if(p_scannermanager->BuildLUT() )
1261  {
1262  Cerr("***** castor-proj() -> A problem occurred while generating/reading the LUT !" << endl);
1263  Exit(EXIT_FAILURE);
1264  }
1265 
1266 
1267  // Check the scanner manager parameters and initialize the scanner
1268  if (p_scannermanager->CheckParameters())
1269  {
1270  Cerr("***** castor-proj() -> A problem occured while checking scanner manager parameters !" << endl);
1271  Exit(EXIT_FAILURE);
1272  }
1273  if (p_scannermanager->Initialize())
1274  {
1275  Cerr("***** castor-proj() -> A problem occured while initializing scanner !" << endl);
1276  Exit(EXIT_FAILURE);
1277  }
1278 
1279 
1280 
1281  // ============================================================================================================
1282  // Recover image dimensions informations from the interfile header of the image to project
1283  // ============================================================================================================
1284 
1285  Intf_fields IF;
1286  IntfKeyInitFields(&IF);
1287  if(IntfReadHeader(path_to_initial_img, &IF, verbose_data) )
1288  {
1289  Cerr("***** castor-proj() -> An error occurred while trying to read the interfile header of file reading file " << path_to_initial_img << " !" << endl);
1290  Exit(EXIT_FAILURE);
1291  }
1292 
1293  if(verbose_general >= 3) IntfKeyPrintFields(IF);
1294 
1295  nb_VoxX = nb_VoxX<0 ? IF.mtx_size[0] : nb_VoxX ;
1296  nb_VoxY = nb_VoxY<0 ? IF.mtx_size[1] : nb_VoxY ;
1297  nb_VoxZ = nb_VoxZ<0 ? IF.mtx_size[2] : nb_VoxZ ;
1298 
1299  // TODO : nb_beds ?
1300  // If voxel sizes not provided/found, set them to the default value (1mm)
1301  vox_SizeX = IF.vox_size[0]<0 ? 1. : IF.vox_size[0];
1302  vox_SizeY = IF.vox_size[1]<0 ? 1. : IF.vox_size[1];
1303  vox_SizeZ = IF.vox_size[2]<0 ? 1. : IF.vox_size[2];
1304 
1305 
1306  // Get the frame duration from the image in the case it is not provided by the user
1307  if (frame_list == "")
1308  {
1309  vector <string> image_filenames;
1310  if ( IntfIsMHD(path_to_initial_img, image_filenames) )
1311  {
1312  Cerr("***** castor-proj() -> Error : while trying to read the interfile header '"<< path_to_initial_img << "' !" << endl);
1313  Exit(EXIT_FAILURE);
1314  }
1315 
1316  if(IF.image_duration.size()>0) // Image durations have been provided
1317  {
1318  if(image_filenames.size() > 1)
1319  {
1320  // If image files are splitted over a serie of 3D image files, we will have image durations for each image, including gates
1321  // As we need to recover the duration for each gate once, we skip gates if the data contains any
1322  for(int fr=0 ; fr<IF.nb_time_frames ; fr++)
1323  {
1324  // Check nb of frames and parse framing
1325  if(IF.image_duration.size() != (uint32_t)IF.nb_time_frames *
1326  IF.nb_resp_gates *
1327  IF.nb_card_gates )
1328  {
1329  Cerr("***** castor-proj() -> Interfile reading error of the input image :"<<endl);
1330  Cerr(" The number of provided image duration: '"<< IF.image_duration.size()
1331  << "' does not match the actual number of time frames * respiratory gates * cardiac gates ) '"<< IF.nb_time_frames *
1332  IF.nb_resp_gates *
1333  IF.nb_card_gates <<"' !" << endl);
1334  Exit(EXIT_FAILURE);
1335  }
1336 
1337  // Recover image duration and dead frames (if any) in the frame list
1338  stringstream ss_id, ss_df;
1339  // Go to the next frame (skip gate if any)
1340  ss_id << IF.image_duration.at(fr * IF.nb_resp_gates* IF.nb_card_gates);
1341  string frame_duration_str = ss_id.str();
1342 
1343  // Check if there is any dead frame (pause between frame groups)
1344  if(IF.frame_group_pause.size()>0)
1345  {
1346  FLTNB dead_frame_sec = IF.frame_group_pause.at(fr * IF.nb_resp_gates* IF.nb_card_gates);
1347  if(dead_frame_sec > 0.)
1348  {
1349  ss_df << dead_frame_sec;
1350  frame_duration_str.append(",").append(ss_df.str()).append(";0");
1351  }
1352  }
1353 
1354  frame_list.append(frame_duration_str);
1355 
1356  // add comma if not last frame
1357  if(fr+1 < IF.nb_time_frames)
1358  frame_list.append(",");
1359  }
1360  }
1361  else
1362  for(int fr=0 ; fr<IF.nb_time_frames ; fr++)
1363  {
1364  if(IF.image_duration.size() != IF.nb_time_frames)
1365  {
1366  Cerr("***** castor-proj() -> Interfile reading error of the input image :"<<endl);
1367  Cerr(" The number of provided image duration: '"<< IF.image_duration.size()
1368  << "' does not match the actual number of time frames * respiratory gates * cardiac gates ) '"<< IF.nb_time_frames <<"' !" << endl);
1369  Exit(EXIT_FAILURE);
1370  }
1371 
1372  stringstream ss;
1373  ss << IF.image_duration.at(fr);
1374  frame_list.append(ss.str());
1375  // add comma if not last frame
1376  if(fr+1 < IF.nb_time_frames)
1377  frame_list.append(",");
1378  }
1379  }
1380  }
1381 
1382  // Check nb gating
1383  (IF.nb_resp_gates >1) ? nb_resp_gates = IF.nb_resp_gates : 1 ;
1384  (IF.nb_card_gates >1) ? nb_card_gates = IF.nb_card_gates : 1 ;
1385 
1386  //TODO : offsets ?
1387 
1388  // Default initialization of max ring difference with the max number of rings
1389  // TODO : How to deal with potential geometry containing several layers of crystals with different number of rings ?
1390  // Only for PET
1391  if(p_scannermanager->GetScannerType()==SCANNER_PET && max_ring_diff<0) max_ring_diff = p_scannermanager->GetScannerLayerNbRings(0);
1392 
1393  if(nb_VoxX<0 || nb_VoxY<0 || nb_VoxZ<0)
1394  {
1395  ReadDataASCIIFile(p_scannermanager->GetPathToScannerFile(), "voxels number transaxial", &nb_VoxX, 1, KEYWORD_MANDATORY);
1396  ReadDataASCIIFile(p_scannermanager->GetPathToScannerFile(), "voxels number transaxial", &nb_VoxY, 1, KEYWORD_MANDATORY);
1397  ReadDataASCIIFile(p_scannermanager->GetPathToScannerFile(), "voxels number axial", &nb_VoxZ, 1, KEYWORD_MANDATORY);
1398  }
1399 
1400  if ( (fov_SizeX<0 || fov_SizeY<0 || fov_SizeZ<0) && (vox_SizeX<0 || vox_SizeY<0 || vox_SizeZ<0) )
1401  {
1402  ReadDataASCIIFile(p_scannermanager->GetPathToScannerFile(), "field of view transaxial", &fov_SizeX, 1, KEYWORD_MANDATORY);
1403  ReadDataASCIIFile(p_scannermanager->GetPathToScannerFile(), "field of view transaxial", &fov_SizeY, 1, KEYWORD_MANDATORY);
1404  ReadDataASCIIFile(p_scannermanager->GetPathToScannerFile(), "field of view axial", &fov_SizeZ, 1, KEYWORD_MANDATORY);
1405  }
1406 /*
1407  else if ( (fov_SizeX<0 || fov_SizeY<0 || fov_SizeZ<0) && (vox_SizeX>=0 || vox_SizeY>=0 || vox_SizeZ>=0))
1408  {
1409  fov_SizeX = vox_SizeX*(FLTNB)nb_VoxX;
1410  fov_SizeY = vox_SizeY*(FLTNB)nb_VoxY;
1411  fov_SizeZ = vox_SizeZ*(FLTNB)nb_VoxZ;
1412  cout << vox_SizeX << " " << fov_SizeX << endl;
1413  }
1414  else
1415  {
1416  vox_SizeX = fov_SizeX/(FLTNB)nb_VoxX;
1417  vox_SizeY = fov_SizeY/(FLTNB)nb_VoxY;
1418  vox_SizeZ = fov_SizeZ/(FLTNB)nb_VoxZ;
1419  }
1420 */
1421  if (verbose_general>=2)
1422  {
1423  Cout(" Image dimensions for analytical projections: " << endl);
1424  Cout(" voxels number (X,Y,Z): " << nb_VoxX << "," << nb_VoxY << "," << nb_VoxZ << endl);
1425  Cout(" voxels size (X,Y,Z): " << vox_SizeX << "," << vox_SizeY << "," << vox_SizeZ << endl);
1426  Cout(" FOV (X,Y,Z): " << fov_SizeX << "," << fov_SizeY << "," << fov_SizeZ << endl);
1427  }
1428 
1429  // ============================================================================================================
1430  // oImageDimensionsAndQuantification settings
1431  // ============================================================================================================
1432 
1433  // Create oImageDimensionsAndQuantification
1434  oImageDimensionsAndQuantification* p_ImageDimensionsAndQuantification = new oImageDimensionsAndQuantification();
1435  p_ImageDimensionsAndQuantification->SetNbVoxX(nb_VoxX);
1436  p_ImageDimensionsAndQuantification->SetNbVoxY(nb_VoxY);
1437  p_ImageDimensionsAndQuantification->SetNbVoxZ(nb_VoxZ);
1438  p_ImageDimensionsAndQuantification->SetNbThreads(nb_threads);
1439  p_ImageDimensionsAndQuantification->SetNbBeds(nb_beds);
1440  p_ImageDimensionsAndQuantification->SetVoxSizeX(vox_SizeX);
1441  p_ImageDimensionsAndQuantification->SetVoxSizeY(vox_SizeY);
1442  p_ImageDimensionsAndQuantification->SetVoxSizeZ(vox_SizeZ);
1443  p_ImageDimensionsAndQuantification->SetFOVSizeX(fov_SizeX);
1444  p_ImageDimensionsAndQuantification->SetFOVSizeY(fov_SizeY);
1445  p_ImageDimensionsAndQuantification->SetFOVSizeZ(fov_SizeZ);
1446  p_ImageDimensionsAndQuantification->SetFOVOutMasking(0., 0);
1447  p_ImageDimensionsAndQuantification->SetOffsetX(offsetX);
1448  p_ImageDimensionsAndQuantification->SetOffsetY(offsetY);
1449  p_ImageDimensionsAndQuantification->SetOffsetZ(offsetZ);
1450  p_ImageDimensionsAndQuantification->SetVerbose(verbose_data);
1451  p_ImageDimensionsAndQuantification->SetNbRespGates(nb_resp_gates);
1452  p_ImageDimensionsAndQuantification->SetNbCardGates(nb_card_gates);
1453  p_ImageDimensionsAndQuantification->SetFrames(frame_list);
1454  if (p_ImageDimensionsAndQuantification->CheckParameters())
1455  {
1456  Cerr("***** castor-proj() -> A problem occured while checking image dimensions parameters !" << endl);
1457  Exit(EXIT_FAILURE);
1458  }
1459  if (p_ImageDimensionsAndQuantification->Initialize())
1460  {
1461  Cerr("***** castor-proj() -> A problem occured while initializing image dimensions !" << endl);
1462  Exit(EXIT_FAILURE);
1463  }
1464 
1465  // Initialization of DynamicDataManager class, related 4D data splitting management
1466  if (p_ImageDimensionsAndQuantification->InitDynamicData( "", 0, 0, 0, 0, nb_resp_gates, nb_card_gates ) )
1467  {
1468  Cerr("***** castor-proj() -> A problem occured while initializing Dynamic data manager's class !" << endl);
1469  Exit(EXIT_FAILURE);
1470  }
1471 
1472 
1473  // ----------------------------------------------------------------------------------------
1474  // Random Number Generator initialization: (we first required to know the number of threads to use from p_ImageDimensionsAndQuantification)
1475  // ----------------------------------------------------------------------------------------
1476 
1478  p_RNG->SetVerbose(verbose_general);
1479  // Use a user-provided seed to initialize the RNG if one has been provided. Use random number otherwise.
1480  (seed_RNG>=0 )?
1481  p_RNG->Initialize(seed_RNG, p_ImageDimensionsAndQuantification->GetNbThreadsForProjection()):
1482  p_RNG->Initialize(p_ImageDimensionsAndQuantification->GetNbThreadsForProjection());
1483 
1484 
1485  // ============================================================================================================
1486  // Datafile initialization
1487  // ============================================================================================================
1488 
1489  vDataFile** p_DataFile = new vDataFile*[nb_beds];
1490 
1491  if (p_scannermanager->GetScannerType() < 0)
1492  {
1493  Cerr("***** castor-proj() -> Unknown scanner type !" << endl);
1494  Exit(EXIT_FAILURE);
1495  }
1496  else if (p_scannermanager->GetScannerType() == SCANNER_PET)
1497  {
1498  if (datafile_type!=TYPE_PET)
1499  {
1500  Cerr("***** castor-proj() -> Only PET type events are allowed while using a PET scanner !" << endl);
1501  Cerr("***** castor-proj() -> (use '-fileType' to define modality (0=PET, 1=SPECT, 2=CT) )" << endl);
1502  Exit(EXIT_FAILURE);
1503  }
1504  for (int i=0 ; i<nb_beds ; i++)
1505  {
1506  p_DataFile[i] = new iDataFilePET();
1507  if(path_to_attenuation_img != "") // Enable correction factor for attenuation
1508  (dynamic_cast<iDataFilePET*>(p_DataFile[i]))->SetAtnCorrectionFlagOn();
1509  }
1510  }
1511  else if (p_scannermanager->GetScannerType() == SCANNER_SPECT_CONVERGENT)
1512  {
1513  if (datafile_type!=TYPE_SPECT)
1514  {
1515  Cerr("***** castor-proj() -> Only SPECT type events are allowed while using a SPECT scanner !" << endl);
1516  Cerr("***** castor-proj() -> (use '-fileType' to define modality (0=PET, 1=SPECT, 2=CT) )" << endl);
1517  Exit(EXIT_FAILURE);
1518  }
1519  for (int i=0 ; i<nb_beds ; i++)
1520  {
1521  p_DataFile[i] = new iDataFileSPECT();
1522  }
1523  }
1524  else
1525  {
1526  Cerr("***** castor-proj() -> Unsupported scanner type !" << endl);
1527  Exit(EXIT_FAILURE);
1528  }
1529 
1530  // Load raw data in memory and do other stuff if needed.
1531  // TODO TM : Even if the loop is there, no support for multi-bed projection as for now (case we have to project an image as in a multi-bad acquisition)
1532  // It would first require some informations about the FOV of each bed (get this from interfile header of image or for command line options ?)
1533  // In the loop, we need to provide the bed in argument of PROJ_InitFile() in order to correctly initialize the beds
1534  for (int bed=0 ; bed<nb_beds ; bed++)
1535  {
1536  p_DataFile[bed]->SetBedIndex(bed);
1537  p_DataFile[bed]->SetPercentageLoad(data_file_percentage_load);
1538  p_DataFile[bed]->SetVerbose(verbose_data);
1539  p_DataFile[bed]->SetDataMode(datafile_mode);
1540  p_DataFile[bed]->SetDataType(datafile_type);
1541  p_DataFile[bed]->SetImageDimensionsAndQuantification(p_ImageDimensionsAndQuantification);
1542  p_DataFile[bed]->PROJ_GetScannerSpecificParameters();
1543  p_DataFile[bed]->PROJ_InitFile();
1544  p_DataFile[bed]->ComputeSizeEvent();
1545 
1546  // The Check parameters function is dedicated to reconstruction, as it checks things such as consistenties between size of the datafile and user-provided information
1547  if(p_DataFile[bed]->CheckParameters())
1548  {
1549  Cerr("***** castor-proj() -> A problem occurred while checking datafile parameters ! Abort." << endl);
1550  Exit(EXIT_FAILURE);
1551  }
1552  p_DataFile[bed]->PrepareDataFile();
1553 
1554  /*
1555  if(p_DataFile[bed]->SetScannerSpecificParameters())
1556  {
1557  Cerr("***** A problem occured while setting scanner parameters from the datafile !" << endl);
1558  Exit(EXIT_FAILURE);
1559  }
1560  */
1561 
1562  //p_DataFile[bed]->InitDynamicData(p_ImageDimensionsAndQuantification->GetNbTimeFrames(), nb_resp_gates, nb_card_gates, "", 0, 0, bed);
1563 
1564  /*
1565  if (p_DataFile[bed]->InitDynamicData(nb_resp_gates, nb_card_gates, "", 0, 0, 0, p_DataFile[bed]) )
1566  {
1567  Cerr("***** A problem occured while initializing Dynamic data !" << endl);
1568  Exit(EXIT_FAILURE);
1569  }
1570  */
1571  }
1572 
1573 
1574  // ============================================================================================================
1575  // Projector manager initialization
1576  // ============================================================================================================
1577 
1578  // Create object
1579  oProjectorManager* p_ProjectorManager = new oProjectorManager();
1580 
1581  // Set all parameters
1582  p_ProjectorManager->SetScanner(p_scannermanager->GetScannerObject());
1583  p_ProjectorManager->SetImageDimensionsAndQuantification(p_ImageDimensionsAndQuantification);
1584  p_ProjectorManager->SetDataFile(p_DataFile[0]);
1585  p_ProjectorManager->SetComputationStrategy(projector_computation_strategy);
1586  p_ProjectorManager->SetOptionsForward(options_projector);
1587  p_ProjectorManager->SetOptionsBackward(options_projector);
1588  p_ProjectorManager->SetOptionsCommon("");
1589  p_ProjectorManager->SetVerbose(verbose_proj);
1590 
1591  // Check parameters
1592  if (p_ProjectorManager->CheckParameters())
1593  {
1594  Cerr("***** castor-proj() -> A problem occured while checking projector manager's parameters !" << endl);
1595  Exit(EXIT_FAILURE);
1596  }
1597 
1598  // Initialize projector manager
1599  if (p_ProjectorManager->Initialize())
1600  {
1601  Cerr("***** castor-proj() -> A problem occured while initializing projector manager !" << endl);
1602  Exit(EXIT_FAILURE);
1603  }
1604 
1605 
1606 
1607 
1608  // ============================================================================================================
1609  // Image Convolver manager initialization
1610  // ============================================================================================================
1611 
1612  // Create object
1613  oImageConvolverManager* p_ImageConvolverManager = new oImageConvolverManager();
1614  // Set all parameters
1615  p_ImageConvolverManager->SetVerbose(verbose_conv);
1616  p_ImageConvolverManager->SetImageDimensionsAndQuantification(p_ImageDimensionsAndQuantification);
1617  p_ImageConvolverManager->SetOptions(options_image_convolver);
1618 
1619  // Check parameters
1620  if (p_ImageConvolverManager->CheckParameters())
1621  {
1622  Cerr("***** castor-proj() -> A problem occured while checking optimizer manager's parameters !" << endl);
1623  Exit(EXIT_FAILURE);
1624  }
1625  // Initialize optimizer manager
1626  if (p_ImageConvolverManager->Initialize())
1627  {
1628  Cerr("***** castor-proj() -> A problem occured while initializing optimizer manager !" << endl);
1629  Exit(EXIT_FAILURE);
1630  }
1631 
1632 
1633 
1634  // ============================================================================================================
1635  // Image Space initialization
1636  // ============================================================================================================
1637 
1638  oImageSpace* p_ImageSpace = new oImageSpace();
1639  p_ImageSpace->SetImageDimensionsAndQuantification(p_ImageDimensionsAndQuantification);
1640  p_ImageSpace->SetVerbose(verbose_data);
1641 
1642 
1643 
1644  // ============================================================================================================
1645  // Init alogorithm and proceed to analytical projection
1646  // ============================================================================================================
1647 
1648  oAnalyticProjection* p_Projection = new oAnalyticProjection();
1649  p_Projection->InitOptimizer(p_ImageDimensionsAndQuantification);
1650  if (p_Projection->InitNoiseModel(noise_model) )
1651  {
1652  Cerr("***** castor-proj() -> A problem occured while initializing noise model !" << endl);
1653  Exit(EXIT_FAILURE);
1654  }
1655 
1656  p_Projection->SetImageDimensionsAndQuantification(p_ImageDimensionsAndQuantification);
1657  p_Projection->SetProjectorManager(p_ProjectorManager);
1658  p_Projection->SetImageConvolverManager(p_ImageConvolverManager);
1659  p_Projection->SetImageSpace(p_ImageSpace);
1660  p_Projection->SetDataFile(p_DataFile);
1661  p_Projection->SetGPUflag(gpu_flag);
1662  p_Projection->SetVerbose(verbose_algo);
1663  p_Projection->SetNbBeds(nb_beds);
1664  p_Projection->SetPathInitImage(path_to_initial_img);
1665  p_Projection->SetPathAtnImage(path_to_attenuation_img);
1666  p_Projection->SetScanner(p_scannermanager->GetScannerObject());
1667 
1668  // Launch analytical projection
1669  if (p_Projection->Launch())
1670  {
1671  Cerr("***** castor-proj() -> Error occured during Projections" << endl;);
1672  }
1673 
1674 
1675  // ============================================================================================================
1676  // End
1677  // ============================================================================================================
1678 
1679  if(SPECT_projection_angles) delete SPECT_projection_angles;
1680 
1681  if (verbose_general>=5) Cout("----- Deleting CASToR objects ... -----" << endl);
1682  delete p_Projection;
1683  delete p_ImageSpace;
1684  delete p_ImageConvolverManager;
1685  delete p_ProjectorManager;
1686 
1687 
1688  // Delete objects
1689  for (int i=0 ; i<nb_beds ; i++) delete p_DataFile[i];
1690  delete[] p_DataFile;
1691  delete p_ImageDimensionsAndQuantification;
1692 
1693 
1694  if (verbose_general>=5) Cout("----- CASToR objects successfully deleted -----" << endl);
1695 
1696  // Ending
1697  if (verbose_general>=1) Cout(endl);
1698  #ifdef CASTOR_MPI
1699  MPI_Finalize();
1700  #endif
1701  return EXIT_SUCCESS;
1702 }
void SetDataMode(int a_dataMode)
set the data mode
Definition: vDataFile.hh:359
This class is designed to be a mother virtual class for Datafile.
Definition: vDataFile.hh:67
This header file is mainly used to declare some macro definitions and all includes needed from the st...
static sScannerManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
void SetImageDimensionsAndQuantification(oImageDimensionsAndQuantification *ap_ImageDimensionsAndQuantification)
Set the member mp_ImageDimensionsAndQuantification to the provided value.
static sRandomNumberGenerator * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
Declaration of class oImageDimensionsAndQuantification.
void SetNbBeds(int a_nbBeds)
Set number of beds (bed positions)
int InitDynamicData(string a_pathTo4DDataSplittingFile, int a_respMotionCorrectionFlag, int a_cardMotionCorrectionFlag, int a_doubleMotionCorrectionFlag, int a_invMotionCorrectionFlag, int a_nbRespGates, int a_nbCardGates)
Call the eponym function from the oDynamicDataManager object in order to initialize its data...
#define MODE_HISTOGRAM
Definition: vDataFile.hh:36
#define TYPE_PET
Definition: vDataFile.hh:51
void SetVerbose(int a_verbose)
set verbosity
void SetFOVSizeZ(FLTNB a_fovSizeZ)
Set the FOV's size along the Z axis, in mm.
static void ShowCommonHelp()
This function is used to print out some help about the use of options common to all projectors...
Definition: vProjector.cc:127
FLTNB vox_size[3]
#define FLTNB
Definition: gVariables.hh:55
This class manages the analytic projection of an image and the computation of the associated datafile...
void SetComputationStrategy(int a_computationStrategy)
Set the computation strategy for the system matrix elements storage.
int Initialize()
A function used to initialize the manager and the projectors or system matrices it manages...
void SetVerbose(int a_verboseLevel)
Set the member m_verboseLevel to the provided value.
void SetImageConvolverManager(oImageConvolverManager *ap_ImageConvolverManager)
Set the Image Convolver Manager Object.
void SetImageDimensionsAndQuantification(oImageDimensionsAndQuantification *ap_ImageDimensionsAndQuantification)
Set the image dimensions in use.
int FindScannerSystem(string a_scannerName)
Look for a file matching with the scanner name in parameter inside the scanner repository.
void SetFOVSizeY(FLTNB a_fovSizeY)
Set the FOV's size along the Y axis, in mm.
void SetNbVoxZ(INTNB a_nbVoxZ)
Set the number of voxels along the Z axis.
void SetDataFile(vDataFile **a2p_DataFile)
Set the list of DataFile.
int BuildScannerObject()
Instantiate the specific scanner object related to the modality, and set verbosity of scanner object...
void IntfKeyPrintFields(Intf_fields a_IF)
Print all the keys of the Intf_fields structure passed in parameter, as well as their values for debu...
void SetNbRespGates(int a_nbRespGates)
Set the number of respiratory gates.
int CheckParameters()
A function used to check the parameters settings.
void SetPathAtnImage(string a_pathToAtnImg)
Set path to an attenuation image.
void SetOffsetY(FLTNB a_offsetY)
Set the image offset along the Y axis, in mm.
void ShowHelpDynamic()
Display command line options related to the dynamic features for castor-proj.
Definition: castor-proj.cc:268
void SetFOVOutMasking(FLTNB a_fovOutPercent, INTNB a_nbSliceOutMask)
Set the output FOV masking settings: transaxial unmasked FOV percent and number of extrem slices to r...
void SetVerbose(int a_verboseLevel)
set verbosity
Definition: vDataFile.hh:394
#define TYPE_TRANSMISSION
Definition: vDataFile.hh:55
void SetOptionsForward(const string &a_optionsForward)
Set the forward projection options contained in the provided string.
void SetPathInitImage(string a_pathToInitialImage)
Set path to an initial image.
virtual int ComputeSizeEvent()=0
This function is implemented in child classes Computation of the size of each event according to th...
vector< FLTNB > frame_group_pause
string GetPathToScannerFile()
void SetVerbose(int a_verboseLevel)
Set verbosity level.
static sOutputManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
void SetOffsetZ(FLTNB a_offsetZ)
Set the image offset along the Z axis, in mm.
#define TYPE_SPECT
Definition: vDataFile.hh:53
#define SCANNER_PET
void ShowHelpProj()
Display command line options related to the Projector module for castor-proj.
Definition: castor-proj.cc:166
void SetSaveLUTFlag(bool a_flag)
Set to on the flag indicating a LUT generated by a geom file should be written on disk or not...
Declaration of class iDataFilePET.
void SetImageDimensionsAndQuantification(oImageDimensionsAndQuantification *ap_ImageDimensionsAndQuantification)
Set the Image Dimensions and Quantification Object.
void SetOptions(vector< string > a_options)
Set the member m_options to the provided value.
void SetOptionsCommon(const string &a_optionsCommon)
Set the common projection options contained in the provided string.
int CheckParameters()
Check if all parameters have been correctly initialized, and call the CheckParameters function of the...
void SetNbVoxX(INTNB a_nbVoxX)
Set the number of voxels along the X axis.
int SetNbThreads(const string &a_nbThreads)
Set the number of threads.
uint16_t nb_resp_gates
Declaration of class iDataFileSPECT.
void Exit(int code)
void SetOptionsBackward(const string &a_optionsBackward)
Set the backward projection options contained in the provided string.
uint16_t nb_card_gates
Declaration of class iScannerPET.
int IntfReadHeader(const string &a_pathToHeaderFile, Intf_fields *ap_IntfFields, int vb)
Read an Interfile header.
int InstantiateScanner()
Instantiate scanner using the related function in the scanner classes.
int LogCommandLine(int argc, char **argv)
Write log file header with the provided command line options and different informations.
static sAddonManager * GetInstance()
int IntfIsMHD(string a_pathToFile, vector< string > &ap_lPathToImgs)
Check if the string in argument contains the path to a Interfile metaheader.
Declaration of class oImageConvolverManager.
int ConvertFromString(const string &a_str, string *a_result)
Copy the 'a_str' string in the position pointed by 'a_result'.
Definition: gOptions.cc:749
#define FIXED_LIST_COMPUTATION_STRATEGY
int InitNoiseModel(string aNoiseModel)
void SetNbVoxY(INTNB a_nbVoxY)
Set the number of voxels along the Y axis.
int CheckConfigDir(const string &a_path)
Set the path to the CASTOR config directory from the given path if not empty or through the existence...
void SetProjectorManager(oProjectorManager *ap_ProjectorManager)
Set the Projector Manager Object.
int CheckParameters()
A function used to check the parameters settings.
void SetNbCardGates(int a_nbCardGates)
Set the number of cardiac gates.
void ShowHelpSPECT()
Definition: castor-proj.cc:301
#define Cerr(MESSAGE)
#define SCANNER_SPECT_CONVERGENT
int BuildLUT()
Call the eponym function of the scanner class.
int Initialize(int a_nbThreads)
Instanciate a number of RNG according to the number of threads used in openMP.
int Initialize()
A function used to initialize all that is needed.
Singleton class that manages output writing on disk (images, sinograms, etc). It also manages loggi...
void SetFrames(const string &a_frameList)
Set the frame list (a string that will be parsed by the InitializeFramingAndQuantification function) ...
void SetFOVSizeX(FLTNB a_fovSizeX)
Set the FOV's size along the X axis, in mm.
void SetDataType(int a_dataType)
set the data type
Definition: vDataFile.hh:366
int ReadDataASCIIFile(const string &a_file, const string &a_keyword, T *ap_return, int a_nbElts, bool a_mandatoryFlag)
Look for "a_nbElts" elts in the "a_file" file matching the "a_keyword" string passed as parameter a...
Definition: gOptions.cc:111
Singleton class that Instantiate and initialize the scanner object.
Declaration of class sScannerManager.
void SetVoxSizeY(FLTNB a_voxSizeY)
Set the voxel's size along the Y axis, in mm.
void ShowHelpComputation()
Display command line options related to the computation settings for castor-proj. ...
Definition: castor-proj.cc:211
virtual int PrepareDataFile()=0
This function is implemented in child classes Store different kind of information inside arrays (da...
int PROJ_SetPETSpecificParameters(int a_maxRingDiff)
Deliver to the PET scanner object all informations provided from the datafile header.
Inherit from vDataFile. Class that manages the reading of a SPECT input file (header + data)...
int main(int argc, char **argv)
Definition: castor-proj.cc:400
This class is designed to manage the different image convolvers and to apply them.
void ShowHelpInput()
Display command line options related to input settings for castor-proj.
Definition: castor-proj.cc:102
void SetOffsetX(FLTNB a_offsetX)
Set the image offset along the X axis, in mm.
int Launch()
Just call either the LaunchCPU or the LaunchGPU function as asked for.
vScanner * GetScannerObject()
void SetVoxSizeX(FLTNB a_voxSizeX)
Set the voxel's size along the X axis, in mm.
void IntfKeyInitFields(Intf_fields *ap_IF)
Init the file of an Interfile fields structure passed in parameter to their default values...
int CheckParameters()
A function used to check the parameters settings.
vector< FLTNB > image_duration
uint32_t mtx_size[7]
void SetVerbose(int a_verboseLevel)
Set Verbosity.
void SetVerbose(int a_verboseLevel)
Set the verbose level.
void SetDataFile(vDataFile *ap_DataFile)
Set a data file in use to later recover some information from it.
void SetImageDimensionsAndQuantification(oImageDimensionsAndQuantification *ap_ImageDimensionsAndQuantification)
set the pointer to the oImageDimensionsAndQuantification object
Definition: oImageSpace.hh:553
void ShowHelpProjector()
Show help about all implemented projectors.
void ShowHelpOutput()
Display command line options related to output settings for castor-proj.
Definition: castor-proj.cc:132
#define KEYWORD_MANDATORY
Definition: gOptions.hh:25
int PROJ_SetSPECTSpecificParameters(uint16_t *ap_nbOfBins, uint32_t a_nbOfProjections, FLTNB a_firstAngle, FLTNB a_stepAngle, FLTNB *ap_projectionAngles, FLTNB a_CORtoDetectorDistance, string a_rotDirection)
Deliver to the SPECT scanner object all informations provided from the acquisition parameters...
static void ShowCommonHelp()
This function does not take any parameter and is used to display some help about the syntax of the op...
Singleton class that generate a thread-safe random generator number for openMP As singleton...
void SetGPUflag(bool a_flagGPU)
Set the GPU flag.
void ShowHelp(int a_returnCode, int a_mpiRank)
Definition: castor-proj.cc:31
uint16_t nb_time_frames
int Initialize()
Initialization : .
Declaration of class sRandomNumberGenerator.
#define SCANNER_SPECT_PINHOLE
void SetVerbose(int a_verbose)
Set the member m_verboseLevel to the provided value.
#define OS_SEP
void SetVoxSizeZ(FLTNB a_voxSizeZ)
Set the voxel's size along the Z axis, in mm.
void ShowHelpMiscellaneous()
Display command line settings related to SPECT modality for castor-proj.
Definition: castor-proj.cc:364
void SetVerbose(int a_verboseLevel)
set verbosity
Definition: oImageSpace.hh:546
This class is designed to manage the projection part of the reconstruction.
void ShowHelpPET()
Definition: castor-proj.cc:339
Interfile fields. This structure contains all the Interfile keys currently managed by CASToR Decl...
void SetBedIndex(int a_bedIndex)
set the bed index corresponding to this data file
Definition: vDataFile.hh:373
Declaration of class sOutputManager.
This class holds all the matrices in the image domain that can be used in the algorithm: image...
Definition: oImageSpace.hh:41
int InitOutputDirectory(const string &a_pathFout, const string &a_pathDout)
Create the output directory if any, extract the base name and create the log file.
This class is designed to manage all dimensions and quantification related stuff. ...
int GetNbThreadsForProjection()
Get the number of threads used for projections.
Declaration of class oAnalyticProjection.
int GetScannerLayerNbRings(int a_layer)
Ask the number of rings to the scanner object for a specific layer. Returns an error if this inform...
void SetNbBeds(int a_nbBeds)
Set the number of bed positions.
This file is used for all kind of different functions designed for options parsing and ASCII file rea...
void SetVerbose(int a_verboseLevel)
set verbosity
Declaration of class iDataFileTransmission.
#define Cout(MESSAGE)
void SetScanner(vScanner *ap_Scanner)
Set the scanner in use.
void SetImageSpace(oImageSpace *ap_ImageSpace)
Set the Image Space Object.
virtual int PROJ_GetScannerSpecificParameters()=0
This function is implemented in child classes It is used to set several variables of the datafile w...
#define CASTOR_VERSION
Definition: gVariables.hh:44
void SetScanner(vScanner *ap_Scanner)
Set the scanner in use.
void SetImageDimensionsAndQuantification(oImageDimensionsAndQuantification *ap_ImageDimensionsAndQuantification)
set the pointer to the oImageDimensionsAndQuantification object
Definition: vDataFile.hh:408
Inherit from vDataFile. Class that manages the reading of a PET input file (header + data)...
Definition: iDataFilePET.hh:26
void SetPercentageLoad(int a_percentageLoad)
Set the percentage of the data file that will be loaded in memory.
Definition: vDataFile.hh:380
int ReadStringOption(const string &a_input, T *ap_return, int a_nbElts, const string &sep, const string &a_option)
Parse the 'a_input' string corresponding to the 'a_option' into 'a_nbElts' elements, using the 'sep' separator. The results are returned in the templated 'ap_return' dynamic templated array. Call "ConvertFromString()" to perform the correct conversion depending on the type of the data to convert.
Definition: gOptions.cc:50
void InitOptimizer(oImageDimensionsAndQuantification *ap_ImageDimensionsAndQuantification)
virtual int PROJ_InitFile()=0
This function is implemented in child classes Initialize the fstream objets for output writing as w...
void GetUserEndianness()
Check user/host computer endianness and write it to the global variable User_Endianness.
int Initialize()
A function used to initialize the manager and all image convolvers it manages.
void ShowHelpImageConvolver()
Show help about all implemented image convolvers.