]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0Reconstructor.cxx
This commit was generated by cvs2svn to compensate for changes in r15986,
[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 "AliRun.h"
24 #include <AliESD.h>
25 #include "AliLog.h"
26 #include <TClonesArray.h>
27 #include "AliT0RecPoint.h"
28 #include "AliRawReader.h"
29 #include "AliT0RawReader.h"
30 #include "AliT0Loader.h"
31 #include "AliT0digit.h"
32 #include "AliT0Reconstructor.h"
33 #include "AliT0Parameters.h"
34 #include "AliCDBLocal.h"
35 #include "AliCDBStorage.h"
36 #include "AliCDBManager.h"
37 #include "AliCDBEntry.h"
38
39 #include <TArrayI.h>
40 #include <TGraph.h>
41
42 ClassImp(AliT0Reconstructor)
43
44   void  AliT0Reconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
45 {
46   //T0 raw data-> digits conversion
47  // reconstruct time information from raw data
48   AliT0RawReader myrawreader(rawReader,digitsTree);
49   if (!myrawreader.Next())
50     AliDebug(1,Form(" no raw data found!! %i", myrawreader.Next()));
51 }
52 void AliT0Reconstructor::Reconstruct(TTree*digitsTree, TTree*clustersTree) const
53 {
54 // T0 digits reconstruction
55 // T0RecPoint writing 
56
57     //Q->T-> coefficients !!!! should be asked!!!
58   //  Float_t ph2MIP=500;
59   Float_t gain[24], timeDelayCFD[24], timeDelayLED[24];
60   Float_t zdetA,zdetC;
61   TObjArray slewingLED;
62     
63   TArrayI * fADC = new TArrayI(24); 
64   TArrayI * fTimeCFD = new TArrayI(24); 
65   TArrayI * fADCLED = new TArrayI(24); 
66   TArrayI * fTimeLED = new TArrayI(24); 
67
68   AliT0Parameters* param = AliT0Parameters::Instance();
69   param->Init();
70
71   Int_t mV2Mip = param->GetmV2Mip();     
72   //mV2Mip = param->GetmV2Mip();     
73   Int_t channelWidth = param->GetChannelWidth() ;  
74   
75   for (Int_t i=0; i<24; i++){
76     timeDelayCFD[i] = param->GetTimeDelayCFD(i);
77     timeDelayLED[i] = param->GetTimeDelayLED(i);
78     gain[i] = param->GetGain(i);
79     //gain[i] = 1;
80     slewingLED.AddAtAndExpand(param->GetSlew(i),i);
81   }
82   zdetC = param->GetZposition(0);
83   zdetA  = param->GetZposition(1);
84     
85   AliDebug(1,Form("Start DIGITS reconstruction "));
86  
87   TBranch *brDigits=digitsTree->GetBranch("T0");
88   AliT0digit *fDigits = new AliT0digit();
89   if (brDigits) {
90     brDigits->SetAddress(&fDigits);
91   }else{
92     cerr<<"EXEC Branch T0 digits not found"<<endl;
93     return;
94   }
95   brDigits->GetEntry(0);
96   fDigits->GetTime(*fTimeCFD);
97   fDigits->GetADC(*fADC);
98   fDigits->GetTimeAmp(*fTimeLED);
99   fDigits->GetADCAmp(*fADCLED);
100
101   Float_t besttimeright=999999;
102   Float_t besttimeleft=999999;
103   Int_t pmtBestRight=99999;
104   Int_t pmtBestLeft=99999;
105   Float_t timeDiff=999999, meanTime=0;
106   
107   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
108   clustersTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
109
110   Float_t time[24], adc[24];
111   for (Int_t ipmt=0; ipmt<24; ipmt++) {
112       
113     if(fTimeCFD->At(ipmt)>0 ){
114       time[ipmt] = channelWidth *( fTimeCFD->At(ipmt)) - 1000*timeDelayCFD[ipmt];
115       Float_t adc_digPs = channelWidth * Float_t (fADCLED->At(ipmt)) ;
116       adc[ipmt] = TMath::Exp(adc_digPs/1000) /gain[ipmt];
117       AliDebug(1,Form(" time %f ps,  adc %f mv in MIP %i\n ",
118                       time[ipmt], adc[ipmt], Int_t (adc[ipmt]/mV2Mip +0.5)));
119       frecpoints->SetTime(ipmt,time[ipmt]);
120       frecpoints->SetAmp(ipmt,adc[ipmt]);
121     }
122     else {
123       time[ipmt] = 0;
124       adc[ipmt] = 0;
125     }
126   }
127   for (Int_t ipmt=0; ipmt<12; ipmt++){
128     if(time[ipmt] > 1 ) {
129       if(time[ipmt]<besttimeleft){
130         besttimeleft=time[ipmt]; //timeleft
131         pmtBestLeft=ipmt;
132       }
133     }
134   }
135   for ( Int_t ipmt=12; ipmt<24; ipmt++){
136     if(time[ipmt] > 1) {
137       if(time[ipmt]<besttimeright) {
138         besttimeright=time[ipmt]; //timeright
139         pmtBestRight=ipmt;}
140     }
141   }
142   if(besttimeright !=999999)  frecpoints->SetTimeBestRight(Int_t(besttimeright));
143   if( besttimeleft != 999999 ) frecpoints->SetTimeBestLeft(Int_t(besttimeleft));
144   AliDebug(1,Form(" besttimeA %f ps,  besttimeC %f ps",besttimeright, besttimeleft));
145   Float_t c = 0.0299792; // cm/ps
146   Float_t vertex = 0;
147   if(besttimeright !=999999 && besttimeleft != 999999 ){
148     timeDiff = besttimeleft - besttimeright;
149     meanTime = (besttimeright + besttimeleft)/2.;
150     vertex = c*(timeDiff)/2.; //-(lenr-lenl))/2;
151     AliDebug(1,Form("  timeDiff %f ps,  meanTime %f ps, vertex %f cm",timeDiff, meanTime,vertex ));
152     frecpoints->SetVertex(vertex);
153     frecpoints->SetMeanTime(Int_t(meanTime));
154     
155   }
156   clustersTree->Fill();
157 }
158
159
160 void AliT0Reconstructor::FillESD(AliRunLoader* runLoader, AliESD *pESD) const
161 {
162
163   /***************************************************
164   Resonstruct digits to vertex position
165   ****************************************************/
166   
167   
168   if (!runLoader) {
169     AliError("Reconstruct >> No run loader");
170     return;
171   }
172   
173   AliDebug(1,Form("Start FillESD T0"));
174
175   AliT0Loader* pStartLoader = (AliT0Loader*) runLoader->GetLoader("T0Loader");
176  
177   pStartLoader->LoadRecPoints("READ");
178
179   TTree *treeR = pStartLoader->TreeR();
180   
181    AliT0RecPoint* frecpoints= new AliT0RecPoint ();
182     if (!frecpoints) {
183     AliError("Reconstruct Fill ESD >> no recpoints found");
184     return;
185   }
186   
187   AliDebug(1,Form("Start FillESD T0"));
188   TBranch *brRec = treeR->GetBranch("T0");
189   if (brRec) {
190     brRec->SetAddress(&frecpoints);
191   }else{
192     cerr<<"EXEC Branch T0 rec not found"<<endl;
193     // exit(111);
194     return;
195   } 
196     
197     brRec->GetEntry(0);
198     Float_t timeStart, Zposition, amp[24], time[24];
199     Int_t i;
200     Zposition = frecpoints -> GetVertex();
201     timeStart = frecpoints -> GetMeanTime();
202     for ( i=0; i<24; i++) {
203        time[i] = Float_t (frecpoints -> GetTime(i)) / 1000.; // ps to ns
204       amp[i] = frecpoints -> GetAmp(i);
205     }
206     pESD->SetT0zVertex(Zposition); //vertex Z position 
207     pESD->SetT0(timeStart);        // interaction time 
208     pESD->SetT0time(time);         // best TOF on each PMT 
209     pESD->SetT0amplitude(amp);     // number of particles(MIPs) on each PMT
210     cout<<" ESD >> "<<Zposition<<" "<<timeStart<<endl;
211
212     pStartLoader->UnloadRecPoints();
213    
214 } // vertex in 3 sigma
215
216
217
218
219
220