CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
src/algorithm/vAlgorithm.cc
Go to the documentation of this file.
1 
8 #include "gVariables.hh"
9 #include "vAlgorithm.hh"
10 #include "iEventHistoCT.hh"
11 
12 // =====================================================================
13 // ---------------------------------------------------------------------
14 // ---------------------------------------------------------------------
15 // =====================================================================
16 
18 {
19  // Default number of iterations and subsets (1 for each)
20  m_nbIterations = 1;
21  mp_nbSubsets = (int*)malloc(1*sizeof(int));
22  mp_nbSubsets[0] = 1;
23  mp_outputIterations = NULL;
24  // Default number of first iterations to be skipped (0)
26  // set all members to default values
27  m_verbose = -1;
28  m_flagGPU = false;
29  mp_ID = NULL;
30  m2p_DataFile= NULL;
31  mp_ProjectorManager = NULL;
32  mp_OptimizerManager = NULL;
33  mp_DeformationManager = NULL;
35  mp_ImageSpace = NULL;
38  m_nbBeds = -1;
39  m_pathToInitialImg = "";
44 
45 }
46 
47 // =====================================================================
48 // ---------------------------------------------------------------------
49 // ---------------------------------------------------------------------
50 // =====================================================================
51 
53 {
54  // Delete some tables
55  if (mp_nbSubsets!=NULL) free(mp_nbSubsets);
57 }
58 
59 // =====================================================================
60 // ---------------------------------------------------------------------
61 // ---------------------------------------------------------------------
62 // =====================================================================
63 
64 int vAlgorithm::SetNbIterationsAndSubsets(const string& a_nbIterationsSubsets)
65 {
66  // If the string is empty, then just keep the defaults
67  if (a_nbIterationsSubsets=="") return 0;
68  // Otherwise, reset the number of iterations to 0
69  m_nbIterations = 0;
70 
71  // For algorithms based on optimization, there can be a number of subiterations that we get here
72  int nb_sub_iter = 1;
74 
75  // Copy the string
76  string buf = a_nbIterationsSubsets;
77 
78  // Loop to search for commas
79  size_t pos_comma;
80  while ((pos_comma=buf.find_first_of(","))!=string::npos)
81  {
82  // Get the substring before the comma
83  string sub_buf = buf.substr(0,pos_comma);
84  // Search for columns
85  size_t pos_column = sub_buf.find_first_of(":");
86  if (pos_column==string::npos || pos_column==0)
87  {
88  Cerr("***** vAlgorithm::SetNbIterationsAndSubsets() -> Syntax problem in number of iterations and subsets !" << endl);
89  return 1;
90  }
91  int iter = atoi( sub_buf.substr(0,pos_column).c_str() );
92  // Multiply number of user-defined iterations by the number of sub steps
93  iter *= nb_sub_iter;
94  // Continue reading the command line
95  int subs = atoi( sub_buf.substr(pos_column+1).c_str() );
96  mp_nbSubsets = (int*)realloc(mp_nbSubsets,(m_nbIterations+iter)*sizeof(int));
97  for (int it=0; it<iter; it++) mp_nbSubsets[m_nbIterations+it] = subs;
98  // Number of iterations takes into account the number of sub steps
99  m_nbIterations += iter;
100  buf = buf.substr(pos_comma+1);
101  }
102 
103  // Last couple of iterations:subsets
104  size_t pos_column = buf.find_first_of(":");
105  if (pos_column==string::npos || pos_column==0)
106  {
107  Cerr("***** vAlgorithm::SetNbIterationsAndSubsets() -> Syntax problem in number of iterations and subsets !" << endl);
108  return 1;
109  }
110  int iter = atoi( buf.substr(0,pos_column).c_str() );
111  // Multiply number of user-defined iterations by the number of sub steps
112  iter *= nb_sub_iter;
113  // Continue reading the command line
114  int subs = atoi( buf.substr(pos_column+1).c_str() );
115  mp_nbSubsets = (int*)realloc(mp_nbSubsets,(m_nbIterations+iter)*sizeof(int));
116  for (int it=0; it<iter; it++) mp_nbSubsets[m_nbIterations+it] = subs;
117  // Number of iterations takes into account the number of sub steps
118  m_nbIterations += iter;
119 
121  {
122  Cout("vAlgorithm::SetNbIterationsAndSubsets() -> Selected numbers of subsets for each iteration:" << endl);
123  Cout(" Iteration / number of subsets "<< endl);
124  for (int it=0 ; it<m_nbIterations/nb_sub_iter ; it++) Cout(" " << it+1 << " / " << mp_nbSubsets[it*nb_sub_iter] << endl); }
125  // End
126  return 0;
127 }
128 
129 // =====================================================================
130 // ---------------------------------------------------------------------
131 // ---------------------------------------------------------------------
132 // =====================================================================
133 
134 int vAlgorithm::SetOutputIterations(const string& a_outputIterations)
135 {
136  if (m_verbose>=VERBOSE_DETAIL) Cout("vAlgorithm::SetOutputIterations ..." << endl);
137 
138  // For algorithms based on optimization, there can be a number of subiterations that we get here
139  int nb_sub_iter = 1;
141 
142  // Allocate the output iterations boolean table
143  mp_outputIterations = (bool*)malloc(m_nbIterations*sizeof(bool));
144  for (int it=0; it<m_nbIterations; it++) mp_outputIterations[it] = false;
145 
146  // If the list is empty, we save all iterations by default
147  if (a_outputIterations=="")
148  {
149  // Only save image at end of iteration, not for all sub steps
150  for (int it=nb_sub_iter-1; it<m_nbIterations; it=it+nb_sub_iter) mp_outputIterations[it] = true;
151  return 0;
152  }
153 
154  // Copy the string
155  string buf = a_outputIterations;
156 
157  // Loop to search for commas
158  size_t pos_comma;
159  int current_iteration = -1;
160  while ((pos_comma=buf.find_first_of(","))!=string::npos)
161  {
162  // Get the substring before the comma
163  string sub_buf = buf.substr(0,pos_comma);
164  // Search for columns
165  size_t pos_column = sub_buf.find_first_of(":");
166  if (pos_column==string::npos || pos_column==0)
167  {
168  Cerr("***** vAlgorithm::SetOutputIterations() -> Syntax problem in output iteration numbers !" << endl);
169  return 1;
170  }
171  int step_iteration = atoi( sub_buf.substr(0,pos_column).c_str() );
172  step_iteration *= nb_sub_iter;
173  if (step_iteration<=0)
174  {
175  Cerr("***** vAlgorithm::SetOutputIterations() -> Iteration step must be strictly positive (here it is " << step_iteration << ") !" << endl);
176  return 1;
177  }
178  int stop_iteration = atoi( sub_buf.substr(pos_column+1).c_str() );
179  stop_iteration *= nb_sub_iter;
180  if (stop_iteration<=current_iteration)
181  {
182  Cerr("***** oIterationAlgorithm::SetOutputIterations() -> Output iteration number must be increasing through provided couples !" << endl);
183  return 1;
184  }
185  // Output iterations taken into account the number of sub steps
186  for (int it=current_iteration+step_iteration; it<stop_iteration; it+=step_iteration) if (it<m_nbIterations) mp_outputIterations[it] = true;
187  current_iteration = stop_iteration-1;
188  buf = buf.substr(pos_comma+1);
189  }
190 
191  // Last couple of iterations:subsets
192  size_t pos_column = buf.find_first_of(":");
193  if (pos_column==string::npos || pos_column==0)
194  {
195  // Special case if -1 is provided, it means we save the last iteration
196  if (buf=="-1")
197  {
198  mp_outputIterations[m_nbIterations-1] = true;
199  // We directly exist here
200  return 0;
201  }
202  else
203  {
204  Cerr("***** vAlgorithm::SetOutputIterations() -> Syntax problem in output iteration numbers !" << endl);
205  return 1;
206  }
207  }
208  int step_iteration = atoi( buf.substr(0,pos_column).c_str() );
209  step_iteration *= nb_sub_iter;
210  if (step_iteration<=0)
211  {
212  Cerr("***** vAlgorithm::SetOutputIterations() -> Iteration step must be strictly positive (here it is " << step_iteration << ") !" << endl);
213  return 1;
214  }
215  int stop_iteration = atoi( buf.substr(pos_column+1).c_str() );
216  stop_iteration *= nb_sub_iter;
217  if (stop_iteration<=current_iteration)
218  {
219  Cerr("***** oIterationAlgorithm::SetOutputIterations() -> Output iteration number must be increasing through provided couples !" << endl);
220  return 1;
221  }
222  // Output iterations taken into account the number of sub steps
223  for (int it=current_iteration+step_iteration; it<stop_iteration; it+=step_iteration) if (it<m_nbIterations) mp_outputIterations[it] = true;
224 
225  // End
226  return 0;
227 }
228 
229 // =====================================================================
230 // ---------------------------------------------------------------------
231 // ---------------------------------------------------------------------
232 // =====================================================================
233 
234 int vAlgorithm::Run()
235 {
236  if(m_verbose>=2) Cout("vAlgorithm::Iterate ..."<< endl);
237 
238  #ifdef CASTOR_GPU
239  if (m_flagGPU)
240  return RunGPU();
241  else
242  return RunCPU();
243  #else
244  return RunCPU();
245  #endif
246 }
247 
248 // =====================================================================
249 // ---------------------------------------------------------------------
250 // ---------------------------------------------------------------------
251 // =====================================================================
252 
253 int vAlgorithm::RunCPU()
254 {
255 // Here is the scheme of all functions call related to each step of the iterative reconstruction:
256 // StepBeforeIterationLoop()
257 // / Loop on iterations
258 // | StepBeforeSubsetLoop(iteration)
259 // | / Loop on subsets
260 // | | StepPreProcessInsideSubsetLoop(iteration,subset)
261 // | | / Loop on bed positions
262 // | | | StepInnerLoopInsideSubsetLoop(iteration,subset,bed)
263 // | | StepPostProcessInsideSubsetLoop(iteration,subset)
264 // | StepAfterSubsetLoop(iteration)
265 // StepAfterIterationLoop()
266 
267  // Synchronize MPI processes
268  #ifdef CASTOR_MPI
269  MPI_Barrier(MPI_COMM_WORLD);
270  #endif
271 
272  // For algorithms based on optimization, there can be a number of subiterations that we get here
273  int nb_sub_iter = 1;
275 
276  // Verbose
277  if (m_verbose>=1) Cout("vAlgorithm::IterateCPU() -> Start algorithm for " << m_nbIterations / nb_sub_iter << " iterations" << endl);
278 
279  // Initial clock for execution time
280  clock_t clock_start_whole = clock();
281  time_t time_start_whole = time(NULL);
282 
283  // Call before iteration function
284  // Call before iteration function
286  {
287  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepBeforeIterationLoop() function !" << endl);
288  return 1;
289  }
290 
291  // Loop on iterations without doing the iterations that the user asked to skip (times the number of sub steps)
292  for (int iteration=m_nbSkippedIterations*nb_sub_iter ; iteration<m_nbIterations ; iteration++)
293  {
294 
295  // Call before subset function
296  if (StepBeforeSubsetLoop(iteration))
297  {
298  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepBeforeSubsetLoop() function !" << endl);
299  return 1;
300  }
301 
302  // Loop on subsets
303  for (int subset=0 ; subset<mp_nbSubsets[iteration] ; subset++)
304  {
305  // Verbose
306  if (m_verbose>=1)
307  {
308  // Number of displayed iteration is current iteration divided by number of sub steps
309  Cout("vAlgorithm::IterateCPU() -> Start iteration " << iteration / nb_sub_iter + 1 << "/" << m_nbIterations / nb_sub_iter
310  << " sub step " << iteration%nb_sub_iter+1 << "/" << nb_sub_iter
311  << " subset " << subset+1 << "/" << mp_nbSubsets[iteration] << endl);
312  }
313 
314  // Clock start for subset execution time
315  clock_t clock_start_subset = clock();
316  time_t time_start_subset = time(NULL);
317 
318  // Call pre-process inside subset loop function
319  if (StepPreProcessInsideSubsetLoop(iteration,subset))
320  {
321  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepPreProcessInsideSubsetLoop() function !" << endl);
322  return 1;
323  }
324 
325  // Synchronize MPI processes
326  #ifdef CASTOR_MPI
327  MPI_Barrier(MPI_COMM_WORLD);
328  #endif
329 
330  // Loop on bed positions
331  for (int bed=0 ; bed<m_nbBeds ; bed++)
332  {
333  // Call the inner loop on events inside subset loop function
334  if (StepInnerLoopInsideSubsetLoop(iteration,subset,bed))
335  {
336  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepInnerLoopInsideSubsetLoop() function !" << endl);
337  return 1;
338  }
339  } // End of beds loop
340 
341  // Synchronize MPI processes
342  #ifdef CASTOR_MPI
343  MPI_Barrier(MPI_COMM_WORLD);
344  #endif
345 
346  // Call post-process inside subset loop function
347  if (StepPostProcessInsideSubsetLoop(iteration,subset))
348  {
349  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepPostProcessInsideSubsetLoop() function !" << endl);
350  return 1;
351  }
352 
353  // Clock stop for subset execution time
354  clock_t clock_stop_subset = clock();
355  time_t time_stop_subset = time(NULL);
356  if (m_verbose>=2) Cout ("vAlgorithm::IterateCPU() -> Time spent for subset " << subset+1 << " | User: " << time_stop_subset-time_start_subset
357  << " sec | CPU: " << (clock_stop_subset-clock_start_subset)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
358 
359 
360  } // End of subsets loop
361 
362 
363  // Call after subset function
364  if (StepAfterSubsetLoop(iteration))
365  {
366  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepAfterSubsetLoop() function !" << endl);
367  return 1;
368  }
369 
370 
371  } // End of iterations loop
372 
373  // Call after iteration function
375  {
376  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepAfterIterationLoop() function !" << endl);
377  return 1;
378  }
379 
380  // Final clock for execution time
381  clock_t clock_stop_whole = clock();
382  time_t time_stop_whole = time(NULL);
383  if (m_verbose>=1) Cout("vAlgorithm::IterateCPU() -> Total time spent | User: " << time_stop_whole-time_start_whole
384  << " sec | CPU: " << (clock_stop_whole-clock_start_whole)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
385 
386 
387  return 0;
388 }
389 
390 // =====================================================================
391 // ---------------------------------------------------------------------
392 // ---------------------------------------------------------------------
393 // =====================================================================
394 
396 {
397  if (m_verbose>=2) Cout("vAlgorithm::StepBeforeIterationLoop ... " << endl);
398 
399  // Check that the image space is already allocated
400  if (mp_ImageSpace==NULL || !mp_ImageSpace->Checked())
401  {
402  Cerr("***** vAlgorithm::StepBeforeIterationLoop() -> Image space has not been created or was not checked !" << endl);
403  return 1;
404  }
405 
406  // Instantiate all the required images from the oImageSpace
412 
414  {
415  Cerr("***** vAlgorithm::StepBeforeIterationLoop()-> Error during attenuation image initialization !" << endl);
416  return 1;
417  }
418 
420  {
421  Cerr("***** vAlgorithm::StepBeforeIterationLoop() -> An error occurred while initializing the sensitivity image !" << endl);
422  return 1;
423  }
424 
426  {
427  Cerr("***** vAlgorithm::StepBeforeIterationLoop()-> Error during multimodal image initialization !" << endl);
428  return 1;
429  }
430 
431  // End
432  return 0;
433 }
434 
435 // =====================================================================
436 // ---------------------------------------------------------------------
437 // ---------------------------------------------------------------------
438 // =====================================================================
439 
440 int vAlgorithm::StepBeforeSubsetLoop(int a_iteration)
441 {
442  // End
443  return 0;
444 }
445 
446 // =====================================================================
447 // ---------------------------------------------------------------------
448 // ---------------------------------------------------------------------
449 // =====================================================================
450 
451 int vAlgorithm::StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
452 {
453  // End
454  return 0;
455 }
456 
457 // =====================================================================
458 // ---------------------------------------------------------------------
459 // ---------------------------------------------------------------------
460 // =====================================================================
461 
462 int vAlgorithm::StepPostProcessInsideSubsetLoop(int a_iteration, int a_subset)
463 {
464  // End
465  return 0;
466 }
467 
468 // =====================================================================
469 // ---------------------------------------------------------------------
470 // ---------------------------------------------------------------------
471 // =====================================================================
472 
473 int vAlgorithm::StepAfterSubsetLoop(int a_iteration)
474 {
475  // End
476  return 0;
477 }
478 
479 // =====================================================================
480 // ---------------------------------------------------------------------
481 // ---------------------------------------------------------------------
482 // =====================================================================
483 
485 {
486  if (m_verbose>=2) Cout("vAlgorithm::StepAfterIterationLoop ... " << endl);
487 
488  // Deallocate everything
495 
496  // Normal end
497  return 0;
498 }
499 
500 // =====================================================================
501 // ---------------------------------------------------------------------
502 // ---------------------------------------------------------------------
503 // =====================================================================
504 
505 int vAlgorithm::StepImageOutput(int a_iteration, int a_subset)
506 {
507  // Normal end
508  return 0;
509 }
510 
511 // =====================================================================
512 // ---------------------------------------------------------------------
513 // ---------------------------------------------------------------------
514 // =====================================================================
515 
516 int vAlgorithm::InitSpecificOptions(string a_specificOptions)
517 {
518  (void)a_specificOptions; // avoid 'unused parameter' compil. warnings
519  // Normal end
520  return 0;
521 }
virtual int InitSpecificOptions(string a_specificOptions)
int SetNbIterationsAndSubsets(const string &a_nbIterationsSubsets)
oImageDimensionsAndQuantification * mp_ID
#define Cerr(MESSAGE)
vector< string > m_pathToMultiModalImg
void DeallocateImage()
Free memory for the main image matrices.
int InitMultiModalImage(const vector< string > &a_pathToMultiModalImage)
int RunCPU()
Perform the iterative loop of the algorithm. Function designed to be executed on the CPU only...
void DeallocateMultiModalImage()
Free memory for the multimodal image.
virtual int StepBeforeSubsetLoop(int a_iteration)
int Run()
Just call either the RunCPU or the RunGPU function as asked for.
void DeallocateSensitivityImage()
Free memory for the sensitivity image matrices.
oOptimizerManager * mp_OptimizerManager
void InstantiateImage()
Allocate memory for the main image matrices.
void InstantiateSensitivityImage(const string &a_pathToSensitivityImage)
int InitAttenuationImage(const string &a_pathToAtnImage)
virtual int StepImageOutput(int a_iteration, int a_subset=-1)=0
void DeallocateVisitedVoxelsImage()
Free memory for the image matrix containing binary information regarding which 3D voxels have been vi...
virtual int StepInnerLoopInsideSubsetLoop(int a_iteration, int a_subset, int a_bed)=0
oImageProcessingManager * mp_ImageProcessingManager
virtual int StepPostProcessInsideSubsetLoop(int a_iteration, int a_subset)
vAlgorithm()
vAlgorithm constructor. Initialize the member variables to their default values.
void InstantiateVisitedVoxelsImage()
Memory allocation and initialization for the image matrix containing binary information regarding whi...
virtual int StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
virtual int StepBeforeIterationLoop()
This function is called at the beginning of the RunCPU function.
virtual ~vAlgorithm()
vAlgorithm destructor.
oImageConvolverManager * mp_ImageConvolverManager
void InstantiateForwardImage()
Allocate memory for the forward image matrices.
int GetNbSubIterationsInOneIteration()
Get the number of sub iterations in one iteration.
oProjectorManager * mp_ProjectorManager
virtual int StepAfterSubsetLoop(int a_iteration)
oDeformationManager * mp_DeformationManager
virtual int StepAfterIterationLoop()
This function is called at the end of the RunCPU function.
void DeallocateOutputImage()
Free memory for the Image matrices dedicated to output writing on disk.
int SetOutputIterations(const string &a_outputIterations)
bool Checked()
Simply check that the image dimensions and verbosity has been set.
int InitSensitivityImage(const string &a_pathToSensitivityImage)
void DeallocateForwardImage()
Free memory for the forward image matrices.
#define Cout(MESSAGE)
oDynamicModelManager * mp_DynamicModelManager
void InstantiateOutputImage()
Instanciate Image matrices dedicated to output writing on disk.