CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
code/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  mp_nbSubsets = NULL;
49  m_currentIteration = -1;
50  m_currentSubset = -1;
51  mp_Penalty = NULL;
54 }
55 
56 // =====================================================================
57 // ---------------------------------------------------------------------
58 // ---------------------------------------------------------------------
59 // =====================================================================
60 
62 {
63  // Free forward and backward buffers
65  {
67  if (m2p_forwardValues[th]) free(m2p_forwardValues[th]);
68  free(m2p_forwardValues);
69  }
71  {
72  for (int th=0; th<mp_ImageDimensionsAndQuantification->GetNbThreadsMax(); th++)
73  {
74  if (m3p_backwardValues[th])
75  {
76  for (int tof=0; tof<m_nbTOFBins; tof++) if (m3p_backwardValues[th][tof]) free(m3p_backwardValues[th][tof]);
77  free(m3p_backwardValues[th]);
78  }
79  }
80  free(m3p_backwardValues);
81  }
82  // Free pointers to attenuation images for SPECT
84  // Delete all image update statistics tables
86  {
94  }
95  // Delete FOM tables
97  {
98  // Delete the log-likelihood table
100  {
101  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
102  {
103  if (m4p_FOMLogLikelihood[fr])
104  {
105  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
106  {
107  if (m4p_FOMLogLikelihood[fr][rg])
108  {
109  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
110  if (m4p_FOMLogLikelihood[fr][rg][cg]) free(m4p_FOMLogLikelihood[fr][rg][cg]);
111  free(m4p_FOMLogLikelihood[fr][rg]);
112  }
113  }
114  free(m4p_FOMLogLikelihood[fr]);
115  }
116  }
117  free(m4p_FOMLogLikelihood);
118  }
119  // Delete the RMSE table
120  if (m4p_FOMRMSE)
121  {
122  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
123  {
124  if (m4p_FOMRMSE[fr])
125  {
126  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
127  {
128  if (m4p_FOMRMSE[fr][rg])
129  {
130  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
131  if (m4p_FOMRMSE[fr][rg][cg]) free(m4p_FOMRMSE[fr][rg][cg]);
132  free(m4p_FOMRMSE[fr][rg]);
133  }
134  }
135  free(m4p_FOMRMSE[fr]);
136  }
137  }
138  free(m4p_FOMRMSE);
139  }
140  // Delete the nbBins table
141  if (m4p_FOMNbBins)
142  {
143  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
144  {
145  if (m4p_FOMNbBins[fr])
146  {
147  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
148  {
149  if (m4p_FOMNbBins[fr][rg])
150  {
151  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
152  if (m4p_FOMNbBins[fr][rg][cg]) free(m4p_FOMNbBins[fr][rg][cg]);
153  free(m4p_FOMNbBins[fr][rg]);
154  }
155  }
156  free(m4p_FOMNbBins[fr]);
157  }
158  }
159  free(m4p_FOMNbBins);
160  }
161  // Delete the nbData table
162  if (m4p_FOMNbData)
163  {
164  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
165  {
166  if (m4p_FOMNbData[fr])
167  {
168  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
169  {
170  if (m4p_FOMNbData[fr][rg])
171  {
172  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
173  if (m4p_FOMNbData[fr][rg][cg]) free(m4p_FOMNbData[fr][rg][cg]);
174  free(m4p_FOMNbData[fr][rg]);
175  }
176  }
177  free(m4p_FOMNbData[fr]);
178  }
179  }
180  free(m4p_FOMNbData);
181  }
182  // Delete the penalty table
183  if (m4p_FOMPenalty)
184  {
185  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
186  {
187  if (m4p_FOMPenalty[fr])
188  {
189  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
190  {
191  if (m4p_FOMPenalty[fr][rg])
192  {
193  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
194  if (m4p_FOMPenalty[fr][rg][cg]) free(m4p_FOMPenalty[fr][rg][cg]);
195  free(m4p_FOMPenalty[fr][rg]);
196  }
197  }
198  free(m4p_FOMPenalty[fr]);
199  }
200  }
201  free(m4p_FOMPenalty);
202  }
203  }
204 }
205 
206 // =====================================================================
207 // ---------------------------------------------------------------------
208 // ---------------------------------------------------------------------
209 // =====================================================================
210 
212 {
213  // Call the specific help function from the children
215  // Say if the child optimizer is compatible with histogram and/or list-mode data
216  if (m_listmodeCompatibility && m_histogramCompatibility) cout << "This optimizer is compatible with both histogram and list-mode data." << endl;
217  else if (m_listmodeCompatibility) cout << "This optimizer is only compatible with list-mode data." << endl;
218  else if (m_histogramCompatibility) cout << "This optimizer is only compatible with histogram data." << endl;
219  else cout << "This optimizer is a total non-sense, not compatible with histogram nor list-mode data !!!" << endl;
220  // Say if the child optimizer is compatible with emission and/or transmission data
221  if (m_emissionCompatibility && m_transmissionCompatibility) cout << "This optimizer is compatible with both emission and transmission data." << endl;
222  else if (m_emissionCompatibility) cout << "This optimizer is only compatible with emission data." << endl;
223  else if (m_transmissionCompatibility) cout << "This optimizer is only compatible with transmission data." << endl;
224  else cout << "This optimizer is a total non-sense, not compatible with emission nor transmission data !!!" << endl;
225 }
226 
227 // =====================================================================
228 // ---------------------------------------------------------------------
229 // ---------------------------------------------------------------------
230 // =====================================================================
231 
233 {
234  // Check image dimensions
236  {
237  Cerr("***** vOptimizer::CheckParameters() -> oImageDimensionsAndQuantification is null !" << endl);
238  return 1;
239  }
240  // Check image space
241  if (mp_ImageSpace==NULL)
242  {
243  Cerr("***** vOptimizer::CheckParameters() -> oImageSpace is null !" << endl);
244  return 1;
245  }
246  // Check data mode
248  {
249  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data mode provided !" << endl);
250  return 1;
251  }
252  // Check data type
254  {
255  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data type provided !" << endl);
256  return 1;
257  }
258  // Check data spec
260  {
261  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data specificity provided (emission or transmission) !" << endl);
262  return 1;
263  }
264  // Check verbose level
265  if (m_verbose<0)
266  {
267  Cerr("***** vOptimizer::CheckParameters() -> Verbose level is negative !" << endl);
268  return 1;
269  }
270  // Listmode incompatibility
272  {
273  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with listmode data !" << endl);
274  return 1;
275  }
276  // Histogram incompatibility
278  {
279  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with histogram data !" << endl);
280  return 1;
281  }
282  // Emission compatability
284  {
285  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with emission data !" << endl);
286  return 1;
287  }
288  // Transmission compatibility
290  {
291  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with transmission data !" << endl);
292  return 1;
293  }
294  // Send an error if FOM computation is required and verbose is null
295  if (m_optimizerFOMFlag && m_verbose==0)
296  {
297  Cerr("***** vOptimizer::CheckParameters() -> Asked to compute FOMs in the data-space whereas the verbose is not set !" << endl);
298  return 1;
299  }
300  // Send an error if image statistics computation is required and verbose is null
302  {
303  Cerr("***** vOptimizer::CheckParameters() -> Asked to compute image update statistics whereas the verbose is not set !" << endl);
304  return 1;
305  }
306  // Send an error if FOM computation with list-mode data (it has no sense)
308  {
309  Cerr("***** vOptimizer::CheckParameters() -> Computing FOMs in the data-space with list-mode data is not possible !" << endl);
310  return 1;
311  }
312  // Call the CheckSpecificParameters function of the child
314  {
315  Cerr("***** vOptimizer::CheckParameters() -> A problem occurred while checking parameters specific to the optimizer module !" << endl);
316  return 1;
317  }
318  // Normal end
319  return 0;
320 }
321 
322 // =====================================================================
323 // ---------------------------------------------------------------------
324 // ---------------------------------------------------------------------
325 // =====================================================================
326 
328 {
329  // Allocate forward buffers (as many as threads and as many as TOF bins)
332  m2p_forwardValues[th] = (FLTNB*)malloc(m_nbTOFBins*sizeof(FLTNB));
333 
334  // Allocate backward buffers (as many as threads then as many as TOF bins and as many as backward images)
336  for (int th=0; th<mp_ImageDimensionsAndQuantification->GetNbThreadsMax(); th++)
337  {
338  m3p_backwardValues[th] = (FLTNB**)malloc(m_nbTOFBins*sizeof(FLTNB*));
339  for (int tof=0; tof<m_nbTOFBins; tof++)
340  m3p_backwardValues[th][tof] = (FLTNB*)malloc(m_nbBackwardImages*sizeof(FLTNB));
341  }
342 
343  // Allocate pointers to the current attenuation image for SPECT attenuation correction
346 
347  // Allocate the thread-safe tables for image update statistics computation
349  {
357  }
358 
359  // We allocate the tables that hold the figures-of-merit if FOM flag is on
360  if (m_optimizerFOMFlag)
361  {
362  m4p_FOMNbBins = (uint64_t****)malloc(mp_ImageDimensionsAndQuantification->GetNbTimeFrames()*sizeof(uint64_t***));
367  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
368  {
369  m4p_FOMNbBins[fr] = (uint64_t***)malloc(mp_ImageDimensionsAndQuantification->GetNbRespGates()*sizeof(uint64_t**));
374  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
375  {
376  m4p_FOMNbBins[fr][rg] = (uint64_t**)malloc(mp_ImageDimensionsAndQuantification->GetNbCardGates()*sizeof(uint64_t*));
381  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
382  {
383  m4p_FOMNbBins[fr][rg][cg] = (uint64_t*)malloc(mp_ImageDimensionsAndQuantification->GetNbThreadsForProjection()*sizeof(uint64_t));
388  }
389  }
390  }
391  }
392 
393  // Call the specific initialization function of the child
394  if (InitializeSpecific())
395  {
396  Cerr("***** vOptimizer::Initialize() -> A problem occurred while initializing stuff specific to the optimizer module !" << endl);
397  return 1;
398  }
399 
400  // Normal end
401  return 0;
402 }
403 
404 // =====================================================================
405 // ---------------------------------------------------------------------
406 // ---------------------------------------------------------------------
407 // =====================================================================
408 
410 {
411  // Simply reset FOMs if activated
412  if (m_optimizerFOMFlag)
413  {
414  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
415  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
416  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
417  {
419  {
420  m4p_FOMNbBins[fr][rg][cg][th] = 0;
421  m4p_FOMNbData[fr][rg][cg][th] = 0.;
422  m4p_FOMLogLikelihood[fr][rg][cg][th] = 0.;
423  m4p_FOMRMSE[fr][rg][cg][th] = 0.;
424  }
426  {
427  m4p_FOMPenalty[fr][rg][cg][th] = 0.;
428  }
429  }
430  }
431  // Then call the specific pre-update-step function from the child optimizer
433  {
434  Cerr("***** vOptimizer::PreDataUpdateStep() -> A problem occurred while applying the specific pre-data-update step of the optimizer module !" << endl);
435  return 1;
436  }
437  // Normal end
438  return 0;
439 }
440 
441 // =====================================================================
442 // ---------------------------------------------------------------------
443 // ---------------------------------------------------------------------
444 // =====================================================================
445 
447 {
448  // By default, just do nothing here. This function is designed to be overloaded by the specific optimizer if needed
449  return 0;
450 }
451 
452 // =====================================================================
453 // ---------------------------------------------------------------------
454 // ---------------------------------------------------------------------
455 // =====================================================================
456 
458 {
459  // Then call the specific post-update-step function from the child optimizer
461  {
462  Cerr("***** vOptimizer::PreImageUpdateStep() -> A problem occurred while applying the specific post-data-update step of the optimizer module !" << endl);
463  return 1;
464  }
465  // Compute optimizer figures-of-merit (data fidelity term and penalty)
466  if (m_optimizerFOMFlag)
467  {
468  // Compute data FOM related to the data fildelity term
469  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
470  {
471  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
472  {
473  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
474  {
475  // Merge thread results
477  {
478  m4p_FOMNbBins[fr][rg][cg][0] += m4p_FOMNbBins[fr][rg][cg][th];
479  m4p_FOMNbData[fr][rg][cg][0] += m4p_FOMNbData[fr][rg][cg][th];
480  m4p_FOMLogLikelihood[fr][rg][cg][0] += m4p_FOMLogLikelihood[fr][rg][cg][th];
481  m4p_FOMRMSE[fr][rg][cg][0] += m4p_FOMRMSE[fr][rg][cg][th];
482  }
483  // Finish mean number of data counts per channel
484  m4p_FOMNbData[fr][rg][cg][0] /= ((HPFLTNB)(m4p_FOMNbBins[fr][rg][cg][0]));
485  // Finish RMSE computation
486  m4p_FOMRMSE[fr][rg][cg][0] /= ((HPFLTNB)(m4p_FOMNbBins[fr][rg][cg][0]));
487  m4p_FOMRMSE[fr][rg][cg][0] = sqrt(m4p_FOMRMSE[fr][rg][cg][0]);
488  }
489  }
490  }
491  // Compute penalty FOM
492  if (mp_Penalty!=NULL)
493  {
494  // As the penalty will be innocently applied to the current image, it is performed on time basis functions
495  // and not frames/gates. As a result, we loop over time basis functions first, compute the penalty, and
496  // distribute its value with respect to the contribution of each time basis function to the frames/gates.
497  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
498  {
499  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
500  {
501  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
502  {
503  // In case a problem occurs during the multi-threaded loop
504  bool problem = false;
505  // The voxel index
506  INTNB v;
507  // Multi-threading over voxels
508  #pragma omp parallel for private(v) schedule(guided)
510  {
511  // Get the thread index
512  int th = 0;
513  #ifdef CASTOR_OMP
514  th = omp_get_thread_num();
515  #endif
516  // Compute the penalty term
517  if (mp_Penalty->LocalPreProcessingStep(tbf,rbf,cbf,v,th))
518  {
519  Cerr("***** vOptimizer::PreImageUpdateStep() -> A problem occurred while precomputing the local step of the penalty for voxel " << v << " !" << endl);
520  problem = true;
521  }
522  HPFLTNB penalty_value = mp_Penalty->ComputePenaltyValue(tbf,rbf,cbf,v,th);
523  // Distribuete the value to frames/gates with respect to basis functions
524  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
525  {
527  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
528  {
530  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
531  {
533  m4p_FOMPenalty[fr][rg][cg][th] += time_basis_coef * resp_basis_coef * card_basis_coef * penalty_value;
534  }
535  }
536  }
537  }
538  // Potential problems during the loop
539  if (problem)
540  {
541  Cerr("***** vOptimizer::PreImageUpdateStep() -> A problem occurred in the multi-threaded loop for the penalty, stop now !" << endl);
542  return 1;
543  }
544  }
545  }
546  }
547  // Merge thread results
548  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
549  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
550  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
551  {
553  m4p_FOMPenalty[fr][rg][cg][0] += m4p_FOMPenalty[fr][rg][cg][th];
554  // Also divide by the number of subsets for balance with respect to the data fidelity term
556  }
557  }
558  // Verbose
559  Cout("vOptimizer::PreImageUpdateStep() -> Optimizer figures-of-merit related to the objective function" << endl);
560  // Number of spaces for printing
561  string spaces_fr = ""; if (mp_ImageDimensionsAndQuantification->GetNbTimeFrames()>1) spaces_fr += " ";
562  string spaces_rg = ""; if (mp_ImageDimensionsAndQuantification->GetNbRespGates()>1) spaces_rg += " ";
563  string spaces_cg = ""; if (mp_ImageDimensionsAndQuantification->GetNbCardGates()>1) spaces_cg += " ";
564  // Loop on frames and gates
565  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
566  {
567  // Verbose
569  Cout(spaces_fr << "Frame " << fr+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbTimeFrames() << endl);
570  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
571  {
572  // Verbose
574  Cout(spaces_fr << spaces_rg << "Respiratory gate " << rg+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbRespGates() << endl);
575  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
576  {
577  // Verbose
579  Cout(spaces_fr << spaces_rg << spaces_cg << "Cardiac gate " << cg+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbCardGates() << endl);
580  // Print out
581  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Number of data channels used in optimization: " << m4p_FOMNbBins[fr][rg][cg][0] << endl);
582  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Mean number of data counts per channel: " << m4p_FOMNbData[fr][rg][cg][0] << endl);
583  if (mp_Penalty!=NULL)
584  {
585  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Penalty: " << m4p_FOMPenalty[fr][rg][cg][0] << endl);
586  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Log-likelihood: " << m4p_FOMLogLikelihood[fr][rg][cg][0]
587  << " (including penalty: " << m4p_FOMLogLikelihood[fr][rg][cg][0] - m4p_FOMPenalty[fr][rg][cg][0] << ")" << endl);
588  Cout(spaces_fr << spaces_rg << spaces_cg << " --> RMSE: " << m4p_FOMRMSE[fr][rg][cg][0]
589  << " (including penalty: " << m4p_FOMRMSE[fr][rg][cg][0] + m4p_FOMPenalty[fr][rg][cg][0] << ")" << endl);
590  }
591  else
592  {
593  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Log-likelihood: " << m4p_FOMLogLikelihood[fr][rg][cg][0] << endl);
594  Cout(spaces_fr << spaces_rg << spaces_cg << " --> RMSE: " << m4p_FOMRMSE[fr][rg][cg][0] << endl);
595  }
596  }
597  }
598  }
599  }
600  // Normal end
601  return 0;
602 }
603 
604 // =====================================================================
605 // ---------------------------------------------------------------------
606 // ---------------------------------------------------------------------
607 // =====================================================================
608 
610 {
611  // By default, just do nothing here. This function is designed to be overloaded by the specific optimizer if needed
612  return 0;
613 }
614 
615 // =====================================================================
616 // ---------------------------------------------------------------------
617 // ---------------------------------------------------------------------
618 // =====================================================================
619 
621  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
622  int a_th )
623 {
624  // =======================================================================
625  // 4D forward projection including framing, respiratory and cardiac gating
626  // =======================================================================
627 
628  // Clean buffers
629  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
630  {
631  m2p_forwardValues[a_th][b] = 0.;
632  for (int img=0; img<m_nbBackwardImages; img++) m3p_backwardValues[a_th][b][img] = 0.;
633  }
634 
635  // Loop over time basis functions
636  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
637  {
638  // Get frame/basis coefficient and continue if null
639  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(tbf,a_timeFrame);
640  if (time_basis_coef==0.) continue;
641  // Loop over respiratory basis functions
642  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
643  {
644  // Get resp_gate/basis coefficient and continue if null
645  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(rbf,a_respGate);
646  if (resp_basis_coef==0.) continue;
647  // Loop over cardiac basis functions
648  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
649  {
650  // Get card_gate_basis coefficient and continue if null
651  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(cbf,a_cardGate);
652  if (card_basis_coef==0.) continue;
653  // Compute global coefficient
654  FLTNB global_basis_coef = time_basis_coef * resp_basis_coef * card_basis_coef;
655  // Loop over TOF bins
656  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
657  {
658  // Set the current TOF bin
659  ap_Line->SetCurrentTOFBin(b);
660  // Project the image and apply dynamic coefficient converions
661  m2p_forwardValues[a_th][b] += ForwardProject(ap_Line,mp_ImageSpace->m4p_forwardImage[tbf][rbf][cbf]) * global_basis_coef;
662  }
663  }
664  }
665  }
666 
667  // ==============================================
668  // Differentiate emission and transmission models
669  // ==============================================
670 
671  // -----------------
672  // Case for emission
673  // -----------------
675  {
676  // Add additive terms (we multiply by the frame duration because the additive terms are provided as rates)
678  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);
679  else
680  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);
681  }
682  // ---------------------
683  // Case for transmission
684  // ---------------------
685  else if (m_dataSpec==SPEC_TRANSMISSION)
686  {
687  int no_tof = 0;
688  // Ignore TOF as it has no sense; add scatter rate multiplied by the frame duration; take blank into account
689  m2p_forwardValues[a_th][no_tof] = ap_Event->GetBlankValue() * exp(-m2p_forwardValues[a_th][no_tof])
690  + ap_Event->GetAdditiveCorrections(no_tof) * mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame);
691  }
692 
693  // End
694  return 0;
695 }
696 
697 // =====================================================================
698 // ---------------------------------------------------------------------
699 // ---------------------------------------------------------------------
700 // =====================================================================
701 
703  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
704  int a_th )
705 {
706  // Deliberately do nothing!
707  return 0;
708 }
709 
710 // =====================================================================
711 // ---------------------------------------------------------------------
712 // ---------------------------------------------------------------------
713 // =====================================================================
714 
716  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
717  int a_th )
718 {
719  // Loop over TOF bins
720  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
721  {
722  // Set the current TOF bin of the projection-line
723  ap_Line->SetCurrentTOFBin(b);
724  // The weight associated to the line to be backprojected
725  FLTNB weight = 0.;
726  // Call the function to compute specific sensitivity
728  (
729  ap_Event->GetEventValue(b), // the measured data
730  m2p_forwardValues[a_th][b], // the forward model
731  &weight, // the sensitivity weight
732  ap_Event->GetMultiplicativeCorrections(), // the multiplicative corrections
733  ap_Event->GetAdditiveCorrections(b)*mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame), // the additive corrections
734  ap_Event->GetBlankValue(), // the blank value
735  mp_ImageDimensionsAndQuantification->GetQuantificationFactor(a_bed,a_timeFrame, a_respGate, a_cardGate), // the quantification factor
736  ap_Line // the projection line
737  ))
738  {
739  Cerr("***** vOptimizer::DataStep3BackwardProjectSensitivity() -> A problem occurred while performing the sensitivity specific operations !" << endl);
740  return 1;
741  }
742  // Back-project the weight into the sensitivity matrix
743  BackwardProject(ap_Line, mp_ImageSpace->m5p_sensitivity[a_th][a_timeFrame][a_respGate][a_cardGate], weight);
744  }
745 
746  // End
747  return 0;
748 }
749 
750 // =====================================================================
751 // ---------------------------------------------------------------------
752 // ---------------------------------------------------------------------
753 // =====================================================================
754 
756  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
757  int a_th )
758 {
759  // Deliberately do nothing!
760  return 0;
761 }
762 
763 // =====================================================================
764 // ---------------------------------------------------------------------
765 // ---------------------------------------------------------------------
766 // =====================================================================
767 
769  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
770  int a_th )
771 {
772  // Loop on TOF bins
773  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
774  {
775  // Set the current TOF bin of the projection-line
776  ap_Line->SetCurrentTOFBin(b);
777  // Then call the pure virtual function for specific data space operations
779  (
780  ap_Event->GetEventValue(b), // the measured data
781  m2p_forwardValues[a_th][b], // the forward model
782  m3p_backwardValues[a_th][b], // the backward values (the result)
783  ap_Event->GetMultiplicativeCorrections(), // the multiplicative corrections
784  ap_Event->GetAdditiveCorrections(b)*mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame), // the additive corrections
785  ap_Event->GetBlankValue(), // the blank value
786  mp_ImageDimensionsAndQuantification->GetQuantificationFactor(a_bed,a_timeFrame, a_respGate, a_cardGate), // the quantification factor
787  ap_Line // the projection line
788  ))
789  {
790  Cerr("***** vOptimizer::DataStep5ComputeCorrections() -> A problem occurred while performing specific data space operations !" << endl);
791  return 1;
792  }
793  }
794 
795  // End
796  return 0;
797 }
798 
799 // =====================================================================
800 // ---------------------------------------------------------------------
801 // ---------------------------------------------------------------------
802 // =====================================================================
803 
805  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
806  int a_th )
807 {
808  // Deliberately do nothing!
809  return 0;
810 }
811 
812 // =====================================================================
813 // ---------------------------------------------------------------------
814 // ---------------------------------------------------------------------
815 // =====================================================================
816 
818  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
819  int a_th )
820 {
821  // ========================================================================
822  // 4D backward projection including framing, respiratory and cardiac gating
823  // ========================================================================
824 
825  // Loop over time basis functions
826  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
827  {
828  // Get frame/basis coefficient and continue if null
829  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(tbf,a_timeFrame);
830  if (time_basis_coef==0.) continue;
831  // Loop over respiratory basis functions
832  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
833  {
834  // Get resp_gate/basis coefficient and continue if null
835  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(rbf,a_respGate);
836  if (resp_basis_coef==0.) continue;
837  // Loop over cardiac basis functions
838  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
839  {
840  // Get card_gate_basis coefficient and continue if null
841  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(cbf,a_cardGate);
842  if (card_basis_coef==0.) continue;
843  // Compute global coefficient
844  FLTNB global_basis_coef = time_basis_coef * resp_basis_coef * card_basis_coef;
845  // Loop over TOF bins
846  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
847  {
848  // Set the current TOF bin of the projection-line
849  ap_Line->SetCurrentTOFBin(b);
850  // Loop over backprojection images
851  for (int img=0; img<m_nbBackwardImages; img++)
852  {
853  // Continue if backproj is null or infinity
854  if (m3p_backwardValues[a_th][b][img]==0. || !isfinite(m3p_backwardValues[a_th][b][img])) continue;
855  // Backward project into the backward image
856  BackwardProject( ap_Line, mp_ImageSpace->m6p_backwardImage[img][a_th][tbf][rbf][cbf] , m3p_backwardValues[a_th][b][img] * global_basis_coef );
857  }
858  }
859  }
860  }
861  }
862 
863  // End
864  return 0;
865 }
866 
867 // =====================================================================
868 // ---------------------------------------------------------------------
869 // ---------------------------------------------------------------------
870 // =====================================================================
871 
873  int a_timeFrame, int a_respGate, int a_cardGate,
874  int a_thread )
875 {
876  // Loop on TOF bins
877  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
878  {
879  // Get the model and the data
880  HPFLTNB model = m2p_forwardValues[a_thread][b];
881  HPFLTNB data = ap_Event->GetEventValue(b);
882  // Update log likelihood only if model is strictly positive (to avoid numerical issues, we use a threshold at e-10)
883  if (model>1.e-10) m4p_FOMLogLikelihood[a_timeFrame][a_respGate][a_cardGate][a_thread] += data*log(model)-model;
884  // Update RMSE
885  m4p_FOMRMSE[a_timeFrame][a_respGate][a_cardGate][a_thread] += (data-model)*(data-model);
886  // Update number of data channels
887  m4p_FOMNbBins[a_timeFrame][a_respGate][a_cardGate][a_thread]++;
888  // Update the number of data counts
889  m4p_FOMNbData[a_timeFrame][a_respGate][a_cardGate][a_thread] += data;
890  }
891  // Normal end
892  return 0;
893 }
894 
895 // =====================================================================
896 // ---------------------------------------------------------------------
897 // ---------------------------------------------------------------------
898 // =====================================================================
899 
901 {
902  // Verbose
903  if (m_verbose>=3) Cout("vOptimizer::UpdateVisitedVoxels() -> Tick visited voxels based on sensitivity" << endl);
904 
905  // The purpose of this mp_visitedVoxelsImage is to keep track within a complete iteration, of which voxels were
906  // visited by some LORs. Everything is based on sensitivity. With list-mode data, it means that the visited voxels
907  // are the same for all subsets. With histogram data, if the sensitivity is computed on the fly from the data, then
908  // the visited voxels may change for each voxel. When updating the voxel values in the ImageUpdateStep, only voxels
909  // with non-zero sensitivity are updated, to avoid artificially decreasing the signal in voxels that may not have
910  // been visited in this subset but which may be visited in another subset (this is true only for histogram data).
911  // But, at the end of each iteration, the mp_visitedVoxelsImage is used to checked for never-visited voxels, in which
912  // case they will be set to 0 to avoid that they keep their initial value. In other words, these never-visited voxels
913  // are somewhat excluded from the reconstruction.
914 
915  // 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
916  // gates and look if voxels were 'visited' by LORs (i.e. the sensitivity is not null);
917  // the information is saved in one static 3D image and then applied to all frames/gates;
918  int thread0 = 0;
919  int resp_gate0 = 0;
920  int card_gate0 = 0;
921  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
922  for (int v=0; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ(); v++)
923  // If the voxel has been visited, we add 1. to the mp_Image->mp_visitedVoxelsImage
924  if (mp_ImageSpace->m5p_sensitivity[thread0][fr][resp_gate0][card_gate0][v]!=0.) mp_ImageSpace->mp_visitedVoxelsImage[v] += 1.;
925 
926  // End
927  return 0;
928 }
929 
930 // =====================================================================
931 // ---------------------------------------------------------------------
932 // ---------------------------------------------------------------------
933 // =====================================================================
934 
936 {
937  // Verbose
938  if (m_verbose>=2 || m_optimizerImageStatFlag) Cout("vOptimizer::ImageUpdateStep() -> Start image update" << endl);
939  // Number of spaces for nice printing
940  string spaces_tbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions()>1) spaces_tbf += " ";
941  string spaces_rbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions()>1) spaces_rbf += " ";
942  string spaces_cbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions()>1) spaces_cbf += " ";
943  // Loop over time basis functions
944  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
945  {
946  // Verbose
948  {
949  if (mp_ImageDimensionsAndQuantification->GetTimeStaticFlag()) Cout(spaces_tbf << "--> Time frame " << tbf+1 << endl);
950  else Cout(spaces_tbf << "--> Time basis function: " << tbf+1 << endl);
951  }
952  // Loop over respiratory basis functions
953  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
954  {
955  // Verbose
957  {
958  if (mp_ImageDimensionsAndQuantification->GetRespStaticFlag()) Cout(spaces_tbf << spaces_rbf << "--> Respiratory gate " << rbf+1 << endl);
959  else Cout(spaces_tbf << spaces_rbf << "--> Respiratory basis function: " << rbf+1 << endl);
960  }
961  // Loop over cardiac basis functions
962  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
963  {
964  // Verbose
966  {
967  if (mp_ImageDimensionsAndQuantification->GetCardStaticFlag()) Cout(spaces_tbf << spaces_rbf << spaces_cbf << "--> Cardiac gate " << cbf+1 << endl);
968  else Cout(spaces_tbf << spaces_rbf << spaces_cbf << "--> Cardiac basis function: " << cbf+1 << endl);
969  }
970  // Set the number of threads for image computation
971  #ifdef CASTOR_OMP
973  #endif
974  // Reset counter for image stats
976  {
978  {
979  mp_imageStatNbVox[th] = 0;
980  mp_imageStatMin[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][0];
981  mp_imageStatMax[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][0];
982  mp_imageStatMean[th] = 0.;
983  mp_imageStatVariance[th] = 0.;
984  mp_correctionStatMean[th] = 0.;
985  mp_correctionStatVariance[th] = 0.;
986  }
987  }
988  // OMP loop over voxels
989  INTNB v;
990  bool error = false;
991  #pragma omp parallel for private(v) schedule(guided)
992  for (v=0 ; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ() ; v++)
993  {
994  // Get the thread index
995  int th = 0;
996  #ifdef CASTOR_OMP
997  th = omp_get_thread_num();
998  #endif
999  // Compute sensitivity
1000  FLTNB sensitivity = ComputeSensitivity(mp_ImageSpace->m5p_sensitivity[0],tbf,rbf,cbf,v);
1001  // If the sensitivity is negative or null, we skip this voxel
1002  if (sensitivity<=0.) continue;
1003  // Keep the current image value in a buffer for update statistics
1004  FLTNB old_image_value = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
1005  // Fill class-local buffer of backward values (thread-safe)
1006  for (int nb=0; nb<m_nbBackwardImages; nb++)
1007  m3p_backwardValues[th][0][nb] = mp_ImageSpace->m6p_backwardImage[nb][0][tbf][rbf][cbf][v];
1008  // Then call the pure virtual function for specific data space operations
1010  (
1011  mp_ImageSpace->m4p_image[tbf][rbf][cbf][v],
1012  // 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
1013  // to perform MLAA reconstruction but using PET basis functions different from diracs for example... In this case, we would get
1014  // segmentation faults. But actually, this function would be overloaded by MLAA anyway...
1015  //mp_Image->m4p_attenuation != NULL ? mp_Image->m4p_attenuation[tbf][rbf][cbf][v] : 0. ,
1016  &mp_ImageSpace->m4p_image[tbf][rbf][cbf][v],
1017  //mp_Image->m4p_attenuation != NULL ? &mp_Image->m4p_attenuation[tbf][rbf][cbf][v] : NULL ,
1018  sensitivity,
1019  m3p_backwardValues[th][0],
1020  v,
1021  tbf,
1022  rbf,
1023  cbf
1024  ))
1025  {
1026  Cerr("***** vOptimizer::ImageUpdateStep() -> A problem occurred while performing image space specific operations of the optimizer !" << endl);
1027  error = true;
1028  }
1029  // Do some image statistics
1031  {
1032  if (mp_imageStatMin[th]>mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]) mp_imageStatMin[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
1033  if (mp_imageStatMax[th]<mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]) mp_imageStatMax[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
1034  // The one pass algorithm is used to compute the variance here for both the image and correction step.
1035  // Do not change anyhting to the order of the operations !
1036  mp_imageStatNbVox[th]++;
1037  HPFLTNB sample = ((HPFLTNB)(mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]));
1038  HPFLTNB delta = sample - mp_imageStatMean[th];
1039  mp_imageStatMean[th] += delta / ((HPFLTNB)(mp_imageStatNbVox[th]));
1040  mp_imageStatVariance[th] += delta * (sample - mp_imageStatMean[th]);
1041  // Same one-pass variance algorithm for the correction step
1042  sample = ((HPFLTNB)(mp_ImageSpace->m4p_image[tbf][rbf][cbf][v] - old_image_value));
1043  delta = sample - mp_correctionStatMean[th];
1044  mp_correctionStatMean[th] += delta / ((HPFLTNB)(mp_imageStatNbVox[th]));
1045  mp_correctionStatVariance[th] += delta * (sample - mp_correctionStatMean[th]);
1046  }
1047  }
1048  // Check for errors (we cannot stop the multi-threaded loop so we do it here)
1049  if (error) return 1;
1050  // Finish image update statistics computations: these operations are specific to the merging of the
1051  // multi-threaded one-pass variance estimations used above.
1053  {
1055  {
1056  // Image minimum and maximum
1059  // Cast the counts into HPFLTNB
1060  HPFLTNB count1 = ((HPFLTNB)(mp_imageStatNbVox[0]));
1061  HPFLTNB count2 = ((HPFLTNB)(mp_imageStatNbVox[th]));
1062  // Update the count
1064  HPFLTNB count12 = ((HPFLTNB)(mp_imageStatNbVox[0]));
1065  // Compute the delta mean
1066  HPFLTNB delta = mp_imageStatMean[0] - mp_imageStatMean[th];
1067  // Update the variance
1068  mp_imageStatVariance[0] += mp_imageStatVariance[th] + delta * delta * count1 * count2 / count12;
1069  // Update the mean
1070  mp_imageStatMean[0] = (count1*mp_imageStatMean[0] + count2*mp_imageStatMean[th]) / count12;
1071  // Do the same for the correction step
1073  mp_correctionStatVariance[0] += mp_correctionStatVariance[th] + delta * delta * count1 * count2 / count12;
1074  mp_correctionStatMean[0] = (count1*mp_correctionStatMean[0] + count2*mp_correctionStatMean[th]) / count12;
1075  }
1076  // Finish variance computation
1079  // Print out
1080  Cout(spaces_tbf << spaces_rbf << spaces_cbf << " --> Image update step "
1081  << " | mean: " << mp_correctionStatMean[0]
1082  << " | stdv: " << sqrt(mp_correctionStatVariance[0]) << endl);
1083  Cout(spaces_tbf << spaces_rbf << spaces_cbf << " --> New image estimate"
1084  << " | mean: " << mp_imageStatMean[0]
1085  << " | stdv: " << sqrt(mp_imageStatVariance[0])
1086  << " | min/max: " << mp_imageStatMin[0] << "/" << mp_imageStatMax[0] << endl);
1087  }
1088  }
1089  }
1090  }
1091 
1092  // End
1093  return 0;
1094 }
1095 
1096 // =====================================================================
1097 // ---------------------------------------------------------------------
1098 // ---------------------------------------------------------------------
1099 // =====================================================================
1100 
1102 {
1103  // Particular case for SPECT with attenuation correction
1104  if (m_dataType==TYPE_SPECT && m2p_attenuationImage[ap_Line->GetThreadNumber()]!=NULL)
1105  return ap_Line->ForwardProjectWithSPECTAttenuation( m2p_attenuationImage[ap_Line->GetThreadNumber()], ap_image );
1106  // Otherwise call the function without attenuation
1107  else
1108  return ap_Line->ForwardProject( ap_image );
1109 }
1110 
1111 // =====================================================================
1112 // ---------------------------------------------------------------------
1113 // ---------------------------------------------------------------------
1114 // =====================================================================
1115 
1116 void vOptimizer::BackwardProject(oProjectionLine* ap_Line, FLTNB* ap_image, FLTNB a_value)
1117 {
1118  // Particular case for SPECT with attenuation correction
1119  if (m_dataType==TYPE_SPECT && m2p_attenuationImage[ap_Line->GetThreadNumber()]!=NULL)
1120  ap_Line->BackwardProjectWithSPECTAttenuation( m2p_attenuationImage[ap_Line->GetThreadNumber()], ap_image, a_value );
1121  // Otherwise call the function without attenuation
1122  else
1123  ap_Line->BackwardProject( ap_image, a_value );
1124 }
1125 
1126 // =====================================================================
1127 // ---------------------------------------------------------------------
1128 // ---------------------------------------------------------------------
1129 // =====================================================================
1130 
1131 FLTNB vOptimizer::ComputeSensitivity( FLTNB**** a4p_sensitivityMatrix,
1132  int a_timeBasisFunction, int a_respBasisFunction, int a_cardBasisFunction,
1133  int a_voxel )
1134 {
1135  // The sensitivity to be computed
1136  FLTNB sensitivity = 0.;
1137  // Loop over time frames
1138  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
1139  {
1140  // Get frame/basis coefficient and continue if null
1141  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(a_timeBasisFunction,fr);
1142  if (time_basis_coef==0.) continue;
1143  // Loop over respiratory gates
1144  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
1145  {
1146  // Get resp_gate/basis coefficient and continue if null
1147  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(a_respBasisFunction,rg);
1148  if (resp_basis_coef==0.) continue;
1149  // Loop over cardiac gates
1150  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
1151  {
1152  // Get card_gate_basis coefficient and continue if null
1153  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(a_cardBasisFunction,cg);
1154  if (card_basis_coef==0.) continue;
1155  // Add actual weight contribution
1156  sensitivity += a4p_sensitivityMatrix[fr][rg][cg][a_voxel] *
1157  time_basis_coef * resp_basis_coef * card_basis_coef;
1158  }
1159  }
1160  }
1161  // If listmode, then divide by the number of subsets
1163  // Return sensitivity
1164  return sensitivity;
1165 }
1166 
1167 // =====================================================================
1168 // ---------------------------------------------------------------------
1169 // ---------------------------------------------------------------------
1170 // =====================================================================
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)
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