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