CASToR  3.0
Tomographic Reconstruction (PET/SPECT/CT)
iIterativeAlgorithm.cc
Go to the documentation of this file.
1 /*
2 This file is part of CASToR.
3 
4  CASToR is free software: you can redistribute it and/or modify it under the
5  terms of the GNU General Public License as published by the Free Software
6  Foundation, either version 3 of the License, or (at your option) any later
7  version.
8 
9  CASToR is distributed in the hope that it will be useful, but WITHOUT ANY
10  WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
12  details.
13 
14  You should have received a copy of the GNU General Public License along with
15  CASToR (in file GNU_GPL.TXT). If not, see <http://www.gnu.org/licenses/>.
16 
17 Copyright 2017-2019 all CASToR contributors listed below:
18 
19  --> Didier BENOIT, Claude COMTAT, Marina FILIPOVIC, Thibaut MERLIN, Mael MILLARDET, Simon STUTE, Valentin VIELZEUF
20 
21 This is CASToR version 3.0.
22 */
23 
30 #include "gVariables.hh"
31 #include "iIterativeAlgorithm.hh"
32 #include "iEventHistoCT.hh"
33 
34 // =====================================================================
35 // ---------------------------------------------------------------------
36 // ---------------------------------------------------------------------
37 // =====================================================================
38 
40 {
41 }
42 
43 // =====================================================================
44 // ---------------------------------------------------------------------
45 // ---------------------------------------------------------------------
46 // =====================================================================
47 
49 {
50 }
51 
52 // =====================================================================
53 // ---------------------------------------------------------------------
54 // ---------------------------------------------------------------------
55 // =====================================================================
56 
58 {
60  {
61  Cerr("***** iIterativeAlgorithm::StepBeforeIterationLoop() -> A problem occurred while calling StepBeforeIterationLoop() function !" << endl);
62  return 1;
63  }
64  if (m_verbose>=2) Cout("iIterativeAlgorithm::StepBeforeIterationLoop ... " << endl);
65 
68 
69  // Main image initialization
71  {
72  Cerr("***** iIterativeAlgorithm::StepBeforeIterationLoop() -> An error occurred while reading the initialization image !" << endl);
73  return 1;
74  }
75 
76  // Set numbers of iterations and subsets to the optimizer
78 
79  // End
80  return 0;
81 }
82 
83 // =====================================================================
84 // ---------------------------------------------------------------------
85 // ---------------------------------------------------------------------
86 // =====================================================================
87 
89 {
90  if (m_verbose>=3) Cout("iIterativeAlgorithm::StepBeforeSubsetLoop ... " << endl);
91 
92  // Set the current iteration to the optimizer
94 
95  // End
96  return 0;
97 }
98 
99 // =====================================================================
100 // ---------------------------------------------------------------------
101 // ---------------------------------------------------------------------
102 // =====================================================================
103 
104 int iIterativeAlgorithm::StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
105 {
106  if (m_verbose>=3) Cout("iIterativeAlgorithm::StepPreProcessInsideSubsetLoop ... " << endl);
107 
108  // Set the current subset to the optimizer
110 
111  // Get the chrono manager singleton pointer
112  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
113 
114  // Initialize the correction backward image(s)
116 
117  // Copy current image in forward-image buffer
119 
120  // Apply image processing to forward image
122  {
123  Cerr("***** iIterativeAlgorithm::StepPreProcessInsideSubsetLoop() -> A problem occurred while applyin image processing to forward image !" << endl);
124  return 1;
125  }
126 
127  // Apply convolver to forward image
128  p_ChronoManager->StartConvolution();
130  {
131  Cerr("***** iIterativeAlgorithm::StepPreProcessInsideSubsetLoop() -> A problem occurred while applying image convolver to forward image !" << endl);
132  return 1;
133  }
134  p_ChronoManager->StopConvolution();
135 
136  // Call the pre-event loop function from the optimizer manager
138 
139  // Initialisation of the backup images for deformation
140  // (bu_fward image initialized with current forward image)
141  // (bu_bward and bu_sens images set to 0)
142 
144 
145  // End
146  return 0;
147 }
148 
149 // =====================================================================
150 // ---------------------------------------------------------------------
151 // ---------------------------------------------------------------------
152 // =====================================================================
153 
154 int iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop(int a_iteration, int a_subset, int a_bed)
155 {
156  // Verbose
158  {
159  if (m_nbBeds>1) Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Start loop over events for bed " << a_bed+1 << endl << flush);
160  else Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Start loop over events" << endl << flush);
161  }
162 
163  // Reinitialize 4D gate indexes
165 
166  // Apply the bed offset for this bed position
168 
169  // Progression (increments of 2%)
171  {
172  cout << "0 10% 20% 30% 40% 50% 60% 70% 80% 90% 100%" << endl;
173  cout << "|----|----|----|----|----|----|----|----|----|----|" << endl;
174  cout << "|" << flush;
175  }
176  int progression_percentage_old = 0;
177  int progression_nb_bars = 0;
178  uint64_t progression_printing_index = 0;
179 
180  // Compute start and stop indices taking MPI into account (the vDataFile does that)
181  int64_t index_start = 0;
182  int64_t index_stop = 0;
183  m2p_DataFile[a_bed]->GetEventIndexStartAndStop(&index_start, &index_stop, a_subset, mp_nbSubsets[a_iteration]);
184 
185  // Synchronize MPI processes
186  #ifdef CASTOR_MPI
187  MPI_Barrier(MPI_COMM_WORLD);
188  #endif
189 
190  // Set the number of threads for projections (right after this loop, we set back the number of threads for image computations)
191  #ifdef CASTOR_OMP
192  omp_set_num_threads(mp_ID->GetNbThreadsForProjection());
193  #endif
194 
195  // This boolean is used to report any problem inside the parallel loop
196  bool problem = false;
197 
198  // Get the chrono manager singleton pointer
199  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
200 
201  // These flags (one per thread) are used to signal when events are beyond the last frame time stop
202  bool* p_end_loop_flags = (bool*)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(bool));
203  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++) p_end_loop_flags[th] = false;
204 
205  // Launch the loop with precomputed start and stop and using a step equal to the number of subsets
206  int64_t index;
207  // Keep the static scheduling with a chunk size at 1, it is important
208  #pragma omp parallel for private(index) schedule(static, 1)
209  for ( index = index_start ; index < index_stop ; index += mp_nbSubsets[a_iteration] )
210  {
211  // Get the thread index
212  int th = 0;
213  #ifdef CASTOR_OMP
214  th = omp_get_thread_num();
215  #endif
216  // Print progression (do not log out with Cout here)
217  if (m_verbose>=2 && th==0 && mp_ID->GetMPIRank()==0)
218  {
219  if (progression_printing_index%1000==0)
220  {
221  int progression_percentage_new = ((int)( (((float)(index-index_start+1))/((float)(index_stop-index_start)) ) * 100.));
222  if (progression_percentage_new>=progression_percentage_old+2) // Increments of 2%
223  {
224  int nb_steps = (progression_percentage_new-progression_percentage_old)/2;
225  for (int i=0; i<nb_steps; i++)
226  {
227  cout << "-" << flush;
228  progression_nb_bars++;
229  }
230  progression_percentage_old += nb_steps*2;
231  }
232  }
233  progression_printing_index++;
234  }
235  // Skip to the end
236  if (p_end_loop_flags[th]) continue;
237 
238  // Step 1: Get the current event for that thread index
239  p_ChronoManager->StartIterativeDataUpdateStep1(th);
240  #ifdef CASTOR_DEBUG
241  if (m_verbose>=4) Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step1: Get current event for that thread index " << endl);
242  #endif
243  vEvent* event = m2p_DataFile[a_bed]->GetEvent(index, th);
244  if (event==NULL)
245  {
246  Cerr("***** iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> An error occurred while getting the event from index "
247  << index << " (thread " << th << ") !" << endl);
248  // Specify that there was a problem
249  problem = true;
250  // We must continue here because we are inside an OpenMP loop
251  continue;
252  }
253  p_ChronoManager->StopIterativeDataUpdateStep1(th);
254 
255  // Step 2: Call the dynamic switch function that updates the current frame and gate numbers, and also detects involuntary patient motion
256  p_ChronoManager->StartIterativeDataUpdateStep2(th);
257  #ifdef CASTOR_DEBUG
258  if (m_verbose>=4)
259  {
260  Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step2: Check for Dynamic event (frame/gate switch, image-based deformation " << endl);
261  }
262  #endif
263  int dynamic_switch_value = mp_ID->DynamicSwitch(index, event->GetTimeInMs(), a_bed, th);
264  // If the DYNAMIC_SWITCH_CONTINUE is returned, then it means that we are not yet at the first frame
265  if ( dynamic_switch_value == DYNAMIC_SWITCH_CONTINUE )
266  {
267  // Then we just skip this event
268  continue;
269  }
270  // Else, if the DYNAMIC_SWITCH_DEFORMATION is returned, then it means that a change of gate or involuntary motion has occurred
271  // and is associated to a deformation
272  else if ( dynamic_switch_value == DYNAMIC_SWITCH_DEFORMATION )
273  {
274  #ifdef CASTOR_OMP
275  // set deformation requirement for this thread
276  #pragma omp critical (deformation_requirements_first)
277  {
278  mp_DeformationManager->SetDeformationRequirement(th);
279  }
280  // block this thread until the deformation is done
281  bool wait_for_deformation = true;
282  while (wait_for_deformation)
283  {
284  // the first thread performs the deformation on the forward image
285  // as soon as all the threads are ready for deformation
286  if (th==0)
287  {
288  bool all_threads_ready = false;
289  //check whether all the threads require deformation
290  #pragma omp critical (deformation_requirements_first)
291  {
292  all_threads_ready = mp_DeformationManager->AllThreadsRequireDeformation();
293  }
294  if (all_threads_ready)
295  {
296  // Perform the deformation
298  // Free all the threads so that they can continue
299  #pragma omp critical (deformation_requirements_second)
300  {
301  mp_DeformationManager->UnsetDeformationRequirements();
302  }
303  }
304  }
305  // check whether the deformation is done and the thread can continue
306  #pragma omp critical (deformation_requirements_second)
307  {
308  wait_for_deformation = mp_DeformationManager->GetDeformationRequirement(th);
309  }
310  }
311  #else
312  // Perform here any resp related deformations on the forward image
314  #endif
315 
316  // Restore progression printing
317  if (th==0 && m_verbose>=3 && mp_ID->GetMPIRank()==0)
318  {
319  cout << "0 10% 20% 30% 40% 50% 60% 70% 80% 90% 100%" << endl;
320  cout << "|----|----|----|----|----|----|----|----|----|----|" << endl;
321  int progression_percentage_new = ((int)( (((float)(index-index_start+1))/((float)(index_stop-index_start)) ) * 100.));
322  cout << "|";
323 
324  for (int i=0; i< progression_percentage_new / 2 ; i++) cout << "-";
325  }
326 
327  }
328  // 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)
329  else if ( dynamic_switch_value == DYNAMIC_SWITCH_END_LOOP )
330  {
331  p_end_loop_flags[th] = true;
332  continue;
333  }
334  p_ChronoManager->StopIterativeDataUpdateStep2(th);
335 
336  // Step 3: Compute the projection line
337  p_ChronoManager->StartIterativeDataUpdateStep3(th);
338  #ifdef CASTOR_DEBUG
339  if (m_verbose>=4) Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step3: Compute the projection line " << endl);
340  #endif
342  if (line==NULL)
343  {
344  Cerr("***** iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> A problem occurred while computing the projection line !" << endl);
345  // Specify that there was a problem
346  problem = true;
347  // We must continue here because we are inside an OpenMP loop
348  continue;
349  }
350  p_ChronoManager->StopIterativeDataUpdateStep3(th);
351 
352  // Step 4: Optimize in the data space (forward-proj, update, backward-proj)
353  p_ChronoManager->StartIterativeDataUpdateStep4(th);
354  #ifdef CASTOR_DEBUG
355  if (m_verbose>=4) Cout("iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step4: Optimize in the data space " << endl);
356  #endif
357  if (line->NotEmptyLine()) mp_OptimizerManager->DataUpdateStep( line,
358  event,
359  a_bed,
363  th);
364  p_ChronoManager->StopIterativeDataUpdateStep4(th);
365  } // End of indices loop (OpenMP stops here)
366 
367  // Synchronize MPI processes
368  #ifdef CASTOR_MPI
369  MPI_Barrier(MPI_COMM_WORLD);
370  #endif
371 
372  // Set back the number of threads for image computation
373  #ifdef CASTOR_OMP
374  omp_set_num_threads(mp_ID->GetNbThreadsForImageComputation());
375  #endif
376 
377  // End of progression printing (do not log out with Cout here)
378  if (m_verbose>=2 && mp_ID->GetMPIRank()==0)
379  {
380  int progression_total_bars = 49;
381  for (int i=0; i<progression_total_bars-progression_nb_bars; i++) cout << "-";
382  cout << "|" << endl;
383  }
384 
385  // If a problem was encountered, then report it here
386  if (problem)
387  {
388  Cerr("***** iIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> A problem occurred inside the parallel loop over events !" << endl);
389  return 1;
390  }
391 
392  // End
393  return 0;
394 }
395 
396 // =====================================================================
397 // ---------------------------------------------------------------------
398 // ---------------------------------------------------------------------
399 // =====================================================================
400 
402 {
403  if (m_verbose>=3) Cout("iIterativeAlgorithm::StepPostProcessInsideSubsetLoop ... " << endl);
404 
405  // Get the chrono manager singleton pointer
406  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
407 
408  // Merge parallel results
410 
411  // If mask provided, apply mask to sensitivity, the CleanNeverVisitedVoxels will take care of the images
413 
414  // Perform here any image-based deformations related to motion on the backward image.
416 
417  // Apply convolver to backward images and sensitivity-if-needed
418  p_ChronoManager->StartConvolution();
420  {
421  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while applying convolver to backward images !" << endl);
422  return 1;
423  }
424  p_ChronoManager->StopConvolution();
425 
426  // Call the pre image update step function from the optimizer manager
428  {
429  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while performing the pre-image-update step !" << endl);
430  return 1;
431  }
432 
433  // Optimize in the image space (apply corrections, MAP and sensitivity); pass the number of subsets for list-mode sensitivity scaling
435  {
436  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while performing the image update step !" << endl);
437  return 1;
438  }
439 
440  // Apply convolver to current estimed images
441  p_ChronoManager->StartConvolution();
443  {
444  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while applying convolver to current estimate images !" << endl);
445  return 1;
446  }
447  p_ChronoManager->StopConvolution();
448 
449  // Post-update Dynamic Modeling step (if enabled). Manage either linear/non linear dynamic model (physiological or not)
450  if (mp_DynamicModelManager->ApplyDynamicModel(mp_ImageSpace, a_iteration, a_subset))
451  {
452  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while applying dynamic model to current estimate images !" << endl);
453  return 1;
454  }
455 
456  // Apply image processing to current estime images
458  {
459  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while applying image processing to current estimate images !" << endl);
460  return 1;
461  }
462 
463  // Save the sensitivity image in histogram mode, if asked for
465  {
466  // Get output manager to build the file name
467  sOutputManager* p_output_manager = sOutputManager::GetInstance();
468  // Build the file name
469  string temp_sens = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
470  stringstream temp_it; temp_it << a_iteration + 1;
471  stringstream temp_ss; temp_ss << a_subset + 1;
472  temp_sens.append("_it").append(temp_it.str()).append("_ss").append(temp_ss.str()).append("_sensitivity");
473  // Save sensitivity
475  }
476 
477  // Save the current image estimate if asked for
478  if (mp_ID->GetMPIRank()==0 && mp_outputIterations[a_iteration] && m_saveImageAfterSubsets)
479  {
480  // Verbose
481  if (m_verbose>=1) Cout("iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> Save image at iteration " << a_iteration+1 << " for subset " << a_subset+1 << endl);
482  // Save image
483  if (StepImageOutput(a_iteration,a_subset))
484  {
485  Cerr("***** iIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occurred while saving images at iteration " << a_iteration+1 << " !" << endl);
486  return 1;
487  }
488  }
489 
490  // End
491  return 0;
492 }
493 
494 // =====================================================================
495 // ---------------------------------------------------------------------
496 // ---------------------------------------------------------------------
497 // =====================================================================
498 
500 {
501  if (m_verbose>=3) Cout("iIterativeAlgorithm::StepAfterSubsetLoop() -> Clean never visited voxels and save images if needed" << endl);
502  // Clean never visited voxels
504  // Save the main image
505  if (mp_ID->GetMPIRank()==0 && mp_outputIterations[a_iteration])
506  {
507  // Verbose
508  if (m_verbose>=1) Cout("iIterativeAlgorithm::StepAfterSubsetLoop() -> Save image at iteration " << a_iteration+1 << endl);
509  // Save image
510  if (StepImageOutput(a_iteration))
511  {
512  Cerr("***** iIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occurred while saving images at iteration " << a_iteration+1 << " !" << endl);
513  return 1;
514  }
515  }
516  // End
517  return 0;
518 }
519 
520 // =====================================================================
521 // ---------------------------------------------------------------------
522 // ---------------------------------------------------------------------
523 // =====================================================================
524 
526 {
528  {
529  Cerr("***** iIterativeAlgorithm::StepAfterIterationLoop() -> A problem occurred while calling StepAfterIterationLoop() function !" << endl);
530  return 1;
531  }
532  if (m_verbose>=2) Cout("iIterativeAlgorithm::StepAfterIterationLoop ... " << endl);
533 
534  // Deallocate everything
537 
538  // Normal end
539  return 0;
540 }
541 
542 // =====================================================================
543 // ---------------------------------------------------------------------
544 // ---------------------------------------------------------------------
545 // =====================================================================
546 
547 int iIterativeAlgorithm::StepImageOutput(int a_iteration, int a_subset)
548 {
549  // Get the chrono manager singleton pointer
550  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
551  // =================================================================================================
552  // Save dynamic coefficient images from the vDynamicModel
553  // =================================================================================================
554  if (mp_DynamicModelManager->SaveParametricImages(a_iteration, a_subset))
555  {
556  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while saving parametric images related to the dynamic model !" << endl);
557  return 1;
558  }
559  // =================================================================================================
560  // Apply pre-processing steps
561  // =================================================================================================
562  // Copy the current image into the forward image
564  // Apply post-convolution if needed
565  p_ChronoManager->StartConvolution();
567  {
568  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while convolving the output image !" << endl);
569  return 1;
570  }
571  p_ChronoManager->StopConvolution();
572  // Apply post-processing if needed
574  {
575  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while applying image processing the output image !" << endl);
576  return 1;
577  }
578  // Apply output transaxial FOV masking
580  {
581  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while applying output FOV masking !" << endl);
582  return 1;
583  }
584  // Apply output flip
586  {
587  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while applying output flip !" << endl);
588  return 1;
589  }
590  // =================================================================================================
591  // If some basis functions are really in use
592  // =================================================================================================
594  {
595  // Save the basis coefficients images if asked for
597  {
598  if (mp_ImageSpace->SaveOutputBasisCoefficientImage(a_iteration, a_subset))
599  {
600  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while saving output image !" << endl);
601  return 1;
602  }
603  }
604  // Compute the output image from the forward image basis functions
606  }
607  // =================================================================================================
608  // Save frames/gates
609  // =================================================================================================
610 
611  // Save output image (note that if no basis functions at all are in use, the the output image already points to the forward image)
612  if (mp_ImageSpace->SaveOutputImage(a_iteration, a_subset))
613  {
614  Cerr("***** iIterativeAlgorithm::StepImageOutput() -> A problem occurred while saving output image !" << endl);
615  return 1;
616  }
617  // Normal end
618  return 0;
619 }
620 
621 // =====================================================================
622 // ---------------------------------------------------------------------
623 // ---------------------------------------------------------------------
624 // =====================================================================
625 
626 int iIterativeAlgorithm::InitSpecificOptions(string a_specificOptions)
627 {
628  (void)a_specificOptions; // avoid 'unused parameter' compil. warnings
629  if (m_verbose>=2) Cout("iIterativeAlgorithm::InitSpecificOptions ... " << endl);
630 
631  return 0;
632 }
oImageConvolverManager * mp_ImageConvolverManager
Definition: vAlgorithm.hh:359
void StopIterativeDataUpdateStep1(int a_thread)
Stop the timer for duration of iterative data update step 1.
int ApplyProcessingForward(oImageSpace *ap_ImageSpace)
A function used to apply image processing modules onto the forward image of the oImageSpace.
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
oOptimizerManager * mp_OptimizerManager
Definition: vAlgorithm.hh:355
void GetEventIndexStartAndStop(int64_t *ap_indexStart, int64_t *ap_indexStop, int a_subsetNum=0, int a_NbSubsets=1)
Compute the index start and stop of the events loop with respect to the current subset and MPI size a...
Definition: vDataFile.cc:678
#define MODE_HISTOGRAM
Definition: vDataFile.hh:58
int PerformDeformation(oImageSpace *ap_Image)
Apply deformations during reconstruction.
oImageProcessingManager * mp_ImageProcessingManager
Definition: vAlgorithm.hh:360
void DeallocateBackwardImageFromDynamicBasis()
Free memory for the backward image matrices.
Definition: oImageSpace.cc:264
FLTNB GetInitialValue()
Return the initial image value used by the optimizer, explaining why the eponym function of vOptimize...
#define DYNAMIC_SWITCH_CONTINUE
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.
Declaration of class oIterativeAlgorithm.
int ConvolveForward(oImageSpace *ap_ImageSpace)
A function used to apply convolvers onto the forward image of the oImageSpace.
int InitSpecificOptions(string a_specificOptions)
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.
oProjectorManager * mp_ProjectorManager
Definition: vAlgorithm.hh:354
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...
bool m_saveSensitivityHistoFlag
Definition: vAlgorithm.hh:367
void SetCurrentSubset(int a_currentSubset)
Set the current subset to the vOptimizer.
void StopIterativeDataUpdateStep3(int a_thread)
Stop the timer for duration of iterative data update step 3.
#define DYNAMIC_SWITCH_END_LOOP
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:219
int m_nbIterations
Definition: vAlgorithm.hh:347
void SetNumbersOfIterationsAndSubsets(int a_nbIterations, int *ap_nbSubsets)
Set these numbers of iterations and subsets to the vOptimizer.
int ApplyMaskToSensitivity()
Apply the mask to the sensitivity image (only for the first thread, the image must be reduced beforeh...
int ApplyDeformationsToBackwardImage(oImageSpace *ap_Image)
Apply final backward deformations on the backward 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)
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)
A function used to apply convolvers onto the output image of the oImageSpace.
int * mp_nbSubsets
Definition: vAlgorithm.hh:348
int InitImage(const string &a_pathToInitialImage, FLTNB a_value)
Initialize the main image, either using:
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)
A function dedicated to the update step in the data space (for each event inside the loop) ...
int InitImageForDeformation(oImageSpace *ap_Image)
If deformation is enabled, ask the Image Space to initialize the temporary backward image for deforma...
void StartIterativeDataUpdateStep1(int a_thread)
Start the timer for duration of iterative data update step 1.
void StopIterativeDataUpdateStep4(int a_thread)
Stop the timer for duration of iterative data update step 4.
bool GetRespStaticFlag()
Get the respiratory static flag that says if the reconstruction has only one respiratory gate or not...
oImageSpace * mp_ImageSpace
Definition: vAlgorithm.hh:358
void StartIterativeDataUpdateStep2(int a_thread)
Start the timer for duration of iterative data update step 2.
int StepPostProcessInsideSubsetLoop(int a_iteration, int a_subset)
vEvent * GetEvent(int64_t a_eventIndex, int a_th=0)
Definition: vDataFile.cc:605
vDataFile ** m2p_DataFile
Definition: vAlgorithm.hh:353
void PrepareForwardImage()
Copy current image matrix in the forward-image buffer matrix.
int SaveOutputBasisCoefficientImage(int a_iteration, int a_subset=-1)
Call the interfile function to write output basis function coefficient image on disk.
#define Cerr(MESSAGE)
bool m_saveDynamicBasisCoefficients
Definition: vAlgorithm.hh:369
const string & GetPathName()
bool GetCardStaticFlag()
Get the cardiac static flag that says if the reconstruction has only one cardiac gate or not...
int StepBeforeSubsetLoop(int a_iteration)
This function is called before starting the subset loop.
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)
This function is called right after starting the subset loop. Apply any kind of processing on the for...
string m_pathToInitialImg
Definition: vAlgorithm.hh:362
int GetNbBackwardImages()
Return the number of backward images used by the optimizer, explaining why the eponym function of vOp...
oDynamicModelManager * mp_DynamicModelManager
Definition: vAlgorithm.hh:357
int SaveOutputImage(int a_iteration, int a_subset=-1)
Call the interfile function to write output image on disk.
Declaration of class iEventHistoCT.
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.
Definition: vAlgorithm.cc:388
bool m_saveImageAfterSubsets
Definition: vAlgorithm.hh:368
void ComputeOutputImage()
void StartIterativeDataUpdateStep3(int a_thread)
Start the timer for duration of iterative data update step 3.
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.
void StopIterativeDataUpdateStep2(int a_thread)
Stop the timer for duration of iterative data update step 2.
#define VERBOSE_NORMAL
int SaveParametricImages(int a_iteration, int a_subset=-1)
Call SaveParametricImages() function of the dynamic model object is &#39;m_UseModel&#39; is on...
void ApplyBedOffset(int a_bed)
Compute the bed offset from the provided bed index and apply it to all projection lines...
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)
If deformation is enabled, ask the Image Space to Instantiate the temporary backward image for deform...
int ApplyProcessingIntra(oImageSpace *ap_ImageSpace)
A function used to apply image processing modules onto the current image of the oImageSpace.
oDeformationManager * mp_DeformationManager
Definition: vAlgorithm.hh:356
int StepBeforeIterationLoop()
This function is called at the beginning of the RunCPU function.
void StartIterativeDataUpdateStep4(int a_thread)
Start the timer for duration of iterative data update step 4.
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.
oImageDimensionsAndQuantification * mp_ID
Definition: vAlgorithm.hh:352
virtual int StepAfterIterationLoop()
This function is called at the end of the RunCPU function.
Definition: vAlgorithm.cc:483
int StepImageOutput(int a_iteration, int a_subset=-1)
This function deals with everything about saving output images from the reconstruction.
This is the base class for reconstructions, containing a framework with iteration and data subset loo...
Definition: vAlgorithm.hh:51
void StopConvolution()
Stop the timer for duration of convolution.
Mother class for the Event objects.
Definition: vEvent.hh:42
int GetNbThreadsForImageComputation()
Get the number of threads used for image operations.
int GetDataMode()
Definition: vDataFile.hh:299
int 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.
int StepAfterSubsetLoop(int a_iteration)
This function is called after finishing the subset loop.
bool * mp_outputIterations
Definition: vAlgorithm.hh:349
#define Cout(MESSAGE)
void SetCurrentIteration(int a_currentIteration)
Set the current iteration to the vOptimizer.
This class is designed to manage some profiling of the code.
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...
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
static sChronoManager * GetInstance()
Instantiate the singleton if not already done, then return the pointer to its instance.
int GetMPIRank()
Get the MPI instance number (the rank)
int ApplyDynamicModel(oImageSpace *ap_ImageS, int a_iteration, int a_subset)
int GetCurrentTimeFrame(int a_th)
call the eponym function from the oDynamicDataManager object
#define DYNAMIC_SWITCH_DEFORMATION