CASToR  2.0
Tomographic Reconstruction (PET/SPECT/CT)
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
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-2018 all CASToR contributors listed below:
18 
19  --> current contributors: Thibaut MERLIN, Simon STUTE, Didier BENOIT, Claude COMTAT, Marina FILIPOVIC, Mael MILLARDET
20  --> past contributors: Valentin VIELZEUF
21 
22 This is CASToR version 2.0.
23 */
24 
31 #include "vOptimizer.hh"
32 
33 // =====================================================================
34 // ---------------------------------------------------------------------
35 // ---------------------------------------------------------------------
36 // =====================================================================
37 
39 {
40  // Affect default values
42  mp_ImageSpace = NULL;
44  m_nbTOFBins = 0;
45  m_initialValue = 0.;
50  m2p_forwardValues = NULL;
51  m3p_backwardValues = NULL;
54  m_optimizerFOMFlag = false;
56  m4p_FOMLogLikelihood = NULL;
57  m4p_FOMNbBins = NULL;
58  m4p_FOMRMSE = NULL;
59  m4p_FOMNbData = 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_penaltyEnergyFunctionDerivativesOrder = 0;
70 }
71 
72 // =====================================================================
73 // ---------------------------------------------------------------------
74 // ---------------------------------------------------------------------
75 // =====================================================================
76 
78 {
79  // Free forward and backward buffers
81  {
83  if (m2p_forwardValues[th]) free(m2p_forwardValues[th]);
84  free(m2p_forwardValues);
85  }
87  {
88  for (int th=0; th<mp_ImageDimensionsAndQuantification->GetNbThreadsMax(); th++)
89  {
90  if (m3p_backwardValues[th])
91  {
92  for (int tof=0; tof<m_nbTOFBins; tof++) if (m3p_backwardValues[th][tof]) free(m3p_backwardValues[th][tof]);
93  free(m3p_backwardValues[th]);
94  }
95  }
96  free(m3p_backwardValues);
97  }
98  // Free pointers to attenuation images for SPECT
100  // Delete all image update statistics tables
102  {
104  if (mp_imageStatMin) free(mp_imageStatMin);
105  if (mp_imageStatMax) free(mp_imageStatMax);
110  }
111  // Delete FOM tables
112  if (m_optimizerFOMFlag)
113  {
114  // Delete the log-likelihood table
116  {
117  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
118  {
119  if (m4p_FOMLogLikelihood[fr])
120  {
121  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
122  {
123  if (m4p_FOMLogLikelihood[fr][rg])
124  {
125  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
126  if (m4p_FOMLogLikelihood[fr][rg][cg]) free(m4p_FOMLogLikelihood[fr][rg][cg]);
127  free(m4p_FOMLogLikelihood[fr][rg]);
128  }
129  }
130  free(m4p_FOMLogLikelihood[fr]);
131  }
132  }
133  free(m4p_FOMLogLikelihood);
134  }
135  // Delete the RMSE table
136  if (m4p_FOMRMSE)
137  {
138  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
139  {
140  if (m4p_FOMRMSE[fr])
141  {
142  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
143  {
144  if (m4p_FOMRMSE[fr][rg])
145  {
146  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
147  if (m4p_FOMRMSE[fr][rg][cg]) free(m4p_FOMRMSE[fr][rg][cg]);
148  free(m4p_FOMRMSE[fr][rg]);
149  }
150  }
151  free(m4p_FOMRMSE[fr]);
152  }
153  }
154  free(m4p_FOMRMSE);
155  }
156  // Delete the nbBins table
157  if (m4p_FOMNbBins)
158  {
159  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
160  {
161  if (m4p_FOMNbBins[fr])
162  {
163  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
164  {
165  if (m4p_FOMNbBins[fr][rg])
166  {
167  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
168  if (m4p_FOMNbBins[fr][rg][cg]) free(m4p_FOMNbBins[fr][rg][cg]);
169  free(m4p_FOMNbBins[fr][rg]);
170  }
171  }
172  free(m4p_FOMNbBins[fr]);
173  }
174  }
175  free(m4p_FOMNbBins);
176  }
177  // Delete the nbData table
178  if (m4p_FOMNbData)
179  {
180  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
181  {
182  if (m4p_FOMNbData[fr])
183  {
184  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
185  {
186  if (m4p_FOMNbData[fr][rg])
187  {
188  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
189  if (m4p_FOMNbData[fr][rg][cg]) free(m4p_FOMNbData[fr][rg][cg]);
190  free(m4p_FOMNbData[fr][rg]);
191  }
192  }
193  free(m4p_FOMNbData[fr]);
194  }
195  }
196  free(m4p_FOMNbData);
197  }
198  }
199 }
200 
201 // =====================================================================
202 // ---------------------------------------------------------------------
203 // ---------------------------------------------------------------------
204 // =====================================================================
205 
207 {
208  // Call the specific help function from the children
210  // Say if the child optimizer is compatible with histogram and/or list-mode data
211  if (m_listmodeCompatibility && m_histogramCompatibility) cout << "This optimizer is compatible with both histogram and list-mode data." << endl;
212  else if (m_listmodeCompatibility) cout << "This optimizer is only compatible with list-mode data." << endl;
213  else if (m_histogramCompatibility) cout << "This optimizer is only compatible with histogram data." << endl;
214  else cout << "This optimizer is a total non-sense, not compatible with histogram nor list-mode data !!!" << endl;
215  // Say if the child optimizer is compatible with emission and/or transmission data
216  if (m_emissionCompatibility && m_transmissionCompatibility) cout << "This optimizer is compatible with both emission and transmission data." << endl;
217  else if (m_emissionCompatibility) cout << "This optimizer is only compatible with emission data." << endl;
218  else if (m_transmissionCompatibility) cout << "This optimizer is only compatible with transmission data." << endl;
219  else cout << "This optimizer is a total non-sense, not compatible with emission nor transmission data !!!" << endl;
220 }
221 
222 // =====================================================================
223 // ---------------------------------------------------------------------
224 // ---------------------------------------------------------------------
225 // =====================================================================
226 
228 {
229  // Check image dimensions
231  {
232  Cerr("***** vOptimizer::CheckParameters() -> oImageDimensionsAndQuantification is null !" << endl);
233  return 1;
234  }
235  // Check image space
236  if (mp_ImageSpace==NULL)
237  {
238  Cerr("***** vOptimizer::CheckParameters() -> oImageSpace is null !" << endl);
239  return 1;
240  }
241  // Check data mode
243  {
244  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data mode provided !" << endl);
245  return 1;
246  }
247  // Check data type
249  {
250  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data type provided !" << endl);
251  return 1;
252  }
253  // Check data spec
255  {
256  Cerr("***** vOptimizer::CheckParameters() -> No or meaningless data specificity provided (emission or transmission) !" << endl);
257  return 1;
258  }
259  // Check verbose level
260  if (m_verbose<0)
261  {
262  Cerr("***** vOptimizer::CheckParameters() -> Verbose level is negative !" << endl);
263  return 1;
264  }
265  // Listmode incompatibility
267  {
268  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with listmode data !" << endl);
269  return 1;
270  }
271  // Histogram incompatibility
273  {
274  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with histogram data !" << endl);
275  return 1;
276  }
277  // Emission compatability
279  {
280  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with emission data !" << endl);
281  return 1;
282  }
283  // Transmission compatibility
285  {
286  Cerr("***** vOptimizer::CheckParameters() -> The selected optimizer is not compatible with transmission data !" << endl);
287  return 1;
288  }
289  // Send an error if FOM computation is required and verbose is null
290  if (m_optimizerFOMFlag && m_verbose==0)
291  {
292  Cerr("***** vOptimizer::CheckParameters() -> Asked to compute FOMs in the data-space whereas the verbose is not set !" << endl);
293  return 1;
294  }
295  // Send an error if image statistics computation is required and verbose is null
297  {
298  Cerr("***** vOptimizer::CheckParameters() -> Asked to compute image update statistics whereas the verbose is not set !" << endl);
299  return 1;
300  }
301  // Send an error if FOM computation with list-mode data (it has no sense)
303  {
304  Cerr("***** vOptimizer::CheckParameters() -> Computing FOMs in the data-space with list-mode data is not possible !" << endl);
305  return 1;
306  }
307  // Call the CheckSpecificParameters function of the child
309  {
310  Cerr("***** vOptimizer::CheckParameters() -> A problem occured while checking parameters specific to the optimizer module !" << endl);
311  return 1;
312  }
313  // Normal end
314  return 0;
315 }
316 
317 // =====================================================================
318 // ---------------------------------------------------------------------
319 // ---------------------------------------------------------------------
320 // =====================================================================
321 
323 {
324  // Allocate forward buffers (as many as threads and as many as TOF bins)
327  m2p_forwardValues[th] = (FLTNB*)malloc(m_nbTOFBins*sizeof(FLTNB));
328 
329  // Allocate backward buffers (as many as threads then as many as TOF bins and as many as backward images)
331  for (int th=0; th<mp_ImageDimensionsAndQuantification->GetNbThreadsMax(); th++)
332  {
333  m3p_backwardValues[th] = (FLTNB**)malloc(m_nbTOFBins*sizeof(FLTNB*));
334  for (int tof=0; tof<m_nbTOFBins; tof++)
335  m3p_backwardValues[th][tof] = (FLTNB*)malloc(m_nbBackwardImages*sizeof(FLTNB));
336  }
337 
338  // Allocate pointers to the current attenuation image for SPECT attenuation correction
341 
342  // Allocate the thread-safe tables for image update statistics computation
344  {
352  }
353 
354  // We allocate the tables that hold the figures-of-merit if FOM flag is on
355  if (m_optimizerFOMFlag)
356  {
357  m4p_FOMNbBins = (uint64_t****)malloc(mp_ImageDimensionsAndQuantification->GetNbTimeFrames()*sizeof(uint64_t***));
361  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
362  {
363  m4p_FOMNbBins[fr] = (uint64_t***)malloc(mp_ImageDimensionsAndQuantification->GetNbRespGates()*sizeof(uint64_t**));
367  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
368  {
369  m4p_FOMNbBins[fr][rg] = (uint64_t**)malloc(mp_ImageDimensionsAndQuantification->GetNbCardGates()*sizeof(uint64_t*));
373  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
374  {
375  m4p_FOMNbBins[fr][rg][cg] = (uint64_t*)malloc(mp_ImageDimensionsAndQuantification->GetNbThreadsForProjection()*sizeof(uint64_t));
379  }
380  }
381  }
382  }
383 
384  // Call the specific initialization function of the child
385  if (InitializeSpecific())
386  {
387  Cerr("***** vOptimizer::Initialize() -> A problem occured while initializing stuff specific to the optimizer module !" << endl);
388  return 1;
389  }
390 
391  // Normal end
392  return 0;
393 }
394 
395 // =====================================================================
396 // ---------------------------------------------------------------------
397 // ---------------------------------------------------------------------
398 // =====================================================================
399 
400 int vOptimizer::PreDataUpdateStep(int a_iteration, int a_nbIterations, int a_subset, int* ap_nbSubsets)
401 {
402  // Simply reset FOMs if activated
403  if (m_optimizerFOMFlag)
404  {
405  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
406  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
407  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
409  {
410  m4p_FOMNbBins[fr][rg][cg][th] = 0;
411  m4p_FOMNbData[fr][rg][cg][th] = 0.;
412  m4p_FOMLogLikelihood[fr][rg][cg][th] = 0.;
413  m4p_FOMRMSE[fr][rg][cg][th] = 0.;
414  }
415  }
416  // Then call the specific pre-update-step function from the child optimizer
417  if (PreDataUpdateSpecificStep(a_iteration, a_nbIterations, a_subset, ap_nbSubsets))
418  {
419  Cerr("***** vOptimizer::PreDataUpdateStep() -> A problem occured while applying the specific pre-data-update step of the optimizer module !" << endl);
420  return 1;
421  }
422  // Normal end
423  return 0;
424 }
425 
426 // =====================================================================
427 // ---------------------------------------------------------------------
428 // ---------------------------------------------------------------------
429 // =====================================================================
430 
431 int vOptimizer::PreDataUpdateSpecificStep(int a_iteration, int a_nbIterations, int a_subset, int* ap_nbSubsets)
432 {
433  // By default, just do nothing here. This function is designed to be overloaded by the specific optimizer if needed
434  return 0;
435 }
436 
437 // =====================================================================
438 // ---------------------------------------------------------------------
439 // ---------------------------------------------------------------------
440 // =====================================================================
441 
442 int vOptimizer::PostDataUpdateStep(int a_iteration, int a_nbIterations, int a_subset, int* ap_nbSubsets)
443 {
444  // Print out the FOMs if activated
445  if (m_optimizerFOMFlag)
446  {
447  // Verbose
448  Cout("vOptimizer::PostDataUpdateStep() -> Data-space figures-of-merit" << endl);
449  // Number of spaces for printing
450  string spaces_fr = ""; if (mp_ImageDimensionsAndQuantification->GetNbTimeFrames()>1) spaces_fr += " ";
451  string spaces_rg = ""; if (mp_ImageDimensionsAndQuantification->GetNbRespGates()>1) spaces_rg += " ";
452  string spaces_cg = ""; if (mp_ImageDimensionsAndQuantification->GetNbCardGates()>1) spaces_cg += " ";
453  // Loop on frames and gates
454  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
455  {
456  // Verbose
458  Cout(spaces_fr << "Frame " << fr+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbTimeFrames() << endl);
459  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
460  {
461  // Verbose
463  Cout(spaces_fr << spaces_rg << "Respiratory gate " << rg+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbRespGates() << endl);
464  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
465  {
466  // Verbose
468  Cout(spaces_fr << spaces_rg << spaces_cg << "Cardiac gate " << cg+1 << " on " << mp_ImageDimensionsAndQuantification->GetNbCardGates() << endl);
469  // Merge thread results
471  {
472  m4p_FOMNbBins[fr][rg][cg][0] += m4p_FOMNbBins[fr][rg][cg][th];
473  m4p_FOMNbData[fr][rg][cg][0] += m4p_FOMNbData[fr][rg][cg][th];
474  m4p_FOMLogLikelihood[fr][rg][cg][0] += m4p_FOMLogLikelihood[fr][rg][cg][th];
475  m4p_FOMRMSE[fr][rg][cg][0] += m4p_FOMRMSE[fr][rg][cg][th];
476  }
477  // Finish mean number of data counts per channel
478  m4p_FOMNbData[fr][rg][cg][0] /= ((HPFLTNB)(m4p_FOMNbBins[fr][rg][cg][0]));
479  // Finish RMSE computation
480  m4p_FOMRMSE[fr][rg][cg][0] /= ((FLTNB)(m4p_FOMNbBins[fr][rg][cg][0]));
481  m4p_FOMRMSE[fr][rg][cg][0] = sqrt(m4p_FOMRMSE[fr][rg][cg][0]);
482  // Print out
483  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Number of data channels used in optimization: " << m4p_FOMNbBins[fr][rg][cg][0] << endl);
484  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Mean number of data counts per channel: " << m4p_FOMNbData[fr][rg][cg][0] << endl);
485  Cout(spaces_fr << spaces_rg << spaces_cg << " --> Log-likelihood: " << m4p_FOMLogLikelihood[fr][rg][cg][0] << endl);
486  Cout(spaces_fr << spaces_rg << spaces_cg << " --> RMSE: " << m4p_FOMRMSE[fr][rg][cg][0] << endl);
487  }
488  }
489  }
490  }
491  // Then call the specific post-update-step function from the child optimizer
492  if (PostDataUpdateSpecificStep(a_iteration, a_nbIterations, a_subset, ap_nbSubsets))
493  {
494  Cerr("***** vOptimizer::PostDataUpdateStep() -> A problem occured while applying the specific post-data-update step of the optimizer module !" << endl);
495  return 1;
496  }
497  // Normal end
498  return 0;
499 }
500 
501 // =====================================================================
502 // ---------------------------------------------------------------------
503 // ---------------------------------------------------------------------
504 // =====================================================================
505 
506 int vOptimizer::PostDataUpdateSpecificStep(int a_iteration, int a_nbIterations, int a_subset, int* ap_nbSubsets)
507 {
508  // By default, just do nothing here. This function is designed to be overloaded by the specific optimizer if needed
509 
510 
511  return 0;
512 }
513 
514 // =====================================================================
515 // ---------------------------------------------------------------------
516 // ---------------------------------------------------------------------
517 // =====================================================================
518 
520  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
521  int a_th )
522 {
523  // =======================================================================
524  // 4D forward projection including framing, respiratory and cardiac gating
525  // =======================================================================
526 
527  // Clean buffers
528  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
529  {
530  m2p_forwardValues[a_th][b] = 0.;
531  for (int img=0; img<m_nbBackwardImages; img++) m3p_backwardValues[a_th][b][img] = 0.;
532  }
533 
534  // Loop over time basis functions
535  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
536  {
537  // Get frame/basis coefficient and continue if null
538  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(tbf,a_timeFrame);
539  if (time_basis_coef==0.) continue;
540  // Loop over respiratory basis functions
541  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
542  {
543  // Get resp_gate/basis coefficient and continue if null
544  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(rbf,a_respGate);
545  if (resp_basis_coef==0.) continue;
546  // Loop over cardiac basis functions
547  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
548  {
549  // Get card_gate_basis coefficient and continue if null
550  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(cbf,a_cardGate);
551  if (card_basis_coef==0.) continue;
552  // Compute global coefficient
553  FLTNB global_basis_coef = time_basis_coef * resp_basis_coef * card_basis_coef;
554  // Loop over TOF bins
555  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
556  {
557  // Set the current TOF bin
558  ap_Line->SetCurrentTOFBin(b);
559  // Project the image and apply dynamic coefficient converions
560  m2p_forwardValues[a_th][b] += ForwardProject(ap_Line,mp_ImageSpace->m4p_forwardImage[tbf][rbf][cbf]) * global_basis_coef;
561  }
562  }
563  }
564  }
565 
566  // ==============================================
567  // Differentiate emission and transmission models
568  // ==============================================
569 
570  // -----------------
571  // Case for emission
572  // -----------------
574  {
575  // Add additive terms (we multiply by the frame duration because the additive terms are provided as rates)
576  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);
577  }
578  // ---------------------
579  // Case for transmission
580  // ---------------------
581  else if (m_dataSpec==SPEC_TRANSMISSION)
582  {
583  int no_tof = 0;
584  // Ignore TOF as it has no sense; add scatter rate multiplied by the frame duration; take blank into account
585  // Note that the image is supposed to be in cm-1 as it is the standard
586  m2p_forwardValues[a_th][no_tof] = ap_Event->GetBlankValue() * exp(-m2p_forwardValues[a_th][no_tof]*0.1)
587  + ap_Event->GetAdditiveCorrections(no_tof) * mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame);
588  }
589 
590  // End
591  return 0;
592 }
593 
594 // =====================================================================
595 // ---------------------------------------------------------------------
596 // ---------------------------------------------------------------------
597 // =====================================================================
598 
600  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
601  int a_iteration, int a_th )
602 {
603  // Deliberately do nothing!
604  return 0;
605 }
606 
607 // =====================================================================
608 // ---------------------------------------------------------------------
609 // ---------------------------------------------------------------------
610 // =====================================================================
611 
613  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
614  int a_th )
615 {
616  // Loop over TOF bins
617  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
618  {
619  // Set the current TOF bin of the projection-line
620  ap_Line->SetCurrentTOFBin(b);
621  // The weight associated to the line to be backprojected
622  FLTNB weight = 0.;
623  // Call the function to compute specific sensitivity
625  (
626  ap_Event->GetEventValue(b), // the measured data
627  m2p_forwardValues[a_th][b], // the forward model
628  &weight, // the sensitivity weight
629  ap_Event->GetMultiplicativeCorrections(), // the multiplicative corrections
630  ap_Event->GetAdditiveCorrections(b)*mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame), // the additive corrections
631  ap_Event->GetBlankValue(), // the blank value
632  mp_ImageDimensionsAndQuantification->GetQuantificationFactor(a_bed,a_timeFrame, a_respGate, a_cardGate), // the quantification factor
633  ap_Line // the projection line
634  ))
635  {
636  Cerr("***** vOptimizer::DataStep3BackwardProjectSensitivity() -> A problem occured while performing the sensitivity specific operations !" << endl);
637  return 1;
638  }
639  // Back-project the weight into the sensitivity matrix
640  BackwardProject(ap_Line, mp_ImageSpace->m5p_sensitivity[a_th][a_timeFrame][a_respGate][a_cardGate], weight);
641  }
642 
643  // End
644  return 0;
645 }
646 
647 // =====================================================================
648 // ---------------------------------------------------------------------
649 // ---------------------------------------------------------------------
650 // =====================================================================
651 
653  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
654  int a_iteration, int a_th )
655 {
656  // Deliberately do nothing!
657  return 0;
658 }
659 
660 // =====================================================================
661 // ---------------------------------------------------------------------
662 // ---------------------------------------------------------------------
663 // =====================================================================
664 
666  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
667  int a_th )
668 {
669  // Loop on TOF bins
670  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
671  {
672  // Set the current TOF bin of the projection-line
673  ap_Line->SetCurrentTOFBin(b);
674  // Then call the pure virtual function for specific data space operations
676  (
677  ap_Event->GetEventValue(b), // the measured data
678  m2p_forwardValues[a_th][b], // the forward model
679  m3p_backwardValues[a_th][b], // the backward values (the result)
680  ap_Event->GetMultiplicativeCorrections(), // the multiplicative corrections
681  ap_Event->GetAdditiveCorrections(b)*mp_ImageDimensionsAndQuantification->GetFrameDurationInSec(a_bed, a_timeFrame), // the additive corrections
682  ap_Event->GetBlankValue(), // the blank value
683  mp_ImageDimensionsAndQuantification->GetQuantificationFactor(a_bed,a_timeFrame, a_respGate, a_cardGate), // the quantification factor
684  ap_Line // the projection line
685  ))
686  {
687  Cerr("***** vOptimizer::DataStep5ComputeCorrections() -> A problem occured while performing specific data space operations !" << endl);
688  return 1;
689  }
690  }
691 
692  // End
693  return 0;
694 }
695 
696 // =====================================================================
697 // ---------------------------------------------------------------------
698 // ---------------------------------------------------------------------
699 // =====================================================================
700 
702  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
703  int a_iteration, int a_th )
704 {
705  // Deliberately do nothing!
706  return 0;
707 }
708 
709 // =====================================================================
710 // ---------------------------------------------------------------------
711 // ---------------------------------------------------------------------
712 // =====================================================================
713 
715  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
716  int a_th )
717 {
718  // ========================================================================
719  // 4D backward projection including framing, respiratory and cardiac gating
720  // ========================================================================
721 
722  // Loop over time basis functions
723  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
724  {
725  // Get frame/basis coefficient and continue if null
726  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(tbf,a_timeFrame);
727  if (time_basis_coef==0.) continue;
728  // Loop over respiratory basis functions
729  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
730  {
731  // Get resp_gate/basis coefficient and continue if null
732  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(rbf,a_respGate);
733  if (resp_basis_coef==0.) continue;
734  // Loop over cardiac basis functions
735  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
736  {
737  // Get card_gate_basis coefficient and continue if null
738  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(cbf,a_cardGate);
739  if (card_basis_coef==0.) continue;
740  // Compute global coefficient
741  FLTNB global_basis_coef = time_basis_coef * resp_basis_coef * card_basis_coef;
742  // Loop over TOF bins
743  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
744  {
745  // Set the current TOF bin of the projection-line
746  ap_Line->SetCurrentTOFBin(b);
747  // Loop over backprojection images
748  for (int img=0; img<m_nbBackwardImages; img++)
749  {
750  // Continue if backproj is null or infinity
751  if (m3p_backwardValues[a_th][b][img]==0. || !isfinite(m3p_backwardValues[a_th][b][img])) continue;
752  // Backward project into the backward image
753  BackwardProject( ap_Line, mp_ImageSpace->m6p_backwardImage[img][a_th][tbf][rbf][cbf] , m3p_backwardValues[a_th][b][img] * global_basis_coef );
754  }
755  }
756  }
757  }
758  }
759 
760  // End
761  return 0;
762 }
763 
764 // =====================================================================
765 // ---------------------------------------------------------------------
766 // ---------------------------------------------------------------------
767 // =====================================================================
768 
770  int a_timeFrame, int a_respGate, int a_cardGate,
771  int a_thread )
772 {
773  // Loop on TOF bins
774  for (int b=0; b<ap_Line->GetNbTOFBins(); b++)
775  {
776  // Get the model and the data
777  FLTNB model = m2p_forwardValues[a_thread][b];
778  FLTNB data = ap_Event->GetEventValue(b);
779  // Update log likelihood only if model is strictly positive (to avoid numerical issues, we use a threshold at e-10)
780  if (model>1.e-10) m4p_FOMLogLikelihood[a_timeFrame][a_respGate][a_cardGate][a_thread] += data*log(model)-model;
781  // Update RMSE
782  m4p_FOMRMSE[a_timeFrame][a_respGate][a_cardGate][a_thread] += (data-model)*(data-model);
783  // Update number of data channels
784  m4p_FOMNbBins[a_timeFrame][a_respGate][a_cardGate][a_thread]++;
785  // Update the number of data counts
786  m4p_FOMNbData[a_timeFrame][a_respGate][a_cardGate][a_thread] += ((HPFLTNB)data);
787  }
788  // Normal end
789  return 0;
790 }
791 
792 // =====================================================================
793 // ---------------------------------------------------------------------
794 // ---------------------------------------------------------------------
795 // =====================================================================
796 
798 {
799  // Verbose
800  if (m_verbose>=3) Cout("vOptimizer::UpdateVisitedVoxels() -> Tick visited voxels based on sensitivity" << endl);
801 
802  // We loop on the sensitivity for the first frame/gates and look if voxels were 'visited' by LORs
803  int thread0 = 0;
804  //int frame0 = 0;
805  //int resp_gate0 = 0;
806  //int card_gate0 = 0;
807 
808  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
809  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
810  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
811  for (int v=0; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ(); v++)
812  // If the voxel has been visited, we add 1. to the mp_Image->mp_visitedVoxelsImage
813  //if (mp_ImageSpace->m5p_sensitivity[thread0][fr][rg][cg][v]!=0.) mp_ImageSpace->mp_visitedVoxelsImage[v] += 1.;
814  if (mp_ImageSpace->m5p_sensitivity[thread0][fr][rg][cg][v]!=0.) mp_ImageSpace->m4p_visitedVoxelsImage[fr][rg][cg][v] += 1.;
815 
816  // End
817  return 0;
818 }
819 
820 // =====================================================================
821 // ---------------------------------------------------------------------
822 // ---------------------------------------------------------------------
823 // =====================================================================
824 
825 int vOptimizer::ImageUpdateStep( int a_iteration, int a_nbSubsets )
826 {
827  // Verbose
828  if (m_verbose>=2 || m_optimizerImageStatFlag) Cout("vOptimizer::ImageUpdateStep() -> Start image update" << endl);
829  // Number of spaces for nice printing
830  string spaces_tbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions()>1) spaces_tbf += " ";
831  string spaces_rbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions()>1) spaces_rbf += " ";
832  string spaces_cbf = ""; if (mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions()>1) spaces_cbf += " ";
833  // Loop over time basis functions
834  for (int tbf=0; tbf<mp_ImageDimensionsAndQuantification->GetNbTimeBasisFunctions(); tbf++)
835  {
836  // Verbose
838  {
839  if (mp_ImageDimensionsAndQuantification->GetTimeStaticFlag()) Cout(spaces_tbf << "--> Time frame " << tbf+1 << endl);
840  else Cout(spaces_tbf << "--> Time basis function: " << tbf+1 << endl);
841  }
842  // Loop over respiratory basis functions
843  for (int rbf=0; rbf<mp_ImageDimensionsAndQuantification->GetNbRespBasisFunctions(); rbf++)
844  {
845  // Verbose
847  {
848  if (mp_ImageDimensionsAndQuantification->GetRespStaticFlag()) Cout(spaces_tbf << spaces_rbf << "--> Respiratory gate " << rbf+1 << endl);
849  else Cout(spaces_tbf << spaces_rbf << "--> Respiratory basis function: " << rbf+1 << endl);
850  }
851  // Loop over cardiac basis functions
852  for (int cbf=0; cbf<mp_ImageDimensionsAndQuantification->GetNbCardBasisFunctions(); cbf++)
853  {
854  // Verbose
856  {
857  if (mp_ImageDimensionsAndQuantification->GetCardStaticFlag()) Cout(spaces_tbf << spaces_rbf << spaces_cbf << "--> Cardiac gate " << cbf+1 << endl);
858  else Cout(spaces_tbf << spaces_rbf << spaces_cbf << "--> Cardiac basis function: " << cbf+1 << endl);
859  }
860  // Set the number of threads for image computation
861  #ifdef CASTOR_OMP
863  #endif
864  // Reset counter for image stats
866  {
868  {
869  mp_imageStatNbVox[th] = 0;
870  mp_imageStatMin[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][0];
871  mp_imageStatMax[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][0];
872  mp_imageStatMean[th] = 0.;
873  mp_imageStatVariance[th] = 0.;
874  mp_correctionStatMean[th] = 0.;
875  mp_correctionStatVariance[th] = 0.;
876  }
877  }
878  // OMP loop over voxels
879  INTNB v;
880  bool error = false;
881  #pragma omp parallel for private(v) schedule(guided)
882  for (v=0 ; v<mp_ImageDimensionsAndQuantification->GetNbVoxXYZ() ; v++)
883  {
884  // Get the thread index
885  int th = 0;
886  #ifdef CASTOR_OMP
887  th = omp_get_thread_num();
888  #endif
889  // Compute sensitivity
890  FLTNB sensitivity = ComputeSensitivity(mp_ImageSpace->m5p_sensitivity[0],tbf,rbf,cbf,v,a_nbSubsets);
891  // If the sensitivity is negative or null, we skip this voxel
892  if (sensitivity<=0.) continue;
893  // Keep the current image value in a buffer for update statistics
894  FLTNB old_image_value = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
895  // Fill class-local buffer of backward values (thread-safe)
896  for (int nb=0; nb<m_nbBackwardImages; nb++)
897  m3p_backwardValues[th][0][nb] = mp_ImageSpace->m6p_backwardImage[nb][0][tbf][rbf][cbf][v];
898  // Then call the pure virtual function for specific data space operations
900  (
901  mp_ImageSpace->m4p_image[tbf][rbf][cbf][v],
902  // 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
903  // to perform MLAA reconstruction but using PET basis functions different from diracs for example... In this case, we would get
904  // segmentation faults. But actually, this function would be overloaded by MLAA anyway...
905  //mp_Image->m4p_attenuation != NULL ? mp_Image->m4p_attenuation[tbf][rbf][cbf][v] : 0. ,
906  &mp_ImageSpace->m4p_image[tbf][rbf][cbf][v],
907  //mp_Image->m4p_attenuation != NULL ? &mp_Image->m4p_attenuation[tbf][rbf][cbf][v] : NULL ,
908  sensitivity,
909  m3p_backwardValues[th][0],
910  v
911  ))
912  {
913  Cerr("***** vOptimizer::ImageUpdateStep() -> A problem occured in while performing image space specific operations of the optimizer !" << endl);
914  error = true;
915  }
916  // Do some image statistics
918  {
919  if (mp_imageStatMin[th]>mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]) mp_imageStatMin[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
920  if (mp_imageStatMax[th]<mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]) mp_imageStatMax[th] = mp_ImageSpace->m4p_image[tbf][rbf][cbf][v];
921  // The one pass algorithm is used to compute the variance here for both the image and correction step.
922  // Do not change anyhting to the order of the operations !
923  mp_imageStatNbVox[th]++;
924  HPFLTNB sample = ((HPFLTNB)(mp_ImageSpace->m4p_image[tbf][rbf][cbf][v]));
925  HPFLTNB delta = sample - mp_imageStatMean[th];
926  mp_imageStatMean[th] += delta / ((HPFLTNB)(mp_imageStatNbVox[th]));
927  mp_imageStatVariance[th] += delta * (sample - mp_imageStatMean[th]);
928  // Same one-pass variance algorithm for the correction step
929  sample = ((HPFLTNB)(mp_ImageSpace->m4p_image[tbf][rbf][cbf][v] - old_image_value));
930  delta = sample - mp_correctionStatMean[th];
931  mp_correctionStatMean[th] += delta / ((HPFLTNB)(mp_imageStatNbVox[th]));
932  mp_correctionStatVariance[th] += delta * (sample - mp_correctionStatMean[th]);
933  }
934  }
935  // Check for errors (we cannot stop the multi-threaded loop so we do it here)
936  if (error) return 1;
937  // Finish image update statistics computations: these operations are specific to the merging of the
938  // multi-threaded one-pass variance estimations used above.
940  {
942  {
943  // Image minimum and maximum
946  // Cast the counts into HPFLTNB
947  HPFLTNB count1 = ((HPFLTNB)(mp_imageStatNbVox[0]));
948  HPFLTNB count2 = ((HPFLTNB)(mp_imageStatNbVox[th]));
949  // Update the count
951  HPFLTNB count12 = ((HPFLTNB)(mp_imageStatNbVox[0]));
952  // Compute the delta mean
953  HPFLTNB delta = mp_imageStatMean[0] - mp_imageStatMean[th];
954  // Update the variance
955  mp_imageStatVariance[0] += mp_imageStatVariance[th] + delta * delta * count1 * count2 / count12;
956  // Update the mean
957  mp_imageStatMean[0] = (count1*mp_imageStatMean[0] + count2*mp_imageStatMean[th]) / count12;
958  // Do the same for the correction step
960  mp_correctionStatVariance[0] += mp_correctionStatVariance[th] + delta * delta * count1 * count2 / count12;
961  mp_correctionStatMean[0] = (count1*mp_correctionStatMean[0] + count2*mp_correctionStatMean[th]) / count12;
962  }
963  // Finish variance computation
966  // Print out
967  Cout(spaces_tbf << spaces_rbf << spaces_cbf << " --> Image update step "
968  << " | mean: " << mp_correctionStatMean[0]
969  << " | stdv: " << sqrt(mp_correctionStatVariance[0]) << endl);
970  Cout(spaces_tbf << spaces_rbf << spaces_cbf << " --> New image estimate"
971  << " | mean: " << mp_imageStatMean[0]
972  << " | stdv: " << sqrt(mp_imageStatVariance[0])
973  << " | min/max: " << mp_imageStatMin[0] << "/" << mp_imageStatMax[0] << endl);
974  }
975  }
976  }
977  }
978 
979  // End
980  return 0;
981 }
982 
983 // =====================================================================
984 // ---------------------------------------------------------------------
985 // ---------------------------------------------------------------------
986 // =====================================================================
987 
989 {
990  // Particular case for SPECT with attenuation correction
991  if (m_dataType==TYPE_SPECT && m2p_attenuationImage[ap_Line->GetThreadNumber()]!=NULL)
992  return ap_Line->ForwardProjectWithSPECTAttenuation( m2p_attenuationImage[ap_Line->GetThreadNumber()], ap_image );
993  // Otherwise call the function without attenuation
994  else
995  return ap_Line->ForwardProject( ap_image );
996 }
997 
998 // =====================================================================
999 // ---------------------------------------------------------------------
1000 // ---------------------------------------------------------------------
1001 // =====================================================================
1002 
1003 void vOptimizer::BackwardProject(oProjectionLine* ap_Line, FLTNB* ap_image, FLTNB a_value)
1004 {
1005  // Particular case for SPECT with attenuation correction
1006  if (m_dataType==TYPE_SPECT && m2p_attenuationImage[ap_Line->GetThreadNumber()]!=NULL)
1007  ap_Line->BackwardProjectWithSPECTAttenuation( m2p_attenuationImage[ap_Line->GetThreadNumber()], ap_image, a_value );
1008  // Otherwise call the function without attenuation
1009  else
1010  ap_Line->BackwardProject( ap_image, a_value );
1011 }
1012 
1013 // =====================================================================
1014 // ---------------------------------------------------------------------
1015 // ---------------------------------------------------------------------
1016 // =====================================================================
1017 
1018 FLTNB vOptimizer::ComputeSensitivity( FLTNB**** a4p_sensitivityMatrix,
1019  int a_timeBasisFunction, int a_respBasisFunction, int a_cardBasisFunction,
1020  int a_voxel, int a_nbSubsets )
1021 {
1022  // The sensitivity to be computed
1023  FLTNB sensitivity = 0.;
1024  // Loop over time frames
1025  for (int fr=0; fr<mp_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
1026  {
1027  // Get frame/basis coefficient and continue if null
1028  FLTNB time_basis_coef = mp_ImageDimensionsAndQuantification->GetTimeBasisCoefficient(a_timeBasisFunction,fr);
1029  if (time_basis_coef==0.) continue;
1030  // Loop over respiratory gates
1031  for (int rg=0; rg<mp_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
1032  {
1033  // Get resp_gate/basis coefficient and continue if null
1034  FLTNB resp_basis_coef = mp_ImageDimensionsAndQuantification->GetRespBasisCoefficient(a_respBasisFunction,rg);
1035  if (resp_basis_coef==0.) continue;
1036  // Loop over cardiac gates
1037  for (int cg=0; cg<mp_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
1038  {
1039  // Get card_gate_basis coefficient and continue if null
1040  FLTNB card_basis_coef = mp_ImageDimensionsAndQuantification->GetCardBasisCoefficient(a_cardBasisFunction,cg);
1041  if (card_basis_coef==0.) continue;
1042  // Add actual weight contribution
1043  sensitivity += a4p_sensitivityMatrix[fr][rg][cg][a_voxel] *
1044  time_basis_coef * resp_basis_coef * card_basis_coef;
1045  }
1046  }
1047  }
1048  // If listmode, then divide by the number of subsets
1049  if (m_dataMode==MODE_LIST) sensitivity /= ((FLTNB)a_nbSubsets);
1050  // Return sensitivity
1051  return sensitivity;
1052 }
1053 
1054 // =====================================================================
1055 // ---------------------------------------------------------------------
1056 // ---------------------------------------------------------------------
1057 // =====================================================================
FLTNB **** m4p_forwardImage
Definition: oImageSpace.hh:88
FLTNB * mp_imageStatMax
Definition: vOptimizer.hh:641
virtual int ImageUpdateStep(int a_iteration, int a_nbSubsets)
A public function used to perform the image update step of the optimizer.
Definition: vOptimizer.cc:825
#define MODE_HISTOGRAM
Definition: vDataFile.hh:59
#define TYPE_PET
Definition: vDataFile.hh:74
virtual FLTNB GetBlankValue()
This is a pure virtual function implemented in the child classes.
Definition: vEvent.hh:197
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:72
#define MODE_LIST
Definition: vDataFile.hh:57
int GetNbCardBasisFunctions()
Get the number of cardiac basis functions.
#define FLTNB
Definition: gVariables.hh:81
FLTNB m_initialValue
Definition: vOptimizer.hh:619
virtual FLTNB GetAdditiveCorrections(int a_bin)=0
Pure virtual function implemented in the child classes.
bool m_listmodeCompatibility
Definition: vOptimizer.hh:620
virtual int DataStep6Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_iteration, int a_thread)
A public function which does nothing but being virtual.
Definition: vOptimizer.cc:701
virtual int PreDataUpdateSpecificStep(int a_iteration, int a_nbIterations, int a_subset, int *ap_nbSubsets)
A private function used to perform any step required by the child optimizer, before the loop on event...
Definition: vOptimizer.cc:431
FLTNB * mp_imageStatMin
Definition: vOptimizer.hh:640
#define HPFLTNB
Definition: gVariables.hh:83
HPFLTNB **** m4p_FOMNbData
Definition: vOptimizer.hh:636
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:625
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:38
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:988
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.
virtual int DataStep4Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_iteration, int a_thread)
A public function which does nothing but being virtual.
Definition: vOptimizer.cc:652
uint64_t **** m4p_FOMNbBins
Definition: vOptimizer.hh:635
int PostDataUpdateStep(int a_iteration, int a_nbIterations, int a_subset, int *ap_nbSubsets)
A public function used to do stuff that need to be done at the beginning of a subset (before the data...
Definition: vOptimizer.cc:442
#define TYPE_SPECT
Definition: vDataFile.hh:76
HPFLTNB * mp_correctionStatMean
Definition: vOptimizer.hh:644
#define TYPE_CT
Definition: vDataFile.hh:78
FLTNB **** m4p_FOMRMSE
Definition: vOptimizer.hh:634
int m_dataMode
Definition: vOptimizer.hh:627
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:629
bool m_optimizerFOMFlag
Definition: vOptimizer.hh:632
bool GetRespStaticFlag()
Get the respiratory static flag that says if the reconstruction has only one respiratory gate or not...
INTNB * mp_imageStatNbVox
Definition: vOptimizer.hh:639
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)
FLTNB ****** m6p_backwardImage
Definition: oImageSpace.hh:95
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:622
FLTNB **** m4p_image
Definition: oImageSpace.hh:81
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_nbBackwardImages
Definition: vOptimizer.hh:615
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:612
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:1003
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:519
Declaration of class vOptimizer.
FLTNB ** m2p_attenuationImage
Definition: vOptimizer.hh:630
bool m_histogramCompatibility
Definition: vOptimizer.hh:621
bool m_optimizerImageStatFlag
Definition: vOptimizer.hh:638
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:769
#define SPEC_UNKNOWN
Definition: vDataFile.hh:89
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.
FLTNB ComputeSensitivity(FLTNB ****a4p_sensitivityImage, int a_timeBasisFunction, int a_respBasisFunction, int a_cardBasisFunction, int a_voxel, int a_nbSubsets)
A function used to compute the sensitivity of a given voxel and a given set of dynamic basis function...
Definition: vOptimizer.cc:1018
HPFLTNB * mp_imageStatVariance
Definition: vOptimizer.hh:643
#define SPEC_TRANSMISSION
Definition: vDataFile.hh:93
virtual ~vOptimizer()
The destructor of vOptimizer.
Definition: vOptimizer.cc:77
#define INTNB
Definition: gVariables.hh:92
bool m_transmissionCompatibility
Definition: vOptimizer.hh:623
int m_nbTOFBins
Definition: vOptimizer.hh:616
This class is designed to manage and store system matrix elements associated to a vEvent...
int GetNbCardGates()
Get the number of cardiac gates.
FLTNB **** m4p_FOMLogLikelihood
Definition: vOptimizer.hh:633
int Initialize()
A public function used to initialize the optimizer.
Definition: vOptimizer.cc:322
oImageSpace * mp_ImageSpace
Definition: vOptimizer.hh:626
int m_dataType
Definition: vOptimizer.hh:628
Mother class for the Event objects.
Definition: vEvent.hh:43
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:665
#define SPEC_EMISSION
Definition: vDataFile.hh:91
int GetNbTOFBins()
This function is used to get the number of TOF bins in use.
int PreDataUpdateStep(int a_iteration, int a_nbIterations, int a_subset, int *ap_nbSubsets)
A public function used to do stuff that need to be done at the beginning of a subset (before the data...
Definition: vOptimizer.cc:400
HPFLTNB * mp_imageStatMean
Definition: vOptimizer.hh:642
int GetNbThreadsForProjection()
Get the number of threads used for projections.
virtual int DataStep2Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_iteration, int a_thread)
A public function which does nothing but being virtual.
Definition: vOptimizer.cc:599
int UpdateVisitedVoxels()
A public function used to update the 'visited' voxels after each subset.
Definition: vOptimizer.cc:797
virtual int PostDataUpdateSpecificStep(int a_iteration, int a_nbIterations, int a_subset, int *ap_nbSubsets)
A private function used to perform any step required by the child optimizer, after the loop on event ...
Definition: vOptimizer.cc:506
#define MODE_UNKNOWN
Definition: vDataFile.hh:55
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.
virtual int ImageSpaceSpecificOperations(FLTNB a_currentImageValue, FLTNB *ap_newImageValue, FLTNB a_sensitivity, FLTNB *ap_correctionValues, INTNB a_voxel)=0
A private function used to update the image value from the provided data.
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:714
void ShowHelp()
A function used to show help about the optimizer.
Definition: vOptimizer.cc:206
HPFLTNB * mp_correctionStatVariance
Definition: vOptimizer.hh:645
FLTNB ** m2p_forwardValues
Definition: vOptimizer.hh:617
int GetNbRespBasisFunctions()
Get the number of respiratory basis functions.
FLTNB ***** m5p_sensitivity
Definition: oImageSpace.hh:105
FLTNB *** m3p_backwardValues
Definition: vOptimizer.hh:618
int CheckParameters()
A public function used to check the parameters settings.
Definition: vOptimizer.cc:227
FLTNB **** m4p_visitedVoxelsImage
Definition: oImageSpace.hh:128
FLTNB GetCardBasisCoefficient(int a_cardBasisFunction, int a_cardGate)
Get the cardiac basis coefficients for the provided cardiac gate and basis function.