]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSRawDecoder.cxx
Updated histogram limits (PHOS energy)
[u/mrichter/AliRoot.git] / PHOS / AliPHOSRawDecoder.cxx
1 /**************************************************************************
2  * Copyright(c) 2007, ALICE Experiment at CERN, All rights reserved.      *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 // This class decodes the stream of ALTRO samples to extract
19 // the PHOS "digits" of current event.
20 // 
21 // Typical use case:
22 //     AliRawReader* rf = new AliRawReaderDate("2006run2211.raw");
23 //     AliPHOSRawDecoder dc(rf);
24 //     while (rf->NextEvent()) {
25 //       dc.SubtractPedestals(kTRUE);
26 //       while ( dc.NextDigit() ) {
27 //         Int_t module = dc.GetModule();
28 //         Int_t column = dc.GetColumn();
29 //         Int_t row = dc.GetRow();
30 //         Double_t amplitude = dc.GetEnergy();
31 //         Double_t time = dc.GetTime();
32 //         Bool_t IsLowGain = dc.IsLowGain();
33 //            ..........
34 //       }
35 //     }
36
37 // Author: Boris Polichtchouk
38
39 // --- ROOT system ---
40 #include "TArrayI.h"
41 #include "TMath.h"
42
43 // --- AliRoot header files ---
44 #include "AliPHOSRawDecoder.h"
45 #include "AliRawReader.h"
46 #include "AliPHOSCalibData.h"
47 #include "AliLog.h"
48
49 ClassImp(AliPHOSRawDecoder)
50
51 //-----------------------------------------------------------------------------
52 AliPHOSRawDecoder::AliPHOSRawDecoder():
53   fRawReader(0),fCaloStream(0),fPedSubtract(kFALSE),fEnergy(-111),fTime(-111),fQuality(0.),fPedestalRMS(0.),
54   fAmpOffset(0),fModule(-1),fColumn(-1),fRow(-1),fNewModule(-1),fNewColumn(-1),fNewRow(-1),fNewAmp(0),fNewTime(0), 
55   fLowGainFlag(kFALSE),fNewLowGainFlag(kFALSE),fOverflow(kFALSE),fSamples(0),fTimes(0),fCalibData(0)
56 {
57   //Default constructor.
58 }
59
60 //-----------------------------------------------------------------------------
61 AliPHOSRawDecoder::AliPHOSRawDecoder(AliRawReader* rawReader,  AliAltroMapping **mapping):
62   fRawReader(0),fCaloStream(0),fPedSubtract(kFALSE),fEnergy(-111),fTime(-111),fQuality(0.),fPedestalRMS(0.),
63   fAmpOffset(0),fModule(-1),fColumn(-1),fRow(-1),fNewModule(-1),fNewColumn(-1),fNewRow(-1),fNewAmp(0),fNewTime(0),
64   fLowGainFlag(kFALSE),fNewLowGainFlag(kFALSE),fOverflow(kFALSE),fSamples(0),fTimes(0),fCalibData(0)
65 {
66   //Construct a decoder object.
67   //Is is user responsibility to provide next raw event 
68   //using AliRawReader::NextEvent().
69
70   fRawReader =  rawReader;
71   fCaloStream = new AliCaloRawStream(rawReader,"PHOS",mapping);
72   fSamples = new TArrayI(100);
73   fTimes = new TArrayI(100);
74 }
75
76 //-----------------------------------------------------------------------------
77 AliPHOSRawDecoder::~AliPHOSRawDecoder()
78 {
79   //Destructor.
80
81   if(fCaloStream){ delete fCaloStream; fCaloStream=0;}
82   if(fSamples){ delete fSamples; fSamples=0 ;}
83   if(fTimes){ delete fTimes; fTimes=0 ;}
84 }
85
86 //-----------------------------------------------------------------------------
87 AliPHOSRawDecoder::AliPHOSRawDecoder(const AliPHOSRawDecoder &phosDecoder ):
88   fRawReader(phosDecoder.fRawReader),fCaloStream(phosDecoder.fCaloStream),
89   fPedSubtract(phosDecoder.fPedSubtract),
90   fEnergy(phosDecoder.fEnergy),fTime(phosDecoder.fTime),fQuality(phosDecoder.fQuality),fPedestalRMS(phosDecoder.fPedestalRMS),
91   fAmpOffset(phosDecoder.fAmpOffset),fModule(phosDecoder.fModule),fColumn(phosDecoder.fColumn),
92   fRow(phosDecoder.fRow),fNewModule(phosDecoder.fNewModule),fNewColumn(phosDecoder.fNewColumn),
93   fNewRow(phosDecoder.fNewRow),fNewAmp(phosDecoder.fNewAmp),fNewTime(phosDecoder.fNewTime),
94   fLowGainFlag(phosDecoder.fLowGainFlag),fNewLowGainFlag(phosDecoder.fNewLowGainFlag),
95   fOverflow(phosDecoder.fOverflow),fSamples(phosDecoder.fSamples),
96   fTimes(phosDecoder.fTimes),fCalibData(phosDecoder.fCalibData) 
97 {
98   //Copy constructor.
99 }
100
101 //-----------------------------------------------------------------------------
102 AliPHOSRawDecoder& AliPHOSRawDecoder::operator = (const AliPHOSRawDecoder &phosDecode)
103 {
104   //Assignment operator.
105
106   if(this != &phosDecode) {
107     fRawReader = phosDecode.fRawReader;
108
109     if(fCaloStream) delete fCaloStream;
110     fCaloStream = phosDecode.fCaloStream;
111
112     fEnergy = phosDecode.fEnergy;
113     fTime = phosDecode.fTime;
114     fQuality = phosDecode.fQuality ;
115     fPedestalRMS = phosDecode.fPedestalRMS ;
116     fAmpOffset = phosDecode.fAmpOffset ;
117     fModule = phosDecode.fModule;
118     fColumn = phosDecode.fColumn;
119     fRow = phosDecode.fRow;
120     fNewModule = phosDecode.fNewModule;
121     fNewColumn = phosDecode.fNewColumn;
122     fNewRow = phosDecode.fNewRow;
123     fNewAmp = phosDecode.fNewAmp ;
124     fNewTime= phosDecode.fNewTime ;
125     fLowGainFlag = phosDecode.fLowGainFlag;
126     fNewLowGainFlag = phosDecode.fNewLowGainFlag;
127     fOverflow = phosDecode.fOverflow ;
128     
129     if(fSamples) delete fSamples;
130     fSamples = phosDecode.fSamples;
131
132     if(fTimes) delete fTimes;
133     fTimes = phosDecode.fTimes;
134     fCalibData = phosDecode.fCalibData; 
135   }
136
137   return *this;
138 }
139
140 //-----------------------------------------------------------------------------
141
142 Bool_t AliPHOSRawDecoder::NextDigit()
143 {
144   //Extract an energy deposited in the crystal,
145   //crystal' position (module,column,row),
146   //time and gain (high or low).
147   
148   AliCaloRawStream* in = fCaloStream;
149   
150   Int_t    iBin     = 0;
151   Int_t    mxSmps   = fSamples->GetSize();
152   Int_t    tLength  = 0;
153   fEnergy = -111;
154   Float_t pedMean = 0;
155   Float_t pedRMS = 0;
156   Int_t   nPed = 0;
157   Float_t baseLine = 1.0;
158   const Int_t kPreSamples = 10;
159   
160   fSamples->Reset();
161   while ( in->Next() ) { 
162
163      if(!tLength) {
164        tLength = in->GetTimeLength();
165        if(tLength>mxSmps) {
166          fSamples->Set(tLength);
167        }
168      }
169      
170      // Fit the full sample
171      if(in->IsNewHWAddress() && iBin>0) {
172        
173        iBin=0;
174        //First remember new sample
175        fNewLowGainFlag = in->IsLowGain();
176        fNewModule = in->GetModule()+1;
177        fNewRow    = in->GetRow()   +1;
178        fNewColumn = in->GetColumn()+1;
179        fNewAmp = in->GetSignal() ;
180        fNewTime=in->GetTime() ;                                                                                                                               
181        
182        // Temporarily we take the energy as a maximum amplitude
183        // and the pedestal from the 0th point (30 Aug 2006).
184        // Time is not evaluated for the moment (12.01.2007). 
185        // Take is as a first time bin multiplied by the sample tick time
186        
187        if(fPedSubtract) {
188          if (nPed > 0){
189            fPedestalRMS=(pedRMS-pedMean*pedMean/nPed)/nPed ;
190            if(fPedestalRMS > 0.) 
191             fPedestalRMS = TMath::Sqrt(fPedestalRMS) ;
192            fEnergy -= (Double_t)(pedMean/nPed); // pedestal subtraction
193          }
194          else
195            return kFALSE;
196        }
197        else{
198          //take pedestals from DB
199          Double_t pedestal = (Double_t) fAmpOffset ;
200          if(fCalibData){
201            Float_t truePed = fCalibData->GetADCpedestalEmc(fModule, fColumn, fRow) ;
202            Int_t   altroSettings = fCalibData->GetAltroOffsetEmc(fModule, fColumn, fRow) ;
203            pedestal += truePed - altroSettings ;
204          }
205          else{
206 //           printf("AliPHOSRawDecoder::NextDigit() Can not read data from OCDB \n") ;
207          }
208          fEnergy-=pedestal ;
209        }
210        if (fEnergy < baseLine) fEnergy = 0;
211
212        return kTRUE;
213      }
214
215      fLowGainFlag = in->IsLowGain();
216 //     fTime =   in->GetTime();
217      fTime = 1;
218      fModule = in->GetModule()+1;
219      fRow    = in->GetRow()   +1;
220      fColumn = in->GetColumn()+1;
221
222     if(fLowGainFlag==fNewLowGainFlag && fModule==fNewModule &&
223        fRow==fNewRow && fColumn==fNewColumn ){
224        if(fNewAmp>fEnergy)  fEnergy = (Double_t)fNewAmp ;
225        fNewModule=-1 ;
226     } 
227
228      //Calculate pedestal if necessary
229      if(fPedSubtract && in->GetTime() < kPreSamples) {
230        pedMean += in->GetSignal();
231        pedRMS+=in->GetSignal()*in->GetSignal() ;
232        nPed++;
233      }
234      if((Double_t)in->GetSignal() > fEnergy) fEnergy = (Double_t)in->GetSignal();
235      iBin++ ;
236      
237    } // in.Next()
238    
239    return kFALSE;
240 }