ROOT logo
// -*- mode: c++ -*-
#ifndef ALICALORAWANALYZER_H
#define ALICALORAWANALYZER_H
/**************************************************************************
 * This file is property of and copyright by                              *
 * the Relatvistic Heavy Ion Group (RHIG), Yale University, US, 2009      *
 *                                                                        *
 * Primary Author: Per Thomas Hille <p.t.hille@fys.uio.no>                *
 *                                                                        *
 * Contributors are mentioned in the code where appropriate.              *
 * Please report bugs to p.t.hille@fys.uio.no                             *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/

//Base class for extraction 
//of signal amplitude and peak position
//From CALO Calorimeter RAW data

#include "Rtypes.h"
#include "TObject.h"
#include <vector>
#include "TObjArray.h"
#include "AliCaloFitResults.h"
#include "AliCaloConstants.h"
using namespace ALTRO;
using namespace CALO;


class AliCaloBunchInfo;


class  AliCaloRawAnalyzer : public TObject
{
public:
  AliCaloRawAnalyzer(const char *name="AliCaloRawAnalyzer", const char *nameshort="RawAna");
  virtual ~AliCaloRawAnalyzer() { ; }

  virtual AliCaloFitResults Evaluate( const std::vector<AliCaloBunchInfo> &/*bunchvector*/, 
                                      UInt_t /*altrocfg1*/,  UInt_t /*altrocfg2*/ )  = 0;

  static void PrintBunches( const std::vector<AliCaloBunchInfo> &bunchvector );
  static void PrintBunch  ( const AliCaloBunchInfo &bunch );
  
  int PreFitEvaluateSamples( const std::vector<AliCaloBunchInfo>  &bunchvector, 
                             UInt_t altrocfg1, UInt_t altrocfg2, Int_t & index,
                             Float_t & maxf, short & maxamp, short & maxampindex,
                             Float_t & ped, int & first, int & last, int acut);
  
  void SetTimeConstraint  (int min, int max );
  void SetVerbose         (bool verbose = true){ fVerbose     = verbose; }
  void SetIsZeroSuppressed(bool iszs = true)   { fIsZerosupressed = iszs ; }
  void SetAmpCut     (Float_t cut)             { fAmpCut      = cut ; }
  void SetFitArrayCut(Int_t cut)               { fFitArrayCut = cut ; }
  void SetNsampleCut (Int_t cut)               { fNsampleCut  = cut ; }
  void SetOverflowCut(Int_t cut)               { fOverflowCut = cut ; }
  void SetNsamplePed (Int_t i)                 { fNsamplePed  = i   ; }
  void SetL1Phase    (Double_t phase)          { fL1Phase     = phase ; }

  bool    GetIsZeroSuppressed() const { return fIsZerosupressed;}
  Float_t GetAmpCut()      const { return fAmpCut      ; }
  Int_t   GetFitArrayCut() const { return fFitArrayCut ; }
  Int_t   GetNsampleCut()  const { return fNsampleCut  ; }
  Int_t   GetOverflowCut() const { return fOverflowCut ; }
  Int_t   GetNsamplePed()  const { return fNsamplePed  ; }

  // access to array info
  Double_t GetReversed(const int i) const { return fReversed[i]; }
  const char * GetAlgoName()   const { return fName      ; }
  const char * GetAlgoAbbr()   const { return fNameShort ; }
  Algo::fitAlgorithm GetAlgo() const { return fAlgo      ; }

  Double_t CalculateChi2(const Double_t amp, const Double_t time,
                         const Int_t first, const Int_t last,
                         const Double_t adcErr=1,
                         const Double_t tau=2.35) const;
  
  void CalculateMeanAndRMS(const Int_t first, const Int_t last,
                           Double_t & mean, Double_t & rms);
  
  short Max( const AliCaloBunchInfo *const bunch, int * maxindex) const;
  
  UShort_t Max(const UShort_t *data, const int length ) const;
  
  bool CheckBunchEdgesForMax( const AliCaloBunchInfo *const bunch) const;
  
  bool IsInTimeRange( const int maxindex, const int maxtime, const int mintime ) const;
  
  Float_t  ReverseAndSubtractPed( const AliCaloBunchInfo *bunch,
                                  UInt_t altrocfg1,  UInt_t altrocfg2,
                                  double *  outarray ) const;
  
  int  SelectBunch( const std::vector<AliCaloBunchInfo> &bunchvector,
                    short * maxampbin, short * maxamplitude );
  
  void SelectSubarray( const Double_t *date, int length, short maxindex,
                       int * first, int * last, int cut) const;
  
  Float_t EvaluatePedestal(const UShort_t * const data, const int length ) const;
  
  // Used in AliCaloRawAnalyzerFitter
  Float_t GetTau()         const { return fTau    ; }
  void    SetTau   (Float_t tau) { fTau = tau     ; }
  Bool_t  GetFixTau()      const { return fFixTau ; }
  void    SetFixTau(Bool_t b)    { fFixTau = b    ; }

protected:
  Double_t fReversed[ALTROMAXSAMPLES]; //Reversed sequence of samples (pedestalsubtracted)
  int fMinTimeIndex; //The timebin of the max signal value must be between fMinTimeIndex and fMaxTimeIndex
  int fMaxTimeIndex; //The timebin of the max signal value must be between fMinTimeIndex and fMaxTimeIndex
  int fFitArrayCut;  //Cut on ADC value (after ped. subtraction) for signals used for fit
  Float_t fAmpCut;   //Max ADC - pedestal must be higher than this befor attemting to extract the amplitude 
  int fNsampleCut;   //Minimum number of sample require before attemting to extract signal parameters 
  int fOverflowCut; // value when ADC starts to saturate
  int fNsamplePed;   //Number of samples used for pedestal calculation (first in bunch) 
  bool fIsZerosupressed; //Wether or not the data is zeros supressed, by default its assumed that the baseline is also subtracted if set to true
  
  bool fVerbose;     //Print debug information to std out if set to true
  char fName[256]; // Name of the algorithm
  char fNameShort[256]; // Abbrevation for the name
  
  Algo::fitAlgorithm fAlgo; // Which algorithm to use

  Double_t fL1Phase; // Phase of the ADC sampling clock relative to the LHC clock
  Double_t fAmp;     // The amplitude in entities of ADC counts
  Double_t fTof;     // The amplitude in entities of ADC counts
  Float_t  fTau;     // Rise time of the signal (peak position = t0 +tau), by defauly it is 235 ns
  Bool_t   fFixTau;  // Fixed fit parameter or not, used in AliCaloRawAnalyzerFitter
  
  ClassDef(AliCaloRawAnalyzer, 3)

};

#endif
 AliCaloRawAnalyzer.h:1
 AliCaloRawAnalyzer.h:2
 AliCaloRawAnalyzer.h:3
 AliCaloRawAnalyzer.h:4
 AliCaloRawAnalyzer.h:5
 AliCaloRawAnalyzer.h:6
 AliCaloRawAnalyzer.h:7
 AliCaloRawAnalyzer.h:8
 AliCaloRawAnalyzer.h:9
 AliCaloRawAnalyzer.h:10
 AliCaloRawAnalyzer.h:11
 AliCaloRawAnalyzer.h:12
 AliCaloRawAnalyzer.h:13
 AliCaloRawAnalyzer.h:14
 AliCaloRawAnalyzer.h:15
 AliCaloRawAnalyzer.h:16
 AliCaloRawAnalyzer.h:17
 AliCaloRawAnalyzer.h:18
 AliCaloRawAnalyzer.h:19
 AliCaloRawAnalyzer.h:20
 AliCaloRawAnalyzer.h:21
 AliCaloRawAnalyzer.h:22
 AliCaloRawAnalyzer.h:23
 AliCaloRawAnalyzer.h:24
 AliCaloRawAnalyzer.h:25
 AliCaloRawAnalyzer.h:26
 AliCaloRawAnalyzer.h:27
 AliCaloRawAnalyzer.h:28
 AliCaloRawAnalyzer.h:29
 AliCaloRawAnalyzer.h:30
 AliCaloRawAnalyzer.h:31
 AliCaloRawAnalyzer.h:32
 AliCaloRawAnalyzer.h:33
 AliCaloRawAnalyzer.h:34
 AliCaloRawAnalyzer.h:35
 AliCaloRawAnalyzer.h:36
 AliCaloRawAnalyzer.h:37
 AliCaloRawAnalyzer.h:38
 AliCaloRawAnalyzer.h:39
 AliCaloRawAnalyzer.h:40
 AliCaloRawAnalyzer.h:41
 AliCaloRawAnalyzer.h:42
 AliCaloRawAnalyzer.h:43
 AliCaloRawAnalyzer.h:44
 AliCaloRawAnalyzer.h:45
 AliCaloRawAnalyzer.h:46
 AliCaloRawAnalyzer.h:47
 AliCaloRawAnalyzer.h:48
 AliCaloRawAnalyzer.h:49
 AliCaloRawAnalyzer.h:50
 AliCaloRawAnalyzer.h:51
 AliCaloRawAnalyzer.h:52
 AliCaloRawAnalyzer.h:53
 AliCaloRawAnalyzer.h:54
 AliCaloRawAnalyzer.h:55
 AliCaloRawAnalyzer.h:56
 AliCaloRawAnalyzer.h:57
 AliCaloRawAnalyzer.h:58
 AliCaloRawAnalyzer.h:59
 AliCaloRawAnalyzer.h:60
 AliCaloRawAnalyzer.h:61
 AliCaloRawAnalyzer.h:62
 AliCaloRawAnalyzer.h:63
 AliCaloRawAnalyzer.h:64
 AliCaloRawAnalyzer.h:65
 AliCaloRawAnalyzer.h:66
 AliCaloRawAnalyzer.h:67
 AliCaloRawAnalyzer.h:68
 AliCaloRawAnalyzer.h:69
 AliCaloRawAnalyzer.h:70
 AliCaloRawAnalyzer.h:71
 AliCaloRawAnalyzer.h:72
 AliCaloRawAnalyzer.h:73
 AliCaloRawAnalyzer.h:74
 AliCaloRawAnalyzer.h:75
 AliCaloRawAnalyzer.h:76
 AliCaloRawAnalyzer.h:77
 AliCaloRawAnalyzer.h:78
 AliCaloRawAnalyzer.h:79
 AliCaloRawAnalyzer.h:80
 AliCaloRawAnalyzer.h:81
 AliCaloRawAnalyzer.h:82
 AliCaloRawAnalyzer.h:83
 AliCaloRawAnalyzer.h:84
 AliCaloRawAnalyzer.h:85
 AliCaloRawAnalyzer.h:86
 AliCaloRawAnalyzer.h:87
 AliCaloRawAnalyzer.h:88
 AliCaloRawAnalyzer.h:89
 AliCaloRawAnalyzer.h:90
 AliCaloRawAnalyzer.h:91
 AliCaloRawAnalyzer.h:92
 AliCaloRawAnalyzer.h:93
 AliCaloRawAnalyzer.h:94
 AliCaloRawAnalyzer.h:95
 AliCaloRawAnalyzer.h:96
 AliCaloRawAnalyzer.h:97
 AliCaloRawAnalyzer.h:98
 AliCaloRawAnalyzer.h:99
 AliCaloRawAnalyzer.h:100
 AliCaloRawAnalyzer.h:101
 AliCaloRawAnalyzer.h:102
 AliCaloRawAnalyzer.h:103
 AliCaloRawAnalyzer.h:104
 AliCaloRawAnalyzer.h:105
 AliCaloRawAnalyzer.h:106
 AliCaloRawAnalyzer.h:107
 AliCaloRawAnalyzer.h:108
 AliCaloRawAnalyzer.h:109
 AliCaloRawAnalyzer.h:110
 AliCaloRawAnalyzer.h:111
 AliCaloRawAnalyzer.h:112
 AliCaloRawAnalyzer.h:113
 AliCaloRawAnalyzer.h:114
 AliCaloRawAnalyzer.h:115
 AliCaloRawAnalyzer.h:116
 AliCaloRawAnalyzer.h:117
 AliCaloRawAnalyzer.h:118
 AliCaloRawAnalyzer.h:119
 AliCaloRawAnalyzer.h:120
 AliCaloRawAnalyzer.h:121
 AliCaloRawAnalyzer.h:122
 AliCaloRawAnalyzer.h:123
 AliCaloRawAnalyzer.h:124
 AliCaloRawAnalyzer.h:125
 AliCaloRawAnalyzer.h:126
 AliCaloRawAnalyzer.h:127
 AliCaloRawAnalyzer.h:128
 AliCaloRawAnalyzer.h:129
 AliCaloRawAnalyzer.h:130
 AliCaloRawAnalyzer.h:131
 AliCaloRawAnalyzer.h:132
 AliCaloRawAnalyzer.h:133
 AliCaloRawAnalyzer.h:134
 AliCaloRawAnalyzer.h:135
 AliCaloRawAnalyzer.h:136
 AliCaloRawAnalyzer.h:137
 AliCaloRawAnalyzer.h:138
 AliCaloRawAnalyzer.h:139
 AliCaloRawAnalyzer.h:140