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