CASToR  1.1
Tomographic Reconstruction (PET/SPECT)
 All Classes Files Functions Variables Typedefs Macros Groups Pages
castor-datafileConversionEx.cc
Go to the documentation of this file.
1 
16 #include "gVariables.hh"
17 #include "gOptions.hh"
18 #include "iDataFilePET.hh"
19 #include "sOutputManager.hh"
20 #include "sScannerManager.hh"
21 #include <iomanip> //std::setprecision
22 
23 #ifdef CASTOR_ROOT
24  #include "TROOT.h"
25  #include "TApplication.h"
26  #include "TGClient.h"
27  #include "TCanvas.h"
28  #include "TSystem.h"
29  #include "TTree.h"
30  #include "TBranch.h"
31  #include "TFile.h"
32 #endif
33 
38 void ShowHelp(int a_returnCode)
39 {
40  cout << endl;
41  cout << "NOTE : This program is an example code providing guidance regarding how to generate a PET CASToR datafile from any system datafile " << endl;
42  cout << " It does NOT perform any conversion by itself and must be adjusted to the conversion of any system dataset !" << endl;
43  cout << endl;
44  cout << endl;
45  cout << "Usage: castor-datafileConversionEx -ih input_datafile : (if datafile is histogram mode)" << endl;
46  cout << " -il input_datafile : (if datafile is list-mode)" << endl;
47  cout << " -o output_datafile" << endl;
48  cout << " -s scanner_alias "<< endl;
49  cout << " [-nc norm_factors_file" << endl;
50  cout << " [-sc scat_factors_file" << endl;
51  cout << " [-rc rdm_factors_file" << endl;
52  cout << " [-ac atn_factors_file" << endl;
53  cout << " [-cf calibration factor" << endl;
54  cout << " [-ist isotope_alias" << endl;
55  //cout << " [-poi flag for POI]" << endl;
56  //cout << " [-tof flag for TOF]" << endl;
57  cout << " [-vb verbosity]" << endl;
58  cout << endl;
59  cout << endl;
60  cout << "[Main settings]:" << endl;
61  cout << " -ih path_to_histo_datafile : provide an input histogram datafile to convert" << endl;
62  cout << " -il path_to_listm_datafile : provide an input list-mode datafile to convert" << endl;
63  cout << " -o path_to_cstr_file.cdh : provide the path to the output file will be created inside this folder (no default)" << endl;
64  cout << " -s scanner_alias : provide the name of the scanner used for to acquire the original data" << endl;
65  cout << " : Must correspond to a .geom or .hscan file in the config/scanner repository." << endl;
66  cout << endl;
67  cout << "[Optional settings]:" << endl;
68  cout << " -nc norm_factors_file : provide a file containing normalization correction factors" << endl;
69  cout << " -sc scat_factors_file : provide a file containing scatter correction factors" << endl;
70  cout << " -rc rdm_factors_file : provide a file containing random correction factors" << endl;
71  cout << " -ac atn_factors_file : provide a file containing attenuation correction factors" << endl;
72  cout << " -cf calibration factor : provide a calibration factor" << endl;
73  cout << " -ist isotope_alias : provide alias of the isotope used in the input datafile"<< endl;
74  cout << " Supported PET isotopes and their parameters are listed in config/misc/isotopes_pet"<< endl;
75  cout << " Other isotopes can be added in the same file"<< endl;
76  cout << endl;
77  cout << "[Miscellaneous settings]:" << endl;
78  cout << " -vb : give the verbosity level, from 0 (no verbose) to above 5 (at the event level) (default: 1)." << endl;
79  cout << endl;
80  #ifdef CASTOR_VERSION
81  cout << " This program is part of the CASToR release version " << CASTOR_VERSION << "." << endl;
82  cout << endl;
83  #endif
84  Exit(a_returnCode);
85 }
86 
87 
88 /*
89  Main program
90 */
91 int main(int argc, char** argv)
92 {
93  // No argument, then show help
94  if (argc==1) ShowHelp(0);
95 
96  // ============================================================================================================
97  // Parameterized variables with their default values
98  // ============================================================================================================
99 
100  // String containing the path to the input data filename provided by the user
101  string path_to_input_file = "";
102 
103  // Datafile mode (histogram or list (default)
104  int data_mode = MODE_LIST;
105 
106  // Path to user provided data
107  string path_to_out_file = "";
108  string path_to_data_file = "";
109  string path_to_header_file = "";
110  string path_to_norm_file = "";
111  string path_to_scat_file = "";
112  string path_to_rdm_file = "";
113  string path_to_acf_file = "";
114  FLTNB calib_factor = 1.;
115 
116  string scanner_alias = "";
117  string istp_alias = "";
118 
119  // Verbosity
120  int vb = 1;
121 
122  // Input datafile is histogram
123  bool is_histogram_flag = 0;
124 
125 
126  // ============================================================================================================
127  // Read command-line parameters
128  // ============================================================================================================
129 
130  for (int i=1; i<argc; i++)
131  {
132  string option = (string)argv[i];
133 
134  // Provide an histogram datafile to convert
135  if (option=="-ih")
136  {
137  if (!path_to_input_file.empty())
138  {
139  Cerr("***** castor-datafileConversionEx :: the following file names have already been provided (either -ih or -il option must be used): " << endl);
140  Cout(path_to_input_file << endl);
141  Exit(EXIT_FAILURE);
142  }
143  else
144  {
145  if (argv[i+1] == NULL)
146  {
147  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
148  Exit(EXIT_FAILURE);
149  }
150  else
151  path_to_input_file = (string)argv[i+1];
152 
153  data_mode = MODE_HISTOGRAM;
154  i++;
155  }
156  }
157 
158  // Provide a list-mode datafile to convert
159  if (option=="-il")
160  {
161  if (!path_to_input_file.empty())
162  {
163  Cerr("***** castor-datafileConversionEx :: the following file names have already been provided (either -ih or -il option must be used): " << endl);
164  Cout(path_to_input_file << endl);
165  Exit(EXIT_FAILURE);
166  }
167  else
168  {
169  if (argv[i+1] == NULL)
170  {
171  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
172  Exit(EXIT_FAILURE);
173  }
174  else
175  path_to_input_file = (string)argv[i+1];
176 
177  data_mode = MODE_LIST;
178  i++;
179  }
180  }
181 
182  // Output CASToR datafile
183  else if (option=="-o")
184  {
185  if (argv[i+1] == NULL)
186  {
187  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
188  Exit(EXIT_FAILURE);
189  }
190  else
191  path_to_out_file = (string)argv[i+1];
192  i++;
193  }
194 
195  // Scanner alias
196  else if (option=="-s")
197  {
198  if (argv[i+1] == NULL)
199  {
200  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
201  Exit(EXIT_FAILURE);
202  }
203  else
204  scanner_alias = (string)argv[i+1];
205  i++;
206  }
207 
208  // Provide a norm correction file
209  else if (option=="-nc")
210  {
211  if (argv[i+1] == NULL)
212  {
213  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
214  Exit(EXIT_FAILURE);
215  }
216  else
217  path_to_norm_file = (string)argv[i+1];
218 
219  i++;
220  }
221 
222  // Provide a scatter correction file
223  else if (option=="-sc")
224  {
225  if (argv[i+1] == NULL)
226  {
227  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
228  Exit(EXIT_FAILURE);
229  }
230  else
231  path_to_scat_file = (string)argv[i+1];
232 
233  i++;
234  }
235 
236  // Provide a rdm correction file
237  else if (option=="-rc")
238  {
239  if (argv[i+1] == NULL)
240  {
241  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
242  Exit(EXIT_FAILURE);
243  }
244  else
245  path_to_rdm_file = (string)argv[i+1];
246 
247  i++;
248  }
249 
250  // Provide an acf correction file
251  else if (option=="-ac")
252  {
253  if (argv[i+1] == NULL)
254  {
255  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
256  Exit(EXIT_FAILURE);
257  }
258  else
259  path_to_acf_file = (string)argv[i+1];
260 
261  i++;
262  }
263 
264  // Provide a calibration factor
265  else if (option=="-cf")
266  {
267  if (argv[i+1] == NULL)
268  {
269  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
270  Exit(EXIT_FAILURE);
271  }
272  else
273  {
274  string val_str = (string)argv[i+1];
275  if(ConvertFromString(val_str, &calib_factor) )
276  {
277  Cerr("***** castor-datafileConversionEx :: Exception when trying to read'" << val_str << " for option: " << option << endl);
278  Exit(EXIT_FAILURE);
279  }
280  }
281  i++;
282  }
283 
284  // Provide an isotope alias
285  else if (option=="-ist")
286  {
287  if (argv[i+1] == NULL)
288  {
289  Cerr("***** castor-datafileConversionEx :: argument missing for option: " << option << endl);
290  Exit(EXIT_FAILURE);
291  }
292  else
293  istp_alias = (string)argv[i+1];
294 
295  i++;
296  }
297 
298  // Verbosity level
299  else if (option=="-vb")
300  {
301  if (i>=argc-1)
302  {
303  Cerr("***** castor-datafileConversionEx :: Argument missing for option: " << option << endl);
304  Exit(EXIT_FAILURE);
305  }
306  vb = atoi(argv[i+1]);
307  i++;
308  }
309 
310  else
311  {
312  Cerr("***** castor-datafileConversionEx :: Unknown option '" << option << "' !" << endl);
313  Exit(EXIT_FAILURE);
314  }
315  }
316 
317 
318  // ============================================================================================================
319  // Basic initialization checks (minimal initializations mandatory for the next steps)
320  // ============================================================================================================
321 
322  // data files
323  if (path_to_input_file.empty() )
324  {
325  Cerr("***** castor-datafileConversionEx :: Please provide at least one data filename (-ih for histogram, -il for list-mode)" << endl);
326  cout << " -input filename : give an input datafile (no default)." << endl;
327  Exit(EXIT_FAILURE);
328  }
329 
330  // output directory
331  if (path_to_out_file.empty() )
332  {
333  Cerr("***** castor-datafileConversionEx :: Please provide the output file name (-o)" << endl);
334  Exit(EXIT_FAILURE);
335  }
336 
337  // scanner
338  if (scanner_alias.empty())
339  {
340  Cerr("***** castor-datafileConversionEx :: Please provide a scanner alias (-s) :" << endl);
341  Cout(" It must correspond to a .geom or .hscan (user-made LUT) file in the config/scanner directory." << endl);
342  ShowHelp(0);
343  Exit(EXIT_FAILURE);
344  }
345 
346  // Throw error by default.
347  // These two lines must be erased when implementing any datafile converter!
348  //Cerr("castor-GATERootToCastor :: This script is only for demonstration purposes " << endl);
349  //Exit(EXIT_FAILURE);
350 
351 
352  // ============================================================================================================
353  // SOutputManager object initialisation:
354  // (Used for log file)
355  // ============================================================================================================
356 
357  sOutputManager* p_outputManager = sOutputManager::GetInstance();
358  // Set verbose level
359  p_outputManager->SetVerbose(vb);
360  // Set MPI rank
361  p_outputManager->SetMPIRank(0);
362 
363  // Set path to the config directory
364  if (p_outputManager->CheckConfigDir(""))
365  {
366  Cerr("***** castor-datafileConversionEx :: A problem occured while checking for the config directory path !" << endl);
367  Exit(EXIT_FAILURE);
368  }
369  // Initialize output directory and base name
370  if (p_outputManager->InitOutputDirectory(path_to_out_file, ""))
371  {
372  Cerr("***** castor-datafileConversionEx :: A problem occured while initializing output directory !" << endl);
373 
374  Exit(EXIT_FAILURE);
375  }
376  // Log command line
377  if (p_outputManager->LogCommandLine(argc,argv))
378  {
379  Cerr("***** castor-datafileConversionEx :: A problem occured while logging command line arguments !" << endl);
380  Exit(EXIT_FAILURE);
381  }
382 
383  // Output progression options
384  cout << std::fixed;
385  cout << std::setprecision(2);
386 
387  // ============================================================================================================
388  // ScannerManager object initialization
389  // (Basic initialization to locate scanner geometry files)
390  // ============================================================================================================
391 
392  sScannerManager* p_scannermanager = sScannerManager::GetInstance();
393  p_scannermanager->SetVerbose(vb);
394 
395  // Check if the scanner exists and get the name from ScannerManager
396  scanner_alias = (scanner_alias.find(OS_SEP)) ?
397  scanner_alias.substr(scanner_alias.find_last_of(OS_SEP)+1) :
398  scanner_alias;
399 
400  if(p_scannermanager->FindScannerSystem(scanner_alias) )
401  {
402  Cerr("**** castor-datafileConversionEx :: A problem occurred while searching for scanner system !" << endl);
403  Exit(EXIT_FAILURE);
404  }
405 
406  // Build output file names
407  path_to_data_file = path_to_out_file + ".Cdf";
408  path_to_header_file = path_to_out_file + ".Cdh";
409 
410  // ============================================================================================================
411  // Instanciate/Initialize CASToR Datafile
412  // (Enable/Disable optionnal field to be written)
413  // ============================================================================================================
414 
415  // Instantiate & Initialize oImageDimensionsAndQuantification object, required for datafile generation (number of threads)
417 
418  // Instanciate & Initialize iDataFilePET and Event objects
419  iDataFilePET* Out_data_file = NULL;
420  vEvent* Event = NULL;
421 
422  // ===> This variable set the maximum number of LORs contained in one event
423  // ===> in the whole dataset (for compression purposes)
424  // ===> This depends on the dataset to convert
425  uint16_t max_nb_lines_per_event = 1;
426 
427  // Create datafile and set variables
428  Out_data_file = new iDataFilePET();
429 
430  // Set the image oImageDimensionsAndQuantification object (required for initialization)
431  Out_data_file->SetImageDimensionsAndQuantification(p_ID);
432 
433  // Set path to header datafile
434  Out_data_file->SetHeaderDataFileName(path_to_header_file);
435 
436  // Just required for initialization
437  Out_data_file->SetPercentageLoad(0);
438 
439  // Set verbosity
440  Out_data_file->SetVerbose(vb);
441 
442  // Set data mode (list/hitogram)
443  Out_data_file->SetDataMode(data_mode);
444 
445  // Set modality
446  Out_data_file->SetDataType(TYPE_PET);
447 
448  // Set maximum number of LORs contained in one event
449  Out_data_file->SetMaxNumberOfLinesPerEvent(max_nb_lines_per_event);
450 
451  // Set calibration factor
452  Out_data_file->SetCalibrationFactor(calib_factor);
453 
454  // Set isotope
455  if(!istp_alias.empty())
456  {
457  if(p_ID->SetPETIsotope(0, istp_alias) )
458  {
459  Cerr("**** castor-datafileConversionEx :: A problem occurred while checking isotope name !" << endl);
460  Exit(EXIT_FAILURE);
461  }
462  Out_data_file->SetIsotope(istp_alias);
463  }
464 
465  // Allocate the Event structure (containing all information for a listmode
466  // or histogram event, depending on data mode)
467  if(data_mode == MODE_HISTOGRAM)
468  {
469  // Instanciate list-mode Event
470  Event = new iEventHistoPET();
471  Event->SetNbLines(max_nb_lines_per_event);
472  Event->AllocateID();
473  }
474  else
475  {
476  // Instanciate histogram Event
477  Event = new iEventListPET();
478  Event->SetNbLines(max_nb_lines_per_event);
479  Event->AllocateID();
480  }
481 
482 
483  // ============================================================================================================
484  // ===> Read sinogram correction file(s) (if any) and set datafile flag for each correction
485  // ============================================================================================================
486 
487  // ===> Compute sinogram dimensions.
488  // ===> Values must be adjusted to the dataset to convert.
489  uint32_t nbins = 1;
490  uint32_t nangles = 1;
491  uint32_t nsinos = 1;
492  uint32_t nb_cor_factors = nbins
493  * nangles
494  * nsinos;
495 
496  // --- Normalization correction factors ---
497 
498  // FLTNBDATA is either a float (default) or double type.
499  // Defined in /management/gVariables.hh
500  FLTNBDATA *p_norm = NULL;
501  bool norm_flag = false;
502 
503  if( !path_to_norm_file.empty() )
504  {
505  if(vb>=2) Cout("--> Reading the norm correction factors file = " << path_to_norm_file << "..." << endl);
506  ifstream ifile(path_to_norm_file.c_str(), ios::binary | ios::in);
507 
508  if (!ifile)
509  {
510  Cerr("***** castor-datafileConversionEx :: Error reading normalization factor file: "<< path_to_norm_file <<" !" << endl);
511  Exit(EXIT_FAILURE);
512  }
513  else
514  {
515  // Allocate buffer storing the normalisation
516  p_norm = new FLTNBDATA [nb_cor_factors * sizeof( FLTNBDATA )];
517  ifile.read(reinterpret_cast < char* > (p_norm), nb_cor_factors * sizeof( FLTNBDATA));
518  ifile.close();
519  }
520 
521  // ===> Set normalization correction flag in datafile objet
522  // ===> (required to indicate that normalization correction will be provided)
523  Out_data_file->SetNormCorrectionFlagOn();
524  norm_flag = true;
525  }
526 
527  // --- Scatter correction factors file ---
528  FLTNBDATA *p_scat = NULL;
529  bool scat_flag = false;
530 
531  if( !path_to_scat_file.empty() )
532  {
533  if(vb>=2) Cout("--> Reading the scatter correction factors file = " << path_to_scat_file << "..." << endl);
534  ifstream ifile(path_to_scat_file.c_str(), ios::binary | ios::in);
535 
536  if (!ifile)
537  {
538  Cerr("***** castor-datafileConversionEx :: Error reading scatter factor file: "<< path_to_scat_file <<" !" << endl);
539  Exit(EXIT_FAILURE);
540  }
541  else
542  {
543  // Allocating buffer storing the scatter
544  p_scat = new FLTNBDATA [nb_cor_factors * sizeof( FLTNBDATA )];
545  ifile.read(reinterpret_cast < char* > (p_scat), nb_cor_factors * sizeof( FLTNBDATA));
546  ifile.close();
547  }
548 
549  // ===> Set scatter correction flag in datafile objet
550  // ===> (required to indicate that scatter correction will be provided)
551  Out_data_file->SetScatterCorrectionFlagOn();
552  scat_flag = true;
553  }
554 
555  // --- Random correction factors file ---
556  FLTNBDATA *p_rdm = NULL;
557  bool rdm_flag = false;
558 
559  if( !path_to_rdm_file.empty() )
560  {
561  if(vb>=2) Cout("--> Reading the random correction factors file = " << path_to_rdm_file << "..." << endl);
562  ifstream ifile(path_to_rdm_file.c_str(), ios::binary | ios::in);
563 
564  if (!ifile)
565  {
566  Cerr("***** castor-datafileConversionEx :: Error reading random correction factor file: "<< path_to_rdm_file <<" !" << endl);
567  Exit(EXIT_FAILURE);
568  }
569  else
570  {
571  // Allocating buffer storing the random
572  p_rdm = new FLTNBDATA [nb_cor_factors * sizeof( FLTNBDATA )];
573  ifile.read(reinterpret_cast < char* > (p_rdm), nb_cor_factors * sizeof( FLTNBDATA));
574  ifile.close();
575  }
576 
577  // ===> Set random correction flag in datafile objet
578  // ===> (required to indicate that random correction will be provided)
579  Out_data_file->SetRandomCorrectionFlagOn();
580  rdm_flag = true;
581  }
582 
583  // --- Attenuation correction factors file ---
584  FLTNBDATA *p_acf = NULL;
585  bool acf_flag = false;
586 
587  if( !path_to_acf_file.empty() )
588  {
589  if(vb>=2) Cout("--> Reading the attenuation correction factors file = " << path_to_acf_file << "..." << endl);
590  ifstream ifile(path_to_acf_file.c_str(), ios::binary | ios::in);
591 
592  if (!ifile)
593  {
594  Cerr("***** castor-datafileConversionEx :: Error reading attenuation correction factor file: "<< path_to_acf_file <<" !" << endl);
595  Exit(EXIT_FAILURE);
596  }
597  else
598  {
599  // Allocating buffer storing the acf
600  p_acf = new FLTNBDATA [nb_cor_factors * sizeof( FLTNBDATA )];
601  ifile.read(reinterpret_cast < char* > (p_acf), nb_cor_factors* sizeof( FLTNBDATA));
602  ifile.close();
603  }
604 
605  // ===> Set attenuation correction flag in datafile objet
606  // ===> (required to indicate that attenuation correction will be provided)
607  Out_data_file->SetAtnCorrectionFlagOn();
608  acf_flag = true;
609  }
610 
611  // Initialize datafile object
612  Out_data_file->PROJ_InitFile();
613 
614  // Compute the event size relatively to the mandatory
615  // and optionnal (correction factors) elements
616  Out_data_file->ComputeSizeEvent();
617 
618  // Init Datafile
619  Out_data_file->PrepareDataFile();
620 
621 
622  // ============================================================================================================
623  // *********************************************** CONVERSION *************************************************
624  // ============================================================================================================
625  Cout(" --- Start conversion of datafile : " << path_to_input_file << endl
626  << " CASToR output header datafile: " << path_to_header_file << endl
627  << " CASToR output binary datafile: " << path_to_data_file << endl << endl);
628 
629 
630  // ============================================================================================================
631  // Init variables
632  // ============================================================================================================
633 
634  // Index for progression printing
635  uint32_t printing_index = 0;
636 
637  // ===> Recover the total number of events (number of coincidences or
638  // ===> histogram bins) from the dataset metadata
639  uint32_t nEvents = 0;
640 
641 
642  // ===> Acquisition duration and start time in seconds
643  // ===> This should be recovered from the dataset metadata)
644  FLTNB start_time_sec = 0.;
645  FLTNB stop_time_sec = 1.;
646  FLTNB duration_sec = stop_time_sec - start_time_sec;
647 
648  // ===> Scanner variables (blocs, modules, submodules, crystals etc..)
649  // ===> depending on the system.
650  uint32_t nRsectorsPerRing = 70;
651  uint32_t nb_crystals_trs = 9,
652  nb_crystals_axl = 6;
653  uint32_t nb_crystals_per_ring = nRsectorsPerRing
654  * nb_crystals_trs;
655 
656  // ===> Variables to recover indices of the two crystals for each potential
657  // ===> LORs contributing to the event
658  uint32_t *castor_id1 = new uint32_t[max_nb_lines_per_event];
659  uint32_t *castor_id2 = new uint32_t[max_nb_lines_per_event];
660 
661  // --------------------------------
662  // Initialization of root variables
663  // The following lines are just required for the demonstration dataset used in this sample code,
664  // Must be erased when adjusting this code to the conversion of any dataset
665 
666  // ROOT objects/variables
667  int32_t rSectorID1, rSectorID2;
668  int32_t moduleID1, moduleID2;
669  int32_t crystalID1, crystalID2;
670 
671  #ifdef CASTOR_ROOT
672  TApplication *Tapp = new TApplication("tapp",0,0);
673  TTree* Coincidences = new TTree;
674  TFile* Tfile_root = new TFile(path_to_input_file.c_str(),"READ","ROOT file with histograms");
675  Coincidences = (TTree*)Tfile_root->Get("Coincidences");
676  nEvents += Coincidences->GetEntries();
677 
678  //------------- ROOT branches
679  double_t time1=0;
680  Coincidences->SetBranchAddress("time1",&time1);
681  Coincidences->SetBranchAddress("rsectorID1",&rSectorID1);
682  Coincidences->SetBranchAddress("moduleID1",&moduleID1);
683  Coincidences->SetBranchAddress("crystalID1",&crystalID1);
684 
685  double_t time2=0;
686  Coincidences->SetBranchAddress("time2",&time2);
687  Coincidences->SetBranchAddress("rsectorID2",&rSectorID2);
688  Coincidences->SetBranchAddress("moduleID2",&moduleID2);
689  Coincidences->SetBranchAddress("crystalID2",&crystalID2);
690 
691  // Recover acquisition start time from the first event
692  Coincidences->GetEntry(0);
693  #endif
694  // --------------------------------
695 
696  // Just for progression output
697  uint32_t printing_ratio = (nEvents>1000) ? 1000 : nEvents/10;
698 
699  // Loop on the data elements (histogram bins or LORs)
700  for (uint32_t i=0; i<nEvents ; i++)
701  {
702  // ===> Recover event value (list-mode, or dynamic histogram dataset)
703  uint32_t time_ms = 0;
704  stop_time_sec = (FLTNB)time_ms/1000;
705 
706  // ===> Recover the number of lines in the event (if compression).
707  int nb_lines_in_event = 1;
708 
709  // -----------------------------------------------------------------
710  // The following lines are just required for the demonstration dataset used in this sample code.
711  // Must be erased/adapted when adjusting this code to the conversion of any dataset
712  #ifdef CASTOR_ROOT
713  Coincidences->GetEntry(i);
714  time_ms = time1;
715  #endif
716 
717  // ===> Compute CASToR indexation depending on the system layout
718  // ===> This directly depends on the type of scanner file used to
719  // ===> described the system geometry.
720  // ===> See documentation for more information
721  uint32_t crystals_trs_id1 = (uint32_t)(crystalID1/nb_crystals_trs);
722 
723  uint32_t ringID1 = moduleID1*nb_crystals_axl
724  + crystals_trs_id1;
725 
726  uint32_t crystals_trs_id2 = (uint32_t)(crystalID2/nb_crystals_trs);
727 
728  uint32_t ringID2 = moduleID2*nb_crystals_axl
729  + crystals_trs_id2;
730 
731  // Loop on the LORs contained in the event
732  for(int l=0 ; l<nb_lines_in_event ; l++)
733  {
734  // Recover the CASToR ID position of the two crystals according to
735  // the LUT corresponding to this scanner
736  castor_id1[l] = ringID1*nb_crystals_per_ring
737  + rSectorID1*nb_crystals_trs
738  + crystalID1%nb_crystals_trs ;
739 
740  castor_id2[l] = ringID2*nb_crystals_per_ring
741  + rSectorID2*nb_crystals_trs
742  + crystalID2%nb_crystals_trs ;
743  }
744  // -----------------------------------------------------------------
745 
746  // --- Histo-mode event ---
747  if(is_histogram_flag) // histogram mode event
748  {
749  // Set the number of lines in this event
750  Event->SetNbLines(nb_lines_in_event);
751 
752  // Set CASToR IDs in the event structure for each contributing line
753  for(int l=0 ; l<nb_lines_in_event ; l++)
754  {
755  Event->SetID1(l, castor_id1[l]); // 1st argument : line, 2nd argument : index
756  Event->SetID2(l, castor_id2[l]); // 1st argument : line, 2nd argument : index
757  }
758 
759  // Set information related to each TOF bin
760  int nb_tof_bins = 1;
761  for (int b=0; b<nb_tof_bins; b++)
762  {
763  // ===> Recover and set event value (histogram dataset)
764  FLTNB event_value = 1.;
765  Event->SetEventValue(b, event_value); // 1st argument : TOF bin, 2nd argument : event value
766 
767  // ===> Set scatter correction factor (if any)
768  // ===> Require computation of the scatter correction value for
769  // ===> this LOR from the factors loaded in the "p_scat" buffer
770  if(scat_flag)
771  {
772  FLTNB scat_correction_coeff = 0.; // computed from p_scat
773  ((iEventHistoPET*)Event)->SetScatterRate(b, scat_correction_coeff);
774  }
775  }
776 
777  // ===> Set timestamp of the event in ms (global timestamp for histogram)
778  // ===> i.e : 0, or timestamp of a dynamic frame
779  Event->SetTimeInMs(time_ms); // uint32_t
780 
781  // ===> Set random correction factor (if any)
782  // ===> Require computation of the random correction value for
783  // ===> this LOR from the factors loaded in the "p_rdm" buffer
784  if(rdm_flag)
785  {
786  FLTNB rdm_correction_coeff = 0.; // computed from p_rdm
787  ((iEventPET*)Event)->SetRandomRate(rdm_correction_coeff);
788  }
789 
790  // ===> Set normalization correction factor (if any)
791  // ===> Require computation of the norm correction value for
792  // ===> this LOR from the factors loaded in the "p_norm" buffer
793  if(norm_flag)
794  {
795  FLTNB norm_correction_coeff = 0.; // computed from p_norm
796  ((iEventPET*)Event)->SetNormalizationFactor(norm_correction_coeff);
797  }
798 
799  // ===> Set attenuation correction factor (if any)
800  // ===> Require computation of the atn correction value for
801  // ===> this LOR from the factors loaded in the "p_acf" buffer
802  if(acf_flag)
803  {
804  FLTNB atn_correction_coeff = 0.; // computed from p_acf
805  ((iEventPET*)Event)->SetAttenuationCorrectionFactor(atn_correction_coeff);
806  }
807 
808  // Write event in the datafile
809  // ('0' argument is just thread-related. No change required)
810  Out_data_file->PROJ_WriteEvent(Event, 0);
811  }
812 
813  // --- List-mode event ---
814  else
815  {
816  // Set the number of lines in this event
817  Event->SetNbLines(nb_lines_in_event);
818 
819  // Set CASToR ID in the event structure
820  for(int l=0 ; l<nb_lines_in_event ; l++)
821  {
822  Event->SetID1(l, castor_id1[l]); // 1st argument : line, 2nd argument :index
823  Event->SetID2(l, castor_id2[l]); // 1st argument : line, 2nd argument :index
824  }
825 
826  // ===> Set scatter correction factor (if any)
827  // ===> Require computation of the scatter correction value for
828  // ===> this LOR from the factors loaded in the "p_scat" buffer
829  if(scat_flag)
830  {
831  FLTNB scat_correction_coeff = 0.; // computed from p_scat
832  // First argument '0' by default for list-mode (it is used for histograms with TOF bins).
833  ((iEventListPET*)Event)->SetScatterRate(0, scat_correction_coeff);
834  }
835 
836  // ===> Set timestamp of the event in ms
837  Event->SetTimeInMs(time_ms); // uint32_t
838 
839  // ===> Set random correction factor (if any)
840  // ===> Require computation of the random correction value for
841  // ===> this LOR from the factors loaded in the "p_rdm" buffer
842  if(rdm_flag)
843  {
844  FLTNB rdm_correction_coeff = 0.; // computed from p_rdm
845  ((iEventPET*)Event)->SetRandomRate(rdm_correction_coeff);
846  }
847 
848  // ===> Set normalization correction factor (if any)
849  // ===> Require computation of the norm correction value for
850  // ===> this LOR from the factors loaded in the "p_norm" buffer
851  if(norm_flag)
852  {
853  FLTNB norm_correction_coeff = 0.; // computed from p_norm
854  ((iEventPET*)Event)->SetNormalizationFactor(norm_correction_coeff);
855  }
856 
857  // ===> Set attenuation correction factor (if any)
858  // ===> Require computation of the atn correction value for
859  // ===> this LOR from the factors loaded in the "p_acf" buffer
860  if(acf_flag)
861  {
862  FLTNB atn_correction_coeff = 0.; // computed from p_acf
863  ((iEventPET*)Event)->SetAttenuationCorrectionFactor(atn_correction_coeff);
864  }
865 
866  // Write event in the datafile
867  // ('0' argument is just thread-related. No change required)
868  Out_data_file->PROJ_WriteEvent(Event, 0);
869  }
870 
871  // Progression output
872  if (printing_index%printing_ratio==0)
873  {
874  FLTNB percent = ( ((FLTNB)(printing_index+1))/((FLTNB)nEvents) ) * ((FLTNB)100);
875  cout << "\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b "
876  << percent << " \% ";
877  }
878  printing_index++;
879  }
880 
881  Cout("The simulated dataset contained " << nEvents << " coincidences / bins" << endl);
882 
883  // Set duration and nb of events in datafile
884  Out_data_file->SetStartTime(start_time_sec);
885  Out_data_file->SetDuration(duration_sec);
886  Out_data_file->SetNbEvents(nEvents);
887 
888  // Write datafile header
889  Out_data_file->PROJ_WriteHeader();
890 
891  // Write raw datafile
892  Out_data_file->PROJ_WriteData();
893 
894  // Delete tmp writing file (if any)
895  Out_data_file->PROJ_DeleteTmpDatafile();
896 
897  cout << endl;
898 
899  // ============================================================================================================
900  // End
901  // ============================================================================================================
902 
903  // Free Root objects
904  #ifdef CASTOR_ROOT
905  delete Coincidences;
906  delete Tfile_root;
907  delete Tapp;
908  #endif
909 
910  delete[] castor_id1;
911  delete[] castor_id2;
912 
913  if(p_acf != NULL) delete[] p_acf;
914  if(p_rdm != NULL) delete[] p_rdm;
915  if(p_scat != NULL) delete[] p_scat;
916  if(p_norm != NULL) delete[] p_norm;
917 
918  // Delete objects
919  if (Event)
920  {
921  if (is_histogram_flag)
922  delete (iEventHistoPET*)Event;
923  else
924  delete (iEventListPET*)Event;
925  }
926 
927  delete Out_data_file;
928  delete p_ID;
929 
930  return 0;
931 }
void SetDataMode(int a_dataMode)
set the data mode
Definition: vDataFile.hh:359
This header file is mainly used to declare some macro definitions and all includes needed from the st...
int SetPETIsotope(int a_bed, const string &a_isotope)
Set the PET isotope for the provided bed.
static sScannerManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
#define MODE_HISTOGRAM
Definition: vDataFile.hh:36
#define TYPE_PET
Definition: vDataFile.hh:51
void SetNbLines(uint16_t a_value)
Set the number of lines of the Event.
Definition: vEvent.hh:134
void SetVerbose(int a_verbose)
set verbosity
void SetAtnCorrectionFlagOn()
set to true the flag indicating the presence of attenuation correction factors in the datafile ...
#define MODE_LIST
Definition: vDataFile.hh:34
#define FLTNB
Definition: gVariables.hh:55
int FindScannerSystem(string a_scannerName)
Look for a file matching with the scanner name in parameter inside the scanner repository.
void SetScatterCorrectionFlagOn()
set to true the flag indicating the presence of scatter correction factors in the datafile ...
void SetVerbose(int a_verboseLevel)
set verbosity
Definition: vDataFile.hh:394
void SetDuration(FLTNB a_value)
Definition: vDataFile.hh:463
virtual void SetEventValue(int a_bin, FLTNBDATA a_value)=0
This function is implemented by child classes.
Definition: vEvent.hh:120
static sOutputManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
int ComputeSizeEvent()
Computation of the size of each event according to the mandatory/optional correction fields...
Definition: iDataFilePET.cc:88
Declaration of class iDataFilePET.
void Exit(int code)
int LogCommandLine(int argc, char **argv)
Write log file header with the provided command line options and different informations.
void SetIsotope(string a_value)
initialize the isotope string value
#define FLTNBDATA
Definition: gVariables.hh:59
int PROJ_WriteData()
Write/Merge chunk of data in a general data file.
Definition: vDataFile.cc:928
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
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...
Inherit from iEventPET. Class for PET list-mode events.
#define Cerr(MESSAGE)
Inherit from iEventPET. Class for PET histogram mode events.
int PROJ_WriteEvent(vEvent *ap_Event, int a_th)
Write event according to the chosen type of data.
void SetStartTime(FLTNB a_value)
Definition: vDataFile.hh:455
Singleton class that manages output writing on disk (images, sinograms, etc). It also manages loggi...
void SetDataType(int a_dataType)
set the data type
Definition: vDataFile.hh:366
Singleton class that Instantiate and initialize the scanner object.
void SetCalibrationFactor(FLTNB a_value)
initialize the global calibration factor with a FLTNB value
Definition: vDataFile.hh:439
void SetRandomCorrectionFlagOn()
set to true the flag indicating the presence of random correction factors in the datafile ...
Declaration of class sScannerManager.
void ShowHelp(int a_returnCode)
void SetHeaderDataFileName(const string &a_headerFileName)
set the data header file name
Definition: vDataFile.hh:431
void SetMPIRank(int a_mpiRank)
Initialize the machine index for MPI.
int AllocateID()
Instantiate the mp_ID1 and mp_ID2 indices arrays.
Definition: vEvent.cc:71
int PROJ_InitFile()
Initialize the fstream objets for output writing as well as some other variables specific to the Proj...
int PROJ_DeleteTmpDatafile()
Delete temporary datafile used for multithreaded output writing if needed.
Definition: vDataFile.cc:1008
int PROJ_WriteHeader()
Generate a header file according to the projection and data output informations. Used by Projection...
int main(int argc, char **argv)
Inherit from vEvent. Main PET class for the Event objects.
Definition: iEventPET.hh:22
void SetNbEvents(int64_t a_value)
initialize the number of events with a int64_t value
Definition: vDataFile.hh:447
int PrepareDataFile()
Store different kind of information inside arrays (data relative to specific correction as well as ba...
void SetID2(int a_line, uint32_t a_value)
Set the indice associated with the line index for the 2nd ID of the Event.
Definition: vEvent.hh:150
#define OS_SEP
void SetTimeInMs(uint32_t a_value)
Set the timestamp of the Event.
Definition: vEvent.hh:119
Declaration of class sOutputManager.
Mother class for the Event objects.
Definition: vEvent.hh:23
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.
void SetNormCorrectionFlagOn()
set to true the flag indicating the presence of normalization correction factors in the datafile ...
This class is designed to manage all dimensions and quantification related stuff. ...
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
#define Cout(MESSAGE)
#define CASTOR_VERSION
Definition: gVariables.hh:44
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
void SetMaxNumberOfLinesPerEvent(uint16_t a_value)
set the max number of line per event in the datafile
void SetID1(int a_line, uint32_t a_value)
Set the indice associated with the line index for the 1st ID of the Event.
Definition: vEvent.hh:142