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