]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0Reconstructor.cxx
macro to convert Raw to Digits
[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 //____________________________________________________________________
56
57 AliT0Reconstructor::AliT0Reconstructor(const AliT0Reconstructor &r):
58   fZposition(0)
59   
60 {
61   //
62   // AliT0Reconstructor copy constructor
63   //
64
65   ((AliT0Reconstructor &) r).Copy(*this);
66
67 }
68
69 //_____________________________________________________________________________
70 AliT0Reconstructor &AliT0Reconstructor::operator=(const AliT0Reconstructor &r)
71 {
72   //
73   // Assignment operator
74   //
75
76   if (this != &r) ((AliT0Reconstructor &) r).Copy(*this);
77   return *this;
78
79 }
80
81 //_____________________________________________________________________________
82
83 void AliT0Reconstructor::Reconstruct(TTree*digitsTree, TTree*clustersTree) const
84
85 {
86 // T0 digits reconstruction
87 // T0RecPoint writing 
88
89   
90   Float_t  timeDelayLED[24];
91   Float_t zdetA,zdetC;
92   TObjArray slewingLEDrec;
93   TObjArray walk;
94     
95   TArrayI * timeCFD = new TArrayI(24); 
96   TArrayI * timeLED = new TArrayI(24); 
97   TArrayI * chargeQT0 = new TArrayI(24); 
98   TArrayI * chargeQT1 = new TArrayI(24); 
99   
100   AliT0Parameters* param = AliT0Parameters::Instance();
101   param->Init();
102   AliT0Calibrator *calib=new AliT0Calibrator(); 
103
104   Int_t mV2Mip = param->GetmV2Mip();     
105   //mV2Mip = param->GetmV2Mip();     
106   Int_t channelWidth = param->GetChannelWidth() ;  
107   
108   for (Int_t i=0; i<24; i++){
109     TGraph* gr = param ->GetSlewRec(i);
110     slewingLEDrec.AddAtAndExpand(gr,i) ;  
111   }
112   zdetC = param->GetZposition(0);
113   zdetA  = param->GetZposition(1);
114     
115   AliDebug(1,Form("Start DIGITS reconstruction "));
116
117   TBranch *brDigits=digitsTree->GetBranch("T0");
118   AliT0digit *fDigits = new AliT0digit() ;
119   if (brDigits) {
120     brDigits->SetAddress(&fDigits);
121   }else{
122     cerr<<"EXEC Branch T0 digits not found"<<endl;
123     return;
124   }
125   
126   digitsTree->GetEvent(0);
127   digitsTree->GetEntry(0);
128   brDigits->GetEntry(0);
129   fDigits->GetTimeCFD(*timeCFD);
130   fDigits->GetTimeLED(*timeLED);
131   fDigits->GetQT0(*chargeQT0);
132   fDigits->GetQT1(*chargeQT1);
133   
134   Float_t besttimeA=999999;
135   Float_t besttimeC=999999;
136   Int_t pmtBestA=99999;
137   Int_t pmtBestC=99999;
138   Float_t timeDiff=999999, meanTime=0;
139   
140
141
142    AliT0RecPoint* frecpoints= new AliT0RecPoint ();
143    clustersTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
144
145   Float_t time[24], adc[24];
146   for (Int_t ipmt=0; ipmt<24; ipmt++) {
147     if(timeCFD->At(ipmt)>0 ){
148       Int_t qt0= chargeQT0->At(ipmt);
149       Int_t qt1= chargeQT1->At(ipmt);
150       if((qt1-qt0)>0)  adc[ipmt] = TMath::Exp( Double_t (channelWidth*(qt1-qt0)/1000));
151       time[ipmt] = channelWidth * (calib-> WalkCorrection( ipmt,qt1 , timeCFD->At(ipmt) ) ) ;
152       
153       //LED
154       Double_t sl = (timeLED->At(ipmt) - timeCFD->At(ipmt)- (1000.*timeDelayLED[ipmt]/channelWidth))*channelWidth;
155       Double_t qt=((TGraph*)slewingLEDrec.At(ipmt))->Eval(sl/1000.);
156       frecpoints->SetTime(ipmt,time[ipmt]);
157       frecpoints->SetAmp(ipmt,adc[ipmt]);
158       frecpoints->SetAmpLED(ipmt,qt);
159     }
160     else {
161       time[ipmt] = 0;
162       adc[ipmt] = 0;
163     }
164   }
165
166   for (Int_t ipmt=0; ipmt<12; ipmt++){
167     if(time[ipmt] > 1 ) {
168       if(time[ipmt]<besttimeC){
169         besttimeC=time[ipmt]; //timeC
170         pmtBestC=ipmt;
171       }
172     }
173   }
174   for ( Int_t ipmt=12; ipmt<24; ipmt++){
175     if(time[ipmt] > 1) {
176       if(time[ipmt]<besttimeA) {
177         besttimeA=time[ipmt]; //timeA
178         pmtBestA=ipmt;}
179     }
180   }
181   if(besttimeA !=999999)  frecpoints->SetTimeBestA(Int_t(besttimeA));
182   if( besttimeC != 999999 ) frecpoints->SetTimeBestC(Int_t(besttimeC));
183   AliDebug(1,Form(" besttimeA %f ps,  besttimeC %f ps",besttimeA, besttimeC));
184   Float_t c = 0.0299792; // cm/ps
185   Float_t vertex = 0;
186   if(besttimeA !=999999 && besttimeC != 999999 ){
187     timeDiff = besttimeC - besttimeA;
188     meanTime = (besttimeA + besttimeC)/2.;
189     vertex = c*(timeDiff)/2.; //-(lenr-lenl))/2;
190     AliDebug(1,Form("  timeDiff %f ps,  meanTime %f ps, vertex %f cm",timeDiff, meanTime,vertex ));
191     frecpoints->SetVertex(vertex);
192     frecpoints->SetMeanTime(Int_t(meanTime));
193     
194   }
195   clustersTree->Fill();
196
197   delete timeCFD;
198   delete timeLED;
199   delete chargeQT0; 
200   delete chargeQT1; 
201 }
202
203
204 //_______________________________________________________________________
205
206 void AliT0Reconstructor::Reconstruct(AliRawReader* rawReader, TTree*recTree) const
207 {
208 // T0 raw ->
209 // T0RecPoint writing 
210
211     //Q->T-> coefficients !!!! should be asked!!!
212   Float_t  timeDelayLED[24];
213   Float_t zdetA,zdetC;
214   TObjArray slewingLEDrec;
215   TObjArray walk;
216     
217   TArrayI * timeCFD = new TArrayI(24); 
218   TArrayI * timeLED = new TArrayI(24); 
219   TArrayI * chargeQT0 = new TArrayI(24); 
220   TArrayI * chargeQT1 = new TArrayI(24); 
221
222    AliT0RawReader myrawreader(rawReader);
223    if (!myrawreader.Next())
224      AliDebug(1,Form(" no raw data found!! %i", myrawreader.Next()));
225    Int_t allData[110][5];
226    for (Int_t i=0; i<110; i++) {
227      allData[i][0]=myrawreader.GetData(i,0);
228    }
229
230   AliT0Parameters* param = AliT0Parameters::Instance();
231   param->Init();
232   AliT0Calibrator *calib=new AliT0Calibrator(); 
233
234   Int_t mV2Mip = param->GetmV2Mip();     
235   //mV2Mip = param->GetmV2Mip();     
236   Int_t channelWidth = param->GetChannelWidth() ;  
237     
238   for (Int_t i=0; i<24; i++){
239     TGraph* gr = param ->GetSlewRec(i);
240     slewingLEDrec.AddAtAndExpand(gr,i) ;  
241   }
242   
243     zdetC = param->GetZposition(0);
244     zdetA  = param->GetZposition(1);
245     //  zdetC=GetZdet("C");
246     // zdetA=GetZdet("A");
247   cout<<" !!!!!  zdetC "<<zdetC<<"  zdetA "<< zdetA<<endl;    
248    for (Int_t in=0; in<24; in++)
249      {
250        timeLED->AddAt(allData[in+1][0],in);
251        timeCFD->AddAt(allData[in+25][0],in);
252        chargeQT1->AddAt(allData[in+55][0],in);
253        chargeQT0->AddAt(allData[in+79][0],in);
254        AliDebug(10, Form(" readed Raw %i %i %i %i %i", in, timeLED->At(in),timeCFD->At(in),chargeQT0->At(in),chargeQT1->At(in)));
255      }
256
257   Float_t besttimeA=999999;
258   Float_t besttimeC=999999;
259   Int_t pmtBestA=99999;
260   Int_t pmtBestC=99999;
261   Float_t timeDiff=999999, meanTime=0;
262
263   
264   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
265  
266   recTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
267   
268
269   Float_t time[24], adc[24];
270   for (Int_t ipmt=0; ipmt<24; ipmt++) {
271     if(timeCFD->At(ipmt)>0 ){
272       Int_t qt0= chargeQT0->At(ipmt);
273       Int_t qt1= chargeQT1->At(ipmt);
274       if((qt1-qt0)>0)  adc[ipmt] = TMath::Exp( Double_t (channelWidth*(qt1-qt0)/1000));
275       time[ipmt] = channelWidth * (calib-> WalkCorrection( ipmt,qt1 , timeCFD->At(ipmt) ) ) ;
276       Double_t sl = (timeLED->At(ipmt) - timeCFD->At(ipmt)- (1000.*timeDelayLED[ipmt]/channelWidth))*channelWidth;
277       Double_t qt=((TGraph*)slewingLEDrec.At(ipmt))->Eval(sl/1000.);
278       frecpoints->SetTime(ipmt,time[ipmt]);
279       frecpoints->SetAmp(ipmt,adc[ipmt]);
280       frecpoints->SetAmpLED(ipmt,qt);
281     }
282     else {
283       time[ipmt] = 0;
284       adc[ipmt] = 0;
285     }
286   }
287
288   for (Int_t ipmt=0; ipmt<12; ipmt++){
289     if(time[ipmt] > 1 ) {
290       if(time[ipmt]<besttimeC){
291         besttimeC=time[ipmt]; //timeC
292         pmtBestC=ipmt;
293       }
294     }
295   }
296   for ( Int_t ipmt=12; ipmt<24; ipmt++){
297     if(time[ipmt] > 1) {
298       if(time[ipmt]<besttimeA) {
299         besttimeA=time[ipmt]; //timeA
300         pmtBestA=ipmt;}
301     }
302   }
303   if(besttimeA !=999999)  frecpoints->SetTimeBestA(Int_t(besttimeA));
304   if( besttimeC != 999999 ) frecpoints->SetTimeBestC(Int_t(besttimeC));
305   AliDebug(1,Form(" besttimeA %f ps,  besttimeC %f ps",besttimeA, besttimeC));
306   Float_t c = 0.0299792; // cm/ps
307   Float_t vertex = 0;
308   if(besttimeA !=999999 && besttimeC != 999999 ){
309     timeDiff = besttimeC - besttimeA;
310     meanTime = (besttimeA + besttimeC)/2.;
311     vertex = c*(timeDiff)/2.; //-(lenr-lenl))/2;
312     AliDebug(1,Form("  timeDiff %f ps,  meanTime %f ps, vertex %f cm",timeDiff, meanTime,vertex ));
313     frecpoints->SetVertex(vertex);
314     frecpoints->SetMeanTime(Int_t(meanTime));
315     
316   }
317   recTree->Fill();
318
319
320   delete timeCFD;
321   delete timeLED;
322   delete chargeQT0; 
323   delete chargeQT1; 
324   
325 }
326 //____________________________________________________________
327
328 void AliT0Reconstructor::FillESD(AliRunLoader* runLoader, AliESD *pESD) const
329 {
330
331   /***************************************************
332   Resonstruct digits to vertex position
333   ****************************************************/
334   
335   
336   if (!runLoader) {
337     AliError("Reconstruct >> No run loader");
338     return;
339   }
340   
341   AliDebug(1,Form("Start FillESD T0"));
342
343   AliT0Loader* pStartLoader = (AliT0Loader*) runLoader->GetLoader("T0Loader");
344  
345   pStartLoader->LoadRecPoints("READ");
346
347   TTree *treeR = pStartLoader->TreeR();
348   
349    AliT0RecPoint* frecpoints= new AliT0RecPoint ();
350     if (!frecpoints) {
351     AliError("Reconstruct Fill ESD >> no recpoints found");
352     return;
353   }
354   
355   AliDebug(1,Form("Start FillESD T0"));
356   TBranch *brRec = treeR->GetBranch("T0");
357   if (brRec) {
358     brRec->SetAddress(&frecpoints);
359   }else{
360     cerr<<"EXEC Branch T0 rec not found"<<endl;
361     // exit(111);
362     return;
363   } 
364     
365     brRec->GetEntry(0);
366     Float_t timeStart, Zposition, amp[24], time[24];
367     Int_t mean0 = 12450;
368     Int_t i;
369     Zposition = frecpoints -> GetVertex();
370     timeStart = frecpoints -> GetMeanTime() - mean0;
371     for ( i=0; i<24; i++) {
372        time[i] = Float_t (frecpoints -> GetTime(i)) / 1000.; // ps to ns
373       amp[i] = frecpoints -> GetAmp(i);
374     }
375     pESD->SetT0zVertex(Zposition); //vertex Z position 
376     pESD->SetT0(timeStart);        // interaction time 
377     pESD->SetT0time(time);         // best TOF on each PMT 
378     pESD->SetT0amplitude(amp);     // number of particles(MIPs) on each PMT
379     cout<<" ESD >> "<<Zposition<<" "<<timeStart<<endl;
380
381     pStartLoader->UnloadRecPoints();
382    
383 } // vertex in 3 sigma
384
385
386
387
388
389