CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
src/algorithm/oSensitivityGenerator.cc
Go to the documentation of this file.
1 
8 #include "gVariables.hh"
9 #include "oSensitivityGenerator.hh"
10 #include "iDataFilePET.hh"
11 
12 // =====================================================================
13 // ---------------------------------------------------------------------
14 // ---------------------------------------------------------------------
15 // =====================================================================
16 
18 {
19  // Simply default all members to "empty" values
22  mp_ImageSpace = NULL;
23  mp_Scanner = NULL;
24  mp_ProjectorManager = NULL;
26  mp_DeformationManager = NULL;
27  m2p_DataFile = NULL;
35  m_mumapAttenuationFlag = false;
39  mp_lineCounter = NULL;
40  m_flagGPU = false;
41  m_verbose = -1;
42  m_checked = false;
43  m_initialized = false;
44 }
45 
46 // =====================================================================
47 // ---------------------------------------------------------------------
48 // ---------------------------------------------------------------------
49 // =====================================================================
50 
52 {
53  // Delete all normalization data files if any
55  {
56  // Actual number of allocated bed positions
57  int actual_nb_beds = mp_ImageDimensionsAndQuantification->GetNbBeds();
58  if (m_oneNormalizationFileForAllBeds) actual_nb_beds = 1;
59  // Loop on beds
60  for (int bed=0; bed<actual_nb_beds; bed++) if (m3p_NormalizationDataFile[bed])
61  {
62  // Actual number of allocated frames
63  int actual_nb_frames = mp_ImageDimensionsAndQuantification->GetNbTimeFrames();
64  if (m_oneNormalizationFileForAllFrames) actual_nb_frames = 1;
65  for (int fr=0; fr<actual_nb_frames; fr++)
66  if (m3p_NormalizationDataFile[bed][fr]) delete m3p_NormalizationDataFile[bed][fr];
67  free(m3p_NormalizationDataFile[bed]);
68  }
70  }
71  if (mp_lineCounter) free(mp_lineCounter);
72 }
73 
74 // =====================================================================
75 // ---------------------------------------------------------------------
76 // ---------------------------------------------------------------------
77 // =====================================================================
78 
80 {
81  // Check all mandatory parameters
83  {
84  Cerr("***** oSensitivityGenerator::CheckParameters() -> Image dimensions and quantification object is null !" << endl);
85  return 1;
86  }
87  if (mp_ImageSpace==NULL)
88  {
89  Cerr("***** oSensitivityGenerator::CheckParameters() -> Image space object is null !" << endl);
90  return 1;
91  }
92  if (mp_Scanner==NULL)
93  {
94  Cerr("***** oSensitivityGenerator::CheckParameters() -> Scanner object is null !" << endl);
95  return 1;
96  }
97  if (mp_ProjectorManager==NULL)
98  {
99  Cerr("***** oSensitivityGenerator::CheckParameters() -> Projector manager object is null !" << endl);
100  return 1;
101  }
102  if (mp_ImageConvolverManager==NULL)
103  {
104  Cerr("***** oSensitivityGenerator::CheckParameters() -> Convolver manager object is null !" << endl);
105  return 1;
106  }
107  if (mp_DeformationManager==NULL)
108  {
109  Cerr("***** oSensitivityGenerator::CheckParameters() -> Deformation manager object is null !" << endl);
110  return 1;
111  }
112  if (m2p_DataFile==NULL)
113  {
114  Cerr("***** oSensitivityGenerator::CheckParameters() -> Data file array is null !" << endl);
115  return 1;
116  }
117  for (int b=0; b<mp_ImageDimensionsAndQuantification->GetNbBeds(); b++)
118  {
119  if (m2p_DataFile[b]==NULL)
120  {
121  Cerr("***** oSensitivityGenerator::CheckParameters() -> Data file object for bed " << b+1 << " is null !" << endl);
122  return 1;
123  }
124  }
125  if (m_verbose<0)
126  {
127  Cerr("***** oSensitivityGenerator::CheckParameters() -> Verbose level is negative !" << endl);
128  return 1;
129  }
130  // If compute the sensitivity from a histogram datafile, then check if the datafile is actually an histogram
131  if (m_computeFromHistogramFlag && m2p_DataFile[0]->GetDataMode()!=MODE_HISTOGRAM)
132  {
133  Cerr("***** oSensitivityGenerator::CheckParameters() -> It was asked to compute global sensitivity from the provided datafile whereas it is not a histogram !" << endl);
134  return 1;
135  }
136  // Sensitivity computation from histogram do not work with dynamism for the moment
140  {
141  Cerr("***** oSensitivityGenerator::CheckParameters() -> Global sensitivity computation from histogram does not work wth dynamic data yet !" << endl);
142  return 1;
143  }
144  // Now it is checked
145  m_checked = true;
146  // End
147  return 0;
148 }
149 
150 // =====================================================================
151 // ---------------------------------------------------------------------
152 // ---------------------------------------------------------------------
153 // =====================================================================
154 
156 {
157  // First check that the parameters has been checked !
158  if (!m_checked)
159  {
160  Cerr("***** oSensitivityGenerator::Initialize() -> Must call the CheckParameters() function before initialization !" << endl);
161  return 1;
162  }
163  // Verbose
164  if (m_verbose>=2) Cout("oSensitivityGenerator::Initialize() -> Start initialization" << endl);
165  // For the moment, SPECT and CT sensitivity are not implemented
166  if (m2p_DataFile[0]->GetDataType()==TYPE_SPECT || m2p_DataFile[0]->GetDataType()==TYPE_CT)
167  {
168  Cerr("***** oSensitivityGenerator::Initialize() -> Sensitivity for SPECT and CT not yet implemented !" << endl);
169  return 1;
170  }
171 
172  // Allocate the forward image
174 
175  // Initialization of normalization files:
176  // --> must be done before attenuation files, because attenuation can be inside the normalization data file
177  // --> is not done if sensitivity generated from the histogram datafile
179  {
180  Cerr("***** oSensitivityGenerator::Initialize() -> A problem occurred while initializing the normalization data files !" << endl);
181  return 1;
182  }
183  // Initialization of attenuation related stuff: not done if sensitivity generated from the histogram datafile
185  {
186  Cerr("***** oSensitivityGenerator::Initialize() -> A problem occurred while initializing the attenuation files !" << endl);
187  return 1;
188  }
189  // We want to know if we will need to forward project into the attenuation image, so here the full condition
190  if (
191  m_mumapAttenuationFlag && // 1. We must have an attenuation image that is loaded
192  m2p_DataFile[0]->GetDataType()==TYPE_PET && // 2. This is only for PET
193  m2p_DataFile[0]->GetDataMode()==MODE_LIST && // 3. This is only for sensitivity generation for list-mode datafile
194  ( m3p_NormalizationDataFile==NULL || // 4.1 There is no normalization data file (potentially including ACF)
195  !(dynamic_cast<iDataFilePET*>(m3p_NormalizationDataFile[0][0]))->GetAtnCorrectionFlag() ) // 4.2 OR the normalization data file does not contain ACF
196  )
197  {
198  // We use a specific boolean that says we need to forward project into the attenuation image
200  }
201  else
202  {
203  // Otherwise, we don't
205  }
206 
207  // Allocate the backward image for multi-thread
209  // Create sensitivity image file name
210  sOutputManager* p_outputManager = sOutputManager::GetInstance();
211  m_pathToSensitivityImage = p_outputManager->GetPathName() + p_outputManager->GetBaseName() + "_sensitivity";
212  // Allocate line counters
213  mp_lineCounter = (uint64_t*)calloc(mp_ImageDimensionsAndQuantification->GetNbThreadsForProjection(),sizeof(uint64_t));
214  // End
215  return 0;
216 }
217 
218 // =====================================================================
219 // ---------------------------------------------------------------------
220 // ---------------------------------------------------------------------
221 // =====================================================================
222 
224 {
225  // TODO: check the consistency between all normalization files (atnflag, normfactor, etc)
226 
227  // If no normalization file name provided, then we will exit, but we do some checks before
228  if (mp_pathToNormalizationFileName.size()==0)
229  {
230  // If the ignore-normalization-correction-flag is on, then we can exit the function
232  // In PET, throw an error if normalization correction is in the datafile and is taken into account (only if the datafile is a list-mode)
233  if (m2p_DataFile[0]->GetDataType()==TYPE_PET && m2p_DataFile[0]->GetDataMode()==MODE_LIST)
234  {
235  // Cast the vDataFile into iDataFilePET
236  iDataFilePET* p_pet_file = (dynamic_cast<iDataFilePET*>(m2p_DataFile[0]));
237  // Check if the normalization data are in the file and we do not ignore it
238  if (p_pet_file->GetNormCorrectionFlag())
239  {
240  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> Normalization correction is included in the data file while it is not in the sensitivity computation !" << endl);
241  return 1;
242  }
243  }
244  // We can exit now
245  return 0;
246  }
247 
248  // Verbose
249  if (m_verbose>=2) Cout("oSensitivityGenerator::InitializeNormalizationFiles() -> Initialize normalization files" << endl);
250 
251  // Check that the number of normalization files options is the same as the number of beds
253  {
254  // If only one is provided, then we assume that it is the same for all bed positions
256  // Otherwise, throw an error
257  else
258  {
259  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> Number of normalization files options should be one or equal to the number of data files !" << endl);
260  return 1;
261  }
262  }
263 
264  // Allocate the normalization data files array for the bed positions
266 
267  // Get the scanner manager
268  sScannerManager* p_ScannerManager = sScannerManager::GetInstance();
269 
270  // Loop on the number of beds
271  for (int bed=0; bed<mp_ImageDimensionsAndQuantification->GetNbBeds(); bed++)
272  {
273  // --------------------------------------------------------------------------
274  // Step 0: When only one normalization file names is provided for all beds
275  // --------------------------------------------------------------------------
277  {
278  // If not the first bed
279  if (bed!=0)
280  {
281  // We copy the pointer of the first bed into this one
283  // And we continue the loop
284  continue;
285  }
286  }
287  // --------------------------------------------------------------------------
288  // Step 1: Get filenames and check consistency with number of beds and frames
289  // --------------------------------------------------------------------------
290  // We will get the normalization file names for this bed
291  vector<string> normalization_files;
292  // Compute the bed index associated to the file name; in case of reverse order, we change the index only for datafile management here
293  int bed_data_file_name_index = bed;
294  if (m_inverseDataFileOrderFlag) bed_data_file_name_index = mp_ImageDimensionsAndQuantification->GetNbBeds() - 1 - bed;
295  // Then we search for commas separating the file name associated to each frame
296  size_t comma = mp_pathToNormalizationFileName[bed_data_file_name_index].find_first_of(",");
297  // Loop to get all file names
298  while (comma!=string::npos)
299  {
300  // Get the file name before the comma
301  normalization_files.push_back(mp_pathToNormalizationFileName[bed_data_file_name_index].substr(0,comma));
302  // Remove this part from the collection of file names
303  mp_pathToNormalizationFileName[bed_data_file_name_index] = mp_pathToNormalizationFileName[bed_data_file_name_index].substr(comma+1);
304  // Search for the next comma
305  comma = mp_pathToNormalizationFileName[bed_data_file_name_index].find_first_of(",");
306  }
307  // Get the last file name
308  normalization_files.push_back(mp_pathToNormalizationFileName[bed_data_file_name_index]);
309  // If we have only one frame, then set the m_oneNormalizationFileForAllFrames to true right now
311  // Check that the number of file names found is equal to the number of frames
312  if (normalization_files.size()!=((size_t)(mp_ImageDimensionsAndQuantification->GetNbTimeFrames())))
313  {
314  // If only one is provided, then we assume that it is the same for all frames
315  if (normalization_files.size()==1) m_oneNormalizationFileForAllFrames = true;
316  // Otherwise, throw an error
317  else
318  {
319  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> Number of normalization files for bed " << bed+1 << " should be one or equal to the number of frames !" << endl);
320  return 1;
321  }
322  }
323 
324  // --------------------------------------------------------------------------
325  // Step 2: Create the normalization data file objects
326  // --------------------------------------------------------------------------
327  // Allocate the normalization data file array for the frames
329  // Loop on the number of frames
330  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
331  {
332  // If one normalization file for all frames
334  {
335  // Then, if not the first frame
336  if (fr!=0)
337  {
338  // We copy the pointer of the first frame into this one
340  // And we continue the loop
341  continue;
342  }
343  }
344 
345  // Switch on the scanner type and create the specific data files. Check if norm file is "empt" type , which means
346  // it will be ignored. To say that it will ignored, that NormDatafiles will be NULL;
347  if ((p_ScannerManager->GetScannerType() == SCANNER_PET) && normalization_files[fr]!="empt")
348  {
349  m3p_NormalizationDataFile[bed][fr] = new iDataFilePET();
350  m3p_NormalizationDataFile[bed][fr]->SetHeaderDataFileName(normalization_files[fr]);
351  m3p_NormalizationDataFile[bed][fr]->SetBedIndex(bed);
352  m3p_NormalizationDataFile[bed][fr]->SetVerbose(m2p_DataFile[0]->GetVerbose());
354  bool affect_quantification_flag = false;
355  if (m3p_NormalizationDataFile[bed][fr]->ReadInfoInHeader(affect_quantification_flag))
356  {
357  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> A problem occurred during normalization datafile header reading !" << endl);
358  return 1;
359  }
361  {
362  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> A problem occurred while checking normalization datafile parameters !" << endl);
363  return 1;
364  }
365  if (m3p_NormalizationDataFile[bed][fr]->ComputeSizeEvent())
366  {
367  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> A problem occurred in normalization datafile event size computation !" << endl);
368  return 1;
369  }
370  if (m3p_NormalizationDataFile[bed][fr]->InitializeMappedFile())
371  {
372  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> A problem occurred in normalization datafile initialization !" << endl);
373  return 1;
374  }
375  if (m3p_NormalizationDataFile[bed][fr]->PrepareDataFile())
376  {
377  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> A problem occurred in normalization datafile preparation !" << endl);
378  return 1;
379  }
380  // Check that the normalization data mode is histogram or normalization
381  if (m3p_NormalizationDataFile[bed][fr]->GetDataMode()!=MODE_HISTOGRAM && m3p_NormalizationDataFile[bed][fr]->GetDataMode()!=MODE_NORMALIZATION)
382  {
383  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> Normalization data file used for sensitivity computation must be either histogram or normalization mode !" << endl);
384  return 1;
385  }
386  }
387  else if ((p_ScannerManager->GetScannerType() == SCANNER_PET) && normalization_files[fr]=="empt")
388  {
389  m3p_NormalizationDataFile[bed][fr] = NULL;
390  }
391  else
392  {
393  Cerr("***** oSensitivityGenerator::InitializeNormalizationFiles() -> Unknown scanner type (" << p_ScannerManager->GetScannerType() << ") or not yet implemented !" << endl);
394  return 1;
395  }
396  }
397  }
398 
399  // End
400  return 0;
401 }
402 
403 // =====================================================================
404 // ---------------------------------------------------------------------
405 // ---------------------------------------------------------------------
406 // =====================================================================
407 
409 {
410  // If it is asked to ignore attenuation correction, then nothing to do
412 
413  // If attenuation image file name is empty, then we do some checks before
414  if (m_pathToAttenuationImage=="")
415  {
416  // Case for PET: if no normalization data file OR no ACF in the normalization data file
417  if (m2p_DataFile[0]->GetDataType()==TYPE_PET && (m3p_NormalizationDataFile==NULL || !(dynamic_cast<iDataFilePET*>(m3p_NormalizationDataFile[0][0]))->GetAtnCorrectionFlag()))
418  {
419  // Cast the vDataFile into iDataFilePET
420  iDataFilePET* p_pet_file = (dynamic_cast<iDataFilePET*>(m2p_DataFile[0]));
421  // Throw an error if attenuation correction is in the datafile (at this step we know we do not ignore it)
422  if (p_pet_file->GetAtnCorrectionFlag())
423  {
424  Cerr("***** oSensitivityGenerator::InitializeAttenuationFiles() -> Attenuation correction is included in the data file while it is not in the sensitivity computation !" << endl);
425  return 1;
426  }
427  }
428  // We can exit now
429  return 0;
430  }
431 
432  // In PET, if the ACF is already present in the normalization data file, then say it and ignore mumap
433  if (m2p_DataFile[0]->GetDataType()==TYPE_PET && m3p_NormalizationDataFile!=NULL && (dynamic_cast<iDataFilePET*>(m3p_NormalizationDataFile[0][0]))->GetAtnCorrectionFlag())
434  {
435  Cout("oSensitivityGenerator::InitializeAttenuationFiles() -> Ignore provided mumap and consider ACF provided in the normalization data file" << endl);
436  return 0;
437  }
438 
439  // Verbose
440  if (m_verbose>=2) Cout("oSensitivityGenerator::InitializeAttenuationFiles() -> Allocate and read attenuation images (assumed to be in cm-1)" << endl);
441 
442  // Allocate and read the attenuation image (this puts the image into the m4p_attenuation buffer of oImageSpace)
444  {
445  Cerr("***** oSensitivityGenerator::Initialize() -> A problem occurred while initializing the attenuation image into the image space !" << endl);
446  return 1;
447  }
448 
449  // Allocate the foward image, where the attenuation will be copied and deformed if needed
450  //mp_ImageSpace->LMS_InstantiateForwardImage();
451 
452  // Copy attenuation image into forward-image buffer
455 
456  // We do not need the m4p_attenuation buffer of oImageSpace anymore, so free it now
458 
459  // TODO: Think about adding a field in the attenuation header to specify if the attenuation images have the same spatial resolution as the
460  // scanner in use, in order to apply a smoothing on it. For the moment, we assume it has been done beforehand.
461 
462  // TODO: for the moment, we assume that the attenuation image is in cm-1, maybe use a header field to detect that and convert it
463  // automatically.
464 
465  // Put attenuation flag on
466  m_mumapAttenuationFlag = true;
467 
468  // End
469  return 0;
470 }
471 
472 // =====================================================================
473 // ---------------------------------------------------------------------
474 // ---------------------------------------------------------------------
475 // =====================================================================
476 
478 {
479  // Verbose
480  if (m_verbose>=1) Cout("oSensitivityGenerator::Launch() -> Start the sensitivity computation" << endl);
481  // Simply switches between the CPU and GPU versions
482  #ifdef CASTOR_GPU
483  if (m_flagGPU) return LaunchGPU();
484  else return LaunchCPU();
485  #else
486  return LaunchCPU();
487  #endif
488 }
489 
490 // =====================================================================
491 // ---------------------------------------------------------------------
492 // ---------------------------------------------------------------------
493 // =====================================================================
494 
496 {
497  // Loop on the bed positions
498  for (int bed=0; bed<mp_ImageDimensionsAndQuantification->GetNbBeds(); bed++)
499  {
500  // Reset line counters
502  // Apply the bed offset for this bed position
504  // Case where the computation is performed directly from the histogram datafile
506  {
508  {
509  Cerr("***** oSensitivityGenerator::LaunchCPU() -> A problem occurred while computing the sensitivity using the histogram data file for bed " << bed+1 << " !" << endl);
510  return 1;
511  }
512  }
513  // If a normalization file is provided, then use them to iterate on the elements to be taken into account
514  else if (m3p_NormalizationDataFile!=NULL)
515  {
517  {
518  Cerr("***** oSensitivityGenerator::LaunchCPU() -> A problem occurred while computing the sensitivity using the normalization file for bed " << bed+1 << " !" << endl);
519  return 1;
520  }
521  }
522  // Otherwise, use a generic method from the scanner elements
523  else
524  {
526  {
527  Cerr("***** oSensitivityGenerator::LaunchCPU() -> A problem occurred while computing the sensitivity from scanner elements for bed " << bed+1 << " !" << endl);
528  return 1;
529  }
530  }
531 
532  // Sum up the line counter into the first thread
534  // Verbose
535  if (m_verbose>=2) Cout(" --> Number of effectively projected lines: " << mp_lineCounter[0] << endl);
536  }
537 
538 
539  // Set the number of threads for image computations
540  #ifdef CASTOR_OMP
542  #endif
543 
544  // Reduce all results (merging parallel computation)
545  for (int fr=0 ; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames() ; fr++)
546  for (int rg=0 ; rg<mp_ImageDimensionsAndQuantification->GetNb1stMotImgsForLMS(fr) ; rg++)
547  for (int cg=0 ; cg<mp_ImageDimensionsAndQuantification->GetNb2ndMotImgsForLMS() ; cg++)
548  {
549  // Merge parallel results into the backward image
550  int img_0 = 0;
551  mp_ImageSpace->ReduceBackwardImage(img_0, fr, rg, cg);
552  // Perform any required deformations on the backward image
556  fr,
557  rg,
558  cg);
559  }
560 
561  // Allocate sensitivity image
563  // Copy backward image to sensitivity image
565  // Deallocate backward image
567 
568  // Deallocate the forward image only if attenuation was used
569  //if (m_mumapAttenuationFlag) mp_ImageSpace->LMS_DeallocateForwardImage();
570  // Deallocate forward image
572 
573  // Apply PSF and save only for first instance
575  {
576  // Apply PSF if needed
578  // Save sensitivity image
580  }
581  // Deallocate sensitivity image if the sensitivity is not needed after
583  // End
584  return 0;
585 }
586 
587 // =====================================================================
588 // ---------------------------------------------------------------------
589 // ---------------------------------------------------------------------
590 // =====================================================================
591 
593 {
594  // For the moment, we restrict the validity of this function to 3D datasets
596  {
597  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromHistogramDataFile() -> This functionality is only available for non-gated reconstruction without image-based motion correction. For these types of reconstruction, the sensitivity must be estimated from the scanner geometry instead !" << endl);
598  return 1;
599  }
600 
601  // Verbose
602  if (m_verbose>=2)
603  {
604  if (mp_ImageDimensionsAndQuantification->GetNbBeds()>1) Cout("oSensitivityGenerator::ComputeSensitivityFromHistogramDataFile() -> Start computation for bed " << a_bed+1 << endl);
605  else Cout("oSensitivityGenerator::ComputeSensitivityFromHistogramDataFile() -> Start computation" << endl);
606  }
607 
608  // Initial clock
609  clock_t clock_start = clock();
610  time_t time_start = time(NULL);
611 
612  // In case a problem occurs in the parallel loop
613  bool problem = false;
614  // Set the number of threads for projections
615  #ifdef CASTOR_OMP
617  #endif
618 
619  // Get index start and stop
620  int64_t index_start = 0;
621  int64_t index_stop = 0;
622  m2p_DataFile[a_bed]->GetEventIndexStartAndStop(&index_start, &index_stop);
623 
624  // Reinitialize 4D gate indexes
626 
627  // Launch the loop on all events
628  int64_t index = 0, printing_index = 0;
629  #pragma omp parallel for private(index) schedule(static, 1)
630  for ( index = index_start ; index < index_stop ; index++)
631  {
632  // Get the thread index
633  int th = 0;
634  #ifdef CASTOR_OMP
635  th = omp_get_thread_num();
636  #endif
637  // Verbose
639  {
640  if (printing_index%1000==0)
641  {
642  int percent = ((int)( ((FLTNB)(index-index_start)) * 100. / ((FLTNB)(index_stop-index_start)) ));
643  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 "
644  << percent << " % " << flush;
645  }
646  printing_index++;
647  }
648  // Get the current event for that thread index
649  vEvent* event = m2p_DataFile[a_bed]->GetEvent(index, th);
650  if (event==NULL)
651  {
652  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromHistogramDataFile() -> An error occurred while getting the event from index "
653  << index << " (thread " << th << ") !" << endl);
654  // A problem was encountered
655  problem = true;
656  // We must continue here because we are inside an OpenMP loop
657  continue;
658  }
659 
660  // Compute the projection line
662  if (line==NULL)
663  {
664  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromHistogramDataFile() -> A problem occurred while computing the projection line !" << endl);
665  // Specify that there was a problem
666  problem = true;
667  // We must continue here because we are inside an OpenMP loop
668  continue;
669  }
670 
671  // Process this line
672  int no_frame = 0;
673  int no_resp_gate = 0;
674  int no_card_gate = 0;
675 // if (line->NotEmptyLine() && ProcessThisLine(line, event, a_bed, mp_ImageDimensionsAndQuantification->GetCurrentTimeFrame(th), no_resp_gate, no_card_gate, th))
676  if (line->NotEmptyLine() && ProcessThisLine(line, event, a_bed, no_frame, no_resp_gate, no_card_gate, th))
677  {
678  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromHistogramDataFile() -> A problem occurred while processing a line !" << endl);
679  }
680  }
681 
682  // Synchronize MPI processes
683  #ifdef CASTOR_MPI
684  MPI_Barrier(MPI_COMM_WORLD);
685  #endif
686 
687  // End of progression printing (do not log out with Cout here)
689  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"
690  << " --> 100 % " << endl;
691  // If a problem was encountered, then report it here
692  if (problem)
693  {
694  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromHistogramDataFile() -> A problem occurred inside the parallel loop over events !" << endl);
695  return 1;
696  }
697 
698  // Clock total
699  clock_t clock_stop = clock();
700  time_t time_stop = time(NULL);
702  Cout(" --> Time spent for sensitivity generation for bed " << a_bed+1 << " | User: " << time_stop-time_start
703  << " sec | CPU: " << (clock_stop-clock_start)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
704  else
705  Cout(" --> Time spent for sensitivity generation | User: " << time_stop-time_start
706  << " sec | CPU: " << (clock_stop-clock_start)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
707 
708  // End
709  return 0;
710 }
711 
712 // =====================================================================
713 // ---------------------------------------------------------------------
714 // ---------------------------------------------------------------------
715 // =====================================================================
716 
718 {
719  // Verbose
720  if (m_verbose>=2)
721  {
722  if (mp_ImageDimensionsAndQuantification->GetNbBeds()>1) Cout("oSensitivityGenerator::ComputeSensitivityFromNormalizationFile() -> Start computation for bed " << a_bed+1 << endl);
723  else Cout("oSensitivityGenerator::ComputeSensitivityFromNormalizationFile() -> Start computation" << endl);
724  }
725 
726  /* // For the moment, we restrict the validity of this function to non-gated datasets and without involuntary patient motion correction
727  if (mp_ImageDimensionsAndQuantification->GetDynRecoType() != STATIC_RECO
728  && mp_ImageDimensionsAndQuantification->GetDynRecoType() != DYN_RECO_FRAMING )
729  {
730  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromNormalizationFile() -> This functionality is only available for non-gated reconstruction without image-based motion correction. For these types of reconstruction, the sensitivity must be estimated from the scanner geometry instead !" << endl);
731  return 1;
732  } */
733 
734  // Initial clock
735  clock_t clock_start = clock();
736  time_t time_start = time(NULL);
737 
738  // Loop on frames
739  for (int fr=0 ; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames() ; fr++)
740  {
741  // Verbose
743  cout << " --> Processing frame " << fr+1 << endl;
744 
745  // Check if norm file for this frame is NULL, in which case it is ignored in this step.
746  if (m3p_NormalizationDataFile[a_bed][fr]!=NULL)
747  {
748  // Loop on respiratory gates
749  for (int rg=0 ; rg<mp_ImageDimensionsAndQuantification->GetNb1stMotImgsForLMS(fr) ; rg++)
750  {
751  // Verbose
753  cout << " --> Processing respiratory gate " << rg+1 << endl;
754  // Loop on cardiac gates
755  for (int cg=0 ; cg<mp_ImageDimensionsAndQuantification->GetNb2ndMotImgsForLMS() ; cg++)
756  {
757  // Verbose
759  cout << " --> Processing cardiac gate " << cg+1 << endl;
760 
761  // Perform any forward deformations on the forward image
765  fr,
766  rg,
767  cg);
768 
769  // In case a problem occurs in the parallel loop
770  bool problem = false;
771  // Set the number of threads for projections
772  #ifdef CASTOR_OMP
774  #endif
775 
776  // Get index start and stop
777  int64_t index_start = 0;
778  int64_t index_stop = 0;
779  m3p_NormalizationDataFile[a_bed][fr]->GetEventIndexStartAndStop(&index_start, &index_stop);
780 
781  // Index for progression printing
782  uint64_t progression_index_total = (index_stop-index_start)
786 
787  // Launch the loop on all events
788  int64_t index = 0, printing_index = 0;
789  #pragma omp parallel for private(index) schedule(static, 1)
790  for ( index = index_start ; index < index_stop ; index++)
791  {
792  // Get the thread index
793  int th = 0;
794  #ifdef CASTOR_OMP
795  th = omp_get_thread_num();
796  #endif
797 
798  // Verbose
799  /*
800  if (th==0 && m_verbose>=2 && mp_ImageDimensionsAndQuantification->GetMPIRank()==0)
801  {
802  if (printing_index%1000==0)
803  {
804  int percent = ((int)( ((FLTNB)(index-index_start)) * 100. / ((FLTNB)(index_stop-index_start)) ));
805  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 "
806  << percent << " % " << flush;
807  }
808  printing_index++;
809  }
810  */
811 
813  {
814  if (printing_index%1000==0)
815  {
816  int64_t progression_chunk = (int64_t)((FLTNB)(index_stop-index_start));
817  int64_t progression_index_current = fr * mp_ImageDimensionsAndQuantification->GetNb1stMotImgsForLMS(fr) * mp_ImageDimensionsAndQuantification->GetNb2ndMotImgsForLMS() * progression_chunk
818  + rg * mp_ImageDimensionsAndQuantification->GetNb2ndMotImgsForLMS() * progression_chunk
819  + cg * progression_chunk
820  + index;
821 
822  int64_t percent = (int64_t) (( ((FLTNB)progression_index_current)/((FLTNB)progression_index_total) ) * ((FLTNB)100));
823  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 "
824  << percent << " % " << flush;
825  }
826  printing_index++;
827  }
828 
829  // Get the current event for that thread index
830  vEvent* event = m3p_NormalizationDataFile[a_bed][fr]->GetEvent(index, th);
831 
832  if (event==NULL)
833  {
834  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromNormalizationFile() -> An error occurred while getting the event from index "
835  << index << " (thread " << th << ") !" << endl);
836  // A problem was encountered
837  problem = true;
838  // We must continue here because we are inside an OpenMP loop
839  continue;
840  }
841  // Compute the projection line
843  if (line==NULL)
844  {
845  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromNormalizationFile() -> A problem occurred while computing the projection line !" << endl);
846  // Specify that there was a problem
847  problem = true;
848  // We must continue here because we are inside an OpenMP loop
849  continue;
850  }
851  // Process this line
852  if (line->NotEmptyLine() && ProcessThisLine(line, event, a_bed, fr, rg, cg, th))
853  {
854  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromNormalizationFile() -> A problem occurred while processing a line !" << endl);
855  }
856  }
857  // End of progression printing (do not log out with Cout here)
859  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"
860  << " --> 100 % " << endl;
861  // If a problem was encountered, then report it here
862  if (problem)
863  {
864  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromNormalizationFile() -> A problem occurred inside the parallel loop over events !" << endl);
865  return 1;
866  }
867  }
868  }
869  }
870  }
871 
872  // Clock total
873  clock_t clock_stop = clock();
874  time_t time_stop = time(NULL);
876  Cout(" --> Time spent for sensitivity generation for bed " << a_bed+1 << " | User: " << time_stop-time_start
877  << " sec | CPU: " << (clock_stop-clock_start)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
878  else
879  Cout(" --> Time spent for sensitivity generation | User: " << time_stop-time_start
880  << " sec | CPU: " << (clock_stop-clock_start)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
881 
882  // End
883  return 0;
884 }
885 
886 // =====================================================================
887 // ---------------------------------------------------------------------
888 // ---------------------------------------------------------------------
889 // =====================================================================
890 
892 {
893  // Verbose
894  if (m_verbose>=2)
895  {
896  if (mp_ImageDimensionsAndQuantification->GetNbBeds()>1) Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> Start computation for bed " << a_bed+1 << endl);
897  else Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> Start computation" << endl);
898  }
899 
900  // Get the scanner manager
901  sScannerManager* p_scannerManager = sScannerManager::GetInstance();
902 
903  // Total number of elements
904  int nb_total_elts = mp_Scanner->GetSystemNbElts();
905 
906  // Initial clock
907  clock_t clock_start = clock();
908  time_t time_start = time(NULL);
909 
910  // Initialize main loop start and stop values
911  int64_t main_loop_start_index = 0 ;
912  int64_t main_loop_stop_index = 0;
913 
914  // Check beforehand any issue with the loop start/stop values (not possible in the inner multithreaded loop)
915  if (p_scannerManager->PROJ_GetModalityStopValueMainLoop()>0) main_loop_stop_index = p_scannerManager->PROJ_GetModalityStopValueMainLoop();
916  else
917  {
918  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromScanner() -> An error occurred when trying to initialize main loop stop index !" << endl);
919  return 1;
920  }
921  if (p_scannerManager->PROJ_GetModalityStartValueInnerLoop(0)<=0)
922  {
923  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromScanner() -> An error occurred when trying to initialize inner loop start index !" << endl);
924  return 1;
925  }
926 
927  // Prepare pre-computed sums of events to avoid the exponential evolution of the percentage (for PET)
928  // TODO Perhaps replace this with a call to scannerManager for error management
929  int64_t* progression_elts_array = new int64_t[nb_total_elts*sizeof(int64_t)];
930  progression_elts_array[0] = 0;
931 
932  for (int idx_elt=1 ; idx_elt<mp_Scanner->GetSystemNbElts() ; idx_elt++)
933  progression_elts_array[idx_elt] = progression_elts_array[idx_elt-1] + (mp_Scanner->GetSystemNbElts()-idx_elt);
934 
935  // Index for progression printing
936  uint64_t printing_index = 0;
937  uint64_t progression_index_total = p_scannerManager->PROJ_GetProgressionFinalValue();
938 
939  //uint64_t progression_index_total = p_scannerManager->PROJ_GetProgressionFinalValue()
940  // * mp_ImageDimensionsAndQuantification->GetNb1stMotImgsForLMS(0)
941  // * mp_ImageDimensionsAndQuantification->GetNb2ndMotImgsForLMS();
942 
943  // Adapt the progression timer to dynamic acquisition
944 
945  // Attenuation image provided and gating/physiological motion correction enabled
948  progression_index_total *= mp_ImageDimensionsAndQuantification->GetNb1stMotImgsForLMS(0)
950 
951  // Attenuation image provided, first frame, and gating/physiological motion correction enabled
953  {
954  progression_index_total = 0 ;
955  for (int fr=0 ; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames() ; fr++)
956  for (int rg=0 ; rg<mp_ImageDimensionsAndQuantification->GetNb1stMotImgsForLMS(fr) ; rg++)
957  if ( ( fr==0 )
958  || ( rg>0 )
960  progression_index_total += p_scannerManager->PROJ_GetProgressionFinalValue();
961  }
962 
963  // Divide by the number of MPI instance
964  progression_index_total /= mp_ImageDimensionsAndQuantification->GetMPISize();
965 
966 
967 
968  // ------------------ Compute here the part that each MPI instance manages --------------------
969  // Synchronize MPI processes
970  #ifdef CASTOR_MPI
971  MPI_Barrier(MPI_COMM_WORLD);
972  #endif
973 
974  int64_t idx_elt =0;
975 
976  while (idx_elt< nb_total_elts
977  &&(uint64_t)progression_elts_array[idx_elt] < (mp_ImageDimensionsAndQuantification->GetMPIRank()+1) * progression_index_total)
978  {
979  // Get starting index for the MPI instance
980  if (mp_ImageDimensionsAndQuantification->GetMPIRank() != 0 // For the first instance, starting index is 0
981  && main_loop_start_index == 0 // Starting index has not been initialized yet (1st instance is discarded with check above)
982  && ((uint64_t)progression_elts_array[idx_elt] >= mp_ImageDimensionsAndQuantification->GetMPIRank() * progression_index_total) )
983  main_loop_start_index = idx_elt;
984 
985  idx_elt++;
986  }
987 
988  // Set the stop element index for the MPI instance
989  main_loop_stop_index = idx_elt;
990 
991  // Get the id of index element for this mpi instance
992  #ifdef CASTOR_MPI
993  if (m_verbose>=2)
994  Cout( "oSensitivityGenerator::ComputeSensitivityFromScanner() -> MPI ID " << mp_ImageDimensionsAndQuantification->GetMPIRank() <<
995  " start/stop: " << main_loop_start_index << "/" << main_loop_stop_index << endl);
996  #endif
997 
998 
999  // Dynamic list-mode sensitivity images are generated according to the nature of the dynamic acquisition:
1000  // - 4D acquisition with dynamic frames : One sensitivity image is generated, then replicated with appropriate quantification for each frames
1001  // - Gated (respiratory/cardiac) acquisition (within each frame with 5D (frame + gates) acquisitions ) :
1002  // ->Atn image provided (each position) : sensitivity image will be generated for each gate
1003  // ->Otherwise : Generation of one sensitivity image, then duplicated for each gate
1004  // - Gated (respiratory/cardiac) acquisition with (voxel-based) motion correction (within each frame with 5D (frame + gates) acquisitions ) :
1005  // ->Atn image of the ref position : For each gate: forward deformation of the sensitivity image, followed by sensitivity image generation
1006  // ->Atn image of each position : TODO ?
1007  // ->Otherwise : Generation of one sensitivity image, then duplicated for each gate
1008  // Backward deformations to the reference position performed on the sensitivity images of each gate.
1009  // An average sensitivity image if the reference position is then produced from these gated image assuming equivalent weighing for each gate
1010  // - Involuntary patient motion (IPM) correction (time-based deformation)
1011  // (within each frame with 5D (frame + gates) acquisitions. Note that contrary to gated motion correction, each frame could contain different number of motion subset) :
1012  // ->Atn image of the ref position : Sensitivity images will be generated for each position resulting from the time-based motion subset
1013  // ->Otherwise : Generation of one sensitivity image, then duplicated for each gate
1014  // Backward deformations to the reference position performed on the sensitivity images of each motion subset.
1015  // An average sensitivity image if the reference position is then produced from these gated image with weighing based on duration of each subset
1016 
1017  // Just a variable to track how many dynamic images have been processed
1018  // for progression feedback calculation
1019  uint16_t nb_dyn_img_processed=0;
1020 
1021  // Loop on frames
1022  for (int fr=0 ; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames() ; fr++)
1023  {
1024  // Verbose
1026  Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> Processing frame " << fr+1 << endl);
1027 
1028  // Loop on 1st motion images (respiratory gates or IPM subset images)
1029  for (int rg=0 ; rg<mp_ImageDimensionsAndQuantification->GetNb1stMotImgsForLMS(fr) ; rg++)
1030  {
1031  // Verbose
1033  Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> Processing 1st motion image for sensitivity " << rg+1 << endl);
1034 
1035  // Loop on cardiac gates
1036  for (int cg=0 ; cg<mp_ImageDimensionsAndQuantification->GetNb2ndMotImgsForLMS() ; cg++)
1037  {
1038  // Verbose
1040  Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> Processing 2nd motion image for sensitivity " << cg+1 << endl);
1041 
1042  // --- Generate a sensivitity image --- //
1043 
1044  if( ( nb_dyn_img_processed == 0 ) // First sensitivity image
1045  || // Attenuation image provided, first frame, and gating/physiological motion correction enabled
1046  ( m_mumapAttenuationFlag // Attenuation image provided
1047  && fr==0 // --> first frame
1048  && ( mp_ImageDimensionsAndQuantification->GetDynRecoType() == DYN_RECO_GATING // --> Gating/physiological motion correction enabled
1050  || // Attenuation image provided , involuntary patient motion (IPM) correction enabled, IPM index (rg)>0 for this frame - or - first IPM index for this frame (rg==0) is different from the last IPM index of the previous frame
1051  ( m_mumapAttenuationFlag // Attenuation image provided
1052  && mp_ImageDimensionsAndQuantification->GetDynRecoType() == DYN_RECO_IPMC // involuntary patient motion correction (IPM) enabled
1053  // IPM index (rg)>0 for this frame - or - first IPM index (rg==0) for this frame is different from the last IPM index of the previous frame
1055  )
1056  )
1057  {
1058  if (mp_ImageDimensionsAndQuantification->GetMPIRank()==0 && m_verbose>=VERBOSE_NORMAL && nb_dyn_img_processed == 0)
1059  Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> Generate sensitivity image" << endl);
1060  else if (mp_ImageDimensionsAndQuantification->GetMPIRank()==0 && m_verbose>=VERBOSE_DETAIL && nb_dyn_img_processed > 0) // Verbose detail only
1061  Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> Generate sensitivity image" << endl);
1062 
1063  // Perform any forward deformations on the forward image (only if attenuation is used)
1068  fr,
1069  rg,
1070  cg);
1071 
1072  // In case a problem occurs in the parallel loop
1073  bool problem = false;
1074  // Set the number of threads for projections
1075  #ifdef CASTOR_OMP
1077  #endif
1078 
1079  // Start the loop
1080  int64_t idx_elt1;
1081 
1082  #pragma omp parallel for private(idx_elt1) schedule(static, 1)
1083  for(idx_elt1=main_loop_start_index ; idx_elt1<main_loop_stop_index ; idx_elt1++)
1084  {
1085  #ifdef CASTOR_MPI
1086  if(idx_elt1==main_loop_start_index && m_verbose>=VERBOSE_NORMAL)
1087  Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> MPI ID " << mp_ImageDimensionsAndQuantification->GetMPIRank()
1088  << " start/stop: " << main_loop_start_index << "/" << main_loop_stop_index << endl);
1089  #endif
1090 
1091  // Get the thread index
1092  int th = 0;
1093  #ifdef CASTOR_OMP
1094  th = omp_get_thread_num();
1095  #endif
1096 
1097  // Initialize inner loop start and stop values
1098  int64_t inner_loop_start_index = p_scannerManager->PROJ_GetModalityStartValueInnerLoop(idx_elt1) ;
1099  int64_t inner_loop_stop_index = mp_Scanner->GetSystemNbElts();
1100 
1101  // Inner loop on scanner elements (idx_elt2) which are used to form LOR using the first scanner element (idx_elt1)
1102  for (int64_t idx_elt2=inner_loop_start_index ; idx_elt2<inner_loop_stop_index ; idx_elt2++)
1103  {
1104  // Print progression (do not log out with Cout here)
1106  {
1107  if (m_verbose>=2 && printing_index%10000==0)
1108  {
1109  int64_t progression_index_current;
1110 
1111  progression_index_current = p_scannerManager->PROJ_GetCurrentProgression(idx_elt1, idx_elt2, progression_elts_array, nb_dyn_img_processed);
1112 
1113  int64_t percent = (int64_t) (( ((FLTNB)progression_index_current)/((FLTNB)progression_index_total) ) * ((FLTNB)100));
1114 
1116  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 "
1117  << percent << " % ";
1118  }
1119  printing_index++;
1120  }
1121 
1122  // Allocate an event using the iDataFile
1123  vEvent* event = m2p_DataFile[a_bed]->PROJ_GenerateEvent(idx_elt1, idx_elt2, th);
1124 
1125  // Check from the scanner requirements that this LOR is allowed
1126  if (!mp_Scanner->IsAvailableLOR(event->GetID1(0), event->GetID2(0))) continue;
1127 
1128  // Generate the projection event and compute the projection line
1130  if (line==NULL)
1131  {
1132  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromScanner() -> A problem occurred while computing the projection line !" << endl);
1133  // Specify that there was a problem
1134  problem = true;
1135  // We must continue here because we are inside an OpenMP loop
1136  continue;
1137  }
1138 
1139  // Process this line
1140  if (line->NotEmptyLine() && ProcessThisLine(line, event, a_bed, fr, rg, cg, th))
1141  {
1142  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromScanner() -> A problem occurred while processing a line !" << endl);
1143  }
1144  }
1145  }
1146 
1147  // Synchronize MPI processes
1148  #ifdef CASTOR_MPI
1149  MPI_Barrier(MPI_COMM_WORLD);
1150  #endif
1151 
1152  // If a problem was encountered, then report it here
1153  if (problem)
1154  {
1155  Cerr("***** oSensitivityGenerator::ComputeSensitivityFromScanner() -> A problem occurred inside the parallel loop over events !" << endl);
1156  return 1;
1157  }
1158 
1159  nb_dyn_img_processed++;
1160  }
1161 
1162  // --- Just duplicate the sensitivity image which has already been generated, with appropriate quantification --- //
1163  else
1164  {
1165  if (m_verbose>=VERBOSE_DETAIL )
1166  Cout("oSensitivityGenerator::ComputeSensitivityFromScanner() -> Duplicate sensitivity image" << endl);
1167 
1168  // Framing, not first frame
1169  if (fr > 0)
1170  {
1171  // IPM is not enabled : duplicate all gates for each frame
1173  {
1175  for (int v=0 ; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ() ; v++)
1176  mp_ImageSpace->m6p_backwardImage[0][th][fr][rg][cg][v] = mp_ImageSpace->m6p_backwardImage[0][th][0][rg][cg][v]
1179 
1180  }
1181  // IPM is enabled : get the image from last IPM index of previous frame
1182  else
1183  {
1185  for (int v=0 ; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ() ; v++)
1189  }
1190  }
1191  // Otherwise : Just duplicate the first image for each frame/gate
1192  else
1193  {
1195  for (int v=0 ; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ() ; v++)
1196  mp_ImageSpace->m6p_backwardImage[0][th][fr][rg][cg][v] = mp_ImageSpace->m6p_backwardImage[0][th][0][0][0][v]
1199  }
1200  }
1201 
1202  }
1203  }
1204  }
1205 
1206  // End of progression printing (do not log out with Cout here)
1208  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"
1209  << " --> 100 % " << endl;
1210 
1211  delete[] progression_elts_array;
1212 
1213  // Clock total
1214  clock_t clock_stop = clock();
1215  time_t time_stop = time(NULL);
1217  Cout(" --> Time spent for sensitivity generation for bed " << a_bed+1 << " | User: " << time_stop-time_start
1218  << " sec | CPU: " << (clock_stop-clock_start)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
1219  else
1220  Cout(" --> Time spent for sensitivity generation | User: " << time_stop-time_start
1221  << " sec | CPU: " << (clock_stop-clock_start)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
1222 
1223  // End
1224  return 0;
1225 }
1226 
1227 // =====================================================================
1228 // ---------------------------------------------------------------------
1229 // ---------------------------------------------------------------------
1230 // =====================================================================
1231 
1232 int oSensitivityGenerator::ProcessThisLine(oProjectionLine* ap_Line, vEvent* ap_Event, int a_bed, int a_frame, int a_respGate, int a_cardGate, int a_thread)
1233 {
1234  // ----------------------------------------------------------------------------------------------
1235  // Part 0: get the multiplicative corrections from the vEvent and the quantification factor
1236  // ----------------------------------------------------------------------------------------------
1237 
1238  // The multiplicative corrections from the vEvent can be:
1239  // 1: generated from vDataFile::GenerateEvent(), this will be a default value, so this will be 1.
1240  // 2: taken from the normalization data file, so this will be the normalization correction factor times the attenuation correction factor
1241  FLTNB multiplicative_factor = ap_Event->GetMultiplicativeCorrections();
1242  // If null, then skip this event
1243  if (multiplicative_factor<=0.) return 0;
1244  // Get the quantification factor
1245  FLTNB quantification_factor = mp_ImageDimensionsAndQuantification->GetQuantificationFactor(a_bed, a_frame, a_respGate, a_cardGate);
1246  // If null, then skip this event
1247  if (quantification_factor<=0.) return 0;
1248 
1249  // ----------------------------------------------------------------------------------------------
1250  // Part 1: deal with the attenuation forward projection for PET if any
1251  // ----------------------------------------------------------------------------------------------
1252 
1253  // The ACF
1254  FLTNB attenuation_correction_factor = 1.;
1255 
1256  // Do a forward projection only if the attenuation image was provided and the ACF was not included in the normalization data file
1258  {
1259  // Cumulative mu (in cm-1)
1260  FLTNB cumulative_mu = 0.;
1261  // Loop on TOF bins
1262  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
1263  {
1264  // Projection operation (the attenuation image in the m4p_forwardImage buffer of the oImageSpace)
1265  for (int vl=0 ; vl<ap_Line->GetCurrentNbVoxels(FORWARD, b) ; vl++)
1266  cumulative_mu += ap_Line->GetVoxelWeights(FORWARD, b, vl) * mp_ImageSpace->m4p_forwardImage[a_frame][a_respGate][a_cardGate][ap_Line->GetVoxelIndex(FORWARD, b,vl)];
1267  }
1268  // Update the ACF (converting the cumulative mu in mm-1)
1269  attenuation_correction_factor = max(((FLTNB)1.),exp(cumulative_mu*((FLTNB)0.1)));
1270  }
1271 
1272  // ----------------------------------------------------------------------------------------------
1273  // Part 2: back-projection of the LOR sensitivity
1274  // ----------------------------------------------------------------------------------------------
1275 
1276  // Compute the global sensitivity (calibration / (ACF * norm))
1277  FLTNB lor_sensitivity = 1. / (quantification_factor * attenuation_correction_factor * multiplicative_factor);
1278 
1279  // Loop on TOF bins
1280  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
1281  {
1282  // Backprojection operation into the backward image
1283  for (int vl=0 ; vl<ap_Line->GetCurrentNbVoxels(BACKWARD, b) ; vl++)
1284  mp_ImageSpace->m6p_backwardImage[0][a_thread][a_frame][a_respGate][a_cardGate][ap_Line->GetVoxelIndex(BACKWARD, b, vl)] += ap_Line->GetVoxelWeights(BACKWARD, b, vl)
1285  * lor_sensitivity;
1286  }
1287 
1288  // Increment line counters
1289  mp_lineCounter[a_thread]++;
1290 
1291  // End
1292  return 0;
1293 }
1294 
1295 // =====================================================================
1296 // ---------------------------------------------------------------------
1297 // ---------------------------------------------------------------------
1298 // =====================================================================
1299 
1301 {
1302  // Get the flag that says if we merge the dynamic files or not
1303  bool merge_dynamic_files = sOutputManager::GetInstance()->MergeDynImages();
1304  // If we merge dynamic file, or we are in a static case (no dynamic), we stay with the hdr extension
1308  merge_dynamic_files )
1309  return m_pathToSensitivityImage + ".hdr";
1310  else
1311  return m_pathToSensitivityImage + ".mhd";
1312 }
This class is designed to be a mother virtual class for DataFile.
bool GetNormCorrectionFlag()
Simply return m_normCorrectionFlag.
static sScannerManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
void GetEventIndexStartAndStop(int64_t *ap_indexStart, int64_t *ap_indexStop, int a_subsetNum=0, int a_NbSubsets=1)
#define MODE_HISTOGRAM
int ProcessThisLine(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_frame, int a_respGate, int a_cardGate, int a_thread)
void LMS_CopyAtnToForwardImage(bool a_use1stMotion, bool a_use2ndMotion)
int64_t PROJ_GetModalityStartValueInnerLoop(int64_t a_elt1)
#define MODE_NORMALIZATION
#define Cerr(MESSAGE)
#define FORWARD_DEFORMATION
oSensitivityGenerator()
The constructor of oSensitivityGenerator.
int Initialize()
A public function used to initialize the sensitivity generator.
bool UseDeformationInv()
Indicate if the involuntary patient motion deformation is enabled.
int LMS_SaveSensitivityImage(const string &a_pathToSensitivityImage, oDeformationManager *ap_DeformationManager)
string GetPathToSensitivityImage()
This function return the path to the sensitivity image.
int64_t PROJ_GetProgressionFinalValue()
Get numerator value according to the modality to compute percent progression during the projection pr...
#define BACKWARD_DEFORMATION
int ApplyDeformationForSensitivityGeneration(oImageSpace *ap_Image, int a_defDirection, int idx, int fr, int rg, int cg)
void SetVerbose(int a_verboseLevel)
bool GetIgnoreAttnCorrectionFlag()
Get the boolean that says if the attenuation correction is ignored or not.
~oSensitivityGenerator()
The destructor of oSensitivityGenerator.
int Launch()
A public function used to launch the sensitivity generator (compute the sensitivity image) ...
virtual FLTNB GetMultiplicativeCorrections()=0
This is a pure virtual function implemented in the child classes.
void LMS_InstantiateForwardImage()
Allocate memory for the forward image matrices (for list-mode sensitivity generation) ...
bool UseDeformationCard()
Indicate if the cardiac motion deformation is enabled.
static sOutputManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
int InitAttenuationImage(const string &a_pathToAtnImage)
bool GetIgnoreNormCorrectionFlag()
Get the boolean that says if the normalization correction is ignored or not.
vEvent * GetEvent(int64_t a_eventIndex, int a_th=0)
INTNB GetVoxelIndex(int a_direction, int a_TOFBin, INTNB a_voxelInLine)
Singleton class that manages output writing on disk (images, sinograms, etc). It also manages loggi...
INTNB GetCurrentNbVoxels(int a_direction, int a_TOFBin)
void LMS_DeallocateSensitivityImage()
Free memory for the sensitivity image matrices (for list-mode sensitivity generation) ...
Singleton class that Instantiate and initialize the scanner object.
int GetMPISize()
Get the MPI size (the number of MPI instances)
void ResetCurrentDynamicIndices()
Call the eponym function from the oDynamicDataManager class using the member object.
int64_t PROJ_GetCurrentProgression(int64_t a_elt1, int64_t a_elt2, int64_t *ap_nbEltsArray, uint16_t a_nbDynImgProcessed)
oProjectionLine * ComputeProjectionLine(vEvent *ap_Event, int a_th)
FLTNB GetQuantificationFactor(int a_bed, int a_frame, int a_respGate, int a_cardGate)
int InitializeAttenuationFiles()
Initialize the attenuation images provided for sensitivity computation.
void SetHeaderDataFileName(const string &a_headerFileName)
vEvent * PROJ_GenerateEvent(int idx_elt1, int idx_elt2, int a_th)
bool NotEmptyLine()
This function is used to know if the line contains any voxel contribution.
int LaunchCPU()
Launch the computation of the sensitivity image (CPU version)
virtual int GetSystemNbElts()=0
This is a pure virtual method that must be implemented by children.
void ReduceBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
This class is designed to manage and store system matrix elements associated to a vEvent...
FLTNB GetVoxelWeights(int a_direction, int a_TOFBin, INTNB a_voxelInLine)
void InstantiateBackwardImageFromDynamicBins()
Allocate memory for the backward image matrices and initialize them.
void LMS_DeallocateAttenuationImage()
Free memory for the Attenuation image matrices (for analytical projection or list-mode sensitivity ge...
bool MergeDynImages()
Indicate if a dynamic serie of 3D images should be merged in one file (true) or written on disk as on...
void SetBedIndex(int a_bedIndex)
bool UseDeformationResp()
Indicate if the respiratory motion deformation is enabled.
Mother class for the Event objects.
int GetNbThreadsForImageComputation()
Get the number of threads used for image operations.
int GetNb2ndMotImgsForLMS()
call the eponym function from the oDynamicDataManager object
int InitializeNormalizationFiles()
Initialize the normalization datafiles provided for sensitivity computation.
int CheckParameters()
A public function used to check the parameters settings.
int GetNbTOFBins()
This function is used to get the number of TOF bins in use.
int GetNbThreadsForProjection()
Get the number of threads used for projections.
int64_t PROJ_GetModalityStopValueMainLoop()
Get the stop value for the main loop of analytic projection depending on the modality.
int ConvolveSensitivity(oImageSpace *ap_ImageSpace)
oImageDimensionsAndQuantification * mp_ImageDimensionsAndQuantification
virtual int IsAvailableLOR(int a_elt1, int a_elt2)
bool GetAtnCorrectionFlag()
Simply return m_atnCorrectionFlag.
void SetImageDimensionsAndQuantification(oImageDimensionsAndQuantification *ap_ImageDimensionsAndQuantification)
void DeallocateBackwardImageFromDynamicBins()
Free memory of the backward image matrices.
Inherit from vDataFile. Class that manages the reading of a PET input file (header + data)...
void LMS_InstantiateSensitivityImage()
Allocate memory for the sensitivity image matrices (for list-mode sensitivity generation) ...
#define Cout(MESSAGE)
void LMS_DeallocateForwardImage()
Free memory for the forward image matrices (for list-mode sensitivity generation) ...