CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
src/optimizer/vOptimizer.cc
Go to the documentation of this file.
1 
8 #include "vOptimizer.hh"
9 
10 // =====================================================================
11 // ---------------------------------------------------------------------
12 // ---------------------------------------------------------------------
13 // =====================================================================
14 
16 {
17  // Affect default values
19  mp_ImageSpace = NULL;
21  m_nbTOFBins = 0;
22  m_initialValue = 0.;
27  m2p_forwardValues = NULL;
28  m3p_backwardValues = NULL;
31  m_optimizerFOMFlag = false;
33  m4p_FOMLogLikelihood = NULL;
34  m4p_FOMNbBins = NULL;
35  m4p_FOMRMSE = NULL;
36  m4p_FOMNbData = NULL;
37  m4p_FOMPenalty = NULL;
38  m_verbose = 0;
39  m2p_attenuationImage = NULL;
40  mp_imageStatNbVox = NULL;
41  mp_imageStatMin = NULL;
42  mp_imageStatMax = NULL;
43  mp_imageStatMean = NULL;
44  mp_imageStatVariance = NULL;
45  mp_correctionStatMean = NULL;
47  m_nbIterations = -1;
48  m_nbSubSteps = 1; // We default this to 1 to avoid having to define it for optimizers that do not need it
52  mp_nbSubsets = NULL;
53  mp_outputIterations = NULL;
54  m_currentIteration = -1;
56  m_currentSubStep = -1;
58  m_currentSubset = -1;
59  m_needPreIteration = false;
60  m_isInPreIteration = false;
61  m_needPostIteration = false;
62  m_isInPostIteration = false;
63  mp_Penalty = NULL;
66 }
67 
68 // =====================================================================
69 // ---------------------------------------------------------------------
70 // ---------------------------------------------------------------------
71 // =====================================================================
72 
74 {
75  // Free forward and backward buffers
77  {
79  if (m2p_forwardValues[th]) free(m2p_forwardValues[th]);
80  free(m2p_forwardValues);
81  }
83  {
84  for (int th=0; th<mp_ImageDimensionsAndQuantification->GetNbThreadsMax(); th++)
85  {
86  if (m3p_backwardValues[th])
87  {
88  for (int tof=0; tof<m_nbTOFBins; tof++) if (m3p_backwardValues[th][tof]) free(m3p_backwardValues[th][tof]);
89  free(m3p_backwardValues[th]);
90  }
91  }
92  free(m3p_backwardValues);
93  }
94  // Free tables of sub iterations
96  // Free pointers to attenuation images for SPECT
98  // Delete all image update statistics tables
100  {
102  if (mp_imageStatMin) free(mp_imageStatMin);
103  if (mp_imageStatMax) free(mp_imageStatMax);
108  }
109  // Delete FOM tables
110  if (m_optimizerFOMFlag)
111  {
112  // Delete the log-likelihood table
114  {
115  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
116  {
117  if (m4p_FOMLogLikelihood[fr])
118  {
119  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
120  {
121  if (m4p_FOMLogLikelihood[fr][rg])
122  {
123  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
124  if (m4p_FOMLogLikelihood[fr][rg][cg]) free(m4p_FOMLogLikelihood[fr][rg][cg]);
125  free(m4p_FOMLogLikelihood[fr][rg]);
126  }
127  }
128  free(m4p_FOMLogLikelihood[fr]);
129  }
130  }
131  free(m4p_FOMLogLikelihood);
132  }
133  // Delete the RMSE table
134  if (m4p_FOMRMSE)
135  {
136  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
137  {
138  if (m4p_FOMRMSE[fr])
139  {
140  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
141  {
142  if (m4p_FOMRMSE[fr][rg])
143  {
144  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
145  if (m4p_FOMRMSE[fr][rg][cg]) free(m4p_FOMRMSE[fr][rg][cg]);
146  free(m4p_FOMRMSE[fr][rg]);
147  }
148  }
149  free(m4p_FOMRMSE[fr]);
150  }
151  }
152  free(m4p_FOMRMSE);
153  }
154  // Delete the nbBins table
155  if (m4p_FOMNbBins)
156  {
157  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
158  {
159  if (m4p_FOMNbBins[fr])
160  {
161  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
162  {
163  if (m4p_FOMNbBins[fr][rg])
164  {
165  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
166  if (m4p_FOMNbBins[fr][rg][cg]) free(m4p_FOMNbBins[fr][rg][cg]);
167  free(m4p_FOMNbBins[fr][rg]);
168  }
169  }
170  free(m4p_FOMNbBins[fr]);
171  }
172  }
173  free(m4p_FOMNbBins);
174  }
175  // Delete the nbData table
176  if (m4p_FOMNbData)
177  {
178  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
179  {
180  if (m4p_FOMNbData[fr])
181  {
182  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
183  {
184  if (m4p_FOMNbData[fr][rg])
185  {
186  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
187  if (m4p_FOMNbData[fr][rg][cg]) free(m4p_FOMNbData[fr][rg][cg]);
188  free(m4p_FOMNbData[fr][rg]);
189  }
190  }
191  free(m4p_FOMNbData[fr]);
192  }
193  }
194  free(m4p_FOMNbData);
195  }
196  // Delete the penalty table
197  if (m4p_FOMPenalty)
198  {
199  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
200  {
201  if (m4p_FOMPenalty[fr])
202  {
203  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
204  {
205  if (m4p_FOMPenalty[fr][rg])
206  {
207  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
208  if (m4p_FOMPenalty[fr][rg][cg]) free(m4p_FOMPenalty[fr][rg][cg]);
209  free(m4p_FOMPenalty[fr][rg]);
210  }
211  }
212  free(m4p_FOMPenalty[fr]);
213  }
214  }
215  free(m4p_FOMPenalty);
216  }
217  }
218 }
219 
220 // =====================================================================
221 // ---------------------------------------------------------------------
222 // ---------------------------------------------------------------------
223 // =====================================================================
224 
226 {
227  // Call the specific help function from the children
229  // Say if the child optimizer is compatible with histogram and/or list-mode data
230  if (m_listmodeCompatibility && m_histogramCompatibility) cout << "This optimizer is compatible with both histogram and list-mode data." << endl;
231  else if (m_listmodeCompatibility) cout << "This optimizer is only compatible with list-mode data." << endl;
232  else if (m_histogramCompatibility) cout << "This optimizer is only compatible with histogram data." << endl;
233  else cout << "This optimizer is a total non-sense, not compatible with histogram nor list-mode data !!!" << endl;
234  // Say if the child optimizer is compatible with emission and/or transmission data
235  if (m_emissionCompatibility && m_transmissionCompatibility) cout << "This optimizer is compatible with both emission and transmission data." << endl;
236  else if (m_emissionCompatibility) cout << "This optimizer is only compatible with emission data." << endl;
237  else if (m_transmissionCompatibility) cout << "This optimizer is only compatible with transmission data." << endl;
238  else cout << "This optimizer is a total non-sense, not compatible with emission nor transmission data !!!" << endl;
239 }
240 
241 // =====================================================================
242 // ---------------------------------------------------------------------
243 // ---------------------------------------------------------------------
244 // =====================================================================
245 
247 {
248  // Check image dimensions
250  {
251  Cerr("***** vOptimizer::CheckParameters() -> oImageDimensionsAndQuantification is null !" << endl);
252  return 1;
253  }
254  // Check image space
255  if (mp_ImageSpace==NULL)
256  {
257  Cerr("***** vOptimizer::CheckParameters() -> oImageSpace is null !" << endl);
258  return 1;
259  }
260  // Check data mode
262  {
263  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data mode provided !" << endl);
264  return 1;
265  }
266  // Check data type
268  {
269  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data type provided !" << endl);
270  return 1;
271  }
272  // Check data spec
274  {
275  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data specificity provided (emission or transmission) !" << endl);
276  return 1;
277  }
278  // Check verbose level
279  if (m_verbose<0)
280  {
281  Cerr("***** vOptimizer::CheckParameters() -> Verbose level is negative !" << endl);
282  return 1;
283  }
284  // Listmode incompatibility
286  {
287  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with listmode data !" << endl);
288  return 1;
289  }
290  // Histogram incompatibility
292  {
293  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with histogram data !" << endl);
294  return 1;
295  }
296  // Emission compatability
298  {
299  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with emission data !" << endl);
300  return 1;
301  }
302  // Transmission compatibility
304  {
305  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with transmission data !" << endl);
306  return 1;
307  }
308  // Send an error if FOM computation is required and verbose is null
309  if (m_optimizerFOMFlag && m_verbose==0)
310  {
311  Cerr("***** vOptimizer::CheckParameters() -> Asked to compute FOMs in the data-space whereas the verbose is not set !" << endl);
312  return 1;
313  }
314  // Send an error if image statistics computation is required and verbose is null
316  {
317  Cerr("***** vOptimizer::CheckParameters() -> Asked to compute image update statistics whereas the verbose is not set !" << endl);
318  return 1;
319  }
320  // Send an error if FOM computation with list-mode data (it has no sense)
322  {
323  Cerr("***** vOptimizer::CheckParameters() -> Computing FOMs in the data-space with list-mode data is not possible !" << endl);
324  return 1;
325  }
326  // Send an error if number of sub-steps is lower than 1
327  if (m_nbSubSteps<1)
328  {
329  Cerr("***** vOptimizer::CheckParameters() -> Number of sub-steps cannot be lower than 1 !" << endl);
330  return 1;
331  }
332  // Call the CheckSpecificParameters function of the child
334  {
335  Cerr("***** vOptimizer::CheckParameters() -> A problem occurred while checking parameters specific to the optimizer module !" << endl);
336  return 1;
337  }
338  // Normal end
339  return 0;
340 }
341 
342 // =====================================================================
343 // ---------------------------------------------------------------------
344 // ---------------------------------------------------------------------
345 // =====================================================================
346 
348 {
349  // Allocate forward buffers (as many as threads and as many as TOF bins)
352  m2p_forwardValues[th] = (FLTNB*)malloc(m_nbTOFBins*sizeof(FLTNB));
353 
354  // Allocate backward buffers (as many as threads then as many as TOF bins and as many as backward images)
356  for (int th=0; th<mp_ImageDimensionsAndQuantification->GetNbThreadsMax(); th++)
357  {
358  m3p_backwardValues[th] = (FLTNB**)malloc(m_nbTOFBins*sizeof(FLTNB*));
359  for (int tof=0; tof<m_nbTOFBins; tof++)
360  m3p_backwardValues[th][tof] = (FLTNB*)malloc(m_nbBackwardImages*sizeof(FLTNB));
361  }
362 
363  // Allocate pointers to the current attenuation image for SPECT attenuation correction
366 
367  // Allocate the thread-safe tables for image update statistics computation
369  {
377  }
378 
379  // We allocate the tables that hold the figures-of-merit if FOM flag is on
380  if (m_optimizerFOMFlag)
381  {
382  m4p_FOMNbBins = (uint64_t****)malloc(mp_ImageDimensionsAndQuantification->GetNbTimeFrames()*sizeof(uint64_t***));
387  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
388  {
389  m4p_FOMNbBins[fr] = (uint64_t***)malloc(mp_ImageDimensionsAndQuantification->GetNbRespGates()*sizeof(uint64_t**));
394  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
395  {
396  m4p_FOMNbBins[fr][rg] = (uint64_t**)malloc(mp_ImageDimensionsAndQuantification->GetNbCardGates()*sizeof(uint64_t*));
401  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
402  {
403  m4p_FOMNbBins[fr][rg][cg] = (uint64_t*)malloc(mp_ImageDimensionsAndQuantification->GetNbThreadsForProjection()*sizeof(uint64_t));
408  }
409  }
410  }
411  }
412 
413  // Allocate the tables containing the number of sub-iterations in each sub-step and default it to 1 (because it is not mandatory for specific optimizers to allow for changing these values)
414  mp_nbSubIterationsInSubSteps = (int*)malloc(m_nbSubSteps*sizeof(int));
415  for (int step=0; step<m_nbSubSteps; step++) mp_nbSubIterationsInSubSteps[step] = 1;
416 
417  // Call the specific initialization function of the child
418  if (InitializeSpecific())
419  {
420  Cerr("***** vOptimizer::Initialize() -> A problem occurred while initializing stuff specific to the optimizer module !" << endl);
421  return 1;
422  }
423 
424  // Compute the total number of sub-iterations in one iteration (i.e. combining all sub-steps)
426  for (int step=0; step<m_nbSubSteps; step++) m_nbSubIterationsInOneIteration += mp_nbSubIterationsInSubSteps[step];
427 
428  // Set the number of iterations
430 
431  // Normal end
432  return 0;
433 }
434 
435 // =====================================================================
436 // ---------------------------------------------------------------------
437 // ---------------------------------------------------------------------
438 // =====================================================================
439 
440 void vOptimizer::SetCurrentIteration(int a_currentTotalSubIteration)
441 {
442  // First set the current global sub iteration number
443  m_currentTotalSubIteration = a_currentTotalSubIteration;
444  // Compute the current iteration number
446  // Get the remaining number of sub iterations in the current iteration
448  // Compute the current sub step and current sub iteration
450  {
453  }
454 }
455 
456 // =====================================================================
457 // ---------------------------------------------------------------------
458 // ---------------------------------------------------------------------
459 // =====================================================================
460 
462 {
463  // Simply reset FOMs if activated
464  if (m_optimizerFOMFlag)
465  {
466  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
467  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
468  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
469  {
471  {
472  m4p_FOMNbBins[fr][rg][cg][th] = 0;
473  m4p_FOMNbData[fr][rg][cg][th] = 0.;
474  m4p_FOMLogLikelihood[fr][rg][cg][th] = 0.;
475  m4p_FOMRMSE[fr][rg][cg][th] = 0.;
476  }
478  {
479  m4p_FOMPenalty[fr][rg][cg][th] = 0.;
480  }
481  }
482  }
483  // Then call the specific pre-update-step function from the child optimizer
485  {
486  Cerr("***** vOptimizer::PreDataUpdateStep() -> A problem occurred while applying the specific pre-data-update step of the optimizer module !" << endl);
487  return 1;
488  }
489  // Normal end
490  return 0;
491 }
492 
493 // =====================================================================
494 // ---------------------------------------------------------------------
495 // ---------------------------------------------------------------------
496 // =====================================================================
497 
499 {
500  // By default, just do nothing here. This function is designed to be overloaded by the specific optimizer if needed
501  return 0;
502 }
503 
504 // =====================================================================
505 // ---------------------------------------------------------------------
506 // ---------------------------------------------------------------------
507 // =====================================================================
508 
510 {
511  // Then call the specific post-update-step function from the child optimizer
513  {
514  Cerr("***** vOptimizer::PreImageUpdateStep() -> A problem occurred while applying the specific post-data-update step of the optimizer module !" << endl);
515  return 1;
516  }
517  // Compute optimizer figures-of-merit (data fidelity term and penalty)
518  if (m_optimizerFOMFlag)
519  {
520  // Compute data FOM related to the data fildelity term
521  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
522  {
523  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
524  {
525  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
526  {
527  // Merge thread results
529  {
530  m4p_FOMNbBins[fr][rg][cg][0] += m4p_FOMNbBins[fr][rg][cg][th];
531  m4p_FOMNbData[fr][rg][cg][0] += m4p_FOMNbData[fr][rg][cg][th];
532  m4p_FOMLogLikelihood[fr][rg][cg][0] += m4p_FOMLogLikelihood[fr][rg][cg][th];
533  m4p_FOMRMSE[fr][rg][cg][0] += m4p_FOMRMSE[fr][rg][cg][th];
534  }
535  // Finish mean number of data counts per channel
536  m4p_FOMNbData[fr][rg][cg][0] /= ((HPFLTNB)(m4p_FOMNbBins[fr][rg][cg][0]));
537  // Finish RMSE computation
538  m4p_FOMRMSE[fr][rg][cg][0] /= ((HPFLTNB)(m4p_FOMNbBins[fr][rg][cg][0]));
539  m4p_FOMRMSE[fr][rg][cg][0] = sqrt(m4p_FOMRMSE[fr][rg][cg][0]);
540  }
541  }
542  }
543  // Compute penalty FOM
544  if (mp_Penalty!=NULL)
545  {
546  // As the penalty will be innocently applied to the current image, it is performed on time basis functions
547  // and not frames/gates. As a result, we loop over time basis functions first, compute the penalty, and
548  // distribute its value with respect to the contribution of each time basis function to the frames/gates.
549  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
550  {
551  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
552  {
553  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
554  {
555  // In case a problem occurs during the multi-threaded loop
556  bool problem = false;
557  // The voxel index
558  INTNB v;
559  // Multi-threading over voxels
560  #pragma omp parallel for private(v) schedule(guided)
562  {
563  // Get the thread index
564  int th = 0;
565  #ifdef CASTOR_OMP
566  th = omp_get_thread_num();
567  #endif
568  // Compute the penalty term
569  if (mp_Penalty->LocalPreProcessingStep(tbf,rbf,cbf,v,th))
570  {
571  Cerr("***** vOptimizer::PreImageUpdateStep() -> A problem occurred while precomputing the local step of the penalty for voxel " << v << " !" << endl);
572  problem = true;
573  }
574  HPFLTNB penalty_value = mp_Penalty->ComputePenaltyValue(mp_ImageSpace->m4p_image[tbf][cbf][rbf],v,th);
575  // Distribuete the value to frames/gates with respect to basis functions
576  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
577  {
579  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
580  {
582  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
583  {
585  m4p_FOMPenalty[fr][rg][cg][th] += time_basis_coef * resp_basis_coef * card_basis_coef * penalty_value;
586  }
587  }
588  }
589  }
590  // Potential problems during the loop
591  if (problem)
592  {
593  Cerr("***** vOptimizer::PreImageUpdateStep() -> A problem occurred in the multi-threaded loop for the penalty, stop now !" << endl);
594  return 1;
595  }
596  }
597  }
598  }
599  // Merge thread results
600  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
601  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
602  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
603  {
605  m4p_FOMPenalty[fr][rg][cg][0] += m4p_FOMPenalty[fr][rg][cg][th];
606  // Also divide by the number of subsets for balance with respect to the data fidelity term
608  }
609  }
610  // Verbose
611  Cout("vOptimizer::PreImageUpdateStep() -> Optimizer figures-of-merit related to the objective function" << endl);
612  // Number of spaces for printing
613  string spaces_fr = ""; if (mp_ImageDimensionsAndQuantification->GetNbTimeFrames()>1) spaces_fr += " ";
614  string spaces_rg = ""; if (mp_ImageDimensionsAndQuantification->GetNbRespGates()>1) spaces_rg += " ";
615  string spaces_cg = ""; if (mp_ImageDimensionsAndQuantification->GetNbCardGates()>1) spaces_cg += " ";
616  // Loop on frames and gates
617  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
618  {
619  // Verbose
621  Cout(spaces_fr << "Frame " << fr+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbTimeFrames() << endl);
622  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
623  {
624  // Verbose
626  Cout(spaces_fr << spaces_rg << "Respiratory gate " << rg+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbRespGates() << endl);
627  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
628  {
629  // Verbose
631  Cout(spaces_fr << spaces_rg << spaces_cg << "Cardiac gate " << cg+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbCardGates() << endl);
632  // Print out
633  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Number of data channels used in optimization: " << m4p_FOMNbBins[fr][rg][cg][0] << endl);
634  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Mean number of data counts per channel: " << m4p_FOMNbData[fr][rg][cg][0] << endl);
635  if (mp_Penalty!=NULL)
636  {
637  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Penalty: " << m4p_FOMPenalty[fr][rg][cg][0] << endl);
638  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Log-likelihood: " << m4p_FOMLogLikelihood[fr][rg][cg][0]
639  << " (including penalty: " << m4p_FOMLogLikelihood[fr][rg][cg][0] - m4p_FOMPenalty[fr][rg][cg][0] << ")" << endl);
640  Cout(spaces_fr << spaces_rg << spaces_cg << " --> RMSE: " << m4p_FOMRMSE[fr][rg][cg][0]
641  << " (including penalty: " << m4p_FOMRMSE[fr][rg][cg][0] + m4p_FOMPenalty[fr][rg][cg][0] << ")" << endl);
642  }
643  else
644  {
645  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Log-likelihood: " << m4p_FOMLogLikelihood[fr][rg][cg][0] << endl);
646  Cout(spaces_fr << spaces_rg << spaces_cg << " --> RMSE: " << m4p_FOMRMSE[fr][rg][cg][0] << endl);
647  }
648  }
649  }
650  }
651  }
652  // Normal end
653  return 0;
654 }
655 
656 // =====================================================================
657 // ---------------------------------------------------------------------
658 // ---------------------------------------------------------------------
659 // =====================================================================
660 
662 {
663  // By default, just do nothing here. This function is designed to be overloaded by the specific optimizer if needed
664  return 0;
665 }
666 
667 // =====================================================================
668 // ---------------------------------------------------------------------
669 // ---------------------------------------------------------------------
670 // =====================================================================
671 
673  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
674  int a_th )
675 {
676  // =======================================================================
677  // 4D forward projection including framing, respiratory and cardiac gating
678  // =======================================================================
679 
680  // Clean buffers
681  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
682  {
683  m2p_forwardValues[a_th][b] = 0.;
684  for (int img=0; img<m_nbBackwardImages; img++) m3p_backwardValues[a_th][b][img] = 0.;
685  }
686 
687  // Loop over time basis functions
688  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
689  {
690  // Get frame/basis coefficient and continue if null
691  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(tbf,a_timeFrame);
692  if (time_basis_coef==0.) continue;
693  // Loop over respiratory basis functions
694  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
695  {
696  // Get resp_gate/basis coefficient and continue if null
697  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(rbf,a_respGate);
698  if (resp_basis_coef==0.) continue;
699  // Loop over cardiac basis functions
700  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
701  {
702  // Get card_gate_basis coefficient and continue if null
703  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(cbf,a_cardGate);
704  if (card_basis_coef==0.) continue;
705  // Compute global coefficient
706  FLTNB global_basis_coef = time_basis_coef * resp_basis_coef * card_basis_coef;
707  // Loop over TOF bins
708  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
709  {
710  // Set the current TOF bin
711  ap_Line->SetCurrentTOFBin(b);
712  // Project the image and apply dynamic coefficient converions
713  m2p_forwardValues[a_th][b] += ForwardProject(ap_Line,mp_ImageSpace->m4p_forwardImage[tbf][rbf][cbf]) * global_basis_coef;
714  }
715  }
716  }
717  }
718 
719  // ==============================================
720  // Differentiate emission and transmission models
721  // ==============================================
722 
723  // -----------------
724  // Case for emission
725  // -----------------
727  {
728  // Add additive terms (we multiply by the frame duration because the additive terms are provided as rates)
730  for (int b=0; b<ap_Line->GetNbTOFBins(); b++) m2p_forwardValues[a_th][b] += ap_Event->GetAdditiveCorrections(b) * mp_ImageDimensionsAndQuantification->GetdurationPerGate(a_timeFrame,a_respGate);
731  else
732  for (int b=0; b<ap_Line->GetNbTOFBins(); b++) m2p_forwardValues[a_th][b] += ap_Event->GetAdditiveCorrections(b) * mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame);
733  }
734  // ---------------------
735  // Case for transmission
736  // ---------------------
737  else if (m_dataSpec==SPEC_TRANSMISSION)
738  {
739  int no_tof = 0;
740  // Ignore TOF as it has no sense; add scatter rate multiplied by the frame duration; take blank into account
741  m2p_forwardValues[a_th][no_tof] = ap_Event->GetBlankValue() * exp(-m2p_forwardValues[a_th][no_tof])
742  + ap_Event->GetAdditiveCorrections(no_tof) * mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame);
743  }
744 
745  // End
746  return 0;
747 }
748 
749 // =====================================================================
750 // ---------------------------------------------------------------------
751 // ---------------------------------------------------------------------
752 // =====================================================================
753 
754 int vOptimizer::DataStep2Optional( oProjectionLine* ap_Line, vEvent* ap_Event,
755  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
756  int a_th )
757 {
758  // Deliberately do nothing!
759  return 0;
760 }
761 
762 // =====================================================================
763 // ---------------------------------------------------------------------
764 // ---------------------------------------------------------------------
765 // =====================================================================
766 
768  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
769  int a_th )
770 {
771  // Loop over TOF bins
772  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
773  {
774  // Set the current TOF bin of the projection-line
775  ap_Line->SetCurrentTOFBin(b);
776  // The weight associated to the line to be backprojected
777  FLTNB weight = 0.;
778  // Call the function to compute specific sensitivity
780  (
781  ap_Event->GetEventValue(b), // the measured data
782  m2p_forwardValues[a_th][b], // the forward model
783  &weight, // the sensitivity weight
784  ap_Event->GetMultiplicativeCorrections(), // the multiplicative corrections
785  ap_Event->GetAdditiveCorrections(b)*mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame), // the additive corrections
786  ap_Event->GetBlankValue(), // the blank value
787  mp_ImageDimensionsAndQuantification->GetQuantificationFactor(a_bed,a_timeFrame, a_respGate, a_cardGate), // the quantification factor
788  ap_Line // the projection line
789  ))
790  {
791  Cerr("***** vOptimizer::DataStep3BackwardProjectSensitivity() -> A problem occurred while performing the sensitivity specific operations !" << endl);
792  return 1;
793  }
794  // Back-project the weight into the sensitivity matrix
795  BackwardProject(ap_Line, mp_ImageSpace->m5p_sensitivity[a_th][a_timeFrame][a_respGate][a_cardGate], weight);
796  }
797 
798  // End
799  return 0;
800 }
801 
802 // =====================================================================
803 // ---------------------------------------------------------------------
804 // ---------------------------------------------------------------------
805 // =====================================================================
806 
807 int vOptimizer::DataStep4Optional( oProjectionLine* ap_Line, vEvent* ap_Event,
808  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
809  int a_th )
810 {
811  // Deliberately do nothing!
812  return 0;
813 }
814 
815 // =====================================================================
816 // ---------------------------------------------------------------------
817 // ---------------------------------------------------------------------
818 // =====================================================================
819 
821  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
822  int a_th )
823 {
824  // Loop on TOF bins
825  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
826  {
827  // Set the current TOF bin of the projection-line
828  ap_Line->SetCurrentTOFBin(b);
829  // Then call the pure virtual function for specific data space operations
831  (
832  ap_Event->GetEventValue(b), // the measured data
833  m2p_forwardValues[a_th][b], // the forward model
834  m3p_backwardValues[a_th][b], // the backward values (the result)
835  ap_Event->GetMultiplicativeCorrections(), // the multiplicative corrections
836  ap_Event->GetAdditiveCorrections(b)*mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame), // the additive corrections
837  ap_Event->GetBlankValue(), // the blank value
838  mp_ImageDimensionsAndQuantification->GetQuantificationFactor(a_bed,a_timeFrame, a_respGate, a_cardGate), // the quantification factor
839  ap_Line // the projection line
840  ))
841  {
842  Cerr("***** vOptimizer::DataStep5ComputeCorrections() -> A problem occurred while performing specific data space operations !" << endl);
843  return 1;
844  }
845  }
846 
847  // End
848  return 0;
849 }
850 
851 // =====================================================================
852 // ---------------------------------------------------------------------
853 // ---------------------------------------------------------------------
854 // =====================================================================
855 
856 int vOptimizer::DataStep6Optional( oProjectionLine* ap_Line, vEvent* ap_Event,
857  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
858  int a_th )
859 {
860  // Deliberately do nothing!
861  return 0;
862 }
863 
864 // =====================================================================
865 // ---------------------------------------------------------------------
866 // ---------------------------------------------------------------------
867 // =====================================================================
868 
870  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
871  int a_th )
872 {
873  // ========================================================================
874  // 4D backward projection including framing, respiratory and cardiac gating
875  // ========================================================================
876 
877  // Loop over time basis functions
878  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
879  {
880  // Get frame/basis coefficient and continue if null
881  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(tbf,a_timeFrame);
882  if (time_basis_coef==0.) continue;
883  // Loop over respiratory basis functions
884  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
885  {
886  // Get resp_gate/basis coefficient and continue if null
887  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(rbf,a_respGate);
888  if (resp_basis_coef==0.) continue;
889  // Loop over cardiac basis functions
890  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
891  {
892  // Get card_gate_basis coefficient and continue if null
893  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(cbf,a_cardGate);
894  if (card_basis_coef==0.) continue;
895  // Compute global coefficient
896  FLTNB global_basis_coef = time_basis_coef * resp_basis_coef * card_basis_coef;
897  // Loop over TOF bins
898  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
899  {
900  // Set the current TOF bin of the projection-line
901  ap_Line->SetCurrentTOFBin(b);
902  // Loop over backprojection images
903  for (int img=0; img<m_nbBackwardImages; img++)
904  {
905  // Continue if backproj is null or infinity
906  if (m3p_backwardValues[a_th][b][img]==0. || !isfinite(m3p_backwardValues[a_th][b][img])) continue;
907  // Backward project into the backward image
908  BackwardProject( ap_Line, mp_ImageSpace->m6p_backwardImage[img][a_th][tbf][rbf][cbf] , m3p_backwardValues[a_th][b][img] * global_basis_coef );
909  }
910  }
911  }
912  }
913  }
914 
915  // End
916  return 0;
917 }
918 
919 // =====================================================================
920 // ---------------------------------------------------------------------
921 // ---------------------------------------------------------------------
922 // =====================================================================
923 
925  int a_timeFrame, int a_respGate, int a_cardGate,
926  int a_thread )
927 {
928  // Loop on TOF bins
929  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
930  {
931  // Get the model and the data
932  HPFLTNB model = m2p_forwardValues[a_thread][b];
933  HPFLTNB data = ap_Event->GetEventValue(b);
934  // Update log likelihood only if model is strictly positive (to avoid numerical issues, we use a threshold at e-10)
935  if (model>1.e-10) m4p_FOMLogLikelihood[a_timeFrame][a_respGate][a_cardGate][a_thread] += data*log(model)-model;
936  // Update RMSE
937  m4p_FOMRMSE[a_timeFrame][a_respGate][a_cardGate][a_thread] += (data-model)*(data-model);
938  // Update number of data channels
939  m4p_FOMNbBins[a_timeFrame][a_respGate][a_cardGate][a_thread]++;
940  // Update the number of data counts
941  m4p_FOMNbData[a_timeFrame][a_respGate][a_cardGate][a_thread] += data;
942  }
943  // Normal end
944  return 0;
945 }
946 
947 // =====================================================================
948 // ---------------------------------------------------------------------
949 // ---------------------------------------------------------------------
950 // =====================================================================
951 
953 {
954  // Verbose
955  if (m_verbose>=3) Cout("vOptimizer::UpdateVisitedVoxels() -> Tick visited voxels based on sensitivity" << endl);
956 
957  // The purpose of this mp_visitedVoxelsImage is to keep track within a complete iteration, of which voxels were
958  // visited by some LORs. Everything is based on sensitivity. With list-mode data, it means that the visited voxels
959  // are the same for all subsets. With histogram data, if the sensitivity is computed on the fly from the data, then
960  // the visited voxels may change for each voxel. When updating the voxel values in the ImageUpdateStep, only voxels
961  // with non-zero sensitivity are updated, to avoid artificially decreasing the signal in voxels that may not have
962  // been visited in this subset but which may be visited in another subset (this is true only for histogram data).
963  // But, at the end of each iteration, the mp_visitedVoxelsImage is used to checked for never-visited voxels, in which
964  // case they will be set to 0 to avoid that they keep their initial value. In other words, these never-visited voxels
965  // are somewhat excluded from the reconstruction.
966 
967  // For this matter of visited voxels does not depend on dynamic aspects, we loop on the sensitivity over all frames and only for the first
968  // gates and look if voxels were 'visited' by LORs (i.e. the sensitivity is not null);
969  // the information is saved in one static 3D image and then applied to all frames/gates;
970  int thread0 = 0;
971  int resp_gate0 = 0;
972  int card_gate0 = 0;
973  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
974  for (int v=0; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ(); v++)
975  // If the voxel has been visited, we add 1. to the mp_Image->mp_visitedVoxelsImage
976  if (mp_ImageSpace->m5p_sensitivity[thread0][fr][resp_gate0][card_gate0][v]!=0.) mp_ImageSpace->mp_visitedVoxelsImage[v] += 1.;
977 
978  // End
979  return 0;
980 }
981 
982 // =====================================================================
983 // ---------------------------------------------------------------------
984 // ---------------------------------------------------------------------
985 // =====================================================================
986 
988 {
989  // Verbose
990  if (m_verbose>=2 || m_optimizerImageStatFlag) Cout("vOptimizer::ImageUpdateStep() -> Start image update" << endl);
991  // Number of spaces for nice printing
992  string spaces_tbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions()>1) spaces_tbf += " ";
993  string spaces_rbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions()>1) spaces_rbf += " ";
994  string spaces_cbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions()>1) spaces_cbf += " ";
995  // Loop over time basis functions
996  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
997  {
998  // Verbose
1000  {
1001  if (mp_ImageDimensionsAndQuantification->GetTimeStaticFlag()) Cout(spaces_tbf << "--> Time frame " << tbf+1 << endl);
1002  else Cout(spaces_tbf << "--> Time basis function: " << tbf+1 << endl);
1003  }
1004  // Loop over respiratory basis functions
1005  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
1006  {
1007  // Verbose
1009  {
1010  if (mp_ImageDimensionsAndQuantification->GetRespStaticFlag()) Cout(spaces_tbf << spaces_rbf << "--> Respiratory gate " << rbf+1 << endl);
1011  else Cout(spaces_tbf << spaces_rbf << "--> Respiratory basis function: " << rbf+1 << endl);
1012  }
1013  // Loop over cardiac basis functions
1014  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
1015  {
1016  // Verbose
1018  {
1019  if (mp_ImageDimensionsAndQuantification->GetCardStaticFlag()) Cout(spaces_tbf << spaces_rbf << spaces_cbf << "--> Cardiac gate " << cbf+1 << endl);
1020  else Cout(spaces_tbf << spaces_rbf << spaces_cbf << "--> Cardiac basis function: " << cbf+1 << endl);
1021  }
1022  // Set the number of threads for image computation
1023  #ifdef CASTOR_OMP
1025  #endif
1026  // Reset counter for image stats
1028  {
1030  {
1031  mp_imageStatNbVox[th] = 0;
1032  mp_imageStatMin[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][0];
1033  mp_imageStatMax[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][0];
1034  mp_imageStatMean[th] = 0.;
1035  mp_imageStatVariance[th] = 0.;
1036  mp_correctionStatMean[th] = 0.;
1037  mp_correctionStatVariance[th] = 0.;
1038  }
1039  }
1040  // OMP loop over voxels
1041  INTNB v;
1042  bool error = false;
1043  #pragma omp parallel for private(v) schedule(guided)
1044  for (v=0 ; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ() ; v++)
1045  {
1046  // Get the thread index
1047  int th = 0;
1048  #ifdef CASTOR_OMP
1049  th = omp_get_thread_num();
1050  #endif
1051  // Compute sensitivity
1052  FLTNB sensitivity = ComputeSensitivity(mp_ImageSpace->m5p_sensitivity[0],tbf,rbf,cbf,v);
1053  // If the sensitivity is negative or null, we skip this voxel
1054  if (sensitivity<=0.) continue;
1055  // Keep the current image value in a buffer for update statistics
1056  FLTNB old_image_value = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
1057  // Fill class-local buffer of backward values (thread-safe)
1058  for (int nb=0; nb<m_nbBackwardImages; nb++)
1059  m3p_backwardValues[th][0][nb] = mp_ImageSpace->m6p_backwardImage[nb][0][tbf][rbf][cbf][v];
1060  // Then call the pure virtual function for specific data space operations
1062  (
1063  mp_ImageSpace->m4p_image[tbf][rbf][cbf][v],
1064  // Note: have to be carefull here because we pass rbf and cbf as the actual respgate and cardgate. It may be wrong if one would like
1065  // to perform MLAA reconstruction but using PET basis functions different from diracs for example... In this case, we would get
1066  // segmentation faults. But actually, this function would be overloaded by MLAA anyway...
1067  //mp_Image->m4p_attenuation != NULL ? mp_Image->m4p_attenuation[tbf][rbf][cbf][v] : 0. ,
1068  &mp_ImageSpace->m4p_image[tbf][rbf][cbf][v],
1069  //mp_Image->m4p_attenuation != NULL ? &mp_Image->m4p_attenuation[tbf][rbf][cbf][v] : NULL ,
1070  sensitivity,
1071  m3p_backwardValues[th][0],
1072  v,
1073  tbf,
1074  rbf,
1075  cbf
1076  ))
1077  {
1078  Cerr("***** vOptimizer::ImageUpdateStep() -> A problem occurred while performing image space specific operations of the optimizer !" << endl);
1079  error = true;
1080  }
1081  // Do some image statistics
1083  {
1084  if (mp_imageStatMin[th]>mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]) mp_imageStatMin[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
1085  if (mp_imageStatMax[th]<mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]) mp_imageStatMax[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
1086  // The one pass algorithm is used to compute the variance here for both the image and correction step.
1087  // Do not change anyhting to the order of the operations !
1088  mp_imageStatNbVox[th]++;
1089  HPFLTNB sample = ((HPFLTNB)(mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]));
1090  HPFLTNB delta = sample - mp_imageStatMean[th];
1091  mp_imageStatMean[th] += delta / ((HPFLTNB)(mp_imageStatNbVox[th]));
1092  mp_imageStatVariance[th] += delta * (sample - mp_imageStatMean[th]);
1093  // Same one-pass variance algorithm for the correction step
1094  sample = ((HPFLTNB)(mp_ImageSpace->m4p_image[tbf][rbf][cbf][v] - old_image_value));
1095  delta = sample - mp_correctionStatMean[th];
1096  mp_correctionStatMean[th] += delta / ((HPFLTNB)(mp_imageStatNbVox[th]));
1097  mp_correctionStatVariance[th] += delta * (sample - mp_correctionStatMean[th]);
1098  }
1099  }
1100  // Check for errors (we cannot stop the multi-threaded loop so we do it here)
1101  if (error) return 1;
1102  // Finish image update statistics computations: these operations are specific to the merging of the
1103  // multi-threaded one-pass variance estimations used above.
1105  {
1107  {
1108  // Image minimum and maximum
1111  // Cast the counts into HPFLTNB
1112  HPFLTNB count1 = ((HPFLTNB)(mp_imageStatNbVox[0]));
1113  HPFLTNB count2 = ((HPFLTNB)(mp_imageStatNbVox[th]));
1114  // Update the count
1116  HPFLTNB count12 = ((HPFLTNB)(mp_imageStatNbVox[0]));
1117  // Compute the delta mean
1118  HPFLTNB delta = mp_imageStatMean[0] - mp_imageStatMean[th];
1119  // Update the variance
1120  mp_imageStatVariance[0] += mp_imageStatVariance[th] + delta * delta * count1 * count2 / count12;
1121  // Update the mean
1122  mp_imageStatMean[0] = (count1*mp_imageStatMean[0] + count2*mp_imageStatMean[th]) / count12;
1123  // Do the same for the correction step
1125  mp_correctionStatVariance[0] += mp_correctionStatVariance[th] + delta * delta * count1 * count2 / count12;
1126  mp_correctionStatMean[0] = (count1*mp_correctionStatMean[0] + count2*mp_correctionStatMean[th]) / count12;
1127  }
1128  // Finish variance computation
1131  // Print out
1132  Cout(spaces_tbf << spaces_rbf << spaces_cbf << " --> Image update step "
1133  << " | mean: " << mp_correctionStatMean[0]
1134  << " | stdv: " << sqrt(mp_correctionStatVariance[0]) << endl);
1135  Cout(spaces_tbf << spaces_rbf << spaces_cbf << " --> New image estimate"
1136  << " | mean: " << mp_imageStatMean[0]
1137  << " | stdv: " << sqrt(mp_imageStatVariance[0])
1138  << " | min/max: " << mp_imageStatMin[0] << "/" << mp_imageStatMax[0] << endl);
1139  }
1140  }
1141  }
1142  }
1143 
1144  // End
1145  return 0;
1146 }
1147 
1148 // =====================================================================
1149 // ---------------------------------------------------------------------
1150 // ---------------------------------------------------------------------
1151 // =====================================================================
1152 
1154 {
1155  // Particular case for SPECT with attenuation correction
1156  if (m_dataType==TYPE_SPECT && m2p_attenuationImage[ap_Line->GetThreadNumber()]!=NULL)
1157  return ap_Line->ForwardProjectWithSPECTAttenuation( m2p_attenuationImage[ap_Line->GetThreadNumber()], ap_image );
1158  // Otherwise call the function without attenuation
1159  else
1160  return ap_Line->ForwardProject( ap_image );
1161 }
1162 
1163 // =====================================================================
1164 // ---------------------------------------------------------------------
1165 // ---------------------------------------------------------------------
1166 // =====================================================================
1167 
1168 void vOptimizer::BackwardProject(oProjectionLine* ap_Line, FLTNB* ap_image, FLTNB a_value)
1169 {
1170  // Particular case for SPECT with attenuation correction
1171  if (m_dataType==TYPE_SPECT && m2p_attenuationImage[ap_Line->GetThreadNumber()]!=NULL)
1172  ap_Line->BackwardProjectWithSPECTAttenuation( m2p_attenuationImage[ap_Line->GetThreadNumber()], ap_image, a_value );
1173  // Otherwise call the function without attenuation
1174  else
1175  ap_Line->BackwardProject( ap_image, a_value );
1176 }
1177 
1178 // =====================================================================
1179 // ---------------------------------------------------------------------
1180 // ---------------------------------------------------------------------
1181 // =====================================================================
1182 
1183 FLTNB vOptimizer::ComputeSensitivity( FLTNB**** a4p_sensitivityMatrix,
1184  int a_timeBasisFunction, int a_respBasisFunction, int a_cardBasisFunction,
1185  int a_voxel )
1186 {
1187  // The sensitivity to be computed
1188  FLTNB sensitivity = 0.;
1189  // Loop over time frames
1190  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
1191  {
1192  // Get frame/basis coefficient and continue if null
1193  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(a_timeBasisFunction,fr);
1194  if (time_basis_coef==0.) continue;
1195  // Loop over respiratory gates
1196  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
1197  {
1198  // Get resp_gate/basis coefficient and continue if null
1199  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(a_respBasisFunction,rg);
1200  if (resp_basis_coef==0.) continue;
1201  // Loop over cardiac gates
1202  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
1203  {
1204  // Get card_gate_basis coefficient and continue if null
1205  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(a_cardBasisFunction,cg);
1206  if (card_basis_coef==0.) continue;
1207  // Add actual weight contribution
1208  sensitivity += a4p_sensitivityMatrix[fr][rg][cg][a_voxel] *
1209  time_basis_coef * resp_basis_coef * card_basis_coef;
1210  }
1211  }
1212  }
1213  // If listmode, then divide by the number of subsets
1215  // Return sensitivity
1216  return sensitivity;
1217 }
1218 
1219 // =====================================================================
1220 // ---------------------------------------------------------------------
1221 // ---------------------------------------------------------------------
1222 // =====================================================================
virtual int DataSpaceSpecificOperations(FLTNB a_data, FLTNB a_forwardModel, FLTNB *ap_backwardValues, FLTNB a_multiplicativeCorrections, FLTNB a_additiveCorrections, FLTNB a_blankValue, FLTNB a_quantificationFactor, oProjectionLine *ap_Line)=0
#define MODE_HISTOGRAM
virtual FLTNB GetBlankValue()
This is a pure virtual function implemented in the child classes.
FLTNB ForwardProjectWithSPECTAttenuation(FLTNB *ap_attenuation, FLTNB *ap_image=NULL)
int GetNbCardBasisFunctions()
Get the number of cardiac basis functions.
virtual int SensitivitySpecificOperations(FLTNB a_data, FLTNB a_forwardModel, FLTNB *ap_weight, FLTNB a_multiplicativeCorrections, FLTNB a_additiveCorrections, FLTNB a_blankValue, FLTNB a_quantificationFactor, oProjectionLine *ap_Line)=0
#define Cerr(MESSAGE)
FLTNB GetTimeBasisCoefficient(int a_timeBasisFunction, int a_timeFrame)
virtual int DataStep2Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
oImageDimensionsAndQuantification * mp_ImageDimensionsAndQuantification
FLTNB ComputeSensitivity(FLTNB ****a4p_sensitivityImage, int a_timeBasisFunction, int a_respBasisFunction, int a_cardBasisFunction, int a_voxel)
virtual FLTNB GetAdditiveCorrections(int a_bin)=0
virtual int ImageUpdateStep()
A public function used to perform the image update step of the optimizer.
virtual int LocalPreProcessingStep(int a_tbf, int a_rbf, int a_cbf, INTNB a_voxel, int a_th)
virtual int CheckSpecificParameters()=0
A private function used to check the parameters settings specific to the child optimizer.
int GetNbTimeBasisFunctions()
Get the number of time basis functions.
vOptimizer()
The constructor of vOptimizer.
FLTNB ForwardProject(oProjectionLine *ap_Line, FLTNB *ap_image=NULL)
virtual FLTNB GetMultiplicativeCorrections()=0
This is a pure virtual function implemented in the child classes.
int GetThreadNumber()
This function is used to get the thread number associated to this line.
int PreDataUpdateStep()
A public function used to do stuff that need to be done at the beginning of a subset (before the data...
bool GetRespStaticFlag()
Get the respiratory static flag that says if the reconstruction has only one respiratory gate or not...
virtual int DataStep6Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
void BackwardProject(FLTNB *ap_image, FLTNB a_value)
bool GetCardStaticFlag()
Get the cardiac static flag that says if the reconstruction has only one cardiac gate or not...
virtual int ImageSpaceSpecificOperations(FLTNB a_currentImageValue, FLTNB *ap_newImageValue, FLTNB a_sensitivity, FLTNB *ap_correctionValues, INTNB a_voxel, int tbf=-1, int rbf=-1, int cbf=-1)=0
FLTNB GetRespBasisCoefficient(int a_respBasisFunction, int a_respGate)
void BackwardProjectWithSPECTAttenuation(FLTNB *ap_attenuation, FLTNB *ap_image, FLTNB a_value)
virtual int DataStep3BackwardProjectSensitivity(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
virtual int InitializeSpecific()=0
A private function used to initialize everything specific to the child optimizer. ...
int GetNbThreadsMax()
Get the maximum between the number of threads used for projections and image operations.
void BackwardProject(oProjectionLine *ap_Line, FLTNB *ap_image, FLTNB a_value)
virtual int DataStep1ForwardProjectModel(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
FLTNB GetCardBasisCoefficient(int a_cardBasisFunction, int a_cardGate)
FLTNB GetQuantificationFactor(int a_bed, int a_frame, int a_respGate, int a_cardGate)
virtual int DataStep8ComputeFOM(oProjectionLine *ap_Line, vEvent *ap_Event, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
bool GetTimeStaticFlag()
Get the time static flag that says if the reconstruction has only one frame or not.
#define SPEC_TRANSMISSION
virtual ~vOptimizer()
The destructor of vOptimizer.
virtual int PreImageUpdateSpecificStep()
A private function used to perform any step required by the child optimizer, between the loop on even...
virtual FLTNB GetEventValue(int a_bin)=0
int PreImageUpdateStep()
A public function used to do stuff that need to be done between the loop over events and the image up...
This class is designed to manage and store system matrix elements associated to a vEvent...
int Initialize()
A public function used to initialize the optimizer.
virtual int DataStep4Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
Mother class for the Event objects.
int GetNbThreadsForImageComputation()
Get the number of threads used for image operations.
FLTNB ForwardProject(FLTNB *ap_image=NULL)
virtual int DataStep5ComputeCorrections(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
void SetCurrentIteration(int a_currentIteration)
int GetNbTOFBins()
This function is used to get the number of TOF bins in use.
int GetNbThreadsForProjection()
Get the number of threads used for projections.
int UpdateVisitedVoxels()
A public function used to update the &#39;visited&#39; voxels after each subset.
virtual void ShowHelpSpecific()=0
A function used to show help about the child module.
virtual int DataStep7BackwardProjectCorrections(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
void ShowHelp()
A function used to show help about the optimizer.
#define Cout(MESSAGE)
int GetNbRespBasisFunctions()
Get the number of respiratory basis functions.
int CheckParameters()
A public function used to check the parameters settings.
virtual int PreDataUpdateSpecificStep()
A private function used to perform any step required by the child optimizer, before the loop on event...
virtual FLTNB ComputePenaltyValue(int a_tbf, int a_rbf, int a_cbf, INTNB a_voxel, int a_th)=0