with TOF satisfied request about writng in ESD T0 spread
[u/mrichter/AliRoot.git] / T0 / AliT0Reconstructor.cxx
1
2
3 /**************************************************************************
4  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Author: The ALICE Off-line Project.                                    *
7  * Contributors are mentioned in the code where appropriate.              *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 /* $Id$ */
19 /*********************************************************************
20  *  T0 reconstruction and filling ESD
21  *  - reconstruct mean time (interation time) 
22  *  - vertex position
23  *  -  multiplicity
24  ********************************************************************/
25
26 #include <AliESDEvent.h>
27 #include "AliLog.h"
28 #include "AliT0RecPoint.h"
29 #include "AliRawReader.h"
30 #include "AliT0RawReader.h"
31 #include "AliT0digit.h"
32 #include "AliT0Reconstructor.h"
33 #include "AliT0Parameters.h"
34 #include "AliT0Calibrator.h"
35 #include "AliESDfriend.h"
36 #include "AliESDTZEROfriend.h"
37 #include "AliLog.h"
38 #include "AliCDBEntry.h" 
39 #include "AliCDBManager.h"
40 #include "AliCTPTimeParams.h"
41 #include "AliLHCClockPhase.h"
42 #include "AliT0CalibSeasonTimeShift.h"
43 #include "AliESDRun.h"
44
45 #include <TArrayI.h>
46 #include <TGraph.h>
47 #include <TMath.h>
48 #include <Riostream.h>
49
50 ClassImp(AliT0Reconstructor)
51
52   AliT0Reconstructor:: AliT0Reconstructor(): AliReconstructor(),
53                                              fdZonA(0),
54                                              fdZonC(0),
55                                              fZposition(0),
56                                              fParam(NULL),
57                                              fAmpLEDrec(),
58                                              fQTC(0),
59                                              fAmpLED(0),
60                                              fCalib(),
61                                              fLatencyHPTDC(9000),
62                                              fLatencyL1(0),
63                                              fLatencyL1A(0),
64                                              fLatencyL1C(0),
65                                              fGRPdelays(0),
66                                              fTimeMeanShift(0x0),
67                                              fTimeSigmaShift(0x0),
68                                              fESDTZEROfriend(NULL)
69
70 {
71
72
73   //constructor
74   AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
75   if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
76   AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
77   Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
78
79   AliCDBEntry *entry1 = AliCDBManager::Instance()->Get("GRP/CTP/TimeAlign");
80   if (!entry1) AliFatal("CTP time-alignment is not found in OCDB !");
81   AliCTPTimeParams *ctpTimeAlign = (AliCTPTimeParams*)entry1->GetObject();
82   l1Delay += ((Float_t)ctpTimeAlign->GetDelayL1L0()*25.0);
83  
84   AliCDBEntry *entry4 = AliCDBManager::Instance()->Get("GRP/Calib/LHCClockPhase");
85   if (!entry4) AliFatal("LHC clock-phase shift is not found in OCDB !");
86   AliLHCClockPhase *phase = (AliLHCClockPhase*)entry4->GetObject();
87
88   fGRPdelays = l1Delay - phase->GetMeanPhase();
89
90   AliCDBEntry *entry5 = AliCDBManager::Instance()->Get("T0/Calib/TimeAdjust");
91   if (entry5) {
92     AliT0CalibSeasonTimeShift *timeshift = (AliT0CalibSeasonTimeShift*)entry5->GetObject();
93     fTimeMeanShift = timeshift->GetT0Means();
94     fTimeSigmaShift  = timeshift->GetT0Sigmas();
95     //   timeshift->Print();
96    }
97   else
98     AliWarning("Time Adjust is not found in OCDB !");
99
100   fParam = AliT0Parameters::Instance();
101   fParam->Init();
102  
103   for (Int_t i=0; i<24; i++){
104         TGraph* gr = fParam ->GetAmpLEDRec(i);
105         if (gr) fAmpLEDrec.AddAtAndExpand(gr,i) ; 
106           TGraph* gr1 = fParam ->GetAmpLED(i);
107           if (gr1) fAmpLED.AddAtAndExpand(gr1,i) ; 
108           TGraph* gr2 = fParam ->GetQTC(i);
109           if (gr2) fQTC.AddAtAndExpand(gr2,i) ;         
110   }
111
112   fLatencyL1 = fParam->GetLatencyL1();
113   fLatencyL1A = fParam->GetLatencyL1A();
114   fLatencyL1C = fParam->GetLatencyL1C();
115   fLatencyHPTDC = fParam->GetLatencyHPTDC();
116   AliDebug(2,Form(" LatencyL1 %f latencyL1A %f latencyL1C %f latencyHPTDC %f \n",fLatencyL1, fLatencyL1A, fLatencyL1C, fLatencyHPTDC));
117   
118   // fdZonC = TMath::Abs(fParam->GetZPositionShift("T0/C/PMT1"));
119   //fdZonA = TMath::Abs(fParam->GetZPositionShift("T0/A/PMT15"));
120   //here real Z position
121   fdZonC = TMath::Abs(fParam->GetZPosition("T0/C/PMT1"));
122   fdZonA = TMath::Abs(fParam->GetZPosition("T0/A/PMT15"));
123
124   fCalib = new AliT0Calibrator();
125   fESDTZEROfriend = new AliESDTZEROfriend();
126
127 }
128
129 //_____________________________________________________________________________
130 void AliT0Reconstructor::Reconstruct(TTree*digitsTree, TTree*clustersTree) const
131   
132 {
133   // T0 digits reconstruction
134   Int_t refAmp = Int_t (GetRecoParam()->GetRefAmp());
135   
136   TArrayI * timeCFD = new TArrayI(24); 
137   TArrayI * timeLED = new TArrayI(24); 
138   TArrayI * chargeQT0 = new TArrayI(24); 
139   TArrayI * chargeQT1 = new TArrayI(24); 
140
141  
142   Float_t channelWidth = fParam->GetChannelWidth() ;  
143   Float_t meanVertex = fParam->GetMeanVertex();
144   Float_t c = 0.0299792; // cm/ps
145   Double32_t vertex = 9999999;
146   Double32_t timeDiff=999999, meanTime=999999, timeclock=999999;
147
148   
149   AliDebug(1,Form("Start DIGITS reconstruction "));
150   
151
152   TBranch *brDigits=digitsTree->GetBranch("T0");
153   AliT0digit *fDigits = new AliT0digit() ;
154   if (brDigits) {
155     brDigits->SetAddress(&fDigits);
156   }else{
157     AliError(Form("EXEC Branch T0 digits not found"));
158      return;
159   }
160   
161   digitsTree->GetEvent(0);
162   digitsTree->GetEntry(0);
163   brDigits->GetEntry(0);
164   fDigits->GetTimeCFD(*timeCFD);
165   fDigits->GetTimeLED(*timeLED);
166   fDigits->GetQT0(*chargeQT0);
167   fDigits->GetQT1(*chargeQT1);
168   Int_t onlineMean =  fDigits->MeanTime();
169
170   Bool_t tr[5];
171   for (Int_t i=0; i<5; i++) tr[i]=false; 
172   
173   Double32_t besttimeA=999999;
174   Double32_t besttimeC=999999;
175   Int_t pmtBestA=99999;
176   Int_t pmtBestC=99999;
177   
178   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
179   clustersTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints);
180   
181   Float_t time[24], adc[24];
182   for (Int_t ipmt=0; ipmt<24; ipmt++) {
183     if(timeCFD->At(ipmt)>0 ){
184      if(( chargeQT1->At(ipmt) - chargeQT0->At(ipmt))>0)  
185         adc[ipmt] = chargeQT1->At(ipmt) - chargeQT0->At(ipmt);
186       else
187         adc[ipmt] = 0;
188       
189      time[ipmt] = fCalib-> WalkCorrection(refAmp, ipmt, Int_t(adc[ipmt]),  timeCFD->At(ipmt)) ;
190              
191       Double_t sl = Double_t(timeLED->At(ipmt) - timeCFD->At(ipmt));
192       //    time[ipmt] = fCalib-> WalkCorrection( refAmp,ipmt, Int_t(sl),  timeCFD->At(ipmt) ) ;
193       AliDebug(5,Form(" ipmt %i QTC %i , time in chann %i (led-cfd) %i ",
194                        ipmt, Int_t(adc[ipmt]) ,Int_t(time[ipmt]),Int_t( sl)));
195
196       Double_t ampMip =((TGraph*)fAmpLED.At(ipmt))->Eval(sl);
197       Double_t qtMip = ((TGraph*)fQTC.At(ipmt))->Eval(adc[ipmt]);
198       AliDebug(5,Form("  Amlitude in MIPS LED %f ,  QTC %f in channels %f\n ",ampMip,qtMip, adc[ipmt]));
199       
200       frecpoints->SetTime(ipmt, Float_t(time[ipmt]) );
201       frecpoints->SetAmpLED(ipmt, Float_t( ampMip)); //for cosmic &pp beam 
202       frecpoints->SetAmp(ipmt, Float_t(qtMip));
203       
204     }
205     else {
206       time[ipmt] = 0;
207       adc[ipmt] = 0;
208     }
209   }
210   
211   for (Int_t ipmt=0; ipmt<12; ipmt++){
212     if(time[ipmt] > 1 ) {
213       if(time[ipmt]<besttimeC){
214         besttimeC=time[ipmt]; //timeC
215         pmtBestC=ipmt;
216       }
217     }
218   }
219   for ( Int_t ipmt=12; ipmt<24; ipmt++){
220     if(time[ipmt] > 1) {
221       if(time[ipmt]<besttimeA) {
222         besttimeA=time[ipmt]; //timeA
223         pmtBestA=ipmt;}
224     }
225   }
226   if(besttimeA < 999999) {
227     frecpoints->SetTimeBestA(Int_t(besttimeA *channelWidth - fdZonA/c));
228     tr[1]=true;
229   }
230   if( besttimeC < 999999 ) {
231     frecpoints->SetTimeBestC(Int_t(besttimeC *channelWidth - fdZonA/c));
232     tr[2]=true;
233   }
234   AliDebug(5,Form(" besttimeA %f ch,  besttimeC %f ch",besttimeA, besttimeC));
235   if(besttimeA <999999 && besttimeC < 999999 ){
236     //    timeDiff = (besttimeC - besttimeA)*channelWidth;
237     timeDiff = (besttimeA - besttimeC)*channelWidth;
238     meanTime = (besttimeA + besttimeC)/2;// * channelWidth); 
239     timeclock = meanTime *channelWidth -fdZonA/c ;
240     vertex = meanVertex - c*(timeDiff)/2.;// + (fdZonA - fdZonC)/2;
241     tr[0]=true; 
242   }
243   frecpoints->SetVertex(vertex);
244   frecpoints->SetMeanTime(meanTime);
245   frecpoints->SetT0clock(timeclock);
246   frecpoints->SetT0Trig(tr);
247
248   AliDebug(5,Form("T0 triggers %d %d %d %d %d",tr[0],tr[1],tr[2],tr[3],tr[4]));
249
250   //online mean
251   frecpoints->SetOnlineMean(Int_t(onlineMean));
252   AliDebug(10,Form("  timeDiff %f #channel,  meanTime %f #channel, vertex %f cm online mean %i timeclock %f ps",timeDiff, meanTime,vertex, Int_t(onlineMean), timeclock));
253   
254   
255
256    
257   
258   clustersTree->Fill();
259
260   delete timeCFD;
261   delete timeLED;
262   delete chargeQT0; 
263   delete chargeQT1; 
264 }
265
266
267 //_______________________________________________________________________
268
269 void AliT0Reconstructor::Reconstruct(AliRawReader* rawReader, TTree*recTree) const
270 {
271   // T0 raw ->
272   //
273   // reference amplitude and time ref. point from reco param
274
275   Float_t refAmp = GetRecoParam()->GetRefAmp();
276   Int_t refPoint = 0;
277   //Bad channel
278   Int_t badpmt = GetRecoParam()->GetRefPoint();
279   Int_t low[110], high[110];
280
281   Int_t allData[110][5];
282   
283   Int_t timeCFD[24], timeLED[24], chargeQT0[24], chargeQT1[24];
284   Double32_t timeDiff=999999, meanTime=999999, timeclock=9999999;
285   Float_t c = 29.9792458; // cm/ns
286   Double32_t vertex = 9999999;
287   Int_t onlineMean=0;
288   // Float_t meanVertex = fParam->GetMeanVertex();
289   Float_t meanVertex = 0;
290   for (Int_t i0=0; i0<110; i0++)
291     {
292       for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0; 
293       low[i0] = Int_t (GetRecoParam()->GetLow(i0));     
294       high[i0] = Int_t (GetRecoParam()->GetHigh(i0));
295       }
296    
297   Double32_t besttimeA=9999999;
298   Double32_t besttimeC=9999999;
299   Int_t pmtBestA=99999;
300   Int_t pmtBestC=99999;
301    
302   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
303   
304   recTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints);
305    
306   AliDebug(10," before read data ");
307   AliT0RawReader myrawreader(rawReader);
308
309   UInt_t type =rawReader->GetType();
310
311   if (!myrawreader.Next())
312     AliDebug(1,Form(" no raw data found!!"));
313   else
314     {  
315    for (Int_t i=0; i<24; i++)
316     {
317       timeCFD[i]=0; timeLED[i]=0; chargeQT0[i]=0; chargeQT1[i]=0;
318     }
319      Int_t fBCID=Int_t (rawReader->GetBCID());
320       Int_t trmbunch= myrawreader.GetTRMBunchID();
321       AliDebug(10,Form(" CDH BC ID %i, TRM BC ID %i \n", fBCID, trmbunch ));
322  
323       if(type == 7  ) {  //only physics 
324         for (Int_t i=0; i<107; i++) {
325         for (Int_t iHit=0; iHit<5; iHit++) 
326           {
327             allData[i][iHit] = myrawreader.GetData(i,iHit);
328           }
329         }
330         Int_t ref=0;
331         if (refPoint>0) 
332           ref = allData[refPoint][0]-5000;
333         
334         Float_t channelWidth = fParam->GetChannelWidth() ;  
335         
336         //       Int_t meanT0 = fParam->GetMeanT0();
337         
338         for (Int_t in=0; in<12; in++)  
339           {
340             for (Int_t iHit=0; iHit<5; iHit++) 
341               {
342                 if(allData[in+1][iHit] > low[in+1] && 
343                    allData[in+1][iHit] < high[in+1])
344                   {
345                     timeCFD[in] = allData[in+1][iHit] ; 
346                     break;
347                   }
348               }
349             for (Int_t iHit=0; iHit<5; iHit++) 
350               {
351                 if(allData[in+12+1][iHit] > low[in+12+1] && 
352                    allData[in+1+12][iHit] < high[in+12+1])
353                   {
354                         timeLED[in] = allData[in+12+1][iHit] ;
355                         break;
356                   }
357               }
358             for (Int_t iHit=0; iHit<5; iHit++) 
359               {
360                 if(allData[in+1+56][iHit] > low[in+1+56] && 
361                    allData[in+1+56][iHit] < high[in+1+56])
362                   {
363                     timeCFD[in+12] = allData[in+56+1][iHit] ;
364                     break;
365                   }
366               }
367             
368             for (Int_t iHit=0; iHit<5; iHit++) 
369               {
370                 if(allData[in+1+68][iHit] > low[in+1+68] && 
371                    allData[in+1+68][iHit] < high[in+1+68])
372                   {
373
374                     timeLED[in+12] = allData[in+68+1][iHit] ;
375                     break;
376                   }
377               }
378             AliDebug(5, Form(" readed i %i cfdC %i cfdA %i ledC %i ledA%i ",
379                               in, timeCFD[in],timeCFD[in+12],timeLED[in], 
380                               timeLED[in+12]));   
381             
382           }
383         
384         
385         for (Int_t in=0; in<12;  in++)
386           {
387             chargeQT0[in]=allData[2*in+25][0];
388             chargeQT1[in]=allData[2*in+26][0];
389             AliDebug(10, Form(" readed Raw %i %i %i",
390                               in, chargeQT0[in],chargeQT1[in]));
391           }     
392         for (Int_t in=12; in<24;  in++)
393           {
394             chargeQT0[in]=allData[2*in+57][0];
395             chargeQT1[in]=allData[2*in+58][0];
396             AliDebug(10, Form(" readed Raw %i %i %i",
397                               in, chargeQT0[in],chargeQT1[in]));
398             
399           }
400         
401         for (Int_t iHit=0; iHit<5; iHit++) 
402           {
403               if(allData[49][iHit] > low[49] && 
404                  allData[49][iHit] < high[49]){
405                 onlineMean = allData[49][iHit];
406                 break;
407               }
408           }
409         Double32_t time[24], adc[24],  noncalibtime[24];
410         for (Int_t ipmt=0; ipmt<24; ipmt++) {
411           if(timeCFD[ipmt]>0 && ipmt != badpmt ){
412            //for simulated data
413              //for physics  data
414            if(( chargeQT0[ipmt] - chargeQT1[ipmt])>0)  {
415              adc[ipmt] = chargeQT0[ipmt] - chargeQT1[ipmt];
416            }
417            else
418              adc[ipmt] = 0;
419            time[ipmt] = fCalib-> WalkCorrection(Int_t (refAmp), ipmt, Int_t(adc[ipmt]), timeCFD[ipmt] ) ;
420            
421            Double_t sl = timeLED[ipmt] - timeCFD[ipmt];
422            // time[ipmt] = fCalib-> WalkCorrection( refAmp,ipmt, Int_t(sl), timeCFD[ipmt] ) ;
423            AliDebug(5,Form(" ipmt %i QTC %i , time in chann %i (led-cfd) %i ",
424                             ipmt, Int_t(adc[ipmt]) ,Int_t(time[ipmt]),Int_t( sl)));
425            Double_t ampMip =( (TGraph*)fAmpLED.At(ipmt))->Eval(sl);
426            Double_t qtMip = ((TGraph*)fQTC.At(ipmt))->Eval(adc[ipmt]);
427            AliDebug(10,Form("  Amlitude in MIPS LED %f ; QTC %f;  in channels %f\n ",ampMip,qtMip, adc[ipmt]));
428            //bad peak removing
429              frecpoints->SetTime(ipmt, Float_t(time[ipmt]) );
430              // frecpoints->SetTime(ipmt,Double32_t(timeCFD[ipmt]));
431              frecpoints->SetAmp(ipmt, Double32_t( qtMip)); //for cosmic &pp beam 
432              frecpoints->SetAmpLED(ipmt, Double32_t(ampMip));        
433              noncalibtime[ipmt]= Double32_t (timeCFD[ipmt]);
434          }
435          else {
436            time[ipmt] = 0;
437            adc[ipmt] = 0;
438            noncalibtime[ipmt] = 0;
439          }
440        }
441        fESDTZEROfriend->SetT0timeCorr(noncalibtime) ;     
442        for (Int_t ipmt=0; ipmt<12; ipmt++){
443          if(time[ipmt] > 1 && ipmt != badpmt && adc[ipmt]>0.1 )
444            {
445              if(time[ipmt]<besttimeC){
446                   besttimeC=time[ipmt]; //timeC
447                   pmtBestC=ipmt;
448              }
449            }
450        }
451        for ( Int_t ipmt=12; ipmt<24; ipmt++)
452          {
453            if(time[ipmt] > 1 && ipmt != badpmt && adc[ipmt]>0.1)
454              {
455                if(time[ipmt]<besttimeA) {
456                  besttimeA=time[ipmt]; //timeA
457                  pmtBestA=ipmt;
458                }
459              }
460          }
461        if(besttimeA < 999999) 
462          frecpoints->SetTimeBestA(besttimeA * channelWidth - 1000.*fLatencyHPTDC + 1000.*fLatencyL1A - 1000.*fGRPdelays - fTimeMeanShift[1]);
463
464        if( besttimeC < 999999 ) 
465          frecpoints->SetTimeBestC(besttimeC * channelWidth - 1000.*fLatencyHPTDC +1000.*fLatencyL1C - 1000.*fGRPdelays - fTimeMeanShift[2]);
466        AliDebug(5,Form(" pmtA %i besttimeA %f shift A %f ps, pmtC %i besttimeC %f shiftC %f ps",
467                        pmtBestA,besttimeA, fTimeMeanShift[1],
468
469                        pmtBestC,  besttimeC,fTimeMeanShift[2]));
470         if(besttimeA <999999 && besttimeC < 999999 ){
471          timeDiff = ( besttimeA - besttimeC)* 0.001* channelWidth + fLatencyL1A - fLatencyL1C;
472          timeclock = channelWidth * Float_t( besttimeA+besttimeC)/2. - 1000.*fLatencyHPTDC + 1000.*fLatencyL1 - 1000.*fGRPdelays - fTimeMeanShift[0];  
473          meanTime = (besttimeA+besttimeC-2.*Float_t(ref))/2.;
474          vertex =  meanVertex - c*(timeDiff)/2. ; //+ (fdZonA - fdZonC)/2; 
475         }
476       }  //if phys event       
477       AliDebug(10,Form("  timeDiff %f #channel,  meanTime %f #channel, TOFmean%f  vertex %f cm meanVertex %f online mean %i \n",timeDiff, meanTime,timeclock, vertex,meanVertex, onlineMean));
478       frecpoints->SetT0clock(timeclock);
479       frecpoints->SetVertex(vertex);
480       frecpoints->SetMeanTime(meanTime);
481       frecpoints->SetOnlineMean(Int_t(onlineMean));
482         // Set triggers
483       
484       Bool_t tr[5];
485       Int_t trchan[5]= {50,51,52,55,56};
486       for (Int_t i=0; i<5; i++) tr[i]=false; 
487       for (Int_t itr=0; itr<5; itr++) {
488         for (Int_t iHit=0; iHit<5; iHit++) 
489           {
490             Int_t trr=trchan[itr];
491             if( allData[trr][iHit] > 0) tr[itr]=true;
492           }
493       }
494       frecpoints->SetT0Trig(tr);
495
496       //Set MPD
497       if(allData[53][0]>0 && allData[54][0]) 
498         frecpoints->SetMultA(allData[53][0]-allData[54][0]);
499         if(allData[105][0]>0 && allData[106][0]) 
500           frecpoints->SetMultC(allData[105][0]-allData[106][0]);
501         
502         
503     } // if (else )raw data
504   recTree->Fill();
505   if(frecpoints) delete frecpoints;
506 }
507   
508   
509   //____________________________________________________________
510   
511   void AliT0Reconstructor::FillESD(TTree */*digitsTree*/, TTree *clustersTree, AliESDEvent *pESD) const
512   {
513
514   /***************************************************
515   Resonstruct digits to vertex position
516   ****************************************************/
517   
518   AliDebug(1,Form("Start FillESD T0"));
519   pESD ->SetT0spread(fTimeSigmaShift);
520  
521
522   Float_t channelWidth = fParam->GetChannelWidth() ;  
523   Float_t c = 0.0299792458; // cm/ps
524   Float_t currentVertex=0, shift=0;
525   Int_t ncont=-1;
526   const AliESDVertex* vertex = pESD->GetPrimaryVertex();
527   if (!vertex)        vertex = pESD->GetPrimaryVertexSPD();
528   if (!vertex)        vertex = pESD->GetPrimaryVertexTPC();
529   if (!vertex)        vertex = pESD->GetVertex();
530
531   if (vertex) {
532     AliDebug(2, Form("Got %s (%s) from ESD: %f", 
533                     vertex->GetName(), vertex->GetTitle(), vertex->GetZ()));
534     currentVertex = vertex->GetZ();
535     
536     ncont = vertex->GetNContributors();
537     //  cout<<"@@ spdver "<<spdver<<" ncont "<<ncont<<endl;
538     if(ncont>0 ) {
539       shift = currentVertex/c;
540     }
541   }
542   TTree *treeR = clustersTree;
543   
544   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
545   if (!frecpoints) {
546     AliError("Reconstruct Fill ESD >> no recpoints found");
547     return;
548   }
549   
550   AliDebug(1,Form("Start FillESD T0"));
551   TBranch *brRec = treeR->GetBranch("T0");
552   if (brRec) {
553     brRec->SetAddress(&frecpoints);
554   }else{
555     AliError(Form("EXEC Branch T0 rec not found"));
556     return;
557   } 
558   
559   brRec->GetEntry(0);
560   Double32_t amp[24], time[24], ampQTC[24], timecorr[24];  
561   Double32_t* tcorr;
562   for(Int_t i=0; i<24; i++) 
563     amp[i]=time[i]=ampQTC[i]=timecorr[i]=0;
564
565
566   Double32_t timeClock[3];
567   Double32_t zPosition = frecpoints -> GetVertex();
568   Double32_t timeStart = frecpoints -> GetMeanTime();
569   timeClock[0] = frecpoints -> GetT0clock() ;
570   timeClock[1] = frecpoints -> GetBestTimeA() + shift;
571   timeClock[2] = frecpoints -> GetBestTimeC() - shift;
572   for ( Int_t i=0; i<24; i++) {
573     time[i] =  frecpoints -> GetTime(i); // ps to ns
574     if ( time[i] >1) {
575       ampQTC[i] = frecpoints -> GetAmp(i);
576       amp[i] = frecpoints -> AmpLED(i);
577       AliDebug(1,Form("T0: time %f  ampQTC %f ampLED %f \n", time[i], ampQTC[i], amp[i]));
578    }
579   }
580   Int_t trig= frecpoints ->GetT0Trig();
581   pESD->SetT0Trig(trig);
582   
583   pESD->SetT0zVertex(zPosition); //vertex Z position 
584
585   Double32_t multA=frecpoints ->GetMultA();
586   Double32_t multC=frecpoints ->GetMultC();
587   //  pESD->SetT0MultC(multC);        // multiplicity Cside 
588   //  pESD->SetT0MultA(multA);        // multiplicity Aside 
589   pESD->SetT0(multA); // for backward compatubility
590   pESD->SetT0clock(multC); // for backward compatubility
591
592   for(Int_t i=0; i<3; i++) 
593     pESD->SetT0TOF(i,timeClock[i]);   // interaction time (ns) 
594   pESD->SetT0time(time);         // best TOF on each PMT 
595   pESD->SetT0amplitude(ampQTC);     // number of particles(MIPs) on each PMT
596   
597   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));
598   
599   if (pESD) {
600     
601     AliESDfriend *fr = (AliESDfriend*)pESD->FindListObject("AliESDfriend");
602     if (fr) {
603       AliDebug(1, Form("Writing TZERO friend data to ESD tree"));
604
605       //     if (ncont>2) {
606         tcorr = fESDTZEROfriend->GetT0timeCorr();
607         for ( Int_t i=0; i<24; i++) {
608           if(i<12 && time[i]>1) timecorr[i] = tcorr[i] -  shift/channelWidth;
609           if(i>11 && time[i]>1) timecorr[i] = tcorr[i] +  shift/channelWidth;
610           if(time[i]>1)   AliDebug(10,Form("T0 friend : time %f  ampQTC %f ampLED %f \n", timecorr[i], ampQTC[i], amp[i]));
611         }
612         fESDTZEROfriend->SetT0timeCorr( timecorr) ;
613         fESDTZEROfriend->SetT0ampLEDminCFD(amp);
614         fESDTZEROfriend->SetT0ampQTC(ampQTC);
615         fr->SetTZEROfriend(fESDTZEROfriend);
616         //      }//
617     }
618   }
619      
620
621
622 } // vertex in 3 sigma
623
624
625
626
627
628