]> git.uio.no Git - u/mrichter/AliRoot.git/blame - T0/AliT0Reconstructor.cxx
the same calibration for cosmic and pdc data
[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$ */
72e48d95 17/*********************************************************************
18 * T0 reconstruction and filling ESD
19 * - reconstruct mean time (interation time)
20 * - vertex position
21 * - multiplicity
22 ********************************************************************/
dc7ca31d 23
af885e0f 24#include <AliESDEvent.h>
dc7ca31d 25#include "AliLog.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
34#include <TArrayI.h>
35#include <TGraph.h>
aad72f45 36#include <TMath.h>
b09247a2 37#include <Riostream.h>
dc7ca31d 38
39ClassImp(AliT0Reconstructor)
40
c41ceaac 41 AliT0Reconstructor:: AliT0Reconstructor(): AliReconstructor(),
f16935f7 42 fdZonA(0),
43 fdZonC(0),
44 fZposition(0),
45 fParam(NULL),
2e6a5ee0 46 fAmpLEDrec(),
c883fdf2 47 fQTC(0),
48 fAmpLED(0),
2e6a5ee0 49 fCalib()
e0bba6cc 50{
72e48d95 51 //constructor
52
c41ceaac 53 AliDebug(1,"Start reconstructor ");
74adb36a 54
55 fParam = AliT0Parameters::Instance();
56 fParam->Init();
c883fdf2 57 TString option = GetOption();
58
74adb36a 59 for (Int_t i=0; i<24; i++){
2e6a5ee0 60 TGraph* gr = fParam ->GetAmpLEDRec(i);
29ed1d0e 61 if (gr) fAmpLEDrec.AddAtAndExpand(gr,i) ;
c883fdf2 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 }
29ed1d0e 68
f16935f7 69 fdZonC = TMath::Abs(fParam->GetZPositionShift("T0/C/PMT1"));
70 fdZonA = TMath::Abs(fParam->GetZPositionShift("T0/A/PMT15"));
12e9daf9 71 fCalib = new AliT0Calibrator();
72
c41ceaac 73}
74//____________________________________________________________________
e0bba6cc 75
c41ceaac 76AliT0Reconstructor::AliT0Reconstructor(const AliT0Reconstructor &r):
d76c31f4 77 AliReconstructor(r),
f16935f7 78 fdZonA(0),
79 fdZonC(0),
80 fZposition(0),
81 fParam(NULL),
2e6a5ee0 82 fAmpLEDrec(),
c883fdf2 83 fQTC(0),
84 fAmpLED(0),
2e6a5ee0 85 fCalib()
86
f16935f7 87 {
c41ceaac 88 //
89 // AliT0Reconstructor copy constructor
90 //
e0bba6cc 91
c41ceaac 92 ((AliT0Reconstructor &) r).Copy(*this);
e0bba6cc 93
94}
95
c41ceaac 96//_____________________________________________________________________________
97AliT0Reconstructor &AliT0Reconstructor::operator=(const AliT0Reconstructor &r)
dc7ca31d 98{
c41ceaac 99 //
100 // Assignment operator
101 //
102
103 if (this != &r) ((AliT0Reconstructor &) r).Copy(*this);
104 return *this;
e0bba6cc 105
dc7ca31d 106}
c41ceaac 107
108//_____________________________________________________________________________
109
dc7ca31d 110void AliT0Reconstructor::Reconstruct(TTree*digitsTree, TTree*clustersTree) const
94c27e4f 111
dc7ca31d 112{
94c27e4f 113 // T0 digits reconstruction
114 // T0RecPoint writing
c41ceaac 115
d0bcd1fb 116
c41ceaac 117 TArrayI * timeCFD = new TArrayI(24);
118 TArrayI * timeLED = new TArrayI(24);
119 TArrayI * chargeQT0 = new TArrayI(24);
120 TArrayI * chargeQT1 = new TArrayI(24);
74adb36a 121
d0bcd1fb 122
8955c6b4 123 Float_t channelWidth = fParam->GetChannelWidth() ;
b95e8d87 124 Float_t meanVertex = fParam->GetMeanVertex();
94c27e4f 125
dc7ca31d 126 AliDebug(1,Form("Start DIGITS reconstruction "));
94c27e4f 127
2e6a5ee0 128
d0bcd1fb 129 TBranch *brDigits=digitsTree->GetBranch("T0");
e0bba6cc 130 AliT0digit *fDigits = new AliT0digit() ;
dc7ca31d 131 if (brDigits) {
132 brDigits->SetAddress(&fDigits);
133 }else{
f16935f7 134 AliError(Form("EXEC Branch T0 digits not found"));
135 return;
dc7ca31d 136 }
e0bba6cc 137
c41ceaac 138 digitsTree->GetEvent(0);
139 digitsTree->GetEntry(0);
140 brDigits->GetEntry(0);
141 fDigits->GetTimeCFD(*timeCFD);
142 fDigits->GetTimeLED(*timeLED);
143 fDigits->GetQT0(*chargeQT0);
144 fDigits->GetQT1(*chargeQT1);
446d6ec4 145 Int_t onlineMean = fDigits->MeanTime();
c883fdf2 146
c41ceaac 147
148 Float_t besttimeA=999999;
149 Float_t besttimeC=999999;
150 Int_t pmtBestA=99999;
151 Int_t pmtBestC=99999;
dc7ca31d 152 Float_t timeDiff=999999, meanTime=0;
153
d0bcd1fb 154 // Int_t mv2MIP = fParam-> GetmV2Mip();
c883fdf2 155
e0bba6cc 156
94c27e4f 157 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
158 clustersTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
159
b95e8d87 160 Float_t time[24], adc[24];
dc7ca31d 161 for (Int_t ipmt=0; ipmt<24; ipmt++) {
c41ceaac 162 if(timeCFD->At(ipmt)>0 ){
d0bcd1fb 163 if(( chargeQT1->At(ipmt) - chargeQT0->At(ipmt))>0)
164 adc[ipmt] = chargeQT1->At(ipmt) - chargeQT0->At(ipmt);
a7027400 165 else
d0bcd1fb 166 adc[ipmt] = 0;
167
168 time[ipmt] = fCalib-> WalkCorrection( ipmt, adc[ipmt], Float_t( timeCFD->At(ipmt))) ;
169
170 Double_t sl = Double_t(timeLED->At(ipmt) - timeCFD->At(ipmt));
171 // time[ipmt] = fCalib-> WalkCorrection( ipmt, Int_t(sl), timeCFD[ipmt],"cosmic" ) ;
172 AliDebug(10,Form(" ipmt %i QTC %i , time in chann %i (led-cfd) %i ",
173 ipmt, Int_t(adc[ipmt]) ,Int_t(time[ipmt]),Int_t( sl)));
174 Double_t ampMip =((TGraph*)fAmpLED.At(ipmt))->Eval(sl);
175 Double_t qtMip = ((TGraph*)fQTC.At(ipmt))->Eval(adc[ipmt]);
176 AliDebug(10,Form(" Amlitude in MIPS LED %f , QTC %f in channels %i\n ",ampMip,qtMip, adc[ipmt]));
177
178 frecpoints->SetTime(ipmt, Float_t(time[ipmt]) );
179 frecpoints->SetAmp(ipmt, Float_t( ampMip)); //for cosmic &pp beam
180 frecpoints->SetAmpLED(ipmt, Float_t(qtMip));
181
dc7ca31d 182 }
183 else {
184 time[ipmt] = 0;
185 adc[ipmt] = 0;
186 }
187 }
94c27e4f 188
dc7ca31d 189 for (Int_t ipmt=0; ipmt<12; ipmt++){
190 if(time[ipmt] > 1 ) {
c41ceaac 191 if(time[ipmt]<besttimeC){
192 besttimeC=time[ipmt]; //timeC
193 pmtBestC=ipmt;
dc7ca31d 194 }
195 }
196 }
197 for ( Int_t ipmt=12; ipmt<24; ipmt++){
198 if(time[ipmt] > 1) {
c41ceaac 199 if(time[ipmt]<besttimeA) {
200 besttimeA=time[ipmt]; //timeA
201 pmtBestA=ipmt;}
dc7ca31d 202 }
203 }
c41ceaac 204 if(besttimeA !=999999) frecpoints->SetTimeBestA(Int_t(besttimeA));
205 if( besttimeC != 999999 ) frecpoints->SetTimeBestC(Int_t(besttimeC));
d0bcd1fb 206 AliDebug(10,Form(" besttimeA %f ch, besttimeC %f ch",besttimeA, besttimeC));
dc7ca31d 207 Float_t c = 0.0299792; // cm/ps
208 Float_t vertex = 0;
c41ceaac 209 if(besttimeA !=999999 && besttimeC != 999999 ){
12e9daf9 210 timeDiff = (besttimeC - besttimeA)*channelWidth;
d0bcd1fb 211 meanTime = Float_t((besttimeA + besttimeC)/2);// * channelWidth);
12e9daf9 212 // meanTime = (meanT0 - (besttimeA + besttimeC)/2) * channelWidth;
d0bcd1fb 213 vertex = meanVertex - c*(timeDiff)/2.;// + (fdZonA - fdZonC)/2;
dc7ca31d 214 frecpoints->SetVertex(vertex);
215 frecpoints->SetMeanTime(Int_t(meanTime));
446d6ec4 216 //online mean
d0bcd1fb 217 frecpoints->SetOnlineMean(Int_t(onlineMean));
218 AliDebug(10,Form(" timeDiff %f ps, meanTime %f ps, vertex %f cm online mean %i ps",timeDiff, meanTime,vertex, Int_t(onlineMean)));
dc7ca31d 219
220 }
b95e8d87 221
dc7ca31d 222 clustersTree->Fill();
bd375212 223
224 delete timeCFD;
225 delete timeLED;
226 delete chargeQT0;
227 delete chargeQT1;
dc7ca31d 228}
229
230
c41ceaac 231//_______________________________________________________________________
232
233void AliT0Reconstructor::Reconstruct(AliRawReader* rawReader, TTree*recTree) const
234{
94c27e4f 235 // T0 raw ->
236 // T0RecPoint writing
237
238 //Q->T-> coefficients !!!! should be measured!!!
e8ed1cd0 239 Int_t allData[110][5];
2e6a5ee0 240
e8ed1cd0 241 Int_t timeCFD[24], timeLED[24], chargeQT0[24], chargeQT1[24];
aaa0a98f 242 TString option = GetOption();
d0bcd1fb 243 AliDebug(1,Form("Option: %s\n", option.Data()));
2e6a5ee0 244
2e6a5ee0 245
bce12dc5 246 for (Int_t i0=0; i0<105; i0++)
247 {
248 for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;
249 }
2e6a5ee0 250
bce12dc5 251 Float_t besttimeA=9999999;
252 Float_t besttimeC=9999999;
253 Int_t pmtBestA=99999;
254 Int_t pmtBestC=99999;
255 Float_t timeDiff=9999999, meanTime=0;
256 Double_t qt=0;
d0bcd1fb 257 // Int_t mv2MIP = fParam-> GetmV2Mip();
b95e8d87 258 Float_t meanVertex = fParam->GetMeanVertex();
d0bcd1fb 259 printf("meanVertex %f \n",meanVertex);
b95e8d87 260 // UInt_t type =rawReader->GetType();
bce12dc5 261
262 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
263
264 recTree->Branch( "T0", "AliT0RecPoint" ,&frecpoints, 405,1);
2e6a5ee0 265
266
bce12dc5 267 AliDebug(10," before read data ");
268 AliT0RawReader myrawreader(rawReader);
269 if (!myrawreader.Next())
270 AliDebug(1,Form(" no raw data found!!"));
271 else
272 {
273 for (Int_t i=0; i<105; i++) {
274 for (Int_t iHit=0; iHit<5; iHit++)
275 {
276 allData[i][iHit] = myrawreader.GetData(i,iHit);
277 }
278 }
279
280 Float_t channelWidth = fParam->GetChannelWidth() ;
281
282 // Int_t meanT0 = fParam->GetMeanT0();
d0bcd1fb 283
bce12dc5 284
285 for (Int_t in=0; in<12; in++)
286 {
287 timeCFD[in] = allData[in+1][0] ;
288 timeCFD[in+12] = allData[in+56+1][0] ;
289 timeLED[in] = allData[in+12+1][0] ;
290 timeLED[in+12] = allData[in+68+1][0] ;
d0bcd1fb 291 AliDebug(10, Form(" readed i %i cfdC %i cfdA %i ledC %i ledA%i ",
292 in, timeCFD[in],timeCFD[in+12],timeLED[in],
293 timeLED[in+12]));
bce12dc5 294 }
295
296 for (Int_t in=0; in<12; in++)
297 {
d0bcd1fb 298 chargeQT0[in]=allData[2*in+25][0];
299 chargeQT1[in]=allData[2*in+26][0];
bce12dc5 300 }
301
302 for (Int_t in=12; in<24; in++)
303 {
d0bcd1fb 304 chargeQT0[in]=allData[2*in+57][0];
305 chargeQT1[in]=allData[2*in+58][0];
bce12dc5 306 }
307
d0bcd1fb 308 // } //cosmic with physics event
2e6a5ee0 309 for (Int_t in=0; in<24; in++)
c883fdf2 310 AliDebug(10, Form(" readed Raw %i %i %i %i %i",
311 in, timeLED[in],timeCFD[in],chargeQT0[in],chargeQT1[in]));
446d6ec4 312 Int_t onlineMean = allData[49][0];
2e6a5ee0 313
12e9daf9 314 Float_t time[24], adc[24];
2e6a5ee0 315 for (Int_t ipmt=0; ipmt<24; ipmt++) {
316 if(timeCFD[ipmt]>0 && timeLED[ipmt]>0){
bce12dc5 317 //for simulated data
d0bcd1fb 318 //for physics data
319 if(( chargeQT1[ipmt] - chargeQT0[ipmt])>0)
320 adc[ipmt] = chargeQT1[ipmt] - chargeQT0[ipmt];
321 else
322 adc[ipmt] = 0;
bce12dc5 323
d0bcd1fb 324 time[ipmt] = fCalib-> WalkCorrection( ipmt, adc[ipmt], timeCFD[ipmt],"cosmic" ) ;
bce12dc5 325
d0bcd1fb 326 Double_t sl = timeLED[ipmt] - timeCFD[ipmt];
327 // time[ipmt] = fCalib-> WalkCorrection( ipmt, Int_t(sl), timeCFD[ipmt],"cosmic" ) ;
328 AliDebug(10,Form(" ipmt %i QTC %i , time in chann %i (led-cfd) %i ",
329 ipmt, Int_t(adc[ipmt]) ,Int_t(time[ipmt]),Int_t( sl)));
330 Double_t ampMip =( (TGraph*)fAmpLED.At(ipmt))->Eval(sl);
331 Double_t qtMip = ((TGraph*)fQTC.At(ipmt))->Eval(adc[ipmt]);
332 AliDebug(10,Form(" Amlitude in MIPS LED %f ; QTC %f; in channels %i\n ",ampMip,qtMip, adc[ipmt]));
bce12dc5 333
d0bcd1fb 334 frecpoints->SetTime(ipmt, Float_t(time[ipmt]) );
335 frecpoints->SetAmp(ipmt, Float_t( ampMip)); //for cosmic &pp beam
336 frecpoints->SetAmpLED(ipmt, Float_t(qtMip));
bce12dc5 337
2e6a5ee0 338 }
339 else {
340 time[ipmt] = 0;
341 adc[ipmt] = 0;
342 }
343 }
344
345 for (Int_t ipmt=0; ipmt<12; ipmt++){
346 if(time[ipmt] > 1 ) {
347 if(time[ipmt]<besttimeC){
348 besttimeC=time[ipmt]; //timeC
349 pmtBestC=ipmt;
350 }
351 }
352 }
353 for ( Int_t ipmt=12; ipmt<24; ipmt++){
354 if(time[ipmt] > 1) {
355 if(time[ipmt]<besttimeA) {
356 besttimeA=time[ipmt]; //timeA
357 pmtBestA=ipmt;}
358 }
359 }
360 if(besttimeA !=9999999) frecpoints->SetTimeBestA(Int_t(besttimeA));
361 if( besttimeC != 9999999 ) frecpoints->SetTimeBestC(Int_t(besttimeC));
d0bcd1fb 362 AliDebug(5,Form(" pmtA %i besttimeA %f ps, pmtC %i besttimeC %f ps",
363 pmtBestA,besttimeA, pmtBestC, besttimeC));
364 Float_t c = 0.0299792458; // cm/ps
2e6a5ee0 365 Float_t vertex = 99999;
366 if(besttimeA <9999999 && besttimeC < 9999999 ){
12e9daf9 367 timeDiff = ( besttimeC - besttimeA) *channelWidth;
d0bcd1fb 368 meanTime = Float_t((besttimeA + besttimeC)/2.);
369 onlineMean = onlineMean ;
370 vertex = meanVertex - c*(timeDiff)/2.; //+ (fdZonA - fdZonC)/2;
2e6a5ee0 371 frecpoints->SetVertex(vertex);
372 frecpoints->SetMeanTime(Int_t(meanTime));
1d7681da 373 frecpoints->SetOnlineMean(Int_t(onlineMean));
d0bcd1fb 374 AliDebug(5,Form(" timeDiff %f ps, meanTime %f ps, vertex %f cm meanVertex %f online mean %i ",timeDiff, meanTime,vertex,meanVertex, onlineMean));
2e6a5ee0 375
376 }
bce12dc5 377 } // if (else )raw data
378 recTree->Fill();
379 if(frecpoints) delete frecpoints;
c41ceaac 380}
bce12dc5 381
382
c41ceaac 383//____________________________________________________________
384
d76c31f4 385void AliT0Reconstructor::FillESD(TTree */*digitsTree*/, TTree *clustersTree, AliESDEvent *pESD) const
dc7ca31d 386{
387
388 /***************************************************
389 Resonstruct digits to vertex position
390 ****************************************************/
391
dc7ca31d 392 AliDebug(1,Form("Start FillESD T0"));
393
d76c31f4 394 TTree *treeR = clustersTree;
dc7ca31d 395
396 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
397 if (!frecpoints) {
398 AliError("Reconstruct Fill ESD >> no recpoints found");
399 return;
400 }
401
402 AliDebug(1,Form("Start FillESD T0"));
403 TBranch *brRec = treeR->GetBranch("T0");
404 if (brRec) {
405 brRec->SetAddress(&frecpoints);
406 }else{
f16935f7 407 AliError(Form("EXEC Branch T0 rec not found"));
dc7ca31d 408 return;
409 }
410
411 brRec->GetEntry(0);
f16935f7 412 Float_t amp[24], time[24];
413 Float_t zPosition = frecpoints -> GetVertex();
414 Float_t timeStart = frecpoints -> GetMeanTime() ;
415 for ( Int_t i=0; i<24; i++) {
94c27e4f 416 time[i] = Float_t (frecpoints -> GetTime(i)); // ps to ns
dc7ca31d 417 amp[i] = frecpoints -> GetAmp(i);
418 }
f16935f7 419 pESD->SetT0zVertex(zPosition); //vertex Z position
dc7ca31d 420 pESD->SetT0(timeStart); // interaction time
421 pESD->SetT0time(time); // best TOF on each PMT
422 pESD->SetT0amplitude(amp); // number of particles(MIPs) on each PMT
5325480c 423
f16935f7 424 AliDebug(1,Form(" Z position %f cm, T0 %f ps",zPosition , timeStart));
dc7ca31d 425
dc7ca31d 426} // vertex in 3 sigma
427
428
429
430
431
432