]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0Reconstructor.cxx
reading latency from OCDB
[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 #include "AliESDfriend.h"
34
35 #include <TArrayI.h>
36 #include <TGraph.h>
37 #include <TMath.h>
38 #include <Riostream.h>
39
40 ClassImp(AliT0Reconstructor)
41
42   AliT0Reconstructor:: AliT0Reconstructor(): AliReconstructor(),
43                                              fdZonA(0),
44                                              fdZonC(0),
45                                              fZposition(0),
46                                              fParam(NULL),
47                                              fAmpLEDrec(),
48                                              fQTC(0),
49                                              fAmpLED(0),
50                                              fCalib(),
51                                              fLatencyHPTDC(9000),
52                                              fLatencyL1(0),
53                                              fLatencyL1A(0),
54                                              fLatencyL1C(0)
55
56 {
57   //constructor
58
59   fParam = AliT0Parameters::Instance();
60   fParam->Init();
61  
62   for (Int_t i=0; i<24; i++){
63         TGraph* gr = fParam ->GetAmpLEDRec(i);
64         if (gr) fAmpLEDrec.AddAtAndExpand(gr,i) ; 
65           TGraph* gr1 = fParam ->GetAmpLED(i);
66           if (gr1) fAmpLED.AddAtAndExpand(gr1,i) ; 
67           TGraph* gr2 = fParam ->GetQTC(i);
68           if (gr2) fQTC.AddAtAndExpand(gr2,i) ;         
69   }
70
71   fLatencyL1 = fParam->GetLatencyL1();
72   fLatencyL1A = fParam->GetLatencyL1A();
73   fLatencyL1C = fParam->GetLatencyL1C();
74   fLatencyHPTDC = fParam->GetLatencyHPTDC();
75   AliDebug(10,Form(" LatencyL1 %f latencyL1A %f latencyL1C %f latencyHPTDC %f \n",fLatencyL1, fLatencyL1A, fLatencyL1C, fLatencyHPTDC));
76   
77   // fdZonC = TMath::Abs(fParam->GetZPositionShift("T0/C/PMT1"));
78   //fdZonA = TMath::Abs(fParam->GetZPositionShift("T0/A/PMT15"));
79
80
81   fCalib = new AliT0Calibrator();
82
83 }
84
85 //_____________________________________________________________________________
86 void AliT0Reconstructor::Reconstruct(TTree*digitsTree, TTree*clustersTree) const
87   
88 {
89   // T0 digits reconstruction
90   Int_t refAmp = GetRecoParam()->GetRefAmp();
91   
92   TArrayI * timeCFD = new TArrayI(24); 
93   TArrayI * timeLED = new TArrayI(24); 
94   TArrayI * chargeQT0 = new TArrayI(24); 
95   TArrayI * chargeQT1 = new TArrayI(24); 
96
97  
98   Float_t channelWidth = fParam->GetChannelWidth() ;  
99   Float_t meanVertex = fParam->GetMeanVertex();
100   Float_t c = 0.0299792; // cm/ps
101   Double32_t vertex = 9999999;
102   Double32_t timeDiff=999999, meanTime=999999, timeclock=999999;
103
104   
105   AliDebug(1,Form("Start DIGITS reconstruction "));
106   
107
108   TBranch *brDigits=digitsTree->GetBranch("T0");
109   AliT0digit *fDigits = new AliT0digit() ;
110   if (brDigits) {
111     brDigits->SetAddress(&fDigits);
112   }else{
113     AliError(Form("EXEC Branch T0 digits not found"));
114      return;
115   }
116   
117   digitsTree->GetEvent(0);
118   digitsTree->GetEntry(0);
119   brDigits->GetEntry(0);
120   fDigits->GetTimeCFD(*timeCFD);
121   fDigits->GetTimeLED(*timeLED);
122   fDigits->GetQT0(*chargeQT0);
123   fDigits->GetQT1(*chargeQT1);
124   Int_t onlineMean =  fDigits->MeanTime();
125
126   Bool_t tr[5];
127   for (Int_t i=0; i<5; i++) tr[i]=false; 
128   
129   Double32_t besttimeA=999999;
130   Double32_t besttimeC=999999;
131   Int_t pmtBestA=99999;
132   Int_t pmtBestC=99999;
133   
134   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
135   clustersTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
136   
137   Float_t time[24], adc[24];
138   for (Int_t ipmt=0; ipmt<24; ipmt++) {
139     if(timeCFD->At(ipmt)>0 ){
140      if(( chargeQT1->At(ipmt) - chargeQT0->At(ipmt))>0)  
141         adc[ipmt] = chargeQT1->At(ipmt) - chargeQT0->At(ipmt);
142       else
143         adc[ipmt] = 0;
144       
145      // time[ipmt] = fCalib-> WalkCorrection(refAmp, ipmt, adc[ipmt],  timeCFD->At(ipmt)) ;
146              
147       Double_t sl = Double_t(timeLED->At(ipmt) - timeCFD->At(ipmt));
148       time[ipmt] = fCalib-> WalkCorrection( refAmp,ipmt, Int_t(sl),  timeCFD->At(ipmt) ) ;
149       AliDebug(10,Form(" ipmt %i QTC %i , time in chann %i (led-cfd) %i ",
150                        ipmt, Int_t(adc[ipmt]) ,Int_t(time[ipmt]),Int_t( sl)));
151
152       Double_t ampMip =((TGraph*)fAmpLED.At(ipmt))->Eval(sl);
153       Double_t qtMip = ((TGraph*)fQTC.At(ipmt))->Eval(adc[ipmt]);
154       AliDebug(10,Form("  Amlitude in MIPS LED %f ,  QTC %f in channels %i\n ",ampMip,qtMip, adc[ipmt]));
155       
156       frecpoints->SetTime(ipmt, Float_t(time[ipmt]) );
157       frecpoints->SetAmp(ipmt, Float_t( ampMip)); //for cosmic &pp beam 
158       frecpoints->SetAmpLED(ipmt, Float_t(qtMip));
159       
160     }
161     else {
162       time[ipmt] = 0;
163       adc[ipmt] = 0;
164     }
165   }
166   
167   for (Int_t ipmt=0; ipmt<12; ipmt++){
168     if(time[ipmt] > 1 ) {
169       if(time[ipmt]<besttimeC){
170         besttimeC=time[ipmt]; //timeC
171         pmtBestC=ipmt;
172       }
173     }
174   }
175   for ( Int_t ipmt=12; ipmt<24; ipmt++){
176     if(time[ipmt] > 1) {
177       if(time[ipmt]<besttimeA) {
178         besttimeA=time[ipmt]; //timeA
179         pmtBestA=ipmt;}
180     }
181   }
182   if(besttimeA < 999999) {
183     frecpoints->SetTimeBestA(Int_t(besttimeA *channelWidth));
184     tr[1]=true;
185   }
186   if( besttimeC < 999999 ) {
187     frecpoints->SetTimeBestC(Int_t(besttimeC *channelWidth));
188     tr[2]=true;
189   }
190   AliDebug(10,Form(" besttimeA %f ch,  besttimeC %f ch",besttimeA, besttimeC));
191   if(besttimeA <999999 && besttimeC < 999999 ){
192     //    timeDiff = (besttimeC - besttimeA)*channelWidth;
193     timeDiff = (besttimeA - besttimeC)*channelWidth;
194     meanTime = (besttimeA + besttimeC)/2;// * channelWidth); 
195     timeclock = meanTime *channelWidth ;
196     vertex = meanVertex - c*(timeDiff)/2.;// + (fdZonA - fdZonC)/2;
197     tr[0]=true; 
198   }
199   frecpoints->SetVertex(vertex);
200   frecpoints->SetMeanTime(meanTime);
201   frecpoints->SetT0clock(timeclock);
202   frecpoints->SetT0Trig(tr);
203
204   for (Int_t i=0; i<5; i++) {
205     printf(" T0 trigers %i ",tr[i]);
206   }
207     printf(" \n ");
208
209   //online mean
210   frecpoints->SetOnlineMean(Int_t(onlineMean));
211   AliDebug(10,Form("  timeDiff %i #channel,  meanTime %i #channel, vertex %f cm online mean %i timeclock %i ps",timeDiff, meanTime,vertex, Int_t(onlineMean), timeclock));
212   
213   
214
215    
216   
217   clustersTree->Fill();
218
219   delete timeCFD;
220   delete timeLED;
221   delete chargeQT0; 
222   delete chargeQT1; 
223 }
224
225
226 //_______________________________________________________________________
227
228 void AliT0Reconstructor::Reconstruct(AliRawReader* rawReader, TTree*recTree) const
229 {
230   // T0 raw ->
231   //
232   // reference amplitude and time ref. point from reco param
233
234   Int_t refAmp = GetRecoParam()->GetRefAmp();
235   Int_t refPoint = GetRecoParam()->GetRefPoint();
236
237   Int_t allData[110][5];
238   
239   Int_t timeCFD[24], timeLED[24], chargeQT0[24], chargeQT1[24];
240   Double32_t timeDiff=999999, meanTime=999999, timeclock=999999;
241   Float_t c = 29.9792458; // cm/ns
242   Double32_t vertex = 9999999;
243   Int_t onlineMean=0;
244   // Float_t meanVertex = fParam->GetMeanVertex();
245   Float_t meanVertex = 0;
246
247   cout<<" @@@@ Latency "<<fLatencyL1<<endl;
248   for (Int_t i0=0; i0<105; i0++)
249     {
250       for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;   
251     }
252    
253   Double32_t besttimeA=9999999;
254   Double32_t besttimeC=9999999;
255   Int_t pmtBestA=99999;
256   Int_t pmtBestC=99999;
257    
258   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
259   
260   recTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
261    
262   AliDebug(10," before read data ");
263   AliT0RawReader myrawreader(rawReader);
264
265   UInt_t type =rawReader->GetType();
266
267   if (!myrawreader.Next())
268     AliDebug(1,Form(" no raw data found!!"));
269   else
270     {  
271       if(type == 7) {  //only physics 
272       for (Int_t i=0; i<105; i++) {
273         for (Int_t iHit=0; iHit<5; iHit++) 
274           {
275             allData[i][iHit] = myrawreader.GetData(i,iHit);
276           }
277       }
278       Int_t ref=0;
279       if (refPoint>0) 
280       ref = allData[refPoint][0]-5000;
281
282       Float_t channelWidth = fParam->GetChannelWidth() ;  
283       
284       //       Int_t meanT0 = fParam->GetMeanT0();
285        
286           
287           for (Int_t in=0; in<12; in++)  
288             {
289               timeCFD[in] = allData[in+1][0] ;
290               timeCFD[in+12] = allData[in+56+1][0] ;
291               timeLED[in] = allData[in+12+1][0] ;
292               timeLED[in+12] = allData[in+68+1][0] ;
293               AliDebug(10, Form(" readed i %i cfdC %i cfdA %i ledC %i ledA%i ",
294                                 in, timeCFD[in],timeCFD[in+12],timeLED[in], 
295                                 timeLED[in+12]));   
296             }
297           
298           for (Int_t in=0; in<12;  in++)
299             {
300               chargeQT0[in]=allData[2*in+25][0];
301               chargeQT1[in]=allData[2*in+26][0];
302             }
303           
304            for (Int_t in=12; in<24;  in++)
305              {
306                chargeQT0[in]=allData[2*in+57][0];
307                chargeQT1[in]=allData[2*in+58][0];
308              }
309            
310            //    } //cosmic with physics event
311        for (Int_t in=0; in<24; in++)  
312          AliDebug(10, Form(" readed Raw %i %i %i %i %i",
313                            in, timeLED[in],timeCFD[in],chargeQT0[in],chargeQT1[in]));
314         onlineMean = allData[49][0];       
315        
316        Double32_t time[24], adc[24];
317        for (Int_t ipmt=0; ipmt<24; ipmt++) {
318          if(timeCFD[ipmt]>0 && timeLED[ipmt]>0){
319            //for simulated data
320              //for physics  data
321            if(( chargeQT1[ipmt] - chargeQT0[ipmt])>0)  
322              adc[ipmt] = chargeQT0[ipmt] - chargeQT1[ipmt];
323            else
324              adc[ipmt] = 0;
325            
326
327            //      time[ipmt] = fCalib-> WalkCorrection(refAmp, ipmt, adc[ipmt], timeCFD[ipmt] ) ;
328            
329            Double_t sl = timeLED[ipmt] - timeCFD[ipmt];
330            time[ipmt] = fCalib-> WalkCorrection( refAmp,ipmt, Int_t(sl), timeCFD[ipmt] ) ;
331            AliDebug(10,Form(" ipmt %i QTC %i , time in chann %i (led-cfd) %i ",
332                             ipmt, Int_t(adc[ipmt]) ,Int_t(time[ipmt]),Int_t( sl)));
333            Double_t ampMip =( (TGraph*)fAmpLED.At(ipmt))->Eval(sl);
334            Double_t qtMip = ((TGraph*)fQTC.At(ipmt))->Eval(adc[ipmt]);
335            AliDebug(10,Form("  Amlitude in MIPS LED %f ; QTC %f;  in channels %i\n ",ampMip,qtMip, adc[ipmt]));
336              
337            frecpoints->SetTime(ipmt, Float_t(time[ipmt]) );
338            // frecpoints->SetTime(ipmt,Double32_t(timeCFD[ipmt]));
339            frecpoints->SetAmpLED(ipmt, Double32_t( qtMip)); //for cosmic &pp beam 
340            frecpoints->SetAmp(ipmt, Double32_t(ampMip));
341              
342          }
343          else {
344            time[ipmt] = 0;
345            adc[ipmt] = 0;
346          }
347        }
348        
349        for (Int_t ipmt=0; ipmt<12; ipmt++){
350          if(time[ipmt] > 1 ) {
351            if(time[ipmt]<besttimeC){
352              besttimeC=time[ipmt]; //timeC
353              pmtBestC=ipmt;
354            }
355          }
356        }
357        for ( Int_t ipmt=12; ipmt<24; ipmt++){
358          if(time[ipmt] > 1) {
359            if(time[ipmt]<besttimeA) {
360              besttimeA=time[ipmt]; //timeA
361              pmtBestA=ipmt;}
362          }
363        }
364        if(besttimeA < 999999) 
365          frecpoints->SetTimeBestA(0.001* besttimeA * channelWidth - fLatencyHPTDC + fLatencyL1A);
366        if( besttimeC < 999999 ) 
367          frecpoints->SetTimeBestC( 0.001 *besttimeC * channelWidth - fLatencyHPTDC + fLatencyL1C);
368        AliDebug(10,Form(" pmtA %i besttimeA %f ps, pmtC %i besttimeC %f ps",
369                        pmtBestA,besttimeA, pmtBestC,  besttimeC));
370         if(besttimeA <999999 && besttimeC < 999999 ){
371          timeDiff = ( besttimeA - besttimeC) *0.001 * channelWidth + fLatencyL1A - fLatencyL1C;
372          timeclock = 0.001*channelWidth * Float_t( besttimeA+besttimeC)/2. - fLatencyHPTDC + fLatencyL1;  
373          meanTime = (besttimeA+besttimeC-2.*Float_t(ref))/2.;
374          vertex =  meanVertex - c*(timeDiff)/2. ; //+ (fdZonA - fdZonC)/2; 
375         }
376       }  //if phys event       
377       AliDebug(5,Form("  timeDiff %f #channel,  meanTime %f #channel, TOFmean%f  vertex %f cm meanVertex %f online mean %i \n",timeDiff, meanTime,timeclock, vertex,meanVertex, onlineMean));
378       frecpoints->SetT0clock(timeclock);
379       frecpoints->SetVertex(vertex);
380       frecpoints->SetMeanTime(meanTime);
381       frecpoints->SetOnlineMean(Int_t(onlineMean));
382         // Set triggers
383       
384       Bool_t tr[5];
385       Int_t trchan[5]= {50,51,52,55,56};
386       for (Int_t i=0; i<5; i++) tr[i]=false; 
387       for (Int_t itr=0; itr<5; itr++) {
388         if(allData[trchan[itr]][0]>0) tr[itr]=true;
389         frecpoints->SetT0Trig(tr);
390       }
391       
392     } // if (else )raw data
393   recTree->Fill();
394   if(frecpoints) delete frecpoints;
395 }
396   
397   
398   //____________________________________________________________
399   
400   void AliT0Reconstructor::FillESD(TTree */*digitsTree*/, TTree *clustersTree, AliESDEvent *pESD) const
401   {
402
403   /***************************************************
404   Resonstruct digits to vertex position
405   ****************************************************/
406   
407   AliDebug(1,Form("Start FillESD T0"));
408   Float_t channelWidth = fParam->GetChannelWidth() ;  
409   Float_t c = 29.9792458; // cm/ns
410   Float_t currentVertex=0, shift=0;
411   Int_t ncont=0;
412   const AliESDVertex* vertex = pESD->GetPrimaryVertex();
413   if (!vertex)        vertex = pESD->GetPrimaryVertexSPD();
414   if (!vertex)        vertex = pESD->GetPrimaryVertexTPC();
415   if (!vertex)        vertex = pESD->GetVertex();
416
417   if (vertex) {
418     AliDebug(2, Form("Got %s (%s) from ESD: %f", 
419                     vertex->GetName(), vertex->GetTitle(), vertex->GetZ()));
420     currentVertex = vertex->GetZ();
421     
422     ncont = vertex->GetNContributors();
423     // cout<<" spdver "<<spdver<<" ncont "<<ncont<<endl;
424     if(ncont>2 ) {
425       shift = currentVertex/c;
426       //          cout<<" vertex shif "<<shift<<" vertex "<<spdver<<" IsFromVertexer3D  "<<fverSPD->IsFromVertexer3D()<<endl;
427     }
428   }
429   TTree *treeR = clustersTree;
430   
431    AliT0RecPoint* frecpoints= new AliT0RecPoint ();
432     if (!frecpoints) {
433     AliError("Reconstruct Fill ESD >> no recpoints found");
434     return;
435   }
436   
437   AliDebug(1,Form("Start FillESD T0"));
438   TBranch *brRec = treeR->GetBranch("T0");
439   if (brRec) {
440     brRec->SetAddress(&frecpoints);
441   }else{
442     AliError(Form("EXEC Branch T0 rec not found"));
443     return;
444   } 
445     
446     brRec->GetEntry(0);
447     Double32_t amp[24], time[24], ampQTC[24], timecorr[24];  
448     Double32_t timeClock[3];
449     Double32_t zPosition = frecpoints -> GetVertex();
450     Double32_t timeStart = frecpoints -> GetMeanTime();
451     timeClock[0] = frecpoints -> GetT0clock() ;
452     timeClock[1] = frecpoints -> GetBestTimeA() + shift;
453     timeClock[2] = frecpoints -> GetBestTimeC() - shift;
454     for ( Int_t i=0; i<24; i++) {
455       time[i] =  frecpoints -> GetTime(i); // ps to ns
456       amp[i] = frecpoints -> GetAmp(i);
457       ampQTC[i] = frecpoints -> AmpLED(i);
458     }
459     Int_t trig= frecpoints ->GetT0Trig();
460     pESD->SetT0Trig(trig);
461     
462     pESD->SetT0zVertex(zPosition); //vertex Z position 
463     pESD->SetT0(timeStart);        // interaction time 
464     for(Int_t i=0; i<3; i++) 
465       pESD->SetT0TOF(i,timeClock[i]);   // interaction time (ns) 
466     pESD->SetT0time(time);         // best TOF on each PMT 
467     pESD->SetT0amplitude(amp);     // number of particles(MIPs) on each PMT
468     
469     AliDebug(1,Form("T0: Vertex %f (T0A+T0C)/2 %f #channels T0signal %f ns OrA %f ns OrC %f T0trig %i\n",zPosition, timeStart, timeClock[0], timeClock[1], timeClock[2], trig));
470
471     /* if (pESD) {
472      AliESDfriend *fr = (AliESDfriend*)pESD->FindListObject("AliESDfriend");
473      if (fr) {
474         AliDebug(1, Form("Writing TZERO friend data to ESD tree"));
475
476         if (ncont>2) {
477           for ( Int_t i=0; i<24; i++) {
478             if(i<12) timecorr[i]=time[i]-shift*channelWidth;
479             if(i>11) timecorr[i]=time[i]-shift*channelWidth;
480             fr->SetT0timeCorr(timecorr) ;
481           }
482             fr->SetT0ampLEDminCFD(amp);
483             fr->SetT0ampQTC(ampQTC);
484         }
485     }
486   }
487     */
488     
489     
490 } // vertex in 3 sigma
491
492
493
494
495
496