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