]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0Reconstructor.cxx
fix some coding violations and error calculation
[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   Float_t meanVertex = fParam->GetMeanVertex();
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       Float_t ampLED;    //LED
173       Double_t sl = (timeLED->At(ipmt) - time[ipmt])*channelWidth;
174       if( adc[ipmt]/Float_t(mv2MIP) < 20) 
175         {
176           Double_t qt=((TGraph*)fAmpLEDrec.At(ipmt))->Eval(sl/1000.);
177           ampLED = qt/Float_t(mv2MIP);
178           AliDebug(1,Form(" ipmt %i QTC %f ch QTC in MIP %f, time in chann %f (led-cfd) %f in MIPs %f",
179                           ipmt, adc[ipmt], adc[ipmt]/Float_t(mv2MIP), time[ipmt],sl,ampLED ));
180         }
181       else
182         { AliDebug(1,Form(" amplitude %f MIPs on PMT %i, write QTC in both ", (adc[ipmt]/Float_t(mv2MIP)), ipmt));
183           ampLED=adc[ipmt]/Float_t(mv2MIP);}
184       frecpoints->SetTime(ipmt,time[ipmt]);
185       frecpoints->SetAmp(ipmt,adc[ipmt]/Float_t(mv2MIP));
186       frecpoints->SetAmpLED(ipmt,ampLED);
187     }
188     else {
189       time[ipmt] = 0;
190       adc[ipmt] = 0;
191     }
192   }
193   
194   for (Int_t ipmt=0; ipmt<12; ipmt++){
195     if(time[ipmt] > 1 ) {
196       if(time[ipmt]<besttimeC){
197         besttimeC=time[ipmt]; //timeC
198         pmtBestC=ipmt;
199       }
200     }
201   }
202   for ( Int_t ipmt=12; ipmt<24; ipmt++){
203     if(time[ipmt] > 1) {
204       if(time[ipmt]<besttimeA) {
205         besttimeA=time[ipmt]; //timeA
206         pmtBestA=ipmt;}
207     }
208   }
209   if(besttimeA !=999999)  frecpoints->SetTimeBestA(Int_t(besttimeA));
210   if( besttimeC != 999999 ) frecpoints->SetTimeBestC(Int_t(besttimeC));
211   AliDebug(1,Form(" besttimeA %f ch,  besttimeC %f ch",besttimeA, besttimeC));
212   Float_t c = 0.0299792; // cm/ps
213   Float_t vertex = 0;
214   if(besttimeA !=999999 && besttimeC != 999999 ){
215     timeDiff = (besttimeC - besttimeA)*channelWidth;
216     meanTime = (Float_t((besttimeA + besttimeC -4000)/2) * channelWidth); 
217     //    meanTime = (meanT0 - (besttimeA + besttimeC)/2) * channelWidth;
218     vertex = meanVertex - c*(timeDiff)/2. + (fdZonA - fdZonC)/2; 
219     frecpoints->SetVertex(vertex);
220     frecpoints->SetMeanTime(Int_t(meanTime));
221     //online mean
222     frecpoints->SetOnlineMean(Int_t(onlineMean * channelWidth));
223     AliDebug(1,Form("  timeDiff %f ps,  meanTime %f ps, vertex %f cm online mean %i ps",timeDiff, meanTime,vertex, Int_t(onlineMean * channelWidth )));
224     
225   }
226   
227   clustersTree->Fill();
228
229   delete timeCFD;
230   delete timeLED;
231   delete chargeQT0; 
232   delete chargeQT1; 
233 }
234
235
236 //_______________________________________________________________________
237
238 void AliT0Reconstructor::Reconstruct(AliRawReader* rawReader, TTree*recTree) const
239 {
240   // T0 raw ->
241   // T0RecPoint writing 
242   
243   //Q->T-> coefficients !!!! should be measured!!!
244   Int_t allData[110][5];
245   
246   Int_t timeCFD[24], timeLED[24], chargeQT0[24], chargeQT1[24];
247   TString option = GetOption(); 
248   AliDebug(10,Form("Option: %s\n", option.Data()));
249    
250    
251   for (Int_t i0=0; i0<105; i0++)
252     {
253       for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;   
254     }
255    
256   Float_t besttimeA=9999999;
257   Float_t besttimeC=9999999;
258   Int_t pmtBestA=99999;
259   Int_t pmtBestC=99999;
260   Float_t timeDiff=9999999, meanTime=0;
261   Double_t qt=0;
262   Int_t mv2MIP = fParam-> GetmV2Mip();     
263   Float_t meanVertex = fParam->GetMeanVertex();
264
265   //  UInt_t type =rawReader->GetType();         
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+81][0] ;
296           }
297        }
298       
299       if(option == "cosmic"  ) //&& type == 7 )
300         {
301           
302           for (Int_t in=0; in<12; in++)  
303             {
304               timeCFD[in] = allData[in+1][0] ;
305               timeCFD[in+12] = allData[in+56+1][0] ;
306               timeLED[in] = allData[in+12+1][0] ;
307               timeLED[in+12] = allData[in+68+1][0] ;
308             }
309           
310           for (Int_t in=0; in<12;  in++)
311             {
312               chargeQT1[in]=allData[2*in+25][0];
313               chargeQT0[in]=allData[2*in+26][0];
314             }
315           
316            for (Int_t in=12; in<24;  in++)
317              {
318                chargeQT1[in]=allData[2*in+57][0];
319                chargeQT0[in]=allData[2*in+58][0];
320              }
321            
322          } //cosmic with physics event
323        for (Int_t in=0; in<24; in++)  
324          AliDebug(10, Form(" readed Raw %i %i %i %i %i",
325                            in, timeLED[in],timeCFD[in],chargeQT0[in],chargeQT1[in]));
326        Int_t onlineMean = allData[49][0];       
327        
328        Float_t time[24], adc[24];
329        for (Int_t ipmt=0; ipmt<24; ipmt++) {
330          if(timeCFD[ipmt]>0 && timeLED[ipmt]>0){
331            //for simulated data
332            if(option == "pdc"){
333              Double_t qt0 = Double_t(chargeQT0[ipmt]);
334              Double_t qt1 = Double_t(chargeQT1[ipmt]);
335              if((qt1-qt0)>0)  adc[ipmt] = Int_t(TMath::Exp( Double_t (channelWidth*(qt1-qt0)/1000.)));
336              time[ipmt] = fCalib-> WalkCorrection( ipmt,Int_t(qt1) , timeCFD[ipmt], "pdc" ) ;
337              Double_t sl = (timeLED[ipmt] - time[ipmt])*channelWidth;
338              //pp collison
339              if( adc[ipmt]/Float_t(mv2MIP) < 20) 
340                {
341                  if(fAmpLEDrec.At(ipmt)) 
342                    qt=((TGraph*)fAmpLEDrec.At(ipmt))->Eval(sl/1000.);
343                }
344              else //PbPb
345                { AliDebug(1,Form(" amplitude %f MIPs on PMT %i, write QTC in both ", (adc[ipmt]/Float_t(mv2MIP)), ipmt));
346                  qt=adc[ipmt]/Float_t(mv2MIP);}
347              
348              frecpoints->SetTime(ipmt,time[ipmt]);
349              frecpoints->SetAmp(ipmt,adc[ipmt]/Float_t(mv2MIP));
350              frecpoints->SetAmpLED(ipmt,qt/Float_t(mv2MIP));
351              AliDebug(10,Form(" QTC %f mv,  time in chann %f ampLED %f",adc[ipmt] ,time[ipmt], qt));
352              AliDebug(10,Form("  Amlitude in MIPS LED %f ,  QTC %f \n ", adc[ipmt]/Float_t(mv2MIP),qt/Float_t(mv2MIP)));
353            }
354            
355            
356              //for physics  data
357            if(option == "cosmic") // && type == 7)
358              {
359              if(( chargeQT1[ipmt] - chargeQT0[ipmt])>0)  
360                adc[ipmt] = chargeQT1[ipmt] - chargeQT0[ipmt];
361              else
362                adc[ipmt] = 0;
363              //      time[ipmt] = fCalib-> WalkCorrection( ipmt, adc[ipmt], timeCFD[ipmt],"cosmic" ) ;
364              
365              Double_t sl = timeLED[ipmt] - timeCFD[ipmt];
366              time[ipmt] = fCalib-> WalkCorrection( ipmt, Int_t(sl), timeCFD[ipmt],"cosmic" ) ;
367              AliDebug(10,Form(" ipmt %i QTC %i , time in chann %i (led-cfd) %i ",
368                               ipmt, Int_t(adc[ipmt]) ,Int_t(time[ipmt]),Int_t( sl)));
369              Double_t ampMip =( (TGraph*)fAmpLED.At(ipmt))->Eval(sl);
370              Double_t qtMip = ((TGraph*)fQTC.At(ipmt))->Eval(adc[ipmt]);
371              AliDebug(10,Form("  Amlitude in MIPS LED %f ,  QTC %f \n ",ampMip,qtMip));
372              
373              frecpoints->SetTime(ipmt, Float_t(time[ipmt]) );
374              frecpoints->SetAmp(ipmt, Float_t( ampMip)); //for cosmic &pp beam 
375              frecpoints->SetAmpLED(ipmt, Float_t(qtMip));
376            } //if physic data end
377            
378          }
379          else {
380            time[ipmt] = 0;
381            adc[ipmt] = 0;
382          }
383        }
384        
385        for (Int_t ipmt=0; ipmt<12; ipmt++){
386          if(time[ipmt] > 1 ) {
387            if(time[ipmt]<besttimeC){
388              besttimeC=time[ipmt]; //timeC
389              pmtBestC=ipmt;
390            }
391          }
392        }
393        for ( Int_t ipmt=12; ipmt<24; ipmt++){
394          if(time[ipmt] > 1) {
395            if(time[ipmt]<besttimeA) {
396              besttimeA=time[ipmt]; //timeA
397              pmtBestA=ipmt;}
398          }
399        }
400        if(besttimeA !=9999999)  frecpoints->SetTimeBestA(Int_t(besttimeA));
401        if( besttimeC != 9999999 ) frecpoints->SetTimeBestC(Int_t(besttimeC));
402        AliDebug(1,Form(" besttimeA %f ps,  besttimeC %f ps",besttimeA, besttimeC));
403        Float_t c = 0.0299792; // cm/ps
404        Float_t vertex = 99999;
405        if(besttimeA <9999999 && besttimeC < 9999999 ){
406          timeDiff = ( besttimeC - besttimeA) *channelWidth;
407          if(option == "pdc"){ 
408            //      meanTime = (besttimeA + besttimeC)/2 * channelWidth; 
409            meanTime = (besttimeA + besttimeC-4000.)/2 * channelWidth; 
410            onlineMean = Int_t (onlineMean * channelWidth);
411          }
412          if(option == "cosmic") {
413            meanTime =  Float_t((besttimeA + besttimeC)/2);  
414            onlineMean = onlineMean ;
415          }
416          vertex =  meanVertex -c*(timeDiff)/2.; //+ (fdZonA - fdZonC)/2; 
417          frecpoints->SetVertex(vertex);
418          frecpoints->SetMeanTime(Int_t(meanTime));
419          frecpoints->SetOnlineMean(Int_t(onlineMean));
420          AliDebug(1,Form("  timeDiff %f ps,  meanTime %f ps, vertex %f cm online mean %i ",timeDiff, meanTime,vertex, onlineMean));
421          
422        }
423     } // if (else )raw data
424          recTree->Fill();
425          if(frecpoints) delete frecpoints;
426 }
427        
428        
429 //____________________________________________________________
430
431 void AliT0Reconstructor::FillESD(TTree */*digitsTree*/, TTree *clustersTree, AliESDEvent *pESD) const
432 {
433
434   /***************************************************
435   Resonstruct digits to vertex position
436   ****************************************************/
437   
438   AliDebug(1,Form("Start FillESD T0"));
439
440   TTree *treeR = clustersTree;
441   
442    AliT0RecPoint* frecpoints= new AliT0RecPoint ();
443     if (!frecpoints) {
444     AliError("Reconstruct Fill ESD >> no recpoints found");
445     return;
446   }
447   
448   AliDebug(1,Form("Start FillESD T0"));
449   TBranch *brRec = treeR->GetBranch("T0");
450   if (brRec) {
451     brRec->SetAddress(&frecpoints);
452   }else{
453     AliError(Form("EXEC Branch T0 rec not found"));
454     return;
455   } 
456     
457     brRec->GetEntry(0);
458     Float_t amp[24], time[24];
459     Float_t  zPosition = frecpoints -> GetVertex();
460     Float_t timeStart = frecpoints -> GetMeanTime() ;
461     for ( Int_t i=0; i<24; i++) {
462       time[i] = Float_t (frecpoints -> GetTime(i)); // ps to ns
463       amp[i] = frecpoints -> GetAmp(i);
464     }
465     pESD->SetT0zVertex(zPosition); //vertex Z position 
466     pESD->SetT0(timeStart);        // interaction time 
467     pESD->SetT0time(time);         // best TOF on each PMT 
468     pESD->SetT0amplitude(amp);     // number of particles(MIPs) on each PMT
469  
470     AliDebug(1,Form(" Z position %f cm,  T0  %f ps",zPosition , timeStart));
471
472 } // vertex in 3 sigma
473
474
475
476
477
478