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