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