CASToR  1.1
Tomographic Reconstruction (PET/SPECT)
 All Classes Files Functions Variables Typedefs Macros Groups Pages
oAnalyticProjection.cc
Go to the documentation of this file.
1 
2 /*
3  Implementation of class oAnalyticProjection
4 
5  - separators:
6  - doxygen:
7  - default initialization:
8  - CASTOR_DEBUG:
9  - CASTOR_VERBOSE:
10 */
11 
12 
21 #include "gVariables.hh"
22 #include "oAnalyticProjection.hh"
23 
24 
25 /*
26  \brief oAnalyticProjection constructor.
27  Initialize the member variables to their default values.
28 */
30 {
31  m_verbose = -1;
32  m_flagGPU = false;
33  mp_ID = NULL;
34  m2p_DataFile = NULL;
35  mp_ProjectorManager = NULL;
37  mp_ImageSpace = NULL;
38  mp_ComputeProjection = NULL;
39  m_nbBeds = -1;
40  m_pathToInitialImg = "";
41  m_pathToAtnImg = "";
42  mp_Scanner = NULL;
43 }
44 
45 
46 
47 /*
48  \brief oAnalyticProjection destructor.
49 */
51 {}
52 
53 
54 
55 /*
56  \fn Launch
57  \brief Just call either the LaunchCPU or the LaunchGPU function as asked for
58  \return 0 if success, positive value otherwise
59 */
61 {
62  #ifdef CASTOR_GPU
63  if (m_flagGPU) return LaunchGPU();
64  else return LaunchCPU();
65  #else
66  return LaunchCPU();
67  #endif
68 }
69 
70 
71 
72 
73 // =====================================================================
74 // ---------------------------------------------------------------------
75 // ---------------------------------------------------------------------
76 // =====================================================================
77 /*
78  \fn LaunchGPU
79  \brief Perform the projection by executing each possible LOR of the scanner,
80  call the different object for optimization, projection, update, etc.
81  Function designed to be executed on the GPU only.
82  Returns error by default as it is not implemented
83  \return 0 if success, positive value otherwise
84 */
85 #ifdef CASTOR_GPU
86 int oAnalyticProjection::LaunchGPU()
87 {
88  Cerr("***** oAnalyticProjection::LaunchGPU() -> The GPU specific function is not implemented !" << endl);
89  return 1;
90 }
91 #endif
92 
93 
94 
95 
96 // =====================================================================
97 // ---------------------------------------------------------------------
98 // ---------------------------------------------------------------------
99 // =====================================================================
100 /*
101  \fn LaunchCPU
102  \brief Perform the projection by executing each possible LOR of the scanner,
103  call the different object for optimization, projection, update, etc.
104  Function designed to be executed on the GPU only.
105  Returns error by default as it is not implemented
106  \return 0 if success, positive value otherwise
107 */
109 {
110  // Verbose
111  if (m_verbose>=1) Cout("oAnalyticProjection::LaunchCPU() -> Start analytic projection" << endl);
112 
113  // Spread verbosity
115 
116  // Image Space allocations
119 
120  // Instanciate and initialize projection image for SPECT
123 
124  // Initialize the image to project
126  {
127  Cerr("***** oAnalyticProjection::LaunchCPU()-> Error during image initialization !" << endl);
128  return 1;
129  }
130 
132  {
133  Cerr("***** oAnalyticProjection::LaunchCPU()-> Error during attenuation image initialization !" << endl);
134  return 1;
135  }
136 
137  // Copy current image in forward-image buffer
139 
141  {
142  Cerr("***** oAnalyticProjection::LaunchCPU() -> A problem occured while applying image convolver to forward image !" << endl);
143  return 1;
144  }
145 
146 
147  // Initial clock
148  clock_t clock_start = clock();
149  time_t time_start = time(NULL);
150 
151  sScannerManager* p_scannerManager;
152  p_scannerManager = sScannerManager::GetInstance();
153 
154  for(int bed=0 ; bed<m_nbBeds ; bed++)
155  {
156  // Initialize main loop start and stop values
157  unsigned int main_loop_start_index = 0 ;
158  unsigned int main_loop_stop_index = 0;
159 
160  // Check beforehand any issue with the loop start/stop values
161  // (not possible to do this check in the inner multithreaded loop)
162  if(p_scannerManager->PROJ_GetModalityStopValueMainLoop() > 0)
163  main_loop_stop_index = p_scannerManager->PROJ_GetModalityStopValueMainLoop();
164  else
165  {
166  Cerr("***** oAnalyticProjection::LaunchCPU()-> An error occured when trying to initialize main loop stop index !" << endl);
167  return 1;
168  }
169 
170  if(p_scannerManager->PROJ_GetModalityStartValueInnerLoop(0) < 0)
171  {
172  Cerr("***** oAnalyticProjection::LaunchCPU()-> An error occured when trying to initialize inner loop start index !" << endl);
173  return 1;
174  }
175 
176  // Prepare pre-computed sums of events to avoid the exponential evolution of the percentage (for PET)
177  // TODO Perhaps replace this with a call to scannerManager for potential error management
178  int32_t nb_total_elts = mp_Scanner->GetSystemNbElts();
179  int64_t* elements_sum = (int64_t*)malloc(nb_total_elts*sizeof(int64_t));
180  elements_sum[0] = 0;
181 
182  for (int64_t idx_elt=1 ; idx_elt<nb_total_elts ; idx_elt++)
183  elements_sum[idx_elt] = elements_sum[idx_elt-1] + (uint64_t)(nb_total_elts-idx_elt);
184 
185  uint64_t* total_events = (uint64_t*)calloc(mp_ID->GetNbThreadsForProjection(),sizeof(uint64_t));
186  double* total_prompts = (double*)calloc(mp_ID->GetNbThreadsForProjection(),sizeof(double));
187 
188  // Index for progression printing
189  uint64_t printing_index = 0;
190  uint64_t progression_index_total = p_scannerManager->PROJ_GetProgressionFinalValue()*
194 
195  // Time reference to initialize the timestamp of the events
196  uint32_t timestamp = 0;
197 
198  // look on the dynamic frames
199  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
200  {
201  uint64_t* total_events_in_frame = (uint64_t*)calloc(mp_ID->GetNbThreadsForProjection(),sizeof(uint64_t));
202  double* total_prompts_in_frame = (double*)calloc(mp_ID->GetNbThreadsForProjection(),sizeof(double));
203 
204  for(int rg=0 ; rg<mp_ID->GetNbRespGates() ; rg++)
205  {
206  uint64_t* total_events_in_rgate = (uint64_t*)calloc(mp_ID->GetNbThreadsForProjection(),sizeof(uint64_t));
207  double* total_prompts_in_rgate = (double*)calloc(mp_ID->GetNbThreadsForProjection(),sizeof(double));
208 
209  for(int cg=0 ; cg<mp_ID->GetNbCardGates() ; cg++)
210  {
211  uint64_t* total_events_in_cgate = (uint64_t*)calloc(mp_ID->GetNbThreadsForProjection(),sizeof(uint64_t));
212  double* total_prompts_in_cgate = (double*)calloc(mp_ID->GetNbThreadsForProjection(),sizeof(double));
213 
214  if(mp_ID->GetNbTimeFrames() > 1)
215  if(m_verbose>=0) Cout(endl << "Processing frame #" << fr+1 << " out of " << mp_ID->GetNbTimeFrames() << endl);
216  if(mp_ID->GetNbRespGates() > 1)
217  if(m_verbose>=0) Cout(endl << "Processing resp gate #" << rg+1 << " out of " << mp_ID->GetNbRespGates() << endl);
218  if(mp_ID->GetNbCardGates() > 1)
219  if(m_verbose>=0) Cout(endl << "Processing card gate #" << cg+1 << " out of " << mp_ID->GetNbCardGates() << endl);
220 
221  // Used to catch errors inside the parallel loop
222  bool problem = false;
223 
224  // Multi-threaded loop on the scanner elements
225  // Use (static,1) as dynamic mode actually different output histogram with multithreading.
226  int64_t idx_elt1; //openMP throws error with unsigned
227  #pragma omp parallel for private(idx_elt1) schedule(static, 1)
228  for (idx_elt1=main_loop_start_index ; idx_elt1<(int)main_loop_stop_index ; idx_elt1++)
229  {
230  // Get the thread index
231  int th = 0;
232  #ifdef CASTOR_OMP
233  th = omp_get_thread_num();
234  #endif
235 
236  // Initialize inner loop start and stop values
237  int64_t inner_loop_start_index = p_scannerManager->PROJ_GetModalityStartValueInnerLoop(idx_elt1) ;
238  int64_t inner_loop_stop_index = mp_Scanner->GetSystemNbElts();
239 
240  // Inner loop on scanner elements (idx_elt2) which are used to form LOR using the first scanner element (idx_elt1)
241  for (int64_t idx_elt2=inner_loop_start_index ; idx_elt2<inner_loop_stop_index ; idx_elt2++)
242  {
243  // Print progression (do not log out with Cout here)
244  if (th==0)
245  {
246  if (printing_index%10000==0)
247  {
248  int64_t progression_index_current = p_scannerManager->PROJ_GetCurrentProgression(idx_elt1, idx_elt2, elements_sum,
249  mp_ID->GetNbRespGates(),
250  mp_ID->GetNbCardGates(),
251  fr, rg, cg);
252 
253  FLTNB percent = ( (FLTNB)progression_index_current /((FLTNB)progression_index_total) ) * ((FLTNB)100);
254  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 "
255  << percent << " % ";
256 
257  }
258  printing_index++;
259  }
260 
261  // Allocate an event using the iDataFile
262  vEvent* event = m2p_DataFile[bed]->PROJ_GenerateEvent(idx_elt1, idx_elt2, th);
263 
264  // TODO : Perhaps replace this with a call to scannerManager for potential error management
265  if (!mp_Scanner->IsAvailableLOR(event->GetID1(0), event->GetID2(0))) continue;
266 
267  // Generate the projection event and compute the projection line
269  if (line==NULL)
270  {
271  Cerr("***** oAnalyticProjection::LaunchCPU() -> A problem occured while computing the projection line !" << endl);
272  // Specify that there was a problem
273  problem = true;
274  // We must continue here because we are inside an OpenMP loop
275  continue;
276  }
277 
278  // Optimize in the data space (forward-proj, update, backward-proj) if the line has been performed
279  if (line->NotEmptyLine())
280  {
281  if(mp_ComputeProjection->DataUpdateStep(m2p_DataFile[bed], line, mp_ImageSpace, event, fr, rg, cg, th, timestamp) )
282  {
283  Cerr("***** oAnalyticProjection::LaunchCPU() -> An error occured during the data update step !" << endl);
284  //return 1; TODO : return in open MP structured block (stop all threads)
285  }
286 
287  // Increment number of events
288  total_events[th]++;
289  total_events_in_cgate[th]++;
290  total_events_in_rgate[th]++;
291  total_events_in_frame[th]++;
292  double event_value = event->GetEventValue(0);
293  total_prompts[th] += event_value;
294  total_prompts_in_cgate[th] += event_value;
295  total_prompts_in_rgate[th] += event_value;
296  total_prompts_in_frame[th] += event_value;
297  }
298 
299  }
300 
301  } // End thread loop
302 
303  // If a problem was encountered, then report it here
304  if (problem)
305  {
306  Cerr("***** oAnalyticProjection::LaunchCPU() -> A problem occured inside the parallel loop over events !" << endl);
307  return 1;
308  }
309 
310  // Write data for this frame
311  if (m2p_DataFile[bed]->PROJ_WriteData() )
312  {
313  Cerr("***** oAnalyticProjection::LaunchCPU()-> An error occured during the data writing step" << endl);
314  return 1;
315  }
316 
317  // Merge counters
318  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++) total_events_in_cgate[0] += total_events_in_cgate[th];
319  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++) total_prompts_in_cgate[0] += total_prompts_in_cgate[th];
320 
321  if(m_verbose>=2 && mp_ID->GetNbCardGates() > 1)
322  {
323  Cout(endl << "Total events in cgate #" << cg+1 << ": " << total_events_in_cgate[0] << endl);
324  Cout("Total prompts in cgate #" << cg+1 << ": " << total_prompts_in_cgate[0] << endl);
325  }
326 
327  free(total_prompts_in_cgate);
328  free(total_events_in_cgate);
329  }
330 
331  // Merge counters
332  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++) total_events_in_rgate[0] += total_events_in_rgate[th];
333  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++) total_prompts_in_rgate[0] += total_prompts_in_rgate[th];
334 
335  if(m_verbose>=2 && mp_ID->GetNbRespGates() > 1)
336  {
337  Cout(endl << "Total events in rgate #" << rg+1 << ": " << total_events_in_rgate[0] << endl);
338  Cout("Total prompts in rgate #" << rg+1 << ": " << total_prompts_in_rgate[0] << endl);
339  }
340  free(total_prompts_in_rgate);
341  free(total_events_in_rgate);
342  }
343 
344  // Get time reference for this frame in order to initialize the timestamp (in ms) of the events
345  timestamp += mp_ID->GetFrameDurationInMs(bed,fr);
346 
347  // Merge counters
348  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++) total_events_in_frame[0] += total_events_in_frame[th];
349  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++) total_prompts_in_frame[0] += total_prompts_in_frame[th];
350 
351  if(m_verbose>=2 && mp_ID->GetNbTimeFrames() > 1)
352  {
353  Cout(endl << "Total events in frame #" << fr+1 << ": " << total_events_in_frame[0] << endl);
354  Cout("Total prompts in frame #" << fr+1 << ": " << total_prompts_in_frame[0] << endl);
355  }
356 
357  free(total_prompts_in_frame);
358  free(total_events_in_frame);
359  }
360 
361  // Compute simulated acquisition duration for 4D projection
362  int acquisition_duration = 0;
363  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
364  acquisition_duration += mp_ID->GetFrameDurationInMs(bed,fr);
365 
366  // If 3D, set the acquisition duration to 1 by default
367  if(acquisition_duration <= 0 ) acquisition_duration = 1; // Set default value for acquisition duration
368 
369  // End of progression printing (do not log out with Cout here)
370  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 100 % " << endl;
371 
372  // Save a projection image for output visualization (currently only for SPECT)
374 
375  // Merge counters
376  if (m_verbose>=2)
377  {
378  Cout("Total number of projected events:" << endl);
379  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++) Cout(" --> Thread " << th << " | " << total_events[th] << " events" << endl);
380  }
381  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++) total_events[0] += total_events[th];
382  if (m_verbose>=1) Cout("Final total number of events projected: " << total_events[0] << endl);
383  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++) total_prompts[0] += total_prompts[th];
384  if (m_verbose>=1) Cout("Final total number of prompts: " << total_prompts[0] << endl);
385 
386  // Write header (Cdh)
387  m2p_DataFile[bed]->SetDuration(acquisition_duration);
388  m2p_DataFile[bed]->SetNbEvents(total_events[0]);
391 
392  free(total_prompts);
393  free(total_events);
394  free(elements_sum);
395  }
396 
397  // Free memory
402 
403  // Clock total
404  clock_t clock_stop = clock();
405  time_t time_stop = time(NULL);
406  Cout (endl << endl << " Total Time spent - Analytic Projection | User: " << time_stop-time_start << " sec | CPU: " << (clock_stop-clock_start)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
407 
408  return 0;
409 }
oImageConvolverManager * mp_ImageConvolverManager
void LMS_InstantiateImage()
Allocate memory for the main image matrices (for list-mode sensitivity generation) ...
This header file is mainly used to declare some macro definitions and all includes needed from the st...
static sScannerManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
oAnalyticProjection()
oAnalyticProjection constructor. Initialize the member variables to their default values...
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:55
int PROJ_InitImage(const string &a_pathToInitialImage)
Load the initial image for the analytical projection.
void LMS_PrepareForwardImage()
Copy current image in forward-image buffer (for list-mode sensitivity generation) ...
int ConvolveForward(oImageSpace *ap_ImageSpace)
A function used to apply convolvers onto the forward image of the oImageSpace.
int64_t PROJ_GetProgressionFinalValue()
Get numerator value according to the modality to compute percent progression during the projection pr...
virtual int PROJ_WriteHeader()=0
This function is implemented in child classes. Generate a header file according to the projection a...
virtual uint16_t PROJ_GetSPECTNbProjections()
return the total number of projections for a SPECT acquisition
Definition: vScanner.cc:354
int PROJ_SaveProjectionImage()
Save an image of the projected data (for analytic projection)
void SetDuration(FLTNB a_value)
Definition: vDataFile.hh:463
void LMS_InstantiateForwardImage()
Allocate memory for the forward image matrices (for list-mode sensitivity generation) ...
int InitAttenuationImage(const string &a_pathToAtnImage)
Memory allocation and initialisation for the attenuation image using either :
Definition: oImageSpace.cc:966
virtual uint16_t PROJ_GetSPECTNbPixels()
return the total number of pixels for a SPECT reconstruction
Definition: vScanner.cc:374
~oAnalyticProjection()
oAnalyticProjection destructor.
oComputeProjection * mp_ComputeProjection
void PROJ_DeallocateProjectionImage(int a_nbProjs)
Free memory for the projection image for analytical projection.
#define Cerr(MESSAGE)
#define SCANNER_SPECT_CONVERGENT
Singleton class that Instantiate and initialize the scanner object.
void LMS_DeallocateImage()
Free memory for the main image matrices (for list-mode sensitivity generation)
oProjectionLine * ComputeProjectionLine(vEvent *ap_Event, int a_th)
This function is used to compute system matrix elements from the associated projector or pre-computed...
oImageDimensionsAndQuantification * mp_ID
int Launch()
Just call either the LaunchCPU or the LaunchGPU function as asked for.
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:1056
int PROJ_DeleteTmpDatafile()
Delete temporary datafile used for multithreaded output writing if needed.
Definition: vDataFile.cc:1008
bool NotEmptyLine()
This function is used to know if the line contains any voxel contribution.
void PROJ_InstantiateProjectionImage(int a_nbProjs, int a_nbPixels)
Instanciate and initialize projection image for analytical projection.
int DataUpdateStep(vDataFile *ap_DataFile, oProjectionLine *a2p_Line, oImageSpace *ap_Image, vEvent *ap_Event, int a_fr, int a_rg, int a_cg, int th, uint32_t a_timestamp)
Perform the forward-projection, and send the result to the datafile for writing. Implements random ...
void SetNbEvents(int64_t a_value)
initialize the number of events with a int64_t value
Definition: vDataFile.hh:447
virtual int GetSystemNbElts()=0
This is a pure virtual method that must be implemented by children.
This class is designed to manage and store system matrix elements associated to a vEvent...
void LMS_DeallocateAttenuationImage()
Free memory for the Attenuation image matrices (for analytical projection or list-mode sensitivity ge...
Definition: oImageSpace.cc:928
int GetNbCardGates()
Get the number of cardiac gates.
oProjectorManager * mp_ProjectorManager
Mother class for the Event objects.
Definition: vEvent.hh:23
int GetNbTimeFrames()
Get the number of time frames.
int64_t PROJ_GetCurrentProgression(int64_t a_elt1, int64_t a_elt2, int64_t *ap_nbEltsArray, int a_nbRGates, int a_nbCGates, int a_fr, int a_rg, int a_cg)
Get numerator value according to the modality to compute percent progression during the analytical pr...
int GetNbThreadsForProjection()
Get the number of threads used for projections.
Declaration of class oAnalyticProjection.
uint32_t GetFrameDurationInMs(int a_bed, int a_frame)
Get the frame duration for the given bed, in milliseconds as a uint32_t.
int GetNbRespGates()
Get the number of respiratory gates.
#define Cout(MESSAGE)
int64_t PROJ_GetModalityStopValueMainLoop()
Get the stop value for the main loop of analytic projection depending on the modality.
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:115
void SetVerbose(int a_verboseLevel)
set verbosity
int GetScannerType()
Definition: vScanner.hh:218
void LMS_DeallocateForwardImage()
Free memory for the forward image matrices (for list-mode sensitivity generation) ...