]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0Reconstructor.cxx
violations fixed
[u/mrichter/AliRoot.git] / T0 / AliT0Reconstructor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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  *  T0 reconstruction and filling ESD
19  *  - reconstruct mean time (interation time) 
20  *  - vertex position
21  *  -  multiplicity
22  ********************************************************************/
23
24 #include <AliESDEvent.h>
25 #include "AliLog.h"
26 #include "AliT0RecPoint.h"
27 #include "AliRawReader.h"
28 #include "AliT0RawReader.h"
29 #include "AliT0digit.h"
30 #include "AliT0Reconstructor.h"
31 #include "AliT0Parameters.h"
32 #include "AliT0Calibrator.h"
33
34 #include <TArrayI.h>
35 #include <TGraph.h>
36 #include <TMath.h>
37
38 ClassImp(AliT0Reconstructor)
39
40   AliT0Reconstructor:: AliT0Reconstructor(): AliReconstructor(),
41                                              fdZonA(0),
42                                              fdZonC(0),
43                                              fZposition(0),
44                                              fParam(NULL),
45                                              fAmpLEDrec()
46 {
47   //constructor
48
49  AliDebug(1,"Start reconstructor ");
50   
51   fParam = AliT0Parameters::Instance();
52   fParam->Init();
53   for (Int_t i=0; i<24; i++){
54     TGraph* gr = fParam ->GetAmpLEDRec(i);
55     fAmpLEDrec.AddAtAndExpand(gr,i) ;  
56 //    fTime0vertex[i]= fParam->GetTimeV0(i);
57   }
58   fdZonC = TMath::Abs(fParam->GetZPositionShift("T0/C/PMT1"));
59   fdZonA = TMath::Abs(fParam->GetZPositionShift("T0/A/PMT15"));
60
61 }
62 //____________________________________________________________________
63
64 AliT0Reconstructor::AliT0Reconstructor(const AliT0Reconstructor &r):
65   AliReconstructor(r),
66                                              fdZonA(0),
67                                              fdZonC(0),
68                                              fZposition(0),
69                                              fParam(NULL),
70                                              fAmpLEDrec()
71  {
72   //
73   // AliT0Reconstructor copy constructor
74   //
75
76   ((AliT0Reconstructor &) r).Copy(*this);
77
78 }
79
80 //_____________________________________________________________________________
81 AliT0Reconstructor &AliT0Reconstructor::operator=(const AliT0Reconstructor &r)
82 {
83   //
84   // Assignment operator
85   //
86
87   if (this != &r) ((AliT0Reconstructor &) r).Copy(*this);
88   return *this;
89
90 }
91
92 //_____________________________________________________________________________
93
94 void AliT0Reconstructor::Reconstruct(TTree*digitsTree, TTree*clustersTree) const
95   
96 {
97   // T0 digits reconstruction
98   // T0RecPoint writing 
99   
100    
101   TArrayI * timeCFD = new TArrayI(24); 
102   TArrayI * timeLED = new TArrayI(24); 
103   TArrayI * chargeQT0 = new TArrayI(24); 
104   TArrayI * chargeQT1 = new TArrayI(24); 
105
106   AliT0Calibrator *calib=new AliT0Calibrator(); 
107
108   //  Int_t mV2Mip = param->GetmV2Mip();     
109   //mV2Mip = param->GetmV2Mip();     
110   Float_t channelWidth = fParam->GetChannelWidth() ;  
111   Int_t meanT0 = fParam->GetMeanT0();
112   
113   AliDebug(1,Form("Start DIGITS reconstruction "));
114   
115   TBranch *brDigits=digitsTree->GetBranch("T0");
116   AliT0digit *fDigits = new AliT0digit() ;
117   if (brDigits) {
118     brDigits->SetAddress(&fDigits);
119   }else{
120     AliError(Form("EXEC Branch T0 digits not found"));
121      return;
122   }
123   
124   digitsTree->GetEvent(0);
125   digitsTree->GetEntry(0);
126   brDigits->GetEntry(0);
127   fDigits->GetTimeCFD(*timeCFD);
128   fDigits->GetTimeLED(*timeLED);
129   fDigits->GetQT0(*chargeQT0);
130   fDigits->GetQT1(*chargeQT1);
131
132   
133   Float_t besttimeA=999999;
134   Float_t besttimeC=999999;
135   Int_t pmtBestA=99999;
136   Int_t pmtBestC=99999;
137   Float_t timeDiff=999999, meanTime=0;
138   
139
140
141   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
142   clustersTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
143   
144   Float_t time[24], adc[24];
145   for (Int_t ipmt=0; ipmt<24; ipmt++) {
146     if(timeCFD->At(ipmt)>0 ){
147       Double_t qt0 = Double_t(chargeQT0->At(ipmt));
148       Double_t qt1 = Double_t(chargeQT1->At(ipmt));
149       if((qt1-qt0)>0)  adc[ipmt] = TMath::Exp( Double_t (channelWidth*(qt1-qt0)/1000));
150       time[ipmt] = calib-> WalkCorrection( ipmt,Int_t(qt1) , timeCFD->At(ipmt) ) ;
151       
152       //LED
153       Double_t sl = (timeLED->At(ipmt) - time[ipmt])*channelWidth;
154       Double_t qt=((TGraph*)fAmpLEDrec.At(ipmt))->Eval(sl/1000.);
155       frecpoints->SetTime(ipmt,time[ipmt]);
156       frecpoints->SetAmp(ipmt,adc[ipmt]);
157       frecpoints->SetAmpLED(ipmt,qt);
158     }
159     else {
160       time[ipmt] = 0;
161       adc[ipmt] = 0;
162     }
163   }
164   
165   for (Int_t ipmt=0; ipmt<12; ipmt++){
166     if(time[ipmt] > 1 ) {
167       if(time[ipmt]<besttimeC){
168         besttimeC=time[ipmt]; //timeC
169         pmtBestC=ipmt;
170       }
171     }
172   }
173   for ( Int_t ipmt=12; ipmt<24; ipmt++){
174     if(time[ipmt] > 1) {
175       if(time[ipmt]<besttimeA) {
176         besttimeA=time[ipmt]; //timeA
177         pmtBestA=ipmt;}
178     }
179   }
180   if(besttimeA !=999999)  frecpoints->SetTimeBestA(Int_t(besttimeA));
181   if( besttimeC != 999999 ) frecpoints->SetTimeBestC(Int_t(besttimeC));
182   AliDebug(1,Form(" besttimeA %f ps,  besttimeC %f ps",besttimeA, besttimeC));
183   Float_t c = 0.0299792; // cm/ps
184   Float_t vertex = 0;
185   if(besttimeA !=999999 && besttimeC != 999999 ){
186     timeDiff =(besttimeC - besttimeA)*channelWidth;
187     meanTime = (meanT0 - (besttimeA + besttimeC)/2) * channelWidth;
188     vertex = c*(timeDiff)/2. + (fdZonA - fdZonC)/2; //-(lenr-lenl))/2;
189     AliDebug(1,Form("  timeDiff %f ps,  meanTime %f ps, vertex %f cm",timeDiff, meanTime,vertex ));
190     frecpoints->SetVertex(vertex);
191     frecpoints->SetMeanTime(Int_t(meanTime));
192     
193   }
194   //time in each channel as time[ipmt]-MeanTimeinThisChannel(with vertex=0)
195   for (Int_t ipmt=0; ipmt<24; ipmt++) {
196     if(time[ipmt]>1) {
197 //      time[ipmt] = (time[ipmt] - fTime0vertex[ipmt])*channelWidth;
198       time[ipmt] = time[ipmt] * channelWidth;
199       frecpoints->SetTime(ipmt,time[ipmt]);
200     }
201   }
202   clustersTree->Fill();
203
204   delete timeCFD;
205   delete timeLED;
206   delete chargeQT0; 
207   delete chargeQT1; 
208 }
209
210
211 //_______________________________________________________________________
212
213 void AliT0Reconstructor::Reconstruct(AliRawReader* rawReader, TTree*recTree) const
214 {
215   // T0 raw ->
216   // T0RecPoint writing 
217   
218   //Q->T-> coefficients !!!! should be measured!!!
219   Int_t allData[110][5];
220   
221   TArrayI * timeCFD = new TArrayI(24); 
222   TArrayI * timeLED = new TArrayI(24); 
223   TArrayI * chargeQT0 = new TArrayI(24); 
224   TArrayI * chargeQT1 = new TArrayI(24); 
225
226    for (Int_t i=0; i<110; i++) {
227      allData[i][0]=0;
228    }
229
230    AliT0RawReader myrawreader(rawReader);
231    if (!myrawreader.Next())
232      AliDebug(1,Form(" no raw data found!! %i", myrawreader.Next()));
233    for (Int_t i=0; i<110; i++) {
234      allData[i][0]=myrawreader.GetData(i,0);
235    }
236   AliT0Calibrator *calib = new AliT0Calibrator(); 
237
238   //  Int_t mV2Mip = param->GetmV2Mip();     
239   //mV2Mip = param->GetmV2Mip();     
240   Float_t channelWidth = fParam->GetChannelWidth() ;  
241
242   Int_t meanT0 = fParam->GetMeanT0();
243  
244   for (Int_t in=0; in<24; in++)
245      {
246        timeLED->AddAt(allData[in+1][0],in);
247        timeCFD->AddAt(allData[in+25][0],in);
248        chargeQT1->AddAt(allData[in+57][0],in);
249        chargeQT0->AddAt(allData[in+80][0],in);
250        AliDebug(10, Form(" readed Raw %i %i %i %i %i", in, timeLED->At(in),timeCFD->At(in),chargeQT0->At(in),chargeQT1->At(in)));
251      }
252
253   Float_t besttimeA=999999;
254   Float_t besttimeC=999999;
255   Int_t pmtBestA=99999;
256   Int_t pmtBestC=99999;
257   Float_t timeDiff=999999, meanTime=0;
258
259   
260   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
261  
262   recTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
263   
264
265   Float_t time[24], adc[24];
266   for (Int_t ipmt=0; ipmt<24; ipmt++) {
267     if(timeCFD->At(ipmt)>0 ){
268       Double_t qt0 = Double_t(chargeQT0->At(ipmt));
269       Double_t qt1 = Double_t(chargeQT1->At(ipmt));
270       if((qt1-qt0)>0)  adc[ipmt] = TMath::Exp( Double_t (channelWidth*(qt1-qt0)/1000));
271       //      time[ipmt] = channelWidth * (calib-> WalkCorrection( ipmt,qt1 , timeCFD->At(ipmt) ) ) ;
272       time[ipmt] = calib-> WalkCorrection( ipmt,Int_t(qt1) , timeCFD->At(ipmt) ) ;
273       Double_t sl = (timeLED->At(ipmt) - time[ipmt])*channelWidth;
274       Double_t qt=((TGraph*)fAmpLEDrec.At(ipmt))->Eval(sl/1000.);
275       frecpoints->SetTime(ipmt,time[ipmt]);
276       frecpoints->SetAmp(ipmt,adc[ipmt]);
277       frecpoints->SetAmpLED(ipmt,qt);
278     AliDebug(1,Form(" QTC %f mv,  QTC  %f MIPS time in chann %f time %f ",adc[ipmt], adc[ipmt]/50.,time[ipmt], time[ipmt]*channelWidth));
279      
280     }
281     else {
282       time[ipmt] = 0;
283       adc[ipmt] = 0;
284     }
285   }
286
287   for (Int_t ipmt=0; ipmt<12; ipmt++){
288     if(time[ipmt] > 1 ) {
289       if(time[ipmt]<besttimeC){
290         besttimeC=time[ipmt]; //timeC
291         pmtBestC=ipmt;
292       }
293     }
294   }
295   for ( Int_t ipmt=12; ipmt<24; ipmt++){
296     if(time[ipmt] > 1) {
297       if(time[ipmt]<besttimeA) {
298         besttimeA=time[ipmt]; //timeA
299         pmtBestA=ipmt;}
300     }
301   }
302   if(besttimeA !=999999)  frecpoints->SetTimeBestA(Int_t(besttimeA));
303   if( besttimeC != 999999 ) frecpoints->SetTimeBestC(Int_t(besttimeC));
304   AliDebug(1,Form(" besttimeA %f ps,  besttimeC %f ps",besttimeA, besttimeC));
305   Float_t c = 0.0299792; // cm/ps
306   Float_t vertex = 0;
307   if(besttimeA !=999999 && besttimeC != 999999 ){
308     timeDiff = (besttimeC - besttimeA)*channelWidth;
309     meanTime = (meanT0 - (besttimeA + besttimeC)/2) * channelWidth;   
310     vertex = c*(timeDiff)/2. + (fdZonA - fdZonC)/2; //-(lenr-lenl))/2;
311     AliDebug(1,Form("  timeDiff %f ps,  meanTime %f ps, vertex %f cm",timeDiff, meanTime,vertex ));
312     frecpoints->SetVertex(vertex);
313     frecpoints->SetMeanTime(Int_t(meanTime));
314     
315   }
316   //time in each channel as time[ipmt]-MeanTimeinThisChannel(with vertex=0)
317   /*
318   for (Int_t ipmt=0; ipmt<24; ipmt++) {
319     if(time[ipmt]>1) {
320       time[ipmt] = (time[ipmt] - fTime0vertex[ipmt])*channelWidth;
321       frecpoints->SetTime(ipmt,time[ipmt]);
322     }
323     }
324   */
325   recTree->Fill();
326  
327
328   delete timeCFD;
329   delete timeLED;
330   delete chargeQT0; 
331   delete chargeQT1; 
332   
333 }
334 //____________________________________________________________
335
336 void AliT0Reconstructor::FillESD(TTree */*digitsTree*/, TTree *clustersTree, AliESDEvent *pESD) const
337 {
338
339   /***************************************************
340   Resonstruct digits to vertex position
341   ****************************************************/
342   
343   AliDebug(1,Form("Start FillESD T0"));
344
345   TTree *treeR = clustersTree;
346   
347    AliT0RecPoint* frecpoints= new AliT0RecPoint ();
348     if (!frecpoints) {
349     AliError("Reconstruct Fill ESD >> no recpoints found");
350     return;
351   }
352   
353   AliDebug(1,Form("Start FillESD T0"));
354   TBranch *brRec = treeR->GetBranch("T0");
355   if (brRec) {
356     brRec->SetAddress(&frecpoints);
357   }else{
358     AliError(Form("EXEC Branch T0 rec not found"));
359     return;
360   } 
361     
362     brRec->GetEntry(0);
363     Float_t amp[24], time[24];
364     Float_t  zPosition = frecpoints -> GetVertex();
365     Float_t timeStart = frecpoints -> GetMeanTime() ;
366     for ( Int_t i=0; i<24; i++) {
367       time[i] = Float_t (frecpoints -> GetTime(i)); // ps to ns
368       amp[i] = frecpoints -> GetAmp(i);
369     }
370     pESD->SetT0zVertex(zPosition); //vertex Z position 
371     pESD->SetT0(timeStart);        // interaction time 
372     pESD->SetT0time(time);         // best TOF on each PMT 
373     pESD->SetT0amplitude(amp);     // number of particles(MIPs) on each PMT
374  
375     AliDebug(1,Form(" Z position %f cm,  T0  %f ps",zPosition , timeStart));
376
377 } // vertex in 3 sigma
378
379
380
381
382
383