CASToR  3.2
Tomographic Reconstruction (PET/SPECT/CT)
code/include/optimizer/vOptimizer.hh
Go to the documentation of this file.
1 
8 #ifndef VOPTIMIZER_HH
9 #define VOPTIMIZER_HH 1
10 
11 #include "gVariables.hh"
12 #include "oImageDimensionsAndQuantification.hh"
13 #include "oImageSpace.hh"
14 #include "vDataFile.hh"
15 #include "oProjectionLine.hh"
16 #include "vPenalty.hh"
17 
38 {
39  // -------------------------------------------------------------------
40  // Constructor & Destructor
41  public:
48  vOptimizer();
56  virtual ~vOptimizer();
57 
58 
59  // -------------------------------------------------------------------
60  // Public member functions
61  public:
67  void ShowHelp();
76  int CheckParameters();
85  int Initialize();
96  int UpdateVisitedVoxels();
106  int PreDataUpdateStep();
116  int PreImageUpdateStep();
117 
118 
119  // -------------------------------------------------------------------
120  // Virtual public member functions that may be overloaded by specific child optimizers
121  public:
141  virtual int DataStep1ForwardProjectModel( oProjectionLine* ap_Line, vEvent* ap_Event,
142  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
143  int a_thread );
158  virtual int DataStep2Optional( oProjectionLine* ap_Line, vEvent* ap_Event,
159  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
160  int a_thread );
178  virtual int DataStep3BackwardProjectSensitivity( oProjectionLine* ap_Line, vEvent* ap_Event,
179  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
180  int a_thread );
195  virtual int DataStep4Optional( oProjectionLine* ap_Line, vEvent* ap_Event,
196  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
197  int a_thread );
215  virtual int DataStep5ComputeCorrections( oProjectionLine* ap_Line, vEvent* ap_Event,
216  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
217  int a_thread );
232  virtual int DataStep6Optional( oProjectionLine* ap_Line, vEvent* ap_Event,
233  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
234  int a_thread );
253  virtual int DataStep7BackwardProjectCorrections( oProjectionLine* ap_Line, vEvent* ap_Event,
254  int a_bed, int a_timeFrame, int a_respGate, int a_cardGate,
255  int a_thread );
270  virtual int DataStep8ComputeFOM( oProjectionLine* ap_Line, vEvent* ap_Event,
271  int a_timeFrame, int a_respGate, int a_cardGate,
272  int a_thread );
283  virtual int ImageUpdateStep();
284 
285 
286  // -----------------------------------------------------------------------------------------
287  // Virtual private member functions that may be overloaded by specific child optimizers
288  private:
297  virtual int PreDataUpdateSpecificStep();
306  virtual int PreImageUpdateSpecificStep();
307 
308  // -----------------------------------------------------------------------------------------
309  // Pure virtual public member functions that need to be implemented by child optimizers
310  public:
321  virtual int ReadConfigurationFile( const string& a_configurationFile ) = 0;
332  virtual int ReadOptionsList( const string& a_optionsList ) = 0;
333 
334 
335  // -----------------------------------------------------------------------------------------
336  // Pure virtual private member functions that need to be implemented by child optimizers
337  private:
346  virtual void ShowHelpSpecific() = 0;
355  virtual int CheckSpecificParameters() = 0;
365  virtual int InitializeSpecific() = 0;
384  virtual int SensitivitySpecificOperations( FLTNB a_data, FLTNB a_forwardModel, FLTNB* ap_weight,
385  FLTNB a_multiplicativeCorrections, FLTNB a_additiveCorrections, FLTNB a_blankValue,
386  FLTNB a_quantificationFactor, oProjectionLine* ap_Line ) = 0;
406  virtual int DataSpaceSpecificOperations( FLTNB a_data, FLTNB a_forwardModel, FLTNB* ap_backwardValues,
407  FLTNB a_multiplicativeCorrections, FLTNB a_additiveCorrections, FLTNB a_blankValue,
408  FLTNB a_quantificationFactor, oProjectionLine* ap_Line ) = 0;
427  virtual int ImageSpaceSpecificOperations( FLTNB a_currentImageValue, FLTNB* ap_newImageValue,
428  FLTNB a_sensitivity, FLTNB* ap_correctionValues, INTNB a_voxel,
429  int tbf = -1, int rbf = -1, int cbf = -1 ) = 0;
430 
431 
432  // -------------------------------------------------------------------
433  // Public Get & Set functions
434  public:
440  inline void SetVerbose(int a_verbose)
441  {m_verbose = a_verbose;}
447  inline void SetImageDimensionsAndQuantification(oImageDimensionsAndQuantification* ap_ImageDimensionsAndQuantification)
448  {mp_ImageDimensionsAndQuantification = ap_ImageDimensionsAndQuantification;}
454  inline void SetImageSpace(oImageSpace* ap_ImageSpace)
455  {mp_ImageSpace = ap_ImageSpace;}
461  inline void SetNbTOFBins(int a_nbTOFBins)
462  {m_nbTOFBins = a_nbTOFBins;}
468  inline void SetDataMode(int a_dataMode)
469  {m_dataMode = a_dataMode;}
475  inline void SetDataType(int a_dataType)
476  {m_dataType = a_dataType;}
482  inline void SetDataSpec(int a_dataSpec)
483  {m_dataSpec = a_dataSpec;}
490  inline void SetAttenuationImage(FLTNB* ap_attenuationImage, int a_thread)
491  {m2p_attenuationImage[a_thread] = ap_attenuationImage;}
497  inline void SetFOMFlag(bool a_optimizerFOMFlag)
498  {m_optimizerFOMFlag = a_optimizerFOMFlag;}
504  inline void SetImageStatFlag(bool a_optimizerImageStatFlag)
505  {m_optimizerImageStatFlag = a_optimizerImageStatFlag;}
512  inline void SetNumbersOfIterationsAndSubsets(int a_nbIterations, int* ap_nbSubsets)
513  {m_nbIterations = a_nbIterations; mp_nbSubsets = ap_nbSubsets;}
519  inline void SetCurrentIteration(int a_currentIteration)
520  {m_currentIteration = a_currentIteration;}
526  inline void SetCurrentSubset(int a_currentSubset)
527  {m_currentSubset = a_currentSubset;}
533  inline int GetNbBackwardImages()
534  {return m_nbBackwardImages;}
541  {return m_initialValue;}
547  inline void SetOptimizerID(const string& a_optimizerID)
548  {m_optimizerID = a_optimizerID;}
555  inline const string& GetOptimizerID()
556  {return m_optimizerID;}
569  inline bool GetAcceptPenalty()
576  inline void SetPenalty(vPenalty* ap_penalty)
577  {mp_Penalty = ap_penalty;}
585  {return m_needGlobalSensitivity;}
586 
587 
588  // -------------------------------------------------------------------
589  // Private member functions
590  protected:
606  FLTNB ComputeSensitivity( FLTNB**** a4p_sensitivityImage,
607  int a_timeBasisFunction, int a_respBasisFunction, int a_cardBasisFunction,
608  int a_voxel );
618  FLTNB ForwardProject( oProjectionLine* ap_Line, FLTNB* ap_image = NULL );
629  void BackwardProject( oProjectionLine* ap_Line, FLTNB* ap_image, FLTNB a_value );
630 
631  // -------------------------------------------------------------------
632  // Data members
633  protected:
634 
635  string m_optimizerID;
636  int m_verbose;
653  // Iterations/subsets
658  // Optimizer figures-of-merit computation
662  uint64_t**** m4p_FOMNbBins;
665  // Image update statistics
674  // Penalty
678 };
679 
680 
681 // ----------------------------------------------------------------------
682 // Part of code that manages the auto declaration of children classes
683 // ----------------------------------------------------------------------
684 
685 // Macro for the function that creates the object
686 #define FUNCTION_OPTIMIZER(CLASS) \
687  static vOptimizer *make_optimizer() { return new CLASS(); };
688 
689 // Macro for the class that links the appropriate function to the map of objects
690 #define CLASS_OPTIMIZER(NAME,CLASS) \
691  class NAME##OptimizerCreator \
692  { \
693  public: \
694  NAME##OptimizerCreator() \
695  { sAddonManager::GetInstance()->mp_listOfOptimizers[#NAME] = CLASS::make_optimizer; } \
696  }; \
697  static NAME##OptimizerCreator OptimizerCreator##NAME;
698 
699 #endif
virtual int ReadConfigurationFile(const string &a_configurationFile)=0
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
bool GetNeedGlobalSensitivity()
Get the boolean saying if the sensitivity has to be computed globally for all data channels and not p...
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
void SetImageDimensionsAndQuantification(oImageDimensionsAndQuantification *ap_ImageDimensionsAndQuantification)
virtual int DataStep2Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
oImageDimensionsAndQuantification * mp_ImageDimensionsAndQuantification
FLTNB ComputeSensitivity(FLTNB ****a4p_sensitivityImage, int a_timeBasisFunction, int a_respBasisFunction, int a_cardBasisFunction, int a_voxel)
virtual int ImageUpdateStep()
A public function used to perform the image update step of the optimizer.
virtual int CheckSpecificParameters()=0
A private function used to check the parameters settings specific to the child optimizer.
vOptimizer()
The constructor of vOptimizer.
void SetDataType(int a_dataType)
FLTNB ForwardProject(oProjectionLine *ap_Line, FLTNB *ap_image=NULL)
int PreDataUpdateStep()
A public function used to do stuff that need to be done at the beginning of a subset (before the data...
void SetFOMFlag(bool a_optimizerFOMFlag)
void SetImageSpace(oImageSpace *ap_ImageSpace)
virtual int DataStep6Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
int GetRequiredPenaltyDerivativesOrder()
Get the penalty derivative order needed for this algorithm.
bool GetAcceptPenalty()
Get the boolean saying if the optimizer accepts penalties.
void SetAttenuationImage(FLTNB *ap_attenuationImage, int a_thread)
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
virtual int DataStep3BackwardProjectSensitivity(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
virtual int InitializeSpecific()=0
A private function used to initialize everything specific to the child optimizer. ...
void BackwardProject(oProjectionLine *ap_Line, FLTNB *ap_image, FLTNB a_value)
void SetOptimizerID(const string &a_optimizerID)
virtual int DataStep1ForwardProjectModel(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
void SetDataMode(int a_dataMode)
This class is designed to generically described any penalty applied to MAP algorithms.
virtual int DataStep8ComputeFOM(oProjectionLine *ap_Line, vEvent *ap_Event, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
void SetNbTOFBins(int a_nbTOFBins)
FLTNB GetInitialValue()
Get the initial image value (for initialization)
void SetNumbersOfIterationsAndSubsets(int a_nbIterations, int *ap_nbSubsets)
void SetVerbose(int a_verbose)
virtual ~vOptimizer()
The destructor of vOptimizer.
virtual int PreImageUpdateSpecificStep()
A private function used to perform any step required by the child optimizer, between the loop on even...
int PreImageUpdateStep()
A public function used to do stuff that need to be done between the loop over events and the image up...
This class is designed to generically described any iterative optimizer.
This class is designed to manage and store system matrix elements associated to a vEvent...
void SetImageStatFlag(bool a_optimizerImageStatFlag)
int Initialize()
A public function used to initialize the optimizer.
virtual int DataStep4Optional(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
This class holds all the matrices in the image domain that can be used in the algorithm: image...
Mother class for the Event objects.
virtual int DataStep5ComputeCorrections(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
void SetCurrentIteration(int a_currentIteration)
void SetCurrentSubset(int a_currentSubset)
This class is designed to manage all dimensions and quantification related stuff. ...
int GetNbBackwardImages()
Get the number of backward images used by the specific optimizer.
int UpdateVisitedVoxels()
A public function used to update the 'visited' voxels after each subset.
void SetDataSpec(int a_dataSpec)
virtual void ShowHelpSpecific()=0
A function used to show help about the child module.
virtual int DataStep7BackwardProjectCorrections(oProjectionLine *ap_Line, vEvent *ap_Event, int a_bed, int a_timeFrame, int a_respGate, int a_cardGate, int a_thread)
void ShowHelp()
A function used to show help about the optimizer.
void SetPenalty(vPenalty *ap_penalty)
int CheckParameters()
A public function used to check the parameters settings.
Declaration of class vPenalty.
virtual int PreDataUpdateSpecificStep()
A private function used to perform any step required by the child optimizer, before the loop on event...
virtual int ReadOptionsList(const string &a_optionsList)=0