CASToR  1.1
Tomographic Reconstruction (PET/SPECT)
 All Classes Files Functions Variables Typedefs Macros Groups Pages
oIterativeAlgorithm.cc
Go to the documentation of this file.
1 
2 /*
3  Implementation of class oIterativeAlgorithm
4 
5  - separators: X
6  - doxygen: X
7  - default initialization: X
8  - CASTOR_DEBUG: X
9  - CASTOR_VERBOSE: none
10 */
11 
18 #include "gVariables.hh"
19 #include "oIterativeAlgorithm.hh"
20 
21 bool m_releaseThreads = false;
24 // =====================================================================
25 // ---------------------------------------------------------------------
26 // ---------------------------------------------------------------------
27 // =====================================================================
28 /*
29  \brief oIterativeAlgorithm constructor.
30  Initialize the member variables to their default values.
31 */
33 {
34  // Default number of iterations and subsets (1 for each)
35  m_nbIterations = 1;
36  mp_nbSubsets = (int*)malloc(1*sizeof(int));
37  mp_nbSubsets[0] = 1;
39  mp_outputIterations = NULL;
40  // set all members to default values
41  m_verbose = -1;
42  m_flagGPU = false;
43  mp_ID = NULL;
44  m2p_DataFile= NULL;
45  mp_ProjectorManager = NULL;
46  mp_OptimizerManager = NULL;
47  mp_DeformationManager = NULL;
49  mp_ImageSpace = NULL;
52  m_nbBeds = -1;
53  m_pathToInitialImg = "";
56  m_pathToMaskImg = "";
59 }
60 
61 // =====================================================================
62 // ---------------------------------------------------------------------
63 // ---------------------------------------------------------------------
64 // =====================================================================
65 /*
66  \brief oIterativeAlgorithm destructor.
67 */
69 {
70  // Delete some tables
71  if (mp_nbSubsets!=NULL) free(mp_nbSubsets);
73 }
74 
75 // =====================================================================
76 // ---------------------------------------------------------------------
77 // ---------------------------------------------------------------------
78 // =====================================================================
79 /*
80  \fn SetNbIterationsAndSubsets
81  \param a_nbIterationsSubsets
82  \brief Set the number of iterations and subsets
83  \details The provided string is a list of iteration:subset pairs
84  separated by commas.
85  \return 0 if success, positive value otherwise
86 */
87 int oIterativeAlgorithm::SetNbIterationsAndSubsets(const string& a_nbIterationsSubsets)
88 {
89  // If the string is empty, then just keep the defaults
90  if (a_nbIterationsSubsets=="") return 0;
91  // Otherwise, reset the number of iterations to 0
92  m_nbIterations = 0;
93 
94  // Copy the string
95  string buf = a_nbIterationsSubsets;
96 
97  // Loop to search for commas
98  size_t pos_comma;
99  while ((pos_comma=buf.find_first_of(","))!=string::npos)
100  {
101  // Get the substring before the comma
102  string sub_buf = buf.substr(0,pos_comma);
103  // Search for columns
104  size_t pos_column = sub_buf.find_first_of(":");
105  if (pos_column==string::npos || pos_column==0)
106  {
107  Cerr("***** oIterativeAlgorithm::SetNbIterationsAndSubsets() -> Syntax problem in number of iterations and subsets !" << endl);
108  return 1;
109  }
110  int iter = atoi( sub_buf.substr(0,pos_column).c_str() );
111  int subs = atoi( sub_buf.substr(pos_column+1).c_str() );
112  mp_nbSubsets = (int*)realloc(mp_nbSubsets,(m_nbIterations+iter)*sizeof(int));
113  for (int it=0; it<iter; it++) mp_nbSubsets[m_nbIterations+it] = subs;
114  m_nbIterations += iter;
115  buf = buf.substr(pos_comma+1);
116  }
117 
118  // Last couple of iterations:subsets
119  size_t pos_column = buf.find_first_of(":");
120  if (pos_column==string::npos || pos_column==0)
121  {
122  Cerr("***** oIterativeAlgorithm::SetNbIterationsAndSubsets() -> Syntax problem in number of iterations and subsets !" << endl);
123  return 1;
124  }
125  int iter = atoi( buf.substr(0,pos_column).c_str() );
126  int subs = atoi( buf.substr(pos_column+1).c_str() );
127  mp_nbSubsets = (int*)realloc(mp_nbSubsets,(m_nbIterations+iter)*sizeof(int));
128  for (int it=0; it<iter; it++) mp_nbSubsets[m_nbIterations+it] = subs;
129  m_nbIterations += iter;
130 
131  if (m_verbose>=3)
132  {
133  Cout("oIterativeAlgorithm::SetNbIterationsAndSubsets() -> Selected numbers of subsets for each iteration:" << endl);
134  Cout(" Iteration / number of subsets "<< endl);
135  for (int it=0 ; it<m_nbIterations ; it++) Cout(" " << it+1 << " / " << mp_nbSubsets[it] << endl);
136  }
137  // End
138  return 0;
139 }
140 
141 // =====================================================================
142 // ---------------------------------------------------------------------
143 // ---------------------------------------------------------------------
144 // =====================================================================
145 /*
146  \fn SetOutputIterations
147  \param a_outputIterations
148  \brief Set the selected output iterations
149  \details The provided string is a list of couple a:b separated by commas.
150  It means that we save one iteration over a until b is reached.
151  "b" must be incrementing for each successive couples.
152  If the list is empty, we save all iterations by default.
153  If the list is equal to "-1", then we save only the last iteration.
154  \return 0 if success, positive value otherwise
155 */
156 int oIterativeAlgorithm::SetOutputIterations(const string& a_outputIterations)
157 {
158  if(m_verbose>=2) Cout("oIterativeAlgorithm::SetOutputIterations ..."<< endl);
159 
160  // Allocate the output iterations boolean table
161  mp_outputIterations = (bool*)malloc(m_nbIterations*sizeof(bool));
162  for (int it=0; it<m_nbIterations; it++) mp_outputIterations[it] = false;
163 
164  // If the list is empty, we save all iterations by default
165  if (a_outputIterations=="")
166  {
167  for (int it=0; it<m_nbIterations; it++) mp_outputIterations[it] = true;
168  return 0;
169  }
170 
171  // Copy the string
172  string buf = a_outputIterations;
173 
174  // Loop to search for commas
175  size_t pos_comma;
176  int current_iteration = -1;
177  while ((pos_comma=buf.find_first_of(","))!=string::npos)
178  {
179  // Get the substring before the comma
180  string sub_buf = buf.substr(0,pos_comma);
181  // Search for columns
182  size_t pos_column = sub_buf.find_first_of(":");
183  if (pos_column==string::npos || pos_column==0)
184  {
185  Cerr("***** oIterativeAlgorithm::SetOutputIterations() -> Syntax problem in output iteration numbers !" << endl);
186  return 1;
187  }
188  int step_iteration = atoi( sub_buf.substr(0,pos_column).c_str() );
189  if (step_iteration<=0)
190  {
191  Cerr("***** oIterativeAlgorithm::SetOutputIterations() -> Iteration step must be strictly positive (here it is " << step_iteration << ") !" << endl);
192  return 1;
193  }
194  int stop_iteration = atoi( sub_buf.substr(pos_column+1).c_str() );
195  if (stop_iteration<=current_iteration)
196  {
197  Cerr("***** oIterationAlgorithm::SetOutputIterations() -> Output iteration number must be increasing through provided couples !" << endl);
198  return 1;
199  }
200  for (int it=current_iteration+step_iteration; it<stop_iteration; it+=step_iteration) mp_outputIterations[it] = true;
201  current_iteration = stop_iteration-1;
202  buf = buf.substr(pos_comma+1);
203  }
204 
205  // Last couple of iterations:subsets
206  size_t pos_column = buf.find_first_of(":");
207  if (pos_column==string::npos || pos_column==0)
208  {
209  // Special case if -1 is provided, it means we save the last iteration
210  if (buf=="-1")
211  {
212  mp_outputIterations[m_nbIterations-1] = true;
213  // We directly exist here
214  return 0;
215  }
216  else
217  {
218  Cerr("***** oIterativeAlgorithm::SetOutputIterations() -> Syntax problem in output iteration numbers !" << endl);
219  return 1;
220  }
221  }
222  int step_iteration = atoi( buf.substr(0,pos_column).c_str() );
223  if (step_iteration<=0)
224  {
225  Cerr("***** oIterativeAlgorithm::SetOutputIterations() -> Iteration step must be strictly positive (here it is " << step_iteration << ") !" << endl);
226  return 1;
227  }
228  int stop_iteration = atoi( buf.substr(pos_column+1).c_str() );
229  if (stop_iteration<=current_iteration)
230  {
231  Cerr("***** oIterationAlgorithm::SetOutputIterations() -> Output iteration number must be increasing through provided couples !" << endl);
232  return 1;
233  }
234  for (int it=current_iteration+step_iteration; it<stop_iteration; it+=step_iteration) mp_outputIterations[it] = true;
235 
236  // End
237  return 0;
238 }
239 
240 // =====================================================================
241 // ---------------------------------------------------------------------
242 // ---------------------------------------------------------------------
243 // =====================================================================
244 /*
245  \fn Iterate
246  \brief Just call either the IterateCPU or the IterateGPU function as asked for
247  \return 0 if success, positive value otherwise
248 */
250 {
251  if(m_verbose>=2) Cout("oIterativeAlgorithm::Iterate ..."<< endl);
252 
253  #ifdef CASTOR_GPU
254  if (m_flagGPU)
255  return IterateGPU();
256  else
257  return IterateCPU();
258  #else
259  return IterateCPU();
260  #endif
261 }
262 
263 // =====================================================================
264 // ---------------------------------------------------------------------
265 // ---------------------------------------------------------------------
266 // =====================================================================
267 /*
268  \fn IterateCPU
269  \brief Perform the iterative loop of the algorithm, call the different object
270  for optimization, projection, update, etc.
271  Function designed to be executed on the CPU only.
272  \details Loops over the iterations, subsets, bed position
273  Call functions related to each steps of iterative reconstruction:
274 
275  StepBeforeIterationLoop()
276  / Loop on iterations
277  | StepBeforeSubsetLoop(iteration)
278  | / Loop on subsets
279  | | StepPreProcessInsideSubsetLoop(iteration,subset)
280  | | / Loop on bed positions
281  | | | StepInnerLoopInsideSubsetLoop(iteration,subset,bed)
282  | | StepPostProcessInsideSubsetLoop(iteration,subset)
283  | StepAfterSubsetLoop(iteration)
284  StepAfterIterationLoop()
285 
286  \return 0 if success, positive value otherwise
287 */
289 {
290  // Synchronize MPI processes
291  #ifdef CASTOR_MPI
292  MPI_Barrier(MPI_COMM_WORLD);
293  #endif
294 
295  // Verbose
296  if (m_verbose>=1) Cout("oIterativeAlgorithm::IterateCPU() -> Start algorithm for " << m_nbIterations << " iterations" << endl);
297 
298  // Initial clock for execution time
299  clock_t clock_start_whole = clock();
300  time_t time_start_whole = time(NULL);
301 
302  // Call before iteration function
304  {
305  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepBeforeIterationLoop() function !" << endl);
306  return 1;
307  }
308 
309  // Loop on iterations
310  for (int iteration=0 ; iteration<m_nbIterations ; iteration++)
311  {
312 
313  // Call before subset function
314  if (StepBeforeSubsetLoop(iteration))
315  {
316  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepBeforeSubsetLoop() function !" << endl);
317  return 1;
318  }
319 
320  // Loop on subsets
321  for (int subset=0 ; subset<mp_nbSubsets[iteration] ; subset++)
322  {
323  // Verbose
324  if (m_verbose>=1)
325  {
326  Cout("oIterativeAlgorithm::IterateCPU() -> Start iteration " << iteration+1 << "/" << m_nbIterations
327  << " subset " << subset+1 << "/" << mp_nbSubsets[iteration] << endl);
328  }
329 
330  // Clock start for subset execution time
331  clock_t clock_start_subset = clock();
332  time_t time_start_subset = time(NULL);
333 
334  // Call pre-process inside subset loop function
335  if (StepPreProcessInsideSubsetLoop(iteration,subset))
336  {
337  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepPreProcessInsideSubsetLoop() function !" << endl);
338  return 1;
339  }
340 
341  // Synchronize MPI processes
342  #ifdef CASTOR_MPI
343  MPI_Barrier(MPI_COMM_WORLD);
344  #endif
345 
346  // Loop on bed positions
347  for (int bed=0 ; bed<m_nbBeds ; bed++)
348  {
349  // Call the inner loop on events inside subset loop function
350  if (StepInnerLoopInsideSubsetLoop(iteration,subset,bed))
351  {
352  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepInnerLoopInsideSubsetLoop() function !" << endl);
353  return 1;
354  }
355  } // End of beds loop
356 
357  // Synchronize MPI processes
358  #ifdef CASTOR_MPI
359  MPI_Barrier(MPI_COMM_WORLD);
360  #endif
361 
362  // Call post-process inside subset loop function
363  if (StepPostProcessInsideSubsetLoop(iteration,subset))
364  {
365  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepPostProcessInsideSubsetLoop() function !" << endl);
366  return 1;
367  }
368 
369  // Clock stop for subset execution time
370  clock_t clock_stop_subset = clock();
371  time_t time_stop_subset = time(NULL);
372  if (m_verbose>=2) Cout ("oIterativeAlgorithm::IterateCPU() -> Time spent for subset " << subset+1 << " | User: " << time_stop_subset-time_start_subset
373  << " sec | CPU: " << (clock_stop_subset-clock_start_subset)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
374 
375 
376  } // End of subsets loop
377 
378 
379  // Call after subset function
380  if (StepAfterSubsetLoop(iteration))
381  {
382  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepAfterSubsetLoop() function !" << endl);
383  return 1;
384  }
385 
386 
387  } // End of iterations loop
388 
389  // Call after iteration function
391  {
392  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepAfterIterationLoop() function !" << endl);
393  return 1;
394  }
395 
396  // Final clock for execution time
397  clock_t clock_stop_whole = clock();
398  time_t time_stop_whole = time(NULL);
399  if (m_verbose>=1) Cout("oIterativeAlgorithm::IterateCPU() -> Total time spent | User: " << time_stop_whole-time_start_whole
400  << " sec | CPU: " << (clock_stop_whole-clock_start_whole)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
401 
402 
403  return 0;
404 }
405 
406 // =====================================================================
407 // ---------------------------------------------------------------------
408 // ---------------------------------------------------------------------
409 // =====================================================================
410 /*
411  \fn StepBeforeIterationLoop
412  \brief This function is called at the beginning of the IterateCPU function.
413  \details Initialization and memory allocation for the imageSpace and
414  some managers
415  \return 0 if success, positive value otherwise.
416 */
418 {
419  if (m_verbose>=2) Cout("oIterativeAlgorithm::StepBeforeIterationLoop ... " << endl);
420 
421  // Check that the image space is already allocated
422  if (mp_ImageSpace==NULL || !mp_ImageSpace->Checked())
423  {
424  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop() -> Image space has not been created or was not checked !" << endl);
425  return 1;
426  }
427 
428  // Instantiate all the required images from the oImageSpace
436 
437  // Main image and sensitivity image initialization
439  {
440  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop() -> An error occured while reading the initialization image !" << endl);
441  return 1;
442  }
443 
445  {
446  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop()-> Error during attenuation image initialization !" << endl);
447  return 1;
448  }
449 
451  {
452  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop() -> An error occured while initializing the sensitivity image !" << endl);
453  return 1;
454  }
455 
457  {
458  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop()-> Error during anatomical image initialization !" << endl);
459  return 1;
460  }
461 
463  {
464  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop()-> Error during mask image initialization !" << endl);
465  return 1;
466  }
467 
468  // End
469  return 0;
470 }
471 
472 // =====================================================================
473 // ---------------------------------------------------------------------
474 // ---------------------------------------------------------------------
475 // =====================================================================
476 /*
477  \fn StepBeforeSubsetLoop
478  \param a_iteration : iteration index
479  \brief This function is called before starting the subset loop.
480  \return 0 if success, positive value otherwise.
481 */
483 {
484  if (m_verbose>=3) Cout("oIterativeAlgorithm::StepBeforeSubsetLoop ... " << endl);
485 
486  // End
487  return 0;
488 }
489 
490 // =====================================================================
491 // ---------------------------------------------------------------------
492 // ---------------------------------------------------------------------
493 // =====================================================================
494 /*
495  \fn StepPreProcessInsideSubsetLoop
496  \param a_iteration : iteration index
497  \param a_subset : subset index
498  \brief This function is called right after starting the subset loop.
499  Apply any kind of processing on the forward image before projections
500  \details Copy current main image into forward image matrix
501  Reinitialize backward image and 4D gating indices
502  Apply image processing/convolution on the forward image matrix
503  (image to be projeted)
504  \return 0 if success, positive value otherwise.
505 */
506 int oIterativeAlgorithm::StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
507 {
508  if (m_verbose>=3) Cout("oIterativeAlgorithm::StepPreProcessInsideSubsetLoop ... " << endl);
509 
510  // Reinitialize 4D gate indexes
512 
513  // Initialize the correction backward image(s)
516 
517  // Copy current image in forward-image buffer (apply deformation if needed)
519 
520  // Apply image processing to forward image
522  {
523  Cerr("***** oIterativeAlgorithm::StepPreProcessInsideSubsetLoop() -> A problem occured while applyin image processing to forward image !" << endl);
524  return 1;
525  }
526 
527  // Apply convolver to forward image
529  {
530  Cerr("***** oIterativeAlgorithm::StepPreProcessInsideSubsetLoop() -> A problem occured while applying image convolver to forward image !" << endl);
531  return 1;
532  }
533 
534  // Call the pre-event loop function from the optimizer manager
535  mp_OptimizerManager->PreDataUpdateStep(a_iteration, m_nbIterations, a_subset, mp_nbSubsets[a_iteration]);
536 
537  // End
538  return 0;
539 }
540 
541 // =====================================================================
542 // ---------------------------------------------------------------------
543 // ---------------------------------------------------------------------
544 // =====================================================================
545 /*
546  \fn StepInnerLoopInsideSubsetLoop
547  \param a_iteration : iteration index
548  \param a_subset : subset index
549  \param a_bed : bed position
550  \brief This function is called inside the subset loop and manages the main loop over the events
551  The loop over the events is multithreaded, and involves a thread lock in the case of image-based deformation
552  \details Step 1: Get the current event for that thread index
553  Step 2: Call the dynamic switch function that updates the current frame and gate numbers, and also detects involuntary patient motion
554  Perform image-based deformation if needed (thread lock is required)
555  Step 3: Compute the projection line
556  Step 4: Optimize in the data space (forward-proj, update, backward-proj)
557  \todo Check the correct implementation of thread-locked operation on different architectures
558  \return 0 if success, positive value otherwise.
559 */
560 int oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop(int a_iteration, int a_subset, int a_bed)
561 {
562  // Verbose
564  {
565  if (m_nbBeds>1) Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Start loop over events for bed " << a_bed+1 << endl << flush);
566  else Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Start loop over events" << endl << flush);
567  }
568 
569  // Apply the bed offset for this bed position
571 
572  // Progression (increments of 2%)
574  {
575  cout << "0 10% 20% 30% 40% 50% 60% 70% 80% 90% 100%" << endl;
576  cout << "|----|----|----|----|----|----|----|----|----|----|" << endl;
577  cout << "|" << flush;
578  }
579 
580  int progression_percentage_old = 0;
581  int progression_nb_bars = 0;
582  uint64_t progression_printing_index = 0;
583 
584  // Compute start and stop indices taking MPI into account (the vDataFile does that)
585  int64_t index_start = 0;
586  int64_t index_stop = 0;
587  m2p_DataFile[a_bed]->GetEventIndexStartAndStop(&index_start, &index_stop, a_subset, mp_nbSubsets[a_iteration]);
588 
589  // Multi-threading with OpenMP
590  m_releaseThreads = false;
591  m_nbThreadsWaiting = 0;
592 
593  // Reset the file buffer range before launching the loop over the datafile (this is done only if the percentage load is under 100)
594  m2p_DataFile[a_bed]->ResetBufferRange();
595 
596  // Synchronize MPI processes
597  #ifdef CASTOR_MPI
598  MPI_Barrier(MPI_COMM_WORLD);
599  #endif
600 
601  // Set the number of threads for projections (right after this loop, we set back the number of threads for image computations)
602  #ifdef CASTOR_OMP
603  omp_set_num_threads(mp_ID->GetNbThreadsForProjection());
604  #endif
605 
606  // This boolean is used to report any problem inside the parallel loop
607  bool problem = false;
608 
609  // Compute the last index that each thread will manage here (this is to know when it has finished reading the datafile).
610  // Note that everything here assumes a static scheduling with a chunck size of 1, so do not change that in the parallel loop below!
611  int64_t* p_last_index = (int64_t*)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(int64_t));
612  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
613  {
614  int64_t loop_increment = mp_nbSubsets[a_iteration];
615  int64_t increment_for_this_thread = loop_increment * mp_ID->GetNbThreadsForProjection();
616  int64_t index_start_for_this_thread = index_start + th * loop_increment;
617  int64_t nb_loop_pass_through = (index_stop - 1 - index_start_for_this_thread) / increment_for_this_thread;
618  p_last_index[th] = index_start_for_this_thread + nb_loop_pass_through * increment_for_this_thread;
619  }
620 
621  // Launch the loop with precomputed start and stop and using a step equal to the number of subsets
622  int64_t index;
623  // Keep the static scheduling with a chunk size at 1, it is important
624  #pragma omp parallel for private(index) schedule(static, 1)
625  for ( index = index_start ; index < index_stop ; index += mp_nbSubsets[a_iteration] )
626  {
627  // Get the thread index
628  int th = 0;
629  #ifdef CASTOR_OMP
630  th = omp_get_thread_num();
631  #endif
632  // Print progression (do not log out with Cout here)
633  if (m_verbose>=2 && th==0 && mp_ID->GetMPIRank()==0)
634  {
635  if (progression_printing_index%1000==0)
636  {
637 // int progression_percentage_new = ((int)( (((float)(index-index_start+1))/((float)(index_stop-index_start)) ) *
638 // ( ((float)(a_bed+1))/((float)m_nbBeds) ) * 100. ));
639  int progression_percentage_new = ((int)( (((float)(index-index_start+1))/((float)(index_stop-index_start)) ) * 100.));
640  if (progression_percentage_new>=progression_percentage_old+2) // Increments of 2%
641  {
642  int nb_steps = (progression_percentage_new-progression_percentage_old)/2;
643  for (int i=0; i<nb_steps; i++)
644  {
645  cout << "-" << flush;
646  progression_nb_bars++;
647  }
648  progression_percentage_old += nb_steps*2;
649  }
650  }
651  progression_printing_index++;
652  }
653 
654  // Step 1: Get the current event for that thread index
655  #ifdef CASTOR_DEBUG
656  if (m_verbose>=4) Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step1: Get current event for that thread index " << endl);
657  #endif
658  #ifdef CASTOR_MPI
659  vEvent* event = m2p_DataFile[a_bed]->GetEventWithoutOrderAssumption(index, th);
660  #else
661  // This function is not compatible with MPI, as the first index of the first thread may not be the index of the first event in the portion
662  // of the datafile manage by the MPI instance (due to subsets). Will be investigated in future releases.
663  vEvent* event = m2p_DataFile[a_bed]->GetEventWithAscendingOrderAssumption(index, p_last_index[th], th);
664  #endif
665  if (event==NULL)
666  {
667  Cerr("***** oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> An error occured while getting the event from index "
668  << index << " (thread " << th << ") !" << endl);
669  // Specify that there was a problem
670  problem = true;
671  // We must continue here because we are inside an OpenMP loop
672  continue;
673  }
674 
675  // Step 2: Call the dynamic switch function that updates the current frame and gate numbers, and also detects involuntary patient motion
676  #ifdef CASTOR_DEBUG
677  if (m_verbose>=4)
678  {
679  Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step2: Check for Dynamic event (frame/gate switch, image-based deformation " << endl);
680  }
681  #endif
682  int dynamic_switch_value = mp_ID->DynamicSwitch(index, event->GetTimeInMs(), a_bed, th);
683  // If the DYNAMIC_SWITCH_CONTINUE is returned, then it means that we are not yet at the first frame
684  if ( dynamic_switch_value == DYNAMIC_SWITCH_CONTINUE )
685  {
686  // Then we just skip this event
687  continue;
688  }
689  // Else, if the DYNAMIC_SWITCH_DEFORMATION is returned, then it means that a change of gate or involuntary motion has occured
690  // and is associated to a deformation
691  else if ( dynamic_switch_value == DYNAMIC_SWITCH_DEFORMATION )
692  {
693  // With multi-threads, we should wait all threads to be there to perform the deformation
694  #ifdef CASTOR_OMP
695  // Count the number of threads reaching this point
697  // Loop to lock the threads until the m_releaseThreads condition is enabled
698  while (ThreadBarrierFlag())
699  {
700  // Allow thread 0 to enter once all threads are in the loop
701  if ( (th==0) && ThreadBarrierCheck() )
702  {
703  // Perform here any deformation on the forward image.
705  // Release the threads and set the thread counter to 0
707  }
708  }
709  // Count the number of threads reaching this point
711  // Set the thread lock condition once all threads reached this point
713  #else
714  // Perform here any resp related deformations on the forward image
716  #endif
717  }
718 
719  // Step 3: Compute the projection line
720  #ifdef CASTOR_DEBUG
721  if (m_verbose>=4) Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step3: Compute the projection line " << endl);
722  #endif
724  if (line==NULL)
725  {
726  Cerr("***** oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> A problem occured while computing the projection line !" << endl);
727  // Specify that there was a problem
728  problem = true;
729  // We must continue here because we are inside an OpenMP loop
730  continue;
731  }
732 
733  // Step 4: Optimize in the data space (forward-proj, update, backward-proj)
734  #ifdef CASTOR_DEBUG
735  if (m_verbose>=4) Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step4: Optimize in the data space " << endl);
736  #endif
737  if (line->NotEmptyLine()) mp_OptimizerManager->DataUpdateStep( line,
738  mp_ImageSpace,
739  event,
740  a_bed,
744  a_iteration,
745  th);
746  } // End of indices loop (OpenMP stops here)
747 
748  // Synchronize MPI processes
749  #ifdef CASTOR_MPI
750  MPI_Barrier(MPI_COMM_WORLD);
751  #endif
752 
753  // Set back the number of threads for image computation
754  #ifdef CASTOR_OMP
755  omp_set_num_threads(mp_ID->GetNbThreadsForImageComputation());
756  #endif
757 
758  // Free the last_index table
759  free(p_last_index);
760 
761  // End of progression printing (do not log out with Cout here)
762  if (m_verbose>=2 && mp_ID->GetMPIRank()==0)
763  {
764  int progression_total_bars = 49;
765  for (int i=0; i<progression_total_bars-progression_nb_bars; i++) cout << "-";
766  cout << "|" << endl;
767  }
768 
769  // If a problem was encountered, then report it here
770  if (problem)
771  {
772  Cerr("***** oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> A problem occured inside the parallel loop over events !" << endl);
773  return 1;
774  }
775 
776  // End
777  return 0;
778 }
779 
780 // =====================================================================
781 // ---------------------------------------------------------------------
782 // ---------------------------------------------------------------------
783 // =====================================================================
784 /*
785  \fn StepPreProcessInsideSubsetLoop
786  \param a_iteration : iteration index
787  \param a_subset : subset index
788  \brief This function is called right after starting the subset loop.
789  Apply any kind of image processing on the backward image and
790  main image after backprojections
791  \details Synchronize parallel results
792  Apply image deformation/processing/convolution on the backward image
793  Update the image using the optimizer functions
794  Apply dynamic model/image processing/convolution on the main image
795  \return 0 if success, positive value otherwise.
796 */
798 {
799  if (m_verbose>=3) Cout("oIterativeAlgorithm::StepPostProcessInsideSubsetLoop ... " << endl);
800 
801  // Merge parallel results
803 
804  // Perform here any image-based deformations related to motion on the backward image.
806 
807  // Apply convolver to backward images and sensitivity-if-needed
809  {
810  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying convolver to backward images !" << endl);
811  return 1;
812  }
813 
814  // Call the post data update step function from the optimizer manager
815  mp_OptimizerManager->PostDataUpdateStep(a_iteration, m_nbIterations, a_subset, mp_nbSubsets[a_iteration]);
816 
817  // Optimize in the image space (apply corrections, MAP and sensitivity); pass the number of subsets for list-mode sensitivity scaling
818  mp_OptimizerManager->ImageUpdateStep(mp_ImageSpace, a_iteration, mp_nbSubsets[a_iteration]);
819 
820  // Apply convolver to current estime images
822  {
823  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying convolver to current estimate images !" << endl);
824  return 1;
825  }
826 
827  // Post-update Dynamic Modeling step (if enabled). Manage either linear/non linear dynamic model (physiological or not)
828  if (mp_DynamicModelManager->ApplyDynamicModel(mp_ImageSpace, a_iteration, mp_nbSubsets[a_iteration]))
829  {
830  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying dynamic model to current estimate images !" << endl);
831  return 1;
832  }
833 
834  // Apply image processing to current estime images
836  {
837  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying image processing to current estimate images !" << endl);
838  return 1;
839  }
840 
841  // Save the sensitivity image in histogram mode, if asked for
843  {
844  // Get output manager to build the file name
845  sOutputManager* p_output_manager = sOutputManager::GetInstance();
846  // Build the file name
847  string temp_sens = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
848  stringstream temp_it; temp_it << a_iteration + 1;
849  stringstream temp_ss; temp_ss << a_subset + 1;
850  temp_sens.append("_it").append(temp_it.str()).append("_ss").append(temp_ss.str()).append("_sensitivity");
851  // Save sensitivity
853  }
854 
855  // Save the current image estimate if asked for
856  if (mp_ID->GetMPIRank()==0 && mp_outputIterations[a_iteration] && m_saveImageAfterSubsets)
857  {
858  // Verbose
859  if (m_verbose>=1) Cout("oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> Save image at iteration " << a_iteration+1 << " for subset " << a_subset+1 << endl);
860  // Build image from basis functions
862  // Apply post-convolution if needed
864  {
865  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while convolving the output image !" << endl);
866  return 1;
867  }
868  // Apply post-processing if needed
870  {
871  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying image processing the output image !" << endl);
872  return 1;
873  }
874  // Apply output transaxial FOV masking
876  {
877  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying output FOV masking !" << endl);
878  return 1;
879  }
880  // Apply output flip
882  {
883  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying output flip !" << endl);
884  return 1;
885  }
886  // Save output image
887  if (mp_ImageSpace->SaveOutputImage(a_iteration, a_subset))
888  {
889  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while saving output image !" << endl);
890  return 1;
891  }
892  }
893 
894  // End
895  return 0;
896 }
897 
898 // =====================================================================
899 // ---------------------------------------------------------------------
900 // ---------------------------------------------------------------------
901 // =====================================================================
902 /*
903  \fn StepAfterSubsetLoop
904  \param a_iteration : iteration index
905  \brief This function is called after finishing the subset loop.
906  \details Clean the main images from never visited voxels
907  Apply post-convolution/post-processing if needed
908  Write output images on disk as requested by the user
909  \todo : save parametric images using intrinsic basis functions ?
910  \return 0 if success, positive value otherwise.
911 */
913 {
914  if (m_verbose>=3) Cout("oIterativeAlgorithm::StepAfterSubsetLoop ... " << endl);
915 
916  // Clean never visited voxels
918  // Save the main image
919  if (mp_ID->GetMPIRank()==0 && mp_outputIterations[a_iteration])
920  {
921  if (m_verbose>=1) Cout("oIterativeAlgorithm::StepAfterSubsetLoop() -> Save image at iteration " << a_iteration+1 << endl);
922  // Build image from basis functions
924  // Apply post-convolution if needed
926  {
927  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while convolving the output image !" << endl);
928  return 1;
929  }
930  // Apply post-processing if needed
932  {
933  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while applying image processing the output image !" << endl);
934  return 1;
935  }
936  // Save Dynamic images
938  {
939  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while saving parametric images related to the dynamic model !" << endl);
940  return 1;
941  }
942  // Todo : save parametric images using intrinsic basis functions
943  // Apply output transaxial FOV masking
945  {
946  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while applying output FOV masking !" << endl);
947  return 1;
948  }
949  // Apply output flip
951  {
952  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while applying output flip !" << endl);
953  return 1;
954  }
955  // Save output image
956  if (mp_ImageSpace->SaveOutputImage(a_iteration))
957  {
958  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while saving output image !" << endl);
959  return 1;
960  }
961  }
962 
963  // End
964  return 0;
965 }
966 
967 // =====================================================================
968 // ---------------------------------------------------------------------
969 // ---------------------------------------------------------------------
970 // =====================================================================
971 /*
972  \fn StepAfterIterationLoop
973  \brief This function is called at the end of the IterateCPU function.
974  \details Free Memory for the imageSpace and
975  some managers
976  \return 0 if success, positive value otherwise.
977 */
979 {
980  if (m_verbose>=2) Cout("oIterativeAlgorithm::StepAfterIterationLoop ... " << endl);
981 
982  // Deallocate everything
992 
993  // Normal end
994  return 0;
995 }
996 
997 // =====================================================================
998 // ---------------------------------------------------------------------
999 // ---------------------------------------------------------------------
1000 // =====================================================================
1001 /*
1002  \fn InitSpecificOptions
1003  \param a_specificOptions
1004  \brief Set the selected output iterations
1005  \details Not Implemented yet
1006  \return 0
1007 */
1008 int oIterativeAlgorithm::InitSpecificOptions(string a_specificOptions)
1009 {
1010  if (m_verbose>=2) Cout("oIterativeAlgorithm::InitSpecificOptions ... " << endl);
1011 
1012  return 0;
1013 }
1014 
1015 // =====================================================================
1016 // ---------------------------------------------------------------------
1017 // ---------------------------------------------------------------------
1018 // =====================================================================
1019 // ----- Functions dedicated to barrier for multithreading -----
1020 // TODO we should keep in mind that different compiler may deal differently with OpenMP
1021 // (as the compiler may take liberties and optimize shared variables as register variables,
1022 // leading to non-concurrent observations of the variable.
1023 // And this is not cool as this can completely stop the reconstruction process)
1024 // The omp-flush directives should prevent that, but further tests may be required
1025 
1026 /*
1027  \fn ThreadBarrierIncrement
1028  \brief Increment (thread safe) the m_nbThreadsWaiting variable.
1029 */
1031 {
1032  #ifdef CASTOR_DEBUG
1033  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierIncrement ... " << endl);
1034  #endif
1035 
1036  #pragma omp critical(increment_waiting_threads)
1037  {
1039  }
1040 }
1041 
1042 // =====================================================================
1043 // ---------------------------------------------------------------------
1044 // ---------------------------------------------------------------------
1045 // =====================================================================
1046 /*
1047  \fn ThreadBarrierFlag
1048  \brief Check if the m_releaseThreads boolean is enabled or not
1049  \return true if the threads can be released, false otherwise
1050 */
1052 {
1053  #ifdef CASTOR_DEBUG
1054  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierFlag ... " << endl);
1055  #endif
1056 
1057  // Check any change in m_releaseThreads
1058  #pragma omp flush(m_releaseThreads)
1059  return !m_releaseThreads;
1060 }
1061 
1062 // =====================================================================
1063 // ---------------------------------------------------------------------
1064 // ---------------------------------------------------------------------
1065 // =====================================================================
1066 /*
1067  \fn ThreadBarrierCheck
1068  \brief Check if all the threads are currently in waiting condition
1069  \return true if all the threads are waiting, false otherwise
1070 */
1072 {
1073  #ifdef CASTOR_DEBUG
1074  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierCheck ... " << endl);
1075  #endif
1076 
1077  #pragma omp flush(m_nbThreadsWaiting) // Check any change in m_nbThreadsWaiting
1079 }
1080 
1081 // =====================================================================
1082 // ---------------------------------------------------------------------
1083 // ---------------------------------------------------------------------
1084 // =====================================================================
1085 /*
1086  \fn ThreadBarrierSetOff
1087  \brief Disable the thread locking variable (thread safe), and reset the number of waiting threads to 0
1088 */
1090 {
1091  #ifdef CASTOR_DEBUG
1092  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierSetOff ... " << endl);
1093  #endif
1094 
1095  m_releaseThreads = true;
1096  #pragma omp flush(m_releaseThreads) // Make sure the value of m_releaseThreads is propagated to other threads
1097  m_nbThreadsWaiting = 0;
1098  #pragma omp flush(m_nbThreadsWaiting) // Same for m_nbThreadsWaiting
1099 }
1100 
1101 // =====================================================================
1102 // ---------------------------------------------------------------------
1103 // ---------------------------------------------------------------------
1104 // =====================================================================
1105 /*
1106  \fn ThreadBarrierSetOn
1107  \brief Enable the thread locking variable (thread safe), and reset the number of waiting threads to 0
1108 */
1110 {
1111  #ifdef CASTOR_DEBUG
1112  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierSetOn ... " << endl);
1113  #endif
1114 
1115  m_releaseThreads = false;
1116  #pragma omp flush(m_releaseThreads) // Make sure the value of m_releaseThreads is propagated to other threads
1117  m_nbThreadsWaiting = 0;
1118  #pragma omp flush(m_nbThreadsWaiting) // Same for m_nbThreadsWaiting
1119 }
int ApplyProcessingForward(oImageSpace *ap_ImageSpace)
This header file is mainly used to declare some macro definitions and all includes needed from the st...
int GetCurrentCardImage(int a_th)
call the eponym function from the oDynamicDataManager object
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:873
#define MODE_HISTOGRAM
Definition: vDataFile.hh:36
int PerformDeformation(oImageSpace *ap_Image)
Apply deformations during reconstruction.
int Iterate()
Just call either the IterateCPU or the IterateGPU function as asked for.
int SetNbIterationsAndSubsets(const string &a_nbIterationsSubsets)
Set the number of iterations and subsets.
void DeallocateBackwardImageFromDynamicBasis()
Free memory for the backward image matrices.
Definition: oImageSpace.cc:247
vEvent * GetEventWithoutOrderAssumption(int64_t a_eventIndex, int a_th)
According to the current part of the datafile loaded in memory, either read from this buffer or direc...
Definition: vDataFile.cc:734
virtual int StepInnerLoopInsideSubsetLoop(int a_iteration, int a_subset, int a_bed)
This function is called inside the subset loop and manages the main loop over the events The loop o...
int ApplyDynamicModel(oImageSpace *ap_ImageS, int a_ite, int a_sset)
FLTNB GetInitialValue()
Return the initial image value used by the optimizer, explaining why the eponym function of vOptimize...
#define FLTNB
Definition: gVariables.hh:55
int InitMaskImage(const string &a_pathToImage)
Memory allocation and initialization for the mask image.
Definition: oImageSpace.cc:593
void DeallocateImage()
Free memory for the main image matrices.
Definition: oImageSpace.cc:100
#define DYNAMIC_SWITCH_CONTINUE
void Reduce()
Merge parallel results into the matrix of the backward image matrix of the first thread. Also for MPI.
int ConvolveForward(oImageSpace *ap_ImageSpace)
A function used to apply convolvers onto the forward image of the oImageSpace.
virtual int StepAfterIterationLoop()
This function is called at the end of the IterateCPU function.
int DynamicSwitch(int64_t a_currentEventIndex, uint32_t a_currentTime, int a_bed, int a_th)
Call the eponym function from the oDynamicDataManager class using the member object.
int ApplyOutputFOVMasking()
Mask the outside of the transaxial FOV based on the m_fovOutPercent.
void ThreadBarrierSetOff()
Disable the thread locking variable (thread safe), and reset the number of waiting threads to 0...
void ThreadBarrierSetOn()
Enable the thread locking variable (thread safe), and reset the number of waiting threads to 0...
oIterativeAlgorithm()
oIterativeAlgorithm constructor. Initialize the member variables to their default values...
void DeallocateSensitivityImage()
Free memory for the sensitivity image matrices.
Definition: oImageSpace.cc:444
void InstantiateImage()
Allocate memory for the main image matrices.
Definition: oImageSpace.cc:64
void InstantiateBackwardImageFromDynamicBasis(int a_nbBackwardImages)
Allocate memory for the backward image matrices and set the number of backward images for the whole c...
Definition: oImageSpace.cc:202
int SetOutputIterations(const string &a_outputIterations)
Set the selected output iterations.
int PreDataUpdateStep(int a_iteration, int a_nbIterations, int a_subset, int a_nbSubsets)
A function that simply calls the eponym function from the vOptimizer.
int ApplyDeformationsToBackwardImage(oImageSpace *ap_Image)
Apply final backward deformations on the backward image.
bool ThreadBarrierCheck()
Check if all the threads are currently in waiting condition.
void InstantiateSensitivityImage(const string &a_pathToSensitivityImage)
Allocate the sensitivity image matrices.
Definition: oImageSpace.cc:369
int ImageUpdateStep(oImageSpace *ap_Image, int a_iteration, int a_nbSubsets)
A function dedicated to the update step in the image space (performed after the loop on events) ...
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)
Memory allocation and initialisation for the attenuation image using either :
Definition: oImageSpace.cc:966
int ConvolvePost(oImageSpace *ap_ImageSpace)
A function used to apply convolvers onto the output image of the oImageSpace.
int ConvolveIntra(oImageSpace *ap_ImageSpace)
A function used to apply convolvers onto the current image of the oImageSpace.
int ApplyProcessingPost(oImageSpace *ap_ImageSpace)
int InitImage(const string &a_pathToInitialImage, FLTNB a_value)
Initialize the main image, either using:
void DeallocateVisitedVoxelsImage()
Free memory for the image matrix containing binary information regarding which 3D voxels have been vi...
Definition: oImageSpace.cc:909
void DeallocateAnatomicalImage()
Free memory for the anatomical image.
Definition: oImageSpace.cc:563
Declaration of class oIterativeAlgorithm.
virtual int StepAfterSubsetLoop(int a_iteration)
This function is called after finishing the subset loop.
void DeallocateMaskImage()
Free memory for the mask image.
Definition: oImageSpace.cc:634
oImageProcessingManager * mp_ImageProcessingManager
int m_nbThreadsWaiting
bool ThreadBarrierFlag()
Check if the m_releaseThreads boolean is enabled or not.
oImageConvolverManager * mp_ImageConvolverManager
void PrepareForwardImage()
Copy current image matrix in the forward-image buffer matrix.
#define Cerr(MESSAGE)
const string & GetPathName()
void InstantiateVisitedVoxelsImage()
Memory allocation and initialization for the image matrix containing binary information regarding whi...
Definition: oImageSpace.cc:889
void CleanNeverVisitedVoxels()
Based on the visitedVoxelsImage, clean the never visited voxels in the image. This function must be c...
Singleton class that manages output writing on disk (images, sinograms, etc). It also manages loggi...
oDynamicModelManager * mp_DynamicModelManager
int GetNbBackwardImages()
Return the number of backward images used by the optimizer, explaining why the eponym function of vOp...
int PostDataUpdateStep(int a_iteration, int a_nbIterations, int a_subset, int a_nbSubsets)
A function that simply calls the eponym function from the vOptimizer.
int InitAnatomicalImage(const string &a_pathToAnatomicalImage)
Memory allocation and initialization for the anatomical image matrices.
Definition: oImageSpace.cc:510
virtual int StepPostProcessInsideSubsetLoop(int a_iteration, int a_subset)
int SaveOutputImage(int a_iteration, int a_subset=-1)
Call the interfile function to write output image on disk.
oDeformationManager * mp_DeformationManager
void ResetBufferRange()
Simply set the m_1stIdxArrayEvents and m_lastIdxArrayEvents to -1 only if the percentage load is stri...
Definition: vDataFile.cc:491
int DataUpdateStep(oProjectionLine *ap_Line, oImageSpace *ap_Image, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_iteration, int a_thread)
A function dedicated to the update step in the data space (for each event inside the loop) ...
void ResetCurrentDynamicIndices()
Call the eponym function from the oDynamicDataManager class using the member object.
void ComputeOutputImage()
bool m_releaseThreads
oProjectionLine * ComputeProjectionLine(vEvent *ap_Event, int a_th)
This function is used to compute system matrix elements from the associated projector or pre-computed...
int ConvolveBackward(oImageSpace *ap_ImageSpace)
A function used to apply convolvers onto the backward images of the oImageSpace.
#define VERBOSE_NORMAL
void ApplyBedOffset(int a_bed)
Compute the bed offset from the provided bed index and apply it to all projection lines...
void InstantiateForwardImage()
Allocate memory for the forward image matrices.
Definition: oImageSpace.cc:133
bool NotEmptyLine()
This function is used to know if the line contains any voxel contribution.
vEvent * GetEventWithAscendingOrderAssumption(int64_t a_eventIndex, int64_t a_eventIndexLimit, int a_th)
According to the current part of the datafile loaded in memory, either read from this buffer or direc...
Definition: vDataFile.cc:609
void InstantiateImageForDeformation(oImageSpace *ap_Image)
If deformation is enabled, ask the Image Space to Instantiate the temporary backward image for deform...
int ApplyProcessingIntra(oImageSpace *ap_ImageSpace)
virtual int StepBeforeIterationLoop()
This function is called at the beginning of the IterateCPU function.
This class is designed to manage and store system matrix elements associated to a vEvent...
const string & GetBaseName()
int ApplyOutputFlip()
Just flip the output image.
int IterateCPU()
Perform the iterative loop of the algorithm, call the different object for optimization, projection, update, etc. Function designed to be executed on the CPU only.
void InitImageForDeformation(oImageSpace *ap_Image)
If deformation is enabled, ask the Image Space to initialize the temporary backward image for deforma...
void DeallocateOutputImage()
Free memory for the Image matrices dedicated to output writing on disk.
Definition: oImageSpace.cc:708
virtual int InitSpecificOptions(string a_specificOptions)
Set the selected output iterations.
virtual int StepBeforeSubsetLoop(int a_iteration)
This function is called before starting the subset loop.
Mother class for the Event objects.
Definition: vEvent.hh:23
int GetNbThreadsForImageComputation()
Get the number of threads used for image operations.
int SaveParametricImages(int a_ite)
Call SaveParametricImages() function of the dynamic model object is 'm_UseModel' is on...
int GetDataMode()
Definition: vDataFile.hh:280
virtual int StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
This function is called right after starting the subset loop. Apply any kind of processing on the for...
int GetNbThreadsForProjection()
Get the number of threads used for projections.
oOptimizerManager * mp_OptimizerManager
#define Cout(MESSAGE)
oProjectorManager * mp_ProjectorManager
bool Checked()
Simply check that the image dimensions and verbosity has been set.
Definition: oImageSpace.hh:584
void ThreadBarrierIncrement()
Increment (thread safe) the m_nbThreadsWaiting variable.
int SaveSensitivityImage(const string &a_pathToSensitivityImage)
Call the interfile function to write the sensitivity image on disk.
void InitBackwardImage()
Initialize each voxel of the backward images to 0, also for sensitivity if not loaded (estimated on t...
int InitSensitivityImage(const string &a_pathToSensitivityImage)
Initialization for the sensitivity image matrices.
void DeallocateImageForDeformation(oImageSpace *ap_Image)
If deformation is enabled, ask the Image Space to free memory of the temporary backward image for def...
int GetCurrentRespImage(int a_th)
call the eponym function from the oDynamicDataManager object
void DeallocateForwardImage()
Free memory for the forward image matrices.
Definition: oImageSpace.cc:169
int GetMPIRank()
Get the MPI instance number (the rank)
int GetCurrentTimeFrame(int a_th)
call the eponym function from the oDynamicDataManager object
oImageDimensionsAndQuantification * mp_ID
virtual ~oIterativeAlgorithm()
oIterativeAlgorithm destructor.
void InstantiateOutputImage()
Instanciate Image matrices dedicated to output writing on disk.
Definition: oImageSpace.cc:661
#define DYNAMIC_SWITCH_DEFORMATION