CASToR  2.0
Tomographic Reconstruction (PET/SPECT/CT)
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
oIterativeAlgorithm.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-2018 all CASToR contributors listed below:
18 
19  --> current contributors: Thibaut MERLIN, Simon STUTE, Didier BENOIT, Claude COMTAT, Marina FILIPOVIC, Mael MILLARDET
20  --> past contributors: Valentin VIELZEUF
21 
22 This is CASToR version 2.0.
23 */
24 
31 #include "gVariables.hh"
32 #include "oIterativeAlgorithm.hh"
33 #include "iEventHistoCT.hh"
34 
35 bool m_releaseThreads = false;
38 // =====================================================================
39 // ---------------------------------------------------------------------
40 // ---------------------------------------------------------------------
41 // =====================================================================
42 /*
43  \brief oIterativeAlgorithm constructor.
44  Initialize the member variables to their default values.
45 */
47 {
48  // Default number of iterations and subsets (1 for each)
49  m_nbIterations = 1;
50  mp_nbSubsets = (int*)malloc(1*sizeof(int));
51  mp_nbSubsets[0] = 1;
53  mp_outputIterations = NULL;
54  // set all members to default values
55  m_verbose = -1;
56  m_flagGPU = false;
57  mp_ID = NULL;
58  m2p_DataFile= NULL;
59  mp_ProjectorManager = NULL;
60  mp_OptimizerManager = NULL;
61  mp_DeformationManager = NULL;
63  mp_ImageSpace = NULL;
66  m_nbBeds = -1;
67  m_pathToInitialImg = "";
69  m_pathToMaskImg = "";
72 }
73 
74 // =====================================================================
75 // ---------------------------------------------------------------------
76 // ---------------------------------------------------------------------
77 // =====================================================================
78 /*
79  \brief oIterativeAlgorithm destructor.
80 */
82 {
83  // Delete some tables
84  if (mp_nbSubsets!=NULL) free(mp_nbSubsets);
86 }
87 
88 // =====================================================================
89 // ---------------------------------------------------------------------
90 // ---------------------------------------------------------------------
91 // =====================================================================
92 /*
93  \fn SetNbIterationsAndSubsets
94  \param a_nbIterationsSubsets
95  \brief Set the number of iterations and subsets
96  \details The provided string is a list of iteration:subset pairs
97  separated by commas.
98  \return 0 if success, positive value otherwise
99 */
100 int oIterativeAlgorithm::SetNbIterationsAndSubsets(const string& a_nbIterationsSubsets)
101 {
102  // If the string is empty, then just keep the defaults
103  if (a_nbIterationsSubsets=="") return 0;
104  // Otherwise, reset the number of iterations to 0
105  m_nbIterations = 0;
106 
107  // Copy the string
108  string buf = a_nbIterationsSubsets;
109 
110  // Loop to search for commas
111  size_t pos_comma;
112  while ((pos_comma=buf.find_first_of(","))!=string::npos)
113  {
114  // Get the substring before the comma
115  string sub_buf = buf.substr(0,pos_comma);
116  // Search for columns
117  size_t pos_column = sub_buf.find_first_of(":");
118  if (pos_column==string::npos || pos_column==0)
119  {
120  Cerr("***** oIterativeAlgorithm::SetNbIterationsAndSubsets() -> Syntax problem in number of iterations and subsets !" << endl);
121  return 1;
122  }
123  int iter = atoi( sub_buf.substr(0,pos_column).c_str() );
124  int subs = atoi( sub_buf.substr(pos_column+1).c_str() );
125  mp_nbSubsets = (int*)realloc(mp_nbSubsets,(m_nbIterations+iter)*sizeof(int));
126  for (int it=0; it<iter; it++) mp_nbSubsets[m_nbIterations+it] = subs;
127  m_nbIterations += iter;
128  buf = buf.substr(pos_comma+1);
129  }
130 
131  // Last couple of iterations:subsets
132  size_t pos_column = buf.find_first_of(":");
133  if (pos_column==string::npos || pos_column==0)
134  {
135  Cerr("***** oIterativeAlgorithm::SetNbIterationsAndSubsets() -> Syntax problem in number of iterations and subsets !" << endl);
136  return 1;
137  }
138  int iter = atoi( buf.substr(0,pos_column).c_str() );
139  int subs = atoi( buf.substr(pos_column+1).c_str() );
140  mp_nbSubsets = (int*)realloc(mp_nbSubsets,(m_nbIterations+iter)*sizeof(int));
141  for (int it=0; it<iter; it++) mp_nbSubsets[m_nbIterations+it] = subs;
142  m_nbIterations += iter;
143 
144  if (m_verbose>=3)
145  {
146  Cout("oIterativeAlgorithm::SetNbIterationsAndSubsets() -> Selected numbers of subsets for each iteration:" << endl);
147  Cout(" Iteration / number of subsets "<< endl);
148  for (int it=0 ; it<m_nbIterations ; it++) Cout(" " << it+1 << " / " << mp_nbSubsets[it] << endl);
149  }
150  // End
151  return 0;
152 }
153 
154 // =====================================================================
155 // ---------------------------------------------------------------------
156 // ---------------------------------------------------------------------
157 // =====================================================================
158 /*
159  \fn SetOutputIterations
160  \param a_outputIterations
161  \brief Set the selected output iterations
162  \details The provided string is a list of couple a:b separated by commas.
163  It means that we save one iteration over a until b is reached.
164  "b" must be incrementing for each successive couples.
165  If the list is empty, we save all iterations by default.
166  If the list is equal to "-1", then we save only the last iteration.
167  \return 0 if success, positive value otherwise
168 */
169 int oIterativeAlgorithm::SetOutputIterations(const string& a_outputIterations)
170 {
171  if(m_verbose>=2) Cout("oIterativeAlgorithm::SetOutputIterations ..."<< endl);
172 
173  // Allocate the output iterations boolean table
174  mp_outputIterations = (bool*)malloc(m_nbIterations*sizeof(bool));
175  for (int it=0; it<m_nbIterations; it++) mp_outputIterations[it] = false;
176 
177  // If the list is empty, we save all iterations by default
178  if (a_outputIterations=="")
179  {
180  for (int it=0; it<m_nbIterations; it++) mp_outputIterations[it] = true;
181  return 0;
182  }
183 
184  // Copy the string
185  string buf = a_outputIterations;
186 
187  // Loop to search for commas
188  size_t pos_comma;
189  int current_iteration = -1;
190  while ((pos_comma=buf.find_first_of(","))!=string::npos)
191  {
192  // Get the substring before the comma
193  string sub_buf = buf.substr(0,pos_comma);
194  // Search for columns
195  size_t pos_column = sub_buf.find_first_of(":");
196  if (pos_column==string::npos || pos_column==0)
197  {
198  Cerr("***** oIterativeAlgorithm::SetOutputIterations() -> Syntax problem in output iteration numbers !" << endl);
199  return 1;
200  }
201  int step_iteration = atoi( sub_buf.substr(0,pos_column).c_str() );
202  if (step_iteration<=0)
203  {
204  Cerr("***** oIterativeAlgorithm::SetOutputIterations() -> Iteration step must be strictly positive (here it is " << step_iteration << ") !" << endl);
205  return 1;
206  }
207  int stop_iteration = atoi( sub_buf.substr(pos_column+1).c_str() );
208  if (stop_iteration<=current_iteration)
209  {
210  Cerr("***** oIterationAlgorithm::SetOutputIterations() -> Output iteration number must be increasing through provided couples !" << endl);
211  return 1;
212  }
213  for (int it=current_iteration+step_iteration; it<stop_iteration; it+=step_iteration) mp_outputIterations[it] = true;
214  current_iteration = stop_iteration-1;
215  buf = buf.substr(pos_comma+1);
216  }
217 
218  // Last couple of iterations:subsets
219  size_t pos_column = buf.find_first_of(":");
220  if (pos_column==string::npos || pos_column==0)
221  {
222  // Special case if -1 is provided, it means we save the last iteration
223  if (buf=="-1")
224  {
225  mp_outputIterations[m_nbIterations-1] = true;
226  // We directly exist here
227  return 0;
228  }
229  else
230  {
231  Cerr("***** oIterativeAlgorithm::SetOutputIterations() -> Syntax problem in output iteration numbers !" << endl);
232  return 1;
233  }
234  }
235  int step_iteration = atoi( buf.substr(0,pos_column).c_str() );
236  if (step_iteration<=0)
237  {
238  Cerr("***** oIterativeAlgorithm::SetOutputIterations() -> Iteration step must be strictly positive (here it is " << step_iteration << ") !" << endl);
239  return 1;
240  }
241  int stop_iteration = atoi( buf.substr(pos_column+1).c_str() );
242  if (stop_iteration<=current_iteration)
243  {
244  Cerr("***** oIterationAlgorithm::SetOutputIterations() -> Output iteration number must be increasing through provided couples !" << endl);
245  return 1;
246  }
247  for (int it=current_iteration+step_iteration; it<stop_iteration; it+=step_iteration) mp_outputIterations[it] = true;
248 
249  // End
250  return 0;
251 }
252 
253 // =====================================================================
254 // ---------------------------------------------------------------------
255 // ---------------------------------------------------------------------
256 // =====================================================================
257 /*
258  \fn Iterate
259  \brief Just call either the IterateCPU or the IterateGPU function as asked for
260  \return 0 if success, positive value otherwise
261 */
263 {
264  if(m_verbose>=2) Cout("oIterativeAlgorithm::Iterate ..."<< endl);
265 
266  #ifdef CASTOR_GPU
267  if (m_flagGPU)
268  return IterateGPU();
269  else
270  return IterateCPU();
271  #else
272  return IterateCPU();
273  #endif
274 }
275 
276 // =====================================================================
277 // ---------------------------------------------------------------------
278 // ---------------------------------------------------------------------
279 // =====================================================================
280 /*
281  \fn IterateCPU
282  \brief Perform the iterative loop of the algorithm, call the different object
283  for optimization, projection, update, etc.
284  Function designed to be executed on the CPU only.
285  \details Loops over the iterations, subsets, bed position
286  Call functions related to each steps of iterative reconstruction:
287 
288  StepBeforeIterationLoop()
289  / Loop on iterations
290  | StepBeforeSubsetLoop(iteration)
291  | / Loop on subsets
292  | | StepPreProcessInsideSubsetLoop(iteration,subset)
293  | | / Loop on bed positions
294  | | | StepInnerLoopInsideSubsetLoop(iteration,subset,bed)
295  | | StepPostProcessInsideSubsetLoop(iteration,subset)
296  | StepAfterSubsetLoop(iteration)
297  StepAfterIterationLoop()
298 
299  \return 0 if success, positive value otherwise
300 */
302 {
303  // Synchronize MPI processes
304  #ifdef CASTOR_MPI
305  MPI_Barrier(MPI_COMM_WORLD);
306  #endif
307 
308  // Verbose
309  if (m_verbose>=1) Cout("oIterativeAlgorithm::IterateCPU() -> Start algorithm for " << m_nbIterations << " iterations" << endl);
310 
311  // Initial clock for execution time
312  clock_t clock_start_whole = clock();
313  time_t time_start_whole = time(NULL);
314 
315  // Call before iteration function
317  {
318  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepBeforeIterationLoop() function !" << endl);
319  return 1;
320  }
321 
322  // Loop on iterations
323  for (int iteration=0 ; iteration<m_nbIterations ; iteration++)
324  {
325 
326  // Call before subset function
327  if (StepBeforeSubsetLoop(iteration))
328  {
329  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepBeforeSubsetLoop() function !" << endl);
330  return 1;
331  }
332 
333  // Loop on subsets
334  for (int subset=0 ; subset<mp_nbSubsets[iteration] ; subset++)
335  {
336  // Verbose
337  if (m_verbose>=1)
338  {
339  Cout("oIterativeAlgorithm::IterateCPU() -> Start iteration " << iteration+1 << "/" << m_nbIterations
340  << " subset " << subset+1 << "/" << mp_nbSubsets[iteration] << endl);
341  }
342 
343  // Clock start for subset execution time
344  clock_t clock_start_subset = clock();
345  time_t time_start_subset = time(NULL);
346 
347  // Call pre-process inside subset loop function
348  if (StepPreProcessInsideSubsetLoop(iteration,subset))
349  {
350  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepPreProcessInsideSubsetLoop() function !" << endl);
351  return 1;
352  }
353 
354  // Synchronize MPI processes
355  #ifdef CASTOR_MPI
356  MPI_Barrier(MPI_COMM_WORLD);
357  #endif
358 
359  // Loop on bed positions
360  for (int bed=0 ; bed<m_nbBeds ; bed++)
361  {
362  // Call the inner loop on events inside subset loop function
363  if (StepInnerLoopInsideSubsetLoop(iteration,subset,bed))
364  {
365  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepInnerLoopInsideSubsetLoop() function !" << endl);
366  return 1;
367  }
368  } // End of beds loop
369 
370  // Synchronize MPI processes
371  #ifdef CASTOR_MPI
372  MPI_Barrier(MPI_COMM_WORLD);
373  #endif
374 
375  // Call post-process inside subset loop function
376  if (StepPostProcessInsideSubsetLoop(iteration,subset))
377  {
378  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepPostProcessInsideSubsetLoop() function !" << endl);
379  return 1;
380  }
381 
382  // Clock stop for subset execution time
383  clock_t clock_stop_subset = clock();
384  time_t time_stop_subset = time(NULL);
385  if (m_verbose>=2) Cout ("oIterativeAlgorithm::IterateCPU() -> Time spent for subset " << subset+1 << " | User: " << time_stop_subset-time_start_subset
386  << " sec | CPU: " << (clock_stop_subset-clock_start_subset)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
387 
388 
389  } // End of subsets loop
390 
391 
392  // Call after subset function
393  if (StepAfterSubsetLoop(iteration))
394  {
395  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepAfterSubsetLoop() function !" << endl);
396  return 1;
397  }
398 
399 
400  } // End of iterations loop
401 
402  // Call after iteration function
404  {
405  Cerr("***** oIterativeAlgorithm::IterateCPU() -> A problem occured while calling StepAfterIterationLoop() function !" << endl);
406  return 1;
407  }
408 
409  // Final clock for execution time
410  clock_t clock_stop_whole = clock();
411  time_t time_stop_whole = time(NULL);
412  if (m_verbose>=1) Cout("oIterativeAlgorithm::IterateCPU() -> Total time spent | User: " << time_stop_whole-time_start_whole
413  << " sec | CPU: " << (clock_stop_whole-clock_start_whole)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
414 
415 
416  return 0;
417 }
418 
419 // =====================================================================
420 // ---------------------------------------------------------------------
421 // ---------------------------------------------------------------------
422 // =====================================================================
423 /*
424  \fn StepBeforeIterationLoop
425  \brief This function is called at the beginning of the IterateCPU function.
426  \details Initialization and memory allocation for the imageSpace and
427  some managers
428  \return 0 if success, positive value otherwise.
429 */
431 {
432  if (m_verbose>=2) Cout("oIterativeAlgorithm::StepBeforeIterationLoop ... " << endl);
433 
434  // Check that the image space is already allocated
435  if (mp_ImageSpace==NULL || !mp_ImageSpace->Checked())
436  {
437  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop() -> Image space has not been created or was not checked !" << endl);
438  return 1;
439  }
440 
441  // Instantiate all the required images from the oImageSpace
449 
450  // Main image and sensitivity image initialization
452  {
453  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop() -> An error occured while reading the initialization image !" << endl);
454  return 1;
455  }
456 
458  {
459  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop()-> Error during attenuation image initialization !" << endl);
460  return 1;
461  }
462 
464  {
465  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop() -> An error occured while initializing the sensitivity image !" << endl);
466  return 1;
467  }
468 
470  {
471  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop()-> Error during multimodal image initialization !" << endl);
472  return 1;
473  }
474 
476  {
477  Cerr("***** oIterativeAlgorithm::StepBeforeIterationLoop()-> Error during mask image initialization !" << endl);
478  return 1;
479  }
480 
481  // End
482  return 0;
483 }
484 
485 // =====================================================================
486 // ---------------------------------------------------------------------
487 // ---------------------------------------------------------------------
488 // =====================================================================
489 /*
490  \fn StepBeforeSubsetLoop
491  \param a_iteration : iteration index
492  \brief This function is called before starting the subset loop.
493  \return 0 if success, positive value otherwise.
494 */
496 {
497  if (m_verbose>=3) Cout("oIterativeAlgorithm::StepBeforeSubsetLoop ... " << endl);
498 
499  // End
500  return 0;
501 }
502 
503 // =====================================================================
504 // ---------------------------------------------------------------------
505 // ---------------------------------------------------------------------
506 // =====================================================================
507 /*
508  \fn StepPreProcessInsideSubsetLoop
509  \param a_iteration : iteration index
510  \param a_subset : subset index
511  \brief This function is called right after starting the subset loop.
512  Apply any kind of processing on the forward image before projections
513  \details Copy current main image into forward image matrix
514  Reinitialize backward image and 4D gating indices
515  Apply image processing/convolution on the forward image matrix
516  (image to be projeted)
517  \return 0 if success, positive value otherwise.
518 */
519 int oIterativeAlgorithm::StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
520 {
521  if (m_verbose>=3) Cout("oIterativeAlgorithm::StepPreProcessInsideSubsetLoop ... " << endl);
522 
523  // Get the chrono manager singleton pointer
524  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
525 
526  // Initialize the correction backward image(s)
528 
529  // Copy current image in forward-image buffer (apply deformation if needed)
531 
532  // Apply image processing to forward image
534  {
535  Cerr("***** oIterativeAlgorithm::StepPreProcessInsideSubsetLoop() -> A problem occured while applyin image processing to forward image !" << endl);
536  return 1;
537  }
538 
539  // Apply convolver to forward image
540  p_ChronoManager->StartConvolution();
542  {
543  Cerr("***** oIterativeAlgorithm::StepPreProcessInsideSubsetLoop() -> A problem occured while applying image convolver to forward image !" << endl);
544  return 1;
545  }
546  p_ChronoManager->StopConvolution();
547 
548  // Call the pre-event loop function from the optimizer manager
549  mp_OptimizerManager->PreDataUpdateStep(a_iteration, m_nbIterations, a_subset, mp_nbSubsets);
550 
551  // Initialisation of the backup images for deformation
552  // (bu_fward image initialized with current forward image)
553  // (bu_bward and bu_sens images set to 0)
555 
556  // End
557  return 0;
558 }
559 
560 // =====================================================================
561 // ---------------------------------------------------------------------
562 // ---------------------------------------------------------------------
563 // =====================================================================
564 /*
565  \fn StepInnerLoopInsideSubsetLoop
566  \param a_iteration : iteration index
567  \param a_subset : subset index
568  \param a_bed : bed position
569  \brief This function is called inside the subset loop and manages the main loop over the events
570  The loop over the events is multithreaded, and involves a thread lock in the case of image-based deformation
571  \details Step 1: Get the current event for that thread index
572  Step 2: Call the dynamic switch function that updates the current frame and gate numbers, and also detects involuntary patient motion
573  Perform image-based deformation if needed (thread lock is required)
574  Step 3: Compute the projection line
575  Step 4: Optimize in the data space (forward-proj, update, backward-proj)
576  \todo Check the correct implementation of thread-locked operation on different architectures
577  \return 0 if success, positive value otherwise.
578 */
579 int oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop(int a_iteration, int a_subset, int a_bed)
580 {
581  // Verbose
583  {
584  if (m_nbBeds>1) Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Start loop over events for bed " << a_bed+1 << endl << flush);
585  else Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Start loop over events" << endl << flush);
586  }
587 
588  // Reinitialize 4D gate indexes
590 
591  // Apply the bed offset for this bed position
593 
594  // Progression (increments of 2%)
596  {
597  cout << "0 10% 20% 30% 40% 50% 60% 70% 80% 90% 100%" << endl;
598  cout << "|----|----|----|----|----|----|----|----|----|----|" << endl;
599  cout << "|" << flush;
600  }
601  int progression_percentage_old = 0;
602  int progression_nb_bars = 0;
603  uint64_t progression_printing_index = 0;
604 
605  // Compute start and stop indices taking MPI into account (the vDataFile does that)
606  int64_t index_start = 0;
607  int64_t index_stop = 0;
608  m2p_DataFile[a_bed]->GetEventIndexStartAndStop(&index_start, &index_stop, a_subset, mp_nbSubsets[a_iteration]);
609 
610  // Multi-threading with OpenMP
611  m_releaseThreads = false;
612  m_nbThreadsWaiting = 0;
613 
614  // Synchronize MPI processes
615  #ifdef CASTOR_MPI
616  MPI_Barrier(MPI_COMM_WORLD);
617  #endif
618 
619  // Set the number of threads for projections (right after this loop, we set back the number of threads for image computations)
620  #ifdef CASTOR_OMP
621  omp_set_num_threads(mp_ID->GetNbThreadsForProjection());
622  #endif
623 
624  // This boolean is used to report any problem inside the parallel loop
625  bool problem = false;
626 
627  // Get the chrono manager singleton pointer
628  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
629 
630  // Launch the loop with precomputed start and stop and using a step equal to the number of subsets
631  int64_t index;
632  // Keep the static scheduling with a chunk size at 1, it is important
633  #pragma omp parallel for private(index) schedule(static, 1)
634  for ( index = index_start ; index < index_stop ; index += mp_nbSubsets[a_iteration] )
635  {
636  // Get the thread index
637  int th = 0;
638  #ifdef CASTOR_OMP
639  th = omp_get_thread_num();
640  #endif
641  // Print progression (do not log out with Cout here)
642  if (m_verbose>=2 && th==0 && mp_ID->GetMPIRank()==0)
643  {
644  if (progression_printing_index%1000==0)
645  {
646  int progression_percentage_new = ((int)( (((float)(index-index_start+1))/((float)(index_stop-index_start)) ) * 100.));
647  if (progression_percentage_new>=progression_percentage_old+2) // Increments of 2%
648  {
649  int nb_steps = (progression_percentage_new-progression_percentage_old)/2;
650  for (int i=0; i<nb_steps; i++)
651  {
652  cout << "-" << flush;
653  progression_nb_bars++;
654  }
655  progression_percentage_old += nb_steps*2;
656  }
657  }
658  progression_printing_index++;
659  }
660 
661  // Step 1: Get the current event for that thread index
662  p_ChronoManager->StartIterativeDataUpdateStep1(th);
663  #ifdef CASTOR_DEBUG
664  if (m_verbose>=4) Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step1: Get current event for that thread index " << endl);
665  #endif
666  vEvent* event = m2p_DataFile[a_bed]->GetEvent(index, th);
667  if (event==NULL)
668  {
669  Cerr("***** oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> An error occured while getting the event from index "
670  << index << " (thread " << th << ") !" << endl);
671  // Specify that there was a problem
672  problem = true;
673  // We must continue here because we are inside an OpenMP loop
674  continue;
675  }
676  p_ChronoManager->StopIterativeDataUpdateStep1(th);
677 
678  // Step 2: Call the dynamic switch function that updates the current frame and gate numbers, and also detects involuntary patient motion
679  p_ChronoManager->StartIterativeDataUpdateStep2(th);
680  #ifdef CASTOR_DEBUG
681  if (m_verbose>=4)
682  {
683  Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step2: Check for Dynamic event (frame/gate switch, image-based deformation " << endl);
684  }
685  #endif
686  int dynamic_switch_value = mp_ID->DynamicSwitch(index, event->GetTimeInMs(), a_bed, th);
687  // If the DYNAMIC_SWITCH_CONTINUE is returned, then it means that we are not yet at the first frame
688  if ( dynamic_switch_value == DYNAMIC_SWITCH_CONTINUE )
689  {
690  // Then we just skip this event
691  continue;
692  }
693  // Else, if the DYNAMIC_SWITCH_DEFORMATION is returned, then it means that a change of gate or involuntary motion has occured
694  // and is associated to a deformation
695  else if ( dynamic_switch_value == DYNAMIC_SWITCH_DEFORMATION )
696  {
697  // With multi-threads, we should wait all threads to be there to perform the deformation
698  #ifdef CASTOR_OMP
699  // Count the number of threads reaching this point
701  // Loop to lock the threads until the m_releaseThreads condition is enabled
702  while (ThreadBarrierFlag())
703  {
704  // Allow thread 0 to enter once all threads are in the loop
705  if ( (th==0) && ThreadBarrierCheck() )
706  {
707  // Perform here any deformation on the forward image.
709  // Release the threads and set the thread counter to 0
711  }
712  }
713  // Count the number of threads reaching this point
715  // Set the thread lock condition once all threads reached this point
717  #else
718  // Perform here any resp related deformations on the forward image
720  #endif
721  }
722  p_ChronoManager->StopIterativeDataUpdateStep2(th);
723 
724  // Step 3: Compute the projection line
725  p_ChronoManager->StartIterativeDataUpdateStep3(th);
726  #ifdef CASTOR_DEBUG
727  if (m_verbose>=4) Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step3: Compute the projection line " << endl);
728  #endif
730  if (line==NULL)
731  {
732  Cerr("***** oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> A problem occured while computing the projection line !" << endl);
733  // Specify that there was a problem
734  problem = true;
735  // We must continue here because we are inside an OpenMP loop
736  continue;
737  }
738  p_ChronoManager->StopIterativeDataUpdateStep3(th);
739 
740  // Step 4: Optimize in the data space (forward-proj, update, backward-proj)
741  p_ChronoManager->StartIterativeDataUpdateStep4(th);
742  #ifdef CASTOR_DEBUG
743  if (m_verbose>=4) Cout("oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> Step4: Optimize in the data space " << endl);
744  #endif
745  if (line->NotEmptyLine()) mp_OptimizerManager->DataUpdateStep( line,
746  event,
747  a_bed,
751  a_iteration,
752  th);
753  p_ChronoManager->StopIterativeDataUpdateStep4(th);
754  } // End of indices loop (OpenMP stops here)
755 
756  // Synchronize MPI processes
757  #ifdef CASTOR_MPI
758  MPI_Barrier(MPI_COMM_WORLD);
759  #endif
760 
761  // Set back the number of threads for image computation
762  #ifdef CASTOR_OMP
763  omp_set_num_threads(mp_ID->GetNbThreadsForImageComputation());
764  #endif
765 
766  // End of progression printing (do not log out with Cout here)
767  if (m_verbose>=2 && mp_ID->GetMPIRank()==0)
768  {
769  int progression_total_bars = 49;
770  for (int i=0; i<progression_total_bars-progression_nb_bars; i++) cout << "-";
771  cout << "|" << endl;
772  }
773 
774  // If a problem was encountered, then report it here
775  if (problem)
776  {
777  Cerr("***** oIterativeAlgorithm::StepInnerLoopInsideSubsetLoop() -> A problem occured inside the parallel loop over events !" << endl);
778  return 1;
779  }
780 
781  // End
782  return 0;
783 }
784 
785 // =====================================================================
786 // ---------------------------------------------------------------------
787 // ---------------------------------------------------------------------
788 // =====================================================================
789 /*
790  \fn StepPreProcessInsideSubsetLoop
791  \param a_iteration : iteration index
792  \param a_subset : subset index
793  \brief This function is called right after starting the subset loop.
794  Apply any kind of image processing on the backward image and
795  main image after backprojections
796  \details Synchronize parallel results
797  Apply image deformation/processing/convolution on the backward image
798  Update the image using the optimizer functions
799  Apply dynamic model/image processing/convolution on the main image
800  \return 0 if success, positive value otherwise.
801 */
803 {
804  if (m_verbose>=3) Cout("oIterativeAlgorithm::StepPostProcessInsideSubsetLoop ... " << endl);
805 
806  // Get the chrono manager singleton pointer
807  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
808 
809  // Merge parallel results
811 
812  // If mask provided, apply mask to sensitivity, the CleanNeverVisitedVoxels will take care of the images
814 
815  // Perform here any image-based deformations related to motion on the backward image.
817 
818  // Apply convolver to backward images and sensitivity-if-needed
819  p_ChronoManager->StartConvolution();
821  {
822  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying convolver to backward images !" << endl);
823  return 1;
824  }
825  p_ChronoManager->StopConvolution();
826 
827  // Call the post data update step function from the optimizer manager
828  mp_OptimizerManager->PostDataUpdateStep(a_iteration, m_nbIterations, a_subset, mp_nbSubsets);
829 
830  // Optimize in the image space (apply corrections, MAP and sensitivity); pass the number of subsets for list-mode sensitivity scaling
831  mp_OptimizerManager->ImageUpdateStep(a_iteration, mp_nbSubsets[a_iteration]);
832 
833  // Apply convolver to current estime images
834  p_ChronoManager->StartConvolution();
836  {
837  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying convolver to current estimate images !" << endl);
838  return 1;
839  }
840  p_ChronoManager->StopConvolution();
841 
842  // Post-update Dynamic Modeling step (if enabled). Manage either linear/non linear dynamic model (physiological or not)
843  if (mp_DynamicModelManager->ApplyDynamicModel(mp_ImageSpace, a_iteration, mp_nbSubsets[a_iteration]))
844  {
845  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying dynamic model to current estimate images !" << endl);
846  return 1;
847  }
848 
849  // Apply image processing to current estime images
851  {
852  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying image processing to current estimate images !" << endl);
853  return 1;
854  }
855 
856  // Save the sensitivity image in histogram mode, if asked for
858  {
859  // Get output manager to build the file name
860  sOutputManager* p_output_manager = sOutputManager::GetInstance();
861  // Build the file name
862  string temp_sens = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
863  stringstream temp_it; temp_it << a_iteration + 1;
864  stringstream temp_ss; temp_ss << a_subset + 1;
865  temp_sens.append("_it").append(temp_it.str()).append("_ss").append(temp_ss.str()).append("_sensitivity");
866  // Save sensitivity
868  }
869 
870  // Save the current image estimate if asked for
871  if (mp_ID->GetMPIRank()==0 && mp_outputIterations[a_iteration] && m_saveImageAfterSubsets)
872  {
873  // Verbose
874  if (m_verbose>=1) Cout("oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> Save image at iteration " << a_iteration+1 << " for subset " << a_subset+1 << endl);
875  // Build image from basis functions
877  // Apply post-convolution if needed
878  p_ChronoManager->StartConvolution();
880  {
881  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while convolving the output image !" << endl);
882  return 1;
883  }
884  p_ChronoManager->StopConvolution();
885  // Apply post-processing if needed
887  {
888  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying image processing the output image !" << endl);
889  return 1;
890  }
891  // Save Dynamic images
892  if (mp_DynamicModelManager->SaveParametricImages(a_iteration, a_subset))
893  {
894  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while saving parametric images related to the dynamic model !" << endl);
895  return 1;
896  }
897  // Apply output transaxial FOV masking
899  {
900  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying output FOV masking !" << endl);
901  return 1;
902  }
903  // Apply output flip
905  {
906  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while applying output flip !" << endl);
907  return 1;
908  }
909  // Save output image
910  if (mp_ImageSpace->SaveOutputImage(a_iteration, a_subset))
911  {
912  Cerr("***** oIterativeAlgorithm::StepPostProcessInsideSubsetLoop() -> A problem occured while saving output image !" << endl);
913  return 1;
914  }
915  }
916 
917  // End
918  return 0;
919 }
920 
921 // =====================================================================
922 // ---------------------------------------------------------------------
923 // ---------------------------------------------------------------------
924 // =====================================================================
925 /*
926  \fn StepAfterSubsetLoop
927  \param a_iteration : iteration index
928  \brief This function is called after finishing the subset loop.
929  \details Clean the main images from never visited voxels
930  Apply post-convolution/post-processing if needed
931  Write output images on disk as requested by the user
932  \todo : save parametric images using intrinsic basis functions ?
933  \return 0 if success, positive value otherwise.
934 */
936 {
937  if (m_verbose>=3) Cout("oIterativeAlgorithm::StepAfterSubsetLoop ... " << endl);
938 
939  // Get the chrono manager singleton pointer
940  sChronoManager* p_ChronoManager = sChronoManager::GetInstance();
941 
942  // Clean never visited voxels
944  // Save the main image
945  if (mp_ID->GetMPIRank()==0 && mp_outputIterations[a_iteration])
946  {
947  if (m_verbose>=1) Cout("oIterativeAlgorithm::StepAfterSubsetLoop() -> Save image at iteration " << a_iteration+1 << endl);
948  // Build image from basis functions
950  // Apply post-convolution if needed
951  p_ChronoManager->StartConvolution();
953  {
954  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while convolving the output image !" << endl);
955  return 1;
956  }
957  p_ChronoManager->StopConvolution();
958  // Apply post-processing if needed
960  {
961  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while applying image processing the output image !" << endl);
962  return 1;
963  }
964  // Save Dynamic images
966  {
967  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while saving parametric images related to the dynamic model !" << endl);
968  return 1;
969  }
970  // Todo : save parametric images using intrinsic basis functions
971  // Apply output transaxial FOV masking
973  {
974  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while applying output FOV masking !" << endl);
975  return 1;
976  }
977  // Apply output flip
979  {
980  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while applying output flip !" << endl);
981  return 1;
982  }
983  // Save output image
984  if (mp_ImageSpace->SaveOutputImage(a_iteration))
985  {
986  Cerr("***** oIterativeAlgorithm::StepAfterSubsetLoop() -> A problem occured while saving output image !" << endl);
987  return 1;
988  }
989  }
990 
991  // End
992  return 0;
993 }
994 
995 // =====================================================================
996 // ---------------------------------------------------------------------
997 // ---------------------------------------------------------------------
998 // =====================================================================
999 /*
1000  \fn StepAfterIterationLoop
1001  \brief This function is called at the end of the IterateCPU function.
1002  \details Free Memory for the imageSpace and
1003  some managers
1004  \return 0 if success, positive value otherwise.
1005 */
1007 {
1008  if (m_verbose>=2) Cout("oIterativeAlgorithm::StepAfterIterationLoop ... " << endl);
1009 
1010  // Deallocate everything
1020 
1021  // Normal end
1022  return 0;
1023 }
1024 
1025 // =====================================================================
1026 // ---------------------------------------------------------------------
1027 // ---------------------------------------------------------------------
1028 // =====================================================================
1029 /*
1030  \fn InitSpecificOptions
1031  \param a_specificOptions
1032  \brief Set the selected output iterations
1033  \details Not Implemented yet
1034  \return 0
1035 */
1036 int oIterativeAlgorithm::InitSpecificOptions(string a_specificOptions)
1037 {
1038  if (m_verbose>=2) Cout("oIterativeAlgorithm::InitSpecificOptions ... " << endl);
1039 
1040  return 0;
1041 }
1042 
1043 // =====================================================================
1044 // ---------------------------------------------------------------------
1045 // ---------------------------------------------------------------------
1046 // =====================================================================
1047 // ----- Functions dedicated to barrier for multithreading -----
1048 // TODO we should keep in mind that different compiler may deal differently with OpenMP
1049 // (as the compiler may take liberties and optimize shared variables as register variables,
1050 // leading to non-concurrent observations of the variable.
1051 // And this is not cool as this can completely stop the reconstruction process)
1052 // The omp-flush directives should prevent that, but further tests may be required
1053 
1054 /*
1055  \fn ThreadBarrierIncrement
1056  \brief Increment (thread safe) the m_nbThreadsWaiting variable.
1057 */
1059 {
1060  #ifdef CASTOR_DEBUG
1061  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierIncrement ... " << endl);
1062  #endif
1063 
1064  #pragma omp critical(increment_waiting_threads)
1065  {
1067  }
1068 }
1069 
1070 // =====================================================================
1071 // ---------------------------------------------------------------------
1072 // ---------------------------------------------------------------------
1073 // =====================================================================
1074 /*
1075  \fn ThreadBarrierFlag
1076  \brief Check if the m_releaseThreads boolean is enabled or not
1077  \return true if the threads can be released, false otherwise
1078 */
1080 {
1081  #ifdef CASTOR_DEBUG
1082  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierFlag ... " << endl);
1083  #endif
1084 
1085  // Check any change in m_releaseThreads
1086  #pragma omp flush(m_releaseThreads)
1087  return !m_releaseThreads;
1088 }
1089 
1090 // =====================================================================
1091 // ---------------------------------------------------------------------
1092 // ---------------------------------------------------------------------
1093 // =====================================================================
1094 /*
1095  \fn ThreadBarrierCheck
1096  \brief Check if all the threads are currently in waiting condition
1097  \return true if all the threads are waiting, false otherwise
1098 */
1100 {
1101  #ifdef CASTOR_DEBUG
1102  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierCheck ... " << endl);
1103  #endif
1104 
1105  #pragma omp flush(m_nbThreadsWaiting) // Check any change in m_nbThreadsWaiting
1107 }
1108 
1109 // =====================================================================
1110 // ---------------------------------------------------------------------
1111 // ---------------------------------------------------------------------
1112 // =====================================================================
1113 /*
1114  \fn ThreadBarrierSetOff
1115  \brief Disable the thread locking variable (thread safe), and reset the number of waiting threads to 0
1116 */
1118 {
1119  #ifdef CASTOR_DEBUG
1120  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierSetOff ... " << endl);
1121  #endif
1122 
1123  m_releaseThreads = true;
1124  #pragma omp flush(m_releaseThreads) // Make sure the value of m_releaseThreads is propagated to other threads
1125  m_nbThreadsWaiting = 0;
1126  #pragma omp flush(m_nbThreadsWaiting) // Same for m_nbThreadsWaiting
1127 }
1128 
1129 // =====================================================================
1130 // ---------------------------------------------------------------------
1131 // ---------------------------------------------------------------------
1132 // =====================================================================
1133 /*
1134  \fn ThreadBarrierSetOn
1135  \brief Enable the thread locking variable (thread safe), and reset the number of waiting threads to 0
1136 */
1138 {
1139  #ifdef CASTOR_DEBUG
1140  if (m_verbose>=4) Cout("oIterativeAlgorithm::ThreadBarrierSetOn ... " << endl);
1141  #endif
1142 
1143  m_releaseThreads = false;
1144  #pragma omp flush(m_releaseThreads) // Make sure the value of m_releaseThreads is propagated to other threads
1145  m_nbThreadsWaiting = 0;
1146  #pragma omp flush(m_nbThreadsWaiting) // Same for m_nbThreadsWaiting
1147 }
void StopIterativeDataUpdateStep1(int a_thread)
Stop the timer for duration of iterative data update step 1.
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
vector< string > m_pathToMultiModalImg
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:671
#define MODE_HISTOGRAM
Definition: vDataFile.hh:59
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:266
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...
FLTNB GetInitialValue()
Return the initial image value used by the optimizer, explaining why the eponym function of vOptimize...
#define FLTNB
Definition: gVariables.hh:81
int InitMaskImage(const string &a_pathToImage)
Memory allocation and initialization for the mask image.
Definition: oImageSpace.cc:675
void DeallocateImage()
Free memory for the main image matrices.
Definition: oImageSpace.cc:118
#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 InitMultiModalImage(const vector< string > &a_pathToMultiModalImage)
Memory allocation and initialization for the multimodal image matrices.
Definition: oImageSpace.cc:582
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 DeallocateMultiModalImage()
Free memory for the multimodal image.
Definition: oImageSpace.cc:641
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:463
void StopIterativeDataUpdateStep3(int a_thread)
Stop the timer for duration of iterative data update step 3.
void InstantiateImage()
Allocate memory for the main image matrices.
Definition: oImageSpace.cc:82
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:221
int SetOutputIterations(const string &a_outputIterations)
Set the selected output iterations.
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.
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:388
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 :
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...
Declaration of class oIterativeAlgorithm.
virtual int StepAfterSubsetLoop(int a_iteration)
This function is called after finishing the subset loop.
void StartIterativeDataUpdateStep1(int a_thread)
Start the timer for duration of iterative data update step 1.
void DeallocateMaskImage()
Free memory for the mask image.
Definition: oImageSpace.cc:716
int PostDataUpdateStep(int a_iteration, int a_nbIterations, int a_subset, int *ap_nbSubsets)
A function that simply calls the eponym function from the vOptimizer.
void StopIterativeDataUpdateStep4(int a_thread)
Stop the timer for duration of iterative data update step 4.
oImageProcessingManager * mp_ImageProcessingManager
int m_nbThreadsWaiting
bool ThreadBarrierFlag()
Check if the m_releaseThreads boolean is enabled or not.
void StartIterativeDataUpdateStep2(int a_thread)
Start the timer for duration of iterative data update step 2.
oImageConvolverManager * mp_ImageConvolverManager
vEvent * GetEvent(int64_t a_eventIndex, int a_th=0)
Definition: vDataFile.cc:598
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...
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...
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.
Declaration of class iEventHistoCT.
oDeformationManager * mp_DeformationManager
void ResetCurrentDynamicIndices()
Call the eponym function from the oDynamicDataManager class using the member object.
void ComputeOutputImage()
bool m_releaseThreads
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 'm_UseModel' 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.
void InstantiateForwardImage()
Allocate memory for the forward image matrices.
Definition: oImageSpace.cc:151
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)
int ImageUpdateStep(int a_iteration, int a_nbSubsets)
A function dedicated to the update step in the image space (performed after the loop on events) ...
int DataUpdateStep(oProjectionLine *ap_Line, 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) ...
virtual int StepBeforeIterationLoop()
This function is called at the beginning of the IterateCPU 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.
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:792
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.
void StopConvolution()
Stop the timer for duration of convolution.
Mother class for the Event objects.
Definition: vEvent.hh:43
int GetNbThreadsForImageComputation()
Get the number of threads used for image operations.
int GetDataMode()
Definition: vDataFile.hh:300
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 PreDataUpdateStep(int a_iteration, int a_nbIterations, int a_subset, int *ap_nbSubsets)
A function that simply calls the eponym function from the vOptimizer.
int GetNbThreadsForProjection()
Get the number of threads used for projections.
oOptimizerManager * mp_OptimizerManager
#define Cout(MESSAGE)
This class is designed to manage some profiling of the code.
oProjectorManager * mp_ProjectorManager
bool Checked()
Simply check that the image dimensions and verbosity has been set.
Definition: oImageSpace.hh:654
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:187
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
oImageDimensionsAndQuantification * mp_ID
virtual ~oIterativeAlgorithm()
oIterativeAlgorithm destructor.
void InstantiateOutputImage()
Instanciate Image matrices dedicated to output writing on disk.
Definition: oImageSpace.cc:745
#define DYNAMIC_SWITCH_DEFORMATION