CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
code/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  // set all members to default values
25  m_verbose = -1;
26  m_flagGPU = false;
27  mp_ID = NULL;
28  m2p_DataFile= NULL;
29  mp_ProjectorManager = NULL;
30  mp_OptimizerManager = NULL;
31  mp_DeformationManager = NULL;
33  mp_ImageSpace = NULL;
36  m_nbBeds = -1;
37  m_pathToInitialImg = "";
42 
43 }
44 
45 // =====================================================================
46 // ---------------------------------------------------------------------
47 // ---------------------------------------------------------------------
48 // =====================================================================
49 
51 {
52  // Delete some tables
53  if (mp_nbSubsets!=NULL) free(mp_nbSubsets);
55 }
56 
57 // =====================================================================
58 // ---------------------------------------------------------------------
59 // ---------------------------------------------------------------------
60 // =====================================================================
61 
62 int vAlgorithm::SetNbIterationsAndSubsets(const string& a_nbIterationsSubsets)
63 {
64  // If the string is empty, then just keep the defaults
65  if (a_nbIterationsSubsets=="") return 0;
66  // Otherwise, reset the number of iterations to 0
67  m_nbIterations = 0;
68 
69  // Copy the string
70  string buf = a_nbIterationsSubsets;
71 
72  // Loop to search for commas
73  size_t pos_comma;
74  while ((pos_comma=buf.find_first_of(","))!=string::npos)
75  {
76  // Get the substring before the comma
77  string sub_buf = buf.substr(0,pos_comma);
78  // Search for columns
79  size_t pos_column = sub_buf.find_first_of(":");
80  if (pos_column==string::npos || pos_column==0)
81  {
82  Cerr("***** vAlgorithm::SetNbIterationsAndSubsets() -> Syntax problem in number of iterations and subsets !" << endl);
83  return 1;
84  }
85  int iter = atoi( sub_buf.substr(0,pos_column).c_str() );
86  int subs = atoi( sub_buf.substr(pos_column+1).c_str() );
87  mp_nbSubsets = (int*)realloc(mp_nbSubsets,(m_nbIterations+iter)*sizeof(int));
88  for (int it=0; it<iter; it++) mp_nbSubsets[m_nbIterations+it] = subs;
89  m_nbIterations += iter;
90  buf = buf.substr(pos_comma+1);
91  }
92 
93  // Last couple of iterations:subsets
94  size_t pos_column = buf.find_first_of(":");
95  if (pos_column==string::npos || pos_column==0)
96  {
97  Cerr("***** vAlgorithm::SetNbIterationsAndSubsets() -> Syntax problem in number of iterations and subsets !" << endl);
98  return 1;
99  }
100  int iter = atoi( buf.substr(0,pos_column).c_str() );
101  int subs = atoi( buf.substr(pos_column+1).c_str() );
102  mp_nbSubsets = (int*)realloc(mp_nbSubsets,(m_nbIterations+iter)*sizeof(int));
103  for (int it=0; it<iter; it++) mp_nbSubsets[m_nbIterations+it] = subs;
104  m_nbIterations += iter;
105 
106  if (m_verbose>=3)
107  {
108  Cout("vAlgorithm::SetNbIterationsAndSubsets() -> Selected numbers of subsets for each iteration:" << endl);
109  Cout(" Iteration / number of subsets "<< endl);
110  for (int it=0 ; it<m_nbIterations ; it++) Cout(" " << it+1 << " / " << mp_nbSubsets[it] << endl);
111  }
112  // End
113  return 0;
114 }
115 
116 // =====================================================================
117 // ---------------------------------------------------------------------
118 // ---------------------------------------------------------------------
119 // =====================================================================
120 
121 int vAlgorithm::SetOutputIterations(const string& a_outputIterations)
122 {
123  if (m_verbose>=2) Cout("vAlgorithm::SetOutputIterations ..." << endl);
124 
125  // Allocate the output iterations boolean table
126  mp_outputIterations = (bool*)malloc(m_nbIterations*sizeof(bool));
127  for (int it=0; it<m_nbIterations; it++) mp_outputIterations[it] = false;
128 
129  // If the list is empty, we save all iterations by default
130  if (a_outputIterations=="")
131  {
132  for (int it=0; it<m_nbIterations; it++) mp_outputIterations[it] = true;
133  return 0;
134  }
135 
136  // Copy the string
137  string buf = a_outputIterations;
138 
139  // Loop to search for commas
140  size_t pos_comma;
141  int current_iteration = -1;
142  while ((pos_comma=buf.find_first_of(","))!=string::npos)
143  {
144  // Get the substring before the comma
145  string sub_buf = buf.substr(0,pos_comma);
146  // Search for columns
147  size_t pos_column = sub_buf.find_first_of(":");
148  if (pos_column==string::npos || pos_column==0)
149  {
150  Cerr("***** vAlgorithm::SetOutputIterations() -> Syntax problem in output iteration numbers !" << endl);
151  return 1;
152  }
153  int step_iteration = atoi( sub_buf.substr(0,pos_column).c_str() );
154  if (step_iteration<=0)
155  {
156  Cerr("***** vAlgorithm::SetOutputIterations() -> Iteration step must be strictly positive (here it is " << step_iteration << ") !" << endl);
157  return 1;
158  }
159  int stop_iteration = atoi( sub_buf.substr(pos_column+1).c_str() );
160  if (stop_iteration<=current_iteration)
161  {
162  Cerr("***** oIterationAlgorithm::SetOutputIterations() -> Output iteration number must be increasing through provided couples !" << endl);
163  return 1;
164  }
165  for (int it=current_iteration+step_iteration; it<stop_iteration; it+=step_iteration) if (it<m_nbIterations) mp_outputIterations[it] = true;
166  current_iteration = stop_iteration-1;
167  buf = buf.substr(pos_comma+1);
168  }
169 
170  // Last couple of iterations:subsets
171  size_t pos_column = buf.find_first_of(":");
172  if (pos_column==string::npos || pos_column==0)
173  {
174  // Special case if -1 is provided, it means we save the last iteration
175  if (buf=="-1")
176  {
177  mp_outputIterations[m_nbIterations-1] = true;
178  // We directly exist here
179  return 0;
180  }
181  else
182  {
183  Cerr("***** vAlgorithm::SetOutputIterations() -> Syntax problem in output iteration numbers !" << endl);
184  return 1;
185  }
186  }
187  int step_iteration = atoi( buf.substr(0,pos_column).c_str() );
188  if (step_iteration<=0)
189  {
190  Cerr("***** vAlgorithm::SetOutputIterations() -> Iteration step must be strictly positive (here it is " << step_iteration << ") !" << endl);
191  return 1;
192  }
193  int stop_iteration = atoi( buf.substr(pos_column+1).c_str() );
194  if (stop_iteration<=current_iteration)
195  {
196  Cerr("***** oIterationAlgorithm::SetOutputIterations() -> Output iteration number must be increasing through provided couples !" << endl);
197  return 1;
198  }
199  for (int it=current_iteration+step_iteration; it<stop_iteration; it+=step_iteration) if (it<m_nbIterations) mp_outputIterations[it] = true;
200 
201  // End
202  return 0;
203 }
204 
205 // =====================================================================
206 // ---------------------------------------------------------------------
207 // ---------------------------------------------------------------------
208 // =====================================================================
209 
211 {
212  if(m_verbose>=2) Cout("vAlgorithm::Iterate ..."<< endl);
213 
214  #ifdef CASTOR_GPU
215  if (m_flagGPU)
216  return RunGPU();
217  else
218  return RunCPU();
219  #else
220  return RunCPU();
221  #endif
222 }
223 
224 // =====================================================================
225 // ---------------------------------------------------------------------
226 // ---------------------------------------------------------------------
227 // =====================================================================
228 
230 {
231 // Here is the scheme of all functions call related to each step of the iterative reconstruction:
232 // StepBeforeIterationLoop()
233 // / Loop on iterations
234 // | StepBeforeSubsetLoop(iteration)
235 // | / Loop on subsets
236 // | | StepPreProcessInsideSubsetLoop(iteration,subset)
237 // | | / Loop on bed positions
238 // | | | StepInnerLoopInsideSubsetLoop(iteration,subset,bed)
239 // | | StepPostProcessInsideSubsetLoop(iteration,subset)
240 // | StepAfterSubsetLoop(iteration)
241 // StepAfterIterationLoop()
242 
243  // Synchronize MPI processes
244  #ifdef CASTOR_MPI
245  MPI_Barrier(MPI_COMM_WORLD);
246  #endif
247 
248  // Verbose
249  if (m_verbose>=1) Cout("vAlgorithm::IterateCPU() -> Start algorithm for " << m_nbIterations << " iterations" << endl);
250 
251  // Initial clock for execution time
252  clock_t clock_start_whole = clock();
253  time_t time_start_whole = time(NULL);
254 
255  // Call before iteration function
256  // Call before iteration function
258  {
259  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepBeforeIterationLoop() function !" << endl);
260  return 1;
261  }
262 
263  // Loop on iterations
264  for (int iteration=0 ; iteration<m_nbIterations ; iteration++)
265  {
266 
267  // Call before subset function
268  if (StepBeforeSubsetLoop(iteration))
269  {
270  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepBeforeSubsetLoop() function !" << endl);
271  return 1;
272  }
273 
274  // Loop on subsets
275  for (int subset=0 ; subset<mp_nbSubsets[iteration] ; subset++)
276  {
277  // Verbose
278  if (m_verbose>=1)
279  {
280  Cout("vAlgorithm::IterateCPU() -> Start iteration " << iteration+1 << "/" << m_nbIterations
281  << " subset " << subset+1 << "/" << mp_nbSubsets[iteration] << endl);
282  }
283 
284  // Clock start for subset execution time
285  clock_t clock_start_subset = clock();
286  time_t time_start_subset = time(NULL);
287 
288  // Call pre-process inside subset loop function
289  if (StepPreProcessInsideSubsetLoop(iteration,subset))
290  {
291  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepPreProcessInsideSubsetLoop() function !" << endl);
292  return 1;
293  }
294 
295  // Synchronize MPI processes
296  #ifdef CASTOR_MPI
297  MPI_Barrier(MPI_COMM_WORLD);
298  #endif
299 
300  // Loop on bed positions
301  for (int bed=0 ; bed<m_nbBeds ; bed++)
302  {
303  // Call the inner loop on events inside subset loop function
304  if (StepInnerLoopInsideSubsetLoop(iteration,subset,bed))
305  {
306  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepInnerLoopInsideSubsetLoop() function !" << endl);
307  return 1;
308  }
309  } // End of beds loop
310 
311  // Synchronize MPI processes
312  #ifdef CASTOR_MPI
313  MPI_Barrier(MPI_COMM_WORLD);
314  #endif
315 
316  // Call post-process inside subset loop function
317  if (StepPostProcessInsideSubsetLoop(iteration,subset))
318  {
319  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepPostProcessInsideSubsetLoop() function !" << endl);
320  return 1;
321  }
322 
323  // Clock stop for subset execution time
324  clock_t clock_stop_subset = clock();
325  time_t time_stop_subset = time(NULL);
326  if (m_verbose>=2) Cout ("vAlgorithm::IterateCPU() -> Time spent for subset " << subset+1 << " | User: " << time_stop_subset-time_start_subset
327  << " sec | CPU: " << (clock_stop_subset-clock_start_subset)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
328 
329 
330  } // End of subsets loop
331 
332 
333  // Call after subset function
334  if (StepAfterSubsetLoop(iteration))
335  {
336  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepAfterSubsetLoop() function !" << endl);
337  return 1;
338  }
339 
340 
341  } // End of iterations loop
342 
343  // Call after iteration function
345  {
346  Cerr("***** vAlgorithm::IterateCPU() -> A problem occurred while calling StepAfterIterationLoop() function !" << endl);
347  return 1;
348  }
349 
350  // Final clock for execution time
351  clock_t clock_stop_whole = clock();
352  time_t time_stop_whole = time(NULL);
353  if (m_verbose>=1) Cout("vAlgorithm::IterateCPU() -> Total time spent | User: " << time_stop_whole-time_start_whole
354  << " sec | CPU: " << (clock_stop_whole-clock_start_whole)/((FLTNB)CLOCKS_PER_SEC) << " sec" << endl);
355 
356 
357  return 0;
358 }
359 
360 // =====================================================================
361 // ---------------------------------------------------------------------
362 // ---------------------------------------------------------------------
363 // =====================================================================
364 
366 {
367  if (m_verbose>=2) Cout("vAlgorithm::StepBeforeIterationLoop ... " << endl);
368 
369  // Check that the image space is already allocated
370  if (mp_ImageSpace==NULL || !mp_ImageSpace->Checked())
371  {
372  Cerr("***** vAlgorithm::StepBeforeIterationLoop() -> Image space has not been created or was not checked !" << endl);
373  return 1;
374  }
375 
376  // Instantiate all the required images from the oImageSpace
382 
384  {
385  Cerr("***** vAlgorithm::StepBeforeIterationLoop()-> Error during attenuation image initialization !" << endl);
386  return 1;
387  }
388 
390  {
391  Cerr("***** vAlgorithm::StepBeforeIterationLoop() -> An error occurred while initializing the sensitivity image !" << endl);
392  return 1;
393  }
394 
396  {
397  Cerr("***** vAlgorithm::StepBeforeIterationLoop()-> Error during multimodal image initialization !" << endl);
398  return 1;
399  }
400 
401  // End
402  return 0;
403 }
404 
405 // =====================================================================
406 // ---------------------------------------------------------------------
407 // ---------------------------------------------------------------------
408 // =====================================================================
409 
411 {
412  // End
413  return 0;
414 }
415 
416 // =====================================================================
417 // ---------------------------------------------------------------------
418 // ---------------------------------------------------------------------
419 // =====================================================================
420 
421 int vAlgorithm::StepPreProcessInsideSubsetLoop(int a_iteration, int a_subset)
422 {
423  // End
424  return 0;
425 }
426 
427 // =====================================================================
428 // ---------------------------------------------------------------------
429 // ---------------------------------------------------------------------
430 // =====================================================================
431 
432 int vAlgorithm::StepPostProcessInsideSubsetLoop(int a_iteration, int a_subset)
433 {
434  // End
435  return 0;
436 }
437 
438 // =====================================================================
439 // ---------------------------------------------------------------------
440 // ---------------------------------------------------------------------
441 // =====================================================================
442 
444 {
445  // End
446  return 0;
447 }
448 
449 // =====================================================================
450 // ---------------------------------------------------------------------
451 // ---------------------------------------------------------------------
452 // =====================================================================
453 
455 {
456  if (m_verbose>=2) Cout("vAlgorithm::StepAfterIterationLoop ... " << endl);
457 
458  // Deallocate everything
465 
466  // Normal end
467  return 0;
468 }
469 
470 // =====================================================================
471 // ---------------------------------------------------------------------
472 // ---------------------------------------------------------------------
473 // =====================================================================
474 
475 int vAlgorithm::StepImageOutput(int a_iteration, int a_subset)
476 {
477  // Normal end
478  return 0;
479 }
480 
481 // =====================================================================
482 // ---------------------------------------------------------------------
483 // ---------------------------------------------------------------------
484 // =====================================================================
485 
486 int vAlgorithm::InitSpecificOptions(string a_specificOptions)
487 {
488  (void)a_specificOptions; // avoid 'unused parameter' compil. warnings
489  // Normal end
490  return 0;
491 }
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.
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.