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