]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCCalibRaw.h
Corrected hopefully a bug which causes problem with variable-size HNSPARSE
[u/mrichter/AliRoot.git] / TPC / AliTPCCalibRaw.h
CommitLineData
7442bceb 1#ifndef ALITPCCALIBRAW_H
2#define ALITPCCALIBRAW_H
ddeb9c4f 3/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
5
6////////////////////////////////////////////////////////////////////////////////////////
7// //
8// TPC ALTRO Header analysis //
9// //
10////////////////////////////////////////////////////////////////////////////////////////
11
12#include <TVectorF.h>
13#include <TObjArray.h>
14#include <THnSparse.h>
15
16#include "AliTPCCalibRawBase.h"
17#include "AliTPCCalPad.h"
18#include "AliTPCROC.h"
19
20class TH2C;
be3dbaa0 21class TH1F;
f113dfeb 22class TMap;
1d32f273 23class TGraph;
24class TCanvas;
ddeb9c4f 25
26class AliTPCCalibRaw : public AliTPCCalibRawBase {
27public:
28 AliTPCCalibRaw();
f113dfeb 29 AliTPCCalibRaw(const TMap *config);
ddeb9c4f 30
31 virtual ~AliTPCCalibRaw();
32
be3dbaa0 33 enum {kNRCU=216};
ddeb9c4f 34
35 virtual Int_t Update(const Int_t isector, const Int_t iRow, const Int_t iPad,
36 const Int_t iTimeBin, const Float_t signal);
5312f439 37 virtual void UpdateDDL();
ddeb9c4f 38 virtual void EndEvent();
39 virtual void ResetEvent();
40 virtual void Analyse();
41
42 UInt_t GetNFailL1Phase() const {return fNFailL1Phase;}
be3dbaa0 43 UInt_t GetNFailL1PhaseEvents() const {return fNFailL1PhaseEvent;}
ddeb9c4f 44 Int_t GetPeakDetectionMinus() const {return fPeakDetMinus;}
45 Int_t GetPeakDetectionPlus() const {return fPeakDetPlus;}
46
47 const TVectorF* GetALTROL1PhaseEvents() const {return &fArrALTROL1Phase;}
48
49 const TVectorF *GetALTROL1PhaseEventsRCU(Int_t rcu) const {return (TVectorF*)fArrALTROL1PhaseEvent.At(rcu);}
50 const TVectorF *GetALTROL1PhaseFailEventsRCU(Int_t rcu) const {return (TVectorF*)fArrALTROL1PhaseFailEvent.At(rcu);}
51
1d32f273 52 const TVectorF *GetOccupancyEvent() const {return &fVOccupancyEvent;}
53 const TVectorF *GetOccupancyEventSensitive() const {return &fVOccupancySenEvent;}
54 const TVectorF *GetSignalSumEvent() const {return &fVSignalSumEvent;}
55 const TVectorF *GetSignalSumEventSensitive() const {return &fVSignalSumSenEvent;}
56 const TVectorF *GetFiredPadsSensitive() const {return &fVNfiredPadsSenEvent;}
57 const TVectorF *GetEventTimeStamps() const {return &fVTimeStampEvent;}
58 UInt_t GetFirstTimeStamp() const {return fFirstTimeStamp;}
ddeb9c4f 59
1d32f273 60 void SetRangePeakDetection(Int_t minus, Int_t plus) { fPeakDetMinus=minus; fPeakDetPlus=plus;}
61 //Phase info
ddeb9c4f 62 TH2C *MakeHistL1RCUEvents(Int_t type=0);
63 TH2C *MakeHistL1RCUEventsIROC(Int_t type=0);
64 TH2C *MakeHistL1RCUEventsOROC(Int_t type=0);
be3dbaa0 65 TH1F *MakeHistL1PhaseDist();
66 TVectorF *MakeVectL1PhaseDist();
1d32f273 67 //Occupancy info
68 TGraph* MakeGraphOccupancy(const Int_t type=0, const Int_t xType=0);
7442bceb 69// TGraph* MakeGraphNoiseEvents();
1d32f273 70 TCanvas* MakeCanvasOccupancy(const Int_t xType=1, Bool_t sen=kFALSE);
ddeb9c4f 71
72 const THnSparseI *GetHnDrift() const {return fHnDrift;}
73// AliTPCCalPad *CreateCalPadL1Mean();
74// AliTPCCalPad *CreateCalPadL1RMS();
75
7442bceb 76 void Merge(AliTPCCalibRaw * const sig);
77 virtual Long64_t Merge(TCollection * const list);
78
ddeb9c4f 79private:
80 Int_t fPeakDetMinus; // Consecutive timebins on rising edge to be regarded as a signal
81 Int_t fPeakDetPlus; // Consecutive timebins on falling edge to be regarded as a signal
82 UInt_t fNFailL1Phase; //Number of failures in L1 phase
be3dbaa0 83 UInt_t fNFailL1PhaseEvent; //Number of events with L1 phase failures
ddeb9c4f 84 UInt_t fFirstTimeStamp; //Time Stamp from first event
85 //binning dv hist
86 UInt_t fNSecTime; //Number of seconds per bin in time
87 UInt_t fNBinsTime; //Number of bin in time
88 //processing information
89 Bool_t fPadProcessed; //! if last pead has already been filled for the current pad
90 Int_t fCurrentChannel; //! current channel processed
91 Int_t fCurrentSector; //! current sector processed
92 Int_t fLastSector; //! current sector processed
93 Int_t fCurrentRow; //! current row processed
94 Int_t fCurrentPad; //! current pad processed
95 Int_t fLastTimeBinProc; //! last time bin processed
96 Int_t fPeakTimeBin; //! time bin with local maximum
97 Int_t fLastSignal; //! last signal processed
98 Int_t fNOkPlus; //! number of processed time bins fullfilling peak criteria
99 Int_t fNOkMinus; //! number of processed time bins fullfilling peak criteria
1d32f273 100 Int_t fNanoSec; //! current nano seconds stamp
ddeb9c4f 101//
102 //L1 phase stuff
103 TVectorF fArrCurrentPhaseDist; //!Phase distribution of the current event
be3dbaa0 104 TVectorF fArrCurrentPhase; //!Current phase of all RCUs
105 TVectorF fArrFailEventNumber; //event numbers of failed events;
ddeb9c4f 106 TVectorF fArrALTROL1Phase; //Array of L1 phases on an event bases;
107 TObjArray fArrALTROL1PhaseEvent; //L1 phase for each RCU and event
108 TObjArray fArrALTROL1PhaseFailEvent; //L1 failure for each RCU and event
109 //drift velocity stuff
110 enum {kHnBinsDV=3};
111 THnSparseI *fHnDrift; //Histogram last time bin vs. ROC, Time
1d32f273 112 //occupancy
113 TVectorF fVOccupancyEvent; //occupancy per event (number of samples above threshold)
114 TVectorF fVSignalSumEvent; //occupancy per event (sum of all adc values)
115 TVectorF fVOccupancySenEvent; //occupancy per event (number of samples abouve threshold) in sensitive regions
116 TVectorF fVSignalSumSenEvent; //occupancy per event (sum of all adc values) in sensitive regions
117 TVectorF fVNfiredPadsSenEvent; //number of pads with a signal above threshold in sensitive regions
118 TVectorF fVTimeStampEvent; //timestamp for all events
ddeb9c4f 119
120 TVectorF *MakeArrL1PhaseRCU(Int_t rcu, Bool_t force=kFALSE);
121 TVectorF *MakeArrL1PhaseFailRCU(Int_t rcu, Bool_t force=kFALSE);
122
123 Bool_t IsEdgePad(Int_t sector, Int_t row, Int_t pad) const;
124 void CreateDVhist();
125
7442bceb 126 AliTPCCalibRaw(const AliTPCCalibRaw &calib);
ddeb9c4f 127 AliTPCCalibRaw& operator = (const AliTPCCalibRaw &source);
128
1d32f273 129 ClassDef(AliTPCCalibRaw,3) // Analysis of the Altro header information
ddeb9c4f 130};
131
132//----------------------
133// Inline Functions
134//----------------------
135inline TVectorF *AliTPCCalibRaw::MakeArrL1PhaseRCU(Int_t rcu, Bool_t force)
136{
137 TVectorF *arr=(TVectorF*)fArrALTROL1PhaseEvent.UncheckedAt(rcu);
138 if (!arr && force) {
be3dbaa0 139 arr=new TVectorF(100);
ddeb9c4f 140 fArrALTROL1PhaseEvent.AddAt(arr,rcu);
141 }
142 return arr;
143}
144//
145inline TVectorF *AliTPCCalibRaw::MakeArrL1PhaseFailRCU(Int_t rcu, Bool_t force)
146{
147 TVectorF *arr=(TVectorF*)fArrALTROL1PhaseFailEvent.UncheckedAt(rcu);
148 if (!arr && force) {
1d32f273 149 arr=new TVectorF(100);
ddeb9c4f 150 fArrALTROL1PhaseFailEvent.AddAt(arr,rcu);
151 }
152 return arr;
153}
154//_____________________________________________________________________
155inline Bool_t AliTPCCalibRaw::IsEdgePad(Int_t sector, Int_t row, Int_t pad) const
156{
157 //
158 // return true if pad is on the edge of a row
159 //
160 Int_t edge1 = 0;
161 if ( pad == edge1 ) return kTRUE;
162 Int_t edge2 = fROC->GetNPads(sector,row)-1;
163 if ( pad == edge2 ) return kTRUE;
164
165 return kFALSE;
166}
167
168
169#endif