CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
code/src/algorithm/iIterativeAlgorithm.cc
Go to the documentation of this file.
1 
8 #include "gVariables.hh"
9 #include "iIterativeAlgorithm.hh"
10 #include "iEventHistoCT.hh"
11 
12 // =====================================================================
13 // ---------------------------------------------------------------------
14 // ---------------------------------------------------------------------
15 // =====================================================================
16 
18 {
19 }
20 
21 // =====================================================================
22 // ---------------------------------------------------------------------
23 // ---------------------------------------------------------------------
24 // =====================================================================
25 
27 {
28 }
29 
30 // =====================================================================
31 // ---------------------------------------------------------------------
32 // ---------------------------------------------------------------------
33 // =====================================================================
34 
36 {
38  {
39  Cerr("***** iIterativeAlgorithm::StepBeforeIterationLoop() -> A problem occurred while calling StepBeforeIterationLoop() function !" << endl);
40  return 1;
41  }
42  if (m_verbose>=2) Cout("iIterativeAlgorithm::StepBeforeIterationLoop ... " << endl);
43 
46 
47  // Main image initialization
49  {
50  Cerr("***** iIterativeAlgorithm::StepBeforeIterationLoop() -> An error occurred while reading the initialization image !" << endl);
51  return 1;
52  }
53 
54  // Set numbers of iterations and subsets to the optimizer
56 
57  // End
58  return 0;
59 }
60 
61 // =====================================================================
62 // ---------------------------------------------------------------------
63 // ---------------------------------------------------------------------
64 // =====================================================================
65 
67 {
68  if (m_verbose>=3) Cout("iIterativeAlgorithm::StepBeforeSubsetLoop ... " << endl);
69 
70  // Set the current iteration to the optimizer
72 
73  // End
74  return 0;
75 }
76 
77 // =====================================================================
78 // ---------------------------------------------------------------------
79 // ---------------------------------------------------------------------
80 // =====================================================================
81 
82 int iIterativeAlgorithm::StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
83 {
84  if (m_verbose>=3) Cout("iIterativeAlgorithm::StepPreProcessInsideSubsetLoop ... " << endl);
85 
86  // Set the current subset to the optimizer
88 
89  // Get the chrono manager singleton pointer
90  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
91 
92  // Initialize the correction backward image(s)
94 
95  // Copy current image in forward-image buffer
97 
98  // Apply image processing to forward image
100  {
101  Cerr("***** iIterativeAlgorithm::StepPreProcessInsideSubsetLoop() -> A problem occurred while applyin image processing to forward image !" << endl);
102  return 1;
103  }
104 
105  // Apply convolver to forward image
106  p_ChronoManager->StartConvolution();
108  {
109  Cerr("***** iIterativeAlgorithm::StepPreProcessInsideSubsetLoop() -> A problem occurred while applying image convolver to forward image !" << endl);
110  return 1;
111  }
112  p_ChronoManager->StopConvolution();
113 
114  // Call the pre-event loop function from the optimizer manager
116 
117  // Initialisation of the backup images for deformation
118  // (bu_fward image initialized with current forward image)
119  // (bu_bward and bu_sens images set to 0)
120 
122 
123  // End
124  return 0;
125 }
126 
127 // =====================================================================
128 // ---------------------------------------------------------------------
129 // ---------------------------------------------------------------------
130 // =====================================================================
131 
132 int iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop(int a_iteration, int a_subset, int a_bed)
133 {
134  // Verbose
136  {
137  if (m_nbBeds>1) Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Start loop over events for bed " << a_bed+1 << endl << flush);
138  else Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Start loop over events" << endl << flush);
139  }
140 
141  // Reinitialize 4D gate indexes
143 
144  // Apply the bed offset for this bed position
146 
147  // Progression (increments of 2%)
149  {
150  cout << "0 10% 20% 30% 40% 50% 60% 70% 80% 90% 100%" << endl;
151  cout << "|----|----|----|----|----|----|----|----|----|----|" << endl;
152  cout << "|" << flush;
153  }
154  int progression_percentage_old = 0;
155  int progression_nb_bars = 0;
156  uint64_t progression_printing_index = 0;
157 
158  // Compute start and stop indices taking MPI into account (the vDataFile does that)
159  int64_t index_start = 0;
160  int64_t index_stop = 0;
161  m2p_DataFile[a_bed]->GetEventIndexStartAndStop(&index_start, &index_stop, a_subset, mp_nbSubsets[a_iteration]);
162 
163  // Synchronize MPI processes
164  #ifdef CASTOR_MPI
165  MPI_Barrier(MPI_COMM_WORLD);
166  #endif
167 
168  // Set the number of threads for projections (right after this loop, we set back the number of threads for image computations)
169  #ifdef CASTOR_OMP
170  omp_set_num_threads(mp_ID->GetNbThreadsForProjection());
171  #endif
172 
173  // This boolean is used to report any problem inside the parallel loop
174  bool problem = false;
175 
176  // Get the chrono manager singleton pointer
177  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
178 
179  // These flags (one per thread) are used to signal when events are beyond the last frame time stop
180  bool* p_end_loop_flags = (bool*)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(bool));
181  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++) p_end_loop_flags[th] = false;
182 
183  // Launch the loop with precomputed start and stop and using a step equal to the number of subsets
184  int64_t index;
185  // Keep the static scheduling with a chunk size at 1, it is important
186  #pragma omp parallel for private(index) schedule(static, 1)
187  for ( index = index_start ; index < index_stop ; index += mp_nbSubsets[a_iteration] )
188  {
189  // Get the thread index
190  int th = 0;
191  #ifdef CASTOR_OMP
192  th = omp_get_thread_num();
193  #endif
194  // Print progression (do not log out with Cout here)
195  if (m_verbose>=2 && th==0 && mp_ID->GetMPIRank()==0)
196  {
197  if (progression_printing_index%1000==0)
198  {
199  int progression_percentage_new = ((int)( (((float)(index-index_start+1))/((float)(index_stop-index_start)) ) * 100.));
200  if (progression_percentage_new>=progression_percentage_old+2) // Increments of 2%
201  {
202  int nb_steps = (progression_percentage_new-progression_percentage_old)/2;
203  for (int i=0; i<nb_steps; i++)
204  {
205  cout << "-" << flush;
206  progression_nb_bars++;
207  }
208  progression_percentage_old += nb_steps*2;
209  }
210  }
211  progression_printing_index++;
212  }
213  // Skip to the end
214  if (p_end_loop_flags[th]) continue;
215 
216  // Step 1: Get the current event for that thread index
217  p_ChronoManager->StartIterativeDataUpdateStep1(th);
218  #ifdef CASTOR_DEBUG
219  if (m_verbose>=4) Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step1: Get current event for that thread index " << endl);
220  #endif
221  vEvent* event = m2p_DataFile[a_bed]->GetEvent(index, th);
222  if (event==NULL)
223  {
224  Cerr("***** iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> An error occurred while getting the event from index "
225  << index << " (thread " << th << ") !" << endl);
226  // Specify that there was a problem
227  problem = true;
228  // We must continue here because we are inside an OpenMP loop
229  continue;
230  }
231  p_ChronoManager->StopIterativeDataUpdateStep1(th);
232 
233  // Step 2: Call the dynamic switch function that updates the current frame and gate numbers, and also detects involuntary patient motion
234  p_ChronoManager->StartIterativeDataUpdateStep2(th);
235  #ifdef CASTOR_DEBUG
236  if (m_verbose>=4)
237  {
238  Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step2: Check for Dynamic event (frame/gate switch, image-based deformation " << endl);
239  }
240  #endif
241  int dynamic_switch_value = mp_ID->DynamicSwitch(index, event->GetTimeInMs(), a_bed, th);
242  // If the DYNAMIC_SWITCH_CONTINUE is returned, then it means that we are not yet at the first frame
243  if ( dynamic_switch_value == DYNAMIC_SWITCH_CONTINUE )
244  {
245  // Then we just skip this event
246  continue;
247  }
248  // Else, if the DYNAMIC_SWITCH_DEFORMATION is returned, then it means that a change of gate or involuntary motion has occurred
249  // and is associated to a deformation
250  else if ( dynamic_switch_value == DYNAMIC_SWITCH_DEFORMATION )
251  {
252  #ifdef CASTOR_OMP
253  // set deformation requirement for this thread
254  #pragma omp critical (deformation_requirements_first)
255  {
256  mp_DeformationManager->SetDeformationRequirement(th);
257  }
258  // block this thread until the deformation is done
259  bool wait_for_deformation = true;
260  while (wait_for_deformation)
261  {
262  // the first thread performs the deformation on the forward image
263  // as soon as all the threads are ready for deformation
264  if (th==0)
265  {
266  bool all_threads_ready = false;
267  //check whether all the threads require deformation
268  #pragma omp critical (deformation_requirements_first)
269  {
270  all_threads_ready = mp_DeformationManager->AllThreadsRequireDeformation();
271  }
272  if (all_threads_ready)
273  {
274  // Perform the deformation
276  // Free all the threads so that they can continue
277  #pragma omp critical (deformation_requirements_second)
278  {
279  mp_DeformationManager->UnsetDeformationRequirements();
280  }
281  }
282  }
283  // check whether the deformation is done and the thread can continue
284  #pragma omp critical (deformation_requirements_second)
285  {
286  wait_for_deformation = mp_DeformationManager->GetDeformationRequirement(th);
287  }
288  }
289  #else
290  // Perform here any resp related deformations on the forward image
292  #endif
293 
294  // Restore progression printing
295  if (th==0 && m_verbose>=3 && mp_ID->GetMPIRank()==0)
296  {
297  cout << "0 10% 20% 30% 40% 50% 60% 70% 80% 90% 100%" << endl;
298  cout << "|----|----|----|----|----|----|----|----|----|----|" << endl;
299  int progression_percentage_new = ((int)( (((float)(index-index_start+1))/((float)(index_stop-index_start)) ) * 100.));
300  cout << "|";
301 
302  for (int i=0; i< progression_percentage_new / 2 ; i++) cout << "-";
303  }
304 
305  }
306  // Else, if the DYNAMIC_SWITCH_END_LOOP is returned, set the local end_loop_flag boolean to true (which allows to quickly end the OMP parallel loop)
307  else if ( dynamic_switch_value == DYNAMIC_SWITCH_END_LOOP )
308  {
309  p_end_loop_flags[th] = true;
310  continue;
311  }
312  p_ChronoManager->StopIterativeDataUpdateStep2(th);
313 
314  // Step 3: Compute the projection line
315  p_ChronoManager->StartIterativeDataUpdateStep3(th);
316  #ifdef CASTOR_DEBUG
317  if (m_verbose>=4) Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step3: Compute the projection line " << endl);
318  #endif
320  if (line==NULL)
321  {
322  Cerr("***** iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> A problem occurred while computing the projection line !" << endl);
323  // Specify that there was a problem
324  problem = true;
325  // We must continue here because we are inside an OpenMP loop
326  continue;
327  }
328  p_ChronoManager->StopIterativeDataUpdateStep3(th);
329 
330  // Step 4: Optimize in the data space (forward-proj, update, backward-proj)
331  p_ChronoManager->StartIterativeDataUpdateStep4(th);
332  #ifdef CASTOR_DEBUG
333  if (m_verbose>=4) Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step4: Optimize in the data space " << endl);
334  #endif
335  if (line->NotEmptyLine()) mp_OptimizerManager->DataUpdateStep( line,
336  event,
337  a_bed,
341  th);
342  p_ChronoManager->StopIterativeDataUpdateStep4(th);
343  } // End of indices loop (OpenMP stops here)
344 
345  // Synchronize MPI processes
346  #ifdef CASTOR_MPI
347  MPI_Barrier(MPI_COMM_WORLD);
348  #endif
349 
350  // Set back the number of threads for image computation
351  #ifdef CASTOR_OMP
352  omp_set_num_threads(mp_ID->GetNbThreadsForImageComputation());
353  #endif
354 
355  // End of progression printing (do not log out with Cout here)
356  if (m_verbose>=2 && mp_ID->GetMPIRank()==0)
357  {
358  int progression_total_bars = 49;
359  for (int i=0; i<progression_total_bars-progression_nb_bars; i++) cout << "-";
360  cout << "|" << endl;
361  }
362 
363  // If a problem was encountered, then report it here
364  if (problem)
365  {
366  Cerr("***** iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> A problem occurred inside the parallel loop over events !" << endl);
367  return 1;
368  }
369 
370  // End
371  return 0;
372 }
373 
374 // =====================================================================
375 // ---------------------------------------------------------------------
376 // ---------------------------------------------------------------------
377 // =====================================================================
378 
380 {
381  if (m_verbose>=3) Cout("iIterativeAlgorithm::StepPostProcessInsideSubsetLoop ... " << endl);
382 
383  // Get the chrono manager singleton pointer
384  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
385 
386  // Merge parallel results
388 
389  // Perform here any image-based deformations related to motion on the backward image.
391 
392  // Apply convolver to backward images and sensitivity-if-needed
393  p_ChronoManager->StartConvolution();
395  {
396  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while applying convolver to backward images !" << endl);
397  return 1;
398  }
399  p_ChronoManager->StopConvolution();
400 
401  // Call the pre image update step function from the optimizer manager
403  {
404  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while performing the pre-image-update step !" << endl);
405  return 1;
406  }
407 
408  // Optimize in the image space (apply corrections, MAP and sensitivity); pass the number of subsets for list-mode sensitivity scaling
410  {
411  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while performing the image update step !" << endl);
412  return 1;
413  }
414 
415  // Apply convolver to current estimed images
416  p_ChronoManager->StartConvolution();
418  {
419  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while applying convolver to current estimate images !" << endl);
420  return 1;
421  }
422  p_ChronoManager->StopConvolution();
423 
424  // Post-update Dynamic Modeling step (if enabled). Manage either linear/non linear dynamic model (physiological or not)
425  if (mp_DynamicModelManager->ApplyDynamicModel(mp_ImageSpace, a_iteration, a_subset))
426  {
427  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while applying dynamic model to current estimate images !" << endl);
428  return 1;
429  }
430 
431  // Apply image processing to current estime images
433  {
434  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while applying image processing to current estimate images !" << endl);
435  return 1;
436  }
437 
438  // Save the sensitivity image in histogram mode, if asked for
440  {
441  // Get output manager to build the file name
442  sOutputManager* p_output_manager = sOutputManager::GetInstance();
443  // Build the file name
444  string temp_sens = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
445  stringstream temp_it; temp_it << a_iteration + 1;
446  stringstream temp_ss; temp_ss << a_subset + 1;
447  temp_sens.append("_it").append(temp_it.str()).append("_ss").append(temp_ss.str()).append("_sensitivity");
448  // Save sensitivity
450  }
451 
452  // Save the current image estimate if asked for
453  if (mp_ID->GetMPIRank()==0 && mp_outputIterations[a_iteration] && m_saveImageAfterSubsets)
454  {
455  // Verbose
456  if (m_verbose>=1) Cout("iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> Save image at iteration " << a_iteration+1 << " for subset " << a_subset+1 << endl);
457  // Save image
458  if (StepImageOutput(a_iteration,a_subset))
459  {
460  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while saving images at iteration " << a_iteration+1 << " !" << endl);
461  return 1;
462  }
463  }
464 
465  // End
466  return 0;
467 }
468 
469 // =====================================================================
470 // ---------------------------------------------------------------------
471 // ---------------------------------------------------------------------
472 // =====================================================================
473 
475 {
476  if (m_verbose>=3) Cout("iIterativeAlgorithm::StepAfterSubsetLoop() -> Clean never visited voxels and save images if needed" << endl);
477  // Clean never visited voxels
479  // Save the main image if not already done for each subset
480  if (mp_ID->GetMPIRank()==0 && mp_outputIterations[a_iteration] && !m_saveImageAfterSubsets)
481  {
482  // Verbose
483  if (m_verbose>=1) Cout("iIterativeAlgorithm::StepAfterSubsetLoop() -> Save image at iteration " << a_iteration+1 << endl);
484  // Save image
485  if (StepImageOutput(a_iteration))
486  {
487  Cerr("***** iIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occurred while saving images at iteration " << a_iteration+1 << " !" << endl);
488  return 1;
489  }
490  }
491  // End
492  return 0;
493 }
494 
495 // =====================================================================
496 // ---------------------------------------------------------------------
497 // ---------------------------------------------------------------------
498 // =====================================================================
499 
501 {
503  {
504  Cerr("***** iIterativeAlgorithm::StepAfterIterationLoop() -> A problem occurred while calling StepAfterIterationLoop() function !" << endl);
505  return 1;
506  }
507  if (m_verbose>=2) Cout("iIterativeAlgorithm::StepAfterIterationLoop ... " << endl);
508 
509  // Deallocate everything
512 
513  // Normal end
514  return 0;
515 }
516 
517 // =====================================================================
518 // ---------------------------------------------------------------------
519 // ---------------------------------------------------------------------
520 // =====================================================================
521 
522 int iIterativeAlgorithm::StepImageOutput(int a_iteration, int a_subset)
523 {
524  // Get the chrono manager singleton pointer
525  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
526  // =================================================================================================
527  // Save dynamic coefficient images from the vDynamicModel
528  // =================================================================================================
529  if (mp_DynamicModelManager->SaveParametricImages(a_iteration, a_subset))
530  {
531  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while saving parametric images related to the dynamic model !" << endl);
532  return 1;
533  }
534  // =================================================================================================
535  // Apply pre-processing steps
536  // =================================================================================================
537  // Copy the current image into the forward image
539  // Apply post-convolution if needed
540  p_ChronoManager->StartConvolution();
542  {
543  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while convolving the output image !" << endl);
544  return 1;
545  }
546  p_ChronoManager->StopConvolution();
547  // Apply post-processing if needed
549  {
550  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while applying image processing the output image !" << endl);
551  return 1;
552  }
553  // Apply output transaxial FOV masking
555  {
556  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while applying output FOV masking !" << endl);
557  return 1;
558  }
559  // Apply output flip
561  {
562  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while applying output flip !" << endl);
563  return 1;
564  }
565  // =================================================================================================
566  // If some basis functions are really in use
567  // =================================================================================================
569  {
570  // Save the basis coefficients images if asked for
572  {
573  if (mp_ImageSpace->SaveOutputBasisCoefficientImage(a_iteration, a_subset))
574  {
575  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while saving output image !" << endl);
576  return 1;
577  }
578  }
579  // Compute the output image from the forward image basis functions
581  }
582  // =================================================================================================
583  // Save frames/gates
584  // =================================================================================================
585 
586  // Save output image (note that if no basis functions at all are in use, the the output image already points to the forward image)
587  if (mp_ImageSpace->SaveOutputImage(a_iteration, a_subset))
588  {
589  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while saving output image !" << endl);
590  return 1;
591  }
592  // Normal end
593  return 0;
594 }
595 
596 // =====================================================================
597 // ---------------------------------------------------------------------
598 // ---------------------------------------------------------------------
599 // =====================================================================
600 
601 int iIterativeAlgorithm::InitSpecificOptions(string a_specificOptions)
602 {
603  (void)a_specificOptions; // avoid 'unused parameter' compil. warnings
604  if (m_verbose>=2) Cout("iIterativeAlgorithm::InitSpecificOptions ... " << endl);
605 
606  return 0;
607 }
int ApplyProcessingForward(oImageSpace *ap_ImageSpace)
void GetEventIndexStartAndStop(int64_t *ap_indexStart, int64_t *ap_indexStop, int a_subsetNum=0, int a_NbSubsets=1)
#define MODE_HISTOGRAM
int PerformDeformation(oImageSpace *ap_Image)
void DeallocateBackwardImageFromDynamicBasis()
Free memory for the backward image matrices.
FLTNB GetInitialValue()
Return the initial image value used by the optimizer, explaining why the eponym function of vOptimize...
oImageDimensionsAndQuantification * mp_ID
#define Cerr(MESSAGE)
int PreImageUpdateStep()
A function that simply calls the eponym function from the vOptimizer.
void Reduce()
Merge parallel results into the matrix of the backward image matrix of the first thread. Also for MPI.
void StartIterativeDataUpdateStep4(int a_thread)
int ConvolveForward(oImageSpace *ap_ImageSpace)
int InitSpecificOptions(string a_specificOptions)
int ApplyOutputFOVMasking()
Mask the outside of the transaxial FOV based on the m_fovOutPercent.
int StepInnerLoopInsideSubsetLoop(int a_iteration, int a_subset, int a_bed)
oOptimizerManager * mp_OptimizerManager
void InstantiateBackwardImageFromDynamicBasis(int a_nbBackwardImages)
void StartIterativeDataUpdateStep3(int a_thread)
void StartIterativeDataUpdateStep2(int a_thread)
void StartIterativeDataUpdateStep1(int a_thread)
int ApplyDeformationsToBackwardImage(oImageSpace *ap_Image)
static sOutputManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
int ConvolvePost(oImageSpace *ap_ImageSpace)
int ConvolveIntra(oImageSpace *ap_ImageSpace)
int ApplyProcessingPost(oImageSpace *ap_ImageSpace)
int InitImage(const string &a_pathToInitialImage, FLTNB a_value)
int DynamicSwitch(int64_t a_currentEventIndex, uint32_t a_currentTime, int a_bed, int a_th)
int StepAfterIterationLoop()
This function is called at the end of the RunCPU function.
int DataUpdateStep(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
void StopIterativeDataUpdateStep2(int a_thread)
int InitImageForDeformation(oImageSpace *ap_Image)
oImageProcessingManager * mp_ImageProcessingManager
bool GetRespStaticFlag()
Get the respiratory static flag that says if the reconstruction has only one respiratory gate or not...
int StepPostProcessInsideSubsetLoop(int a_iteration, int a_subset)
vEvent * GetEvent(int64_t a_eventIndex, int a_th=0)
void PrepareForwardImage()
Copy current image matrix in the forward-image buffer matrix.
int SaveOutputBasisCoefficientImage(int a_iteration, int a_subset=-1)
void StopIterativeDataUpdateStep1(int a_thread)
bool GetCardStaticFlag()
Get the cardiac static flag that says if the reconstruction has only one cardiac gate or not...
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...
int StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
int GetNbBackwardImages()
Return the number of backward images used by the optimizer, explaining why the eponym function of vOp...
int SaveOutputImage(int a_iteration, int a_subset=-1)
int PreDataUpdateStep()
A function that simply calls the eponym function from the vOptimizer.
void ResetCurrentDynamicIndices()
Call the eponym function from the oDynamicDataManager class using the member object.
virtual int StepBeforeIterationLoop()
This function is called at the beginning of the RunCPU function.
oProjectionLine * ComputeProjectionLine(vEvent *ap_Event, int a_th)
int ConvolveBackward(oImageSpace *ap_ImageSpace)
oImageConvolverManager * mp_ImageConvolverManager
int SaveParametricImages(int a_iteration, int a_subset=-1)
void SetCurrentIteration(int a_currentIteration)
void StartConvolution()
Start the timer for duration of convolution.
bool GetTimeStaticFlag()
Get the time static flag that says if the reconstruction has only one frame or not.
bool NotEmptyLine()
This function is used to know if the line contains any voxel contribution.
void InstantiateImageForDeformation(oImageSpace *ap_Image)
int ApplyProcessingIntra(oImageSpace *ap_ImageSpace)
oProjectorManager * mp_ProjectorManager
int StepBeforeIterationLoop()
This function is called at the beginning of the RunCPU function.
oDeformationManager * mp_DeformationManager
This class is designed to manage and store system matrix elements associated to a vEvent...
int ApplyOutputFlip()
Just flip the output image.
virtual int StepAfterIterationLoop()
This function is called at the end of the RunCPU function.
int StepImageOutput(int a_iteration, int a_subset=-1)
This is the base class for reconstructions, containing a framework with iteration and data subset loo...
void StopConvolution()
Stop the timer for duration of convolution.
void SetNumbersOfIterationsAndSubsets(int a_nbIterations, int *ap_nbSubsets)
Mother class for the Event objects.
int GetNbThreadsForImageComputation()
Get the number of threads used for image operations.
int ImageUpdateStep()
A function dedicated to the update step in the image space (performed after the loop on events) ...
int GetNbThreadsForProjection()
Get the number of threads used for projections.
This class is designed to manage some profiling of the code.
int SaveSensitivityImage(const string &a_pathToSensitivityImage)
void InitBackwardImage()
Initialize each voxel of the backward images to 0, also for sensitivity if not loaded (estimated on t...
void DeallocateImageForDeformation(oImageSpace *ap_Image)
#define Cout(MESSAGE)
static sChronoManager * GetInstance()
Instantiate the singleton if not already done, then return the pointer to its instance.
int ApplyDynamicModel(oImageSpace *ap_ImageS, int a_iteration, int a_subset)
oDynamicModelManager * mp_DynamicModelManager