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