]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCcalibTime.cxx
bug fix - intorduced with coverity changes
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibTime.cxx
CommitLineData
c74c5f6c 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/*
74235403 17Comments to be written here:
c74c5f6c 18
74235403 191. What do we calibrate.
20
21 Time dependence of gain and drift velocity in order to account for changes in: temperature, pressure, gas composition.
c74c5f6c 22
23 AliTPCcalibTime *calibTime = new AliTPCcalibTime("cosmicTime","cosmicTime",0, 1213.9e+06, 1213.96e+06, 0.04e+04, 0.04e+04);
24
74235403 252. How to interpret results
26
273. Simple example
c74c5f6c 28
74235403 29 a) determine the required time range:
c74c5f6c 30
74235403 31 AliXRDPROOFtoolkit tool;
32 TChain * chain = tool.MakeChain("pass2.txt","esdTree",0,6000);
33 chain->Draw("GetTimeStamp()")
c74c5f6c 34
74235403 35 b) analyse calibration object on Proof in calibration train
c74c5f6c 36
74235403 37 AliTPCcalibTime *calibTime = new AliTPCcalibTime("cosmicTime","cosmicTime", StartTimeStamp, EndTimeStamp, IntegrationTimeVdrift);
c74c5f6c 38
74235403 39 c) plot results
40 .x ~/NimStyle.C
2be25cc9 41 gSystem->Load("libANALYSIS");
42 gSystem->Load("libTPCcalib");
c74c5f6c 43
dde68d36 44 TFile f("CalibObjectsTrain1.root");
45 AliTPCcalibTime *calib = (AliTPCcalibTime *)f->Get("calibTime");
46 calib->GetHistoDrift("all")->Projection(2,0)->Draw()
47 calib->GetFitDrift("all")->Draw("lp")
da6c0bc9 48
74235403 494. Analysis using debug streamers.
da6c0bc9 50
74235403 51 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
52 gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
53 AliXRDPROOFtoolkit tool;
3e55050f 54 TChain * chainTime = tool.MakeChainRandom("time.txt","trackInfo",0,10000);
817766d5 55
3e55050f 56 AliXRDPROOFtoolkit::FilterList("timetpctpc.txt","* tpctpc",1)
57 AliXRDPROOFtoolkit::FilterList("timetoftpc.txt","* toftpc",1)
58 AliXRDPROOFtoolkit::FilterList("timeitstpc.txt","* itstpc",1)
59 AliXRDPROOFtoolkit::FilterList("timelaser.txt","* laserInfo",1)
60 TChain * chainTPCTPC = tool.MakeChainRandom("timetpctpc.txt.Good","tpctpc",0,10000);
817766d5 61 TChain * chainTPCITS = tool.MakeChainRandom("timeitstpc.txt.Good","itstpc",0,10000);
3e55050f 62 TChain * chainTPCTOF = tool.MakeChainRandom("timetoftpc.txt.Good","toftpc",0,10000);
817766d5 63 TChain * chainLaser = tool.MakeChainRandom("timelaser.txt.Good","laserInfo",0,10000);
74235403 64 chainTime->Lookup();
dde68d36 65 chainLaser->Lookup();
c74c5f6c 66*/
67
c74c5f6c 68#include "Riostream.h"
2811495d 69#include "TDatabasePDG.h"
70#include "TGraphErrors.h"
c74c5f6c 71#include "TH1F.h"
c74c5f6c 72#include "THnSparse.h"
73#include "TList.h"
74#include "TMath.h"
2811495d 75#include "TTimeStamp.h"
76#include "TTree.h"
c74c5f6c 77#include "TVectorD.h"
2811495d 78//#include "TChain.h"
79//#include "TFile.h"
80
81#include "AliDCSSensor.h"
82#include "AliDCSSensorArray.h"
c74c5f6c 83#include "AliESDEvent.h"
c74c5f6c 84#include "AliESDInputHandler.h"
2811495d 85#include "AliESDVertex.h"
86#include "AliESDfriend.h"
c74c5f6c 87#include "AliLog.h"
817766d5 88#include "AliRelAlignerKalman.h"
2811495d 89#include "AliTPCCalROC.h"
90#include "AliTPCParam.h"
c74c5f6c 91#include "AliTPCTracklet.h"
da6c0bc9 92#include "AliTPCcalibDB.h"
93#include "AliTPCcalibLaser.h"
2811495d 94#include "AliTPCcalibTime.h"
95#include "AliTPCclusterMI.h"
96#include "AliTPCseed.h"
817766d5 97#include "AliTrackPointArray.h"
2811495d 98#include "AliTracker.h"
70036ef3 99#include "AliKFVertex.h"
660f0543 100#include <AliLog.h>
101
c74c5f6c 102ClassImp(AliTPCcalibTime)
103
104
105AliTPCcalibTime::AliTPCcalibTime()
660f0543 106 :AliTPCcalibBase(),
107 fMemoryMode(1), // 0 -do not fill THnSparse with residuals 1- fill only important QA THn 2 - Fill all THnsparse for calibration
da6c0bc9 108 fLaser(0), // pointer to laser calibration
109 fDz(0), // current delta z
d3ce44cb 110 fCutMaxD(3), // maximal distance in rfi ditection
111 fCutMaxDz(25), // maximal distance in rfi ditection
c74c5f6c 112 fCutTheta(0.03), // maximal distan theta
2be25cc9 113 fCutMinDir(-0.99), // direction vector products
b9ab8e40 114 fCutTracks(2500),
97d3fea0 115 fArrayLaserA(0), //laser fit parameters C
116 fArrayLaserC(0), //laser fit parameters A
dde68d36 117 fArrayDz(0), //NEW! Tmap of V drifts for different triggers
817766d5 118 fAlignITSTPC(0), //alignemnt array ITS TPC match
119 fAlignTRDTPC(0), //alignemnt array TRD TPC match
120 fAlignTOFTPC(0), //alignemnt array TOF TPC match
97d3fea0 121 fTimeKalmanBin(60*15), //time bin width for kalman - 15 minutes default
2be25cc9 122 fTimeBins(0),
123 fTimeStart(0),
124 fTimeEnd(0),
125 fPtBins(0),
126 fPtStart(0),
127 fPtEnd(0),
128 fVdriftBins(0),
129 fVdriftStart(0),
130 fVdriftEnd(0),
131 fRunBins(0),
132 fRunStart(0),
133 fRunEnd(0)
c74c5f6c 134{
6ba68367 135 //
136 // default constructor
137 //
c74c5f6c 138 AliInfo("Default Constructor");
2be25cc9 139 for (Int_t i=0;i<3;i++) {
140 fHistVdriftLaserA[i]=0;
141 fHistVdriftLaserC[i]=0;
142 }
d3ce44cb 143 for (Int_t i=0;i<10;i++) {
74235403 144 fCosmiMatchingHisto[i]=0;
145 }
391ffdb2 146 //
147 for (Int_t i=0;i<5;i++) {
b9908d0b 148 fResHistoTPCCE[i]=0;
391ffdb2 149 fResHistoTPCITS[i]=0;
150 fResHistoTPCTRD[i]=0;
0b736a46 151 fResHistoTPCTOF[i]=0;
391ffdb2 152 fResHistoTPCvertex[i]=0;
70036ef3 153 fTPCVertex[i]=0;
154 }
155 for (Int_t i=0;i<12;i++) {
156 fTPCVertex[i]=0;
391ffdb2 157 }
b1fea14f 158 for (Int_t i=0;i<5;i++) {
159 fTPCVertexCorrelation[i]=0;
160 }
660f0543 161 static Int_t counter=0;
162 if (1) {
163 TTimeStamp s;
164 Int_t time=s;
165 AliInfo(Form("Counter Constructor\t%d\t%d",counter,time));
166 counter++;
167 }
168
c74c5f6c 169}
170
660f0543 171AliTPCcalibTime::AliTPCcalibTime(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeVdrift, Int_t memoryMode)
c74c5f6c 172 :AliTPCcalibBase(),
660f0543 173 fMemoryMode(memoryMode), // 0 -do not fill THnSparse with residuals 1- fill only important QA THn 2 - Fill all THnsparse for calibration
d3ce44cb 174 fLaser(0), // pointer to laser calibration
175 fDz(0), // current delta z
176 fCutMaxD(5*0.5356), // maximal distance in rfi ditection
dde68d36 177 fCutMaxDz(40), // maximal distance in rfi ditection
d3ce44cb 178 fCutTheta(5*0.004644),// maximal distan theta
179 fCutMinDir(-0.99), // direction vector products
b9ab8e40 180 fCutTracks(2500),
97d3fea0 181 fArrayLaserA(new TObjArray(1000)), //laser fit parameters C
182 fArrayLaserC(new TObjArray(1000)), //laser fit parameters A
dde68d36 183 fArrayDz(0), //Tmap of V drifts for different triggers
817766d5 184 fAlignITSTPC(0), //alignemnt array ITS TPC match
185 fAlignTRDTPC(0), //alignemnt array TRD TPC match
186 fAlignTOFTPC(0), //alignemnt array TOF TPC match
97d3fea0 187 fTimeKalmanBin(60*15), //time bin width for kalman - 15 minutes default
2be25cc9 188 fTimeBins(0),
189 fTimeStart(0),
190 fTimeEnd(0),
191 fPtBins(0),
192 fPtStart(0),
193 fPtEnd(0),
194 fVdriftBins(0),
195 fVdriftStart(0),
196 fVdriftEnd(0),
197 fRunBins(0),
198 fRunStart(0),
199 fRunEnd(0)
c74c5f6c 200{
6ba68367 201 //
202 // Non deafaul constructor - to be used in the Calibration setups
203 //
204
c74c5f6c 205 SetName(name);
206 SetTitle(title);
2be25cc9 207 for (Int_t i=0;i<3;i++) {
208 fHistVdriftLaserA[i]=0;
209 fHistVdriftLaserC[i]=0;
210 }
c74c5f6c 211
391ffdb2 212 for (Int_t i=0;i<5;i++) {
b9908d0b 213 fResHistoTPCCE[i]=0;
391ffdb2 214 fResHistoTPCITS[i]=0;
215 fResHistoTPCTRD[i]=0;
0b736a46 216 fResHistoTPCTOF[i]=0;
391ffdb2 217 fResHistoTPCvertex[i]=0;
218 }
219
220
c74c5f6c 221 AliInfo("Non Default Constructor");
74235403 222 fTimeBins =(EndTime-StartTime)/deltaIntegrationTimeVdrift;
2be25cc9 223 fTimeStart =StartTime; //(((TObjString*)(mapGRP->GetValue("fAliceStartTime")))->GetString()).Atoi();
224 fTimeEnd =EndTime; //(((TObjString*)(mapGRP->GetValue("fAliceStopTime")))->GetString()).Atoi();
dde68d36 225 fPtBins = 400;
2be25cc9 226 fPtStart = -0.04;
227 fPtEnd = 0.04;
dde68d36 228 fVdriftBins = 500;
229 fVdriftStart= -0.1;
230 fVdriftEnd = 0.1;
32100b1d 231 fRunBins = 1000001;
dde68d36 232 fRunStart = -1.5;
32100b1d 233 fRunEnd = 999999.5;
2be25cc9 234
235 Int_t binsVdriftLaser[4] = {fTimeBins , fPtBins , fVdriftBins*20, fRunBins };
236 Double_t xminVdriftLaser[4] = {fTimeStart, fPtStart, fVdriftStart , fRunStart};
237 Double_t xmaxVdriftLaser[4] = {fTimeEnd , fPtEnd , fVdriftEnd , fRunEnd };
dde68d36 238 TString axisTitle[4]={
239 "T",
240 "#delta_{P/T}",
241 "value",
242 "run"
243 };
244 TString histoName[3]={
245 "Loffset",
246 "Lcorr",
247 "Lgy"
248 };
2be25cc9 249
dde68d36 250
2be25cc9 251 for (Int_t i=0;i<3;i++) {
252 fHistVdriftLaserA[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
253 fHistVdriftLaserC[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
dde68d36 254 fHistVdriftLaserA[i]->SetName(histoName[i]);
255 fHistVdriftLaserC[i]->SetName(histoName[i]);
256 for (Int_t iaxis=0; iaxis<4;iaxis++){
257 fHistVdriftLaserA[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
258 fHistVdriftLaserC[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
259 }
2be25cc9 260 }
261 fBinsVdrift[0] = fTimeBins;
262 fBinsVdrift[1] = fPtBins;
263 fBinsVdrift[2] = fVdriftBins;
264 fBinsVdrift[3] = fRunBins;
265 fXminVdrift[0] = fTimeStart;
266 fXminVdrift[1] = fPtStart;
267 fXminVdrift[2] = fVdriftStart;
268 fXminVdrift[3] = fRunStart;
269 fXmaxVdrift[0] = fTimeEnd;
270 fXmaxVdrift[1] = fPtEnd;
271 fXmaxVdrift[2] = fVdriftEnd;
272 fXmaxVdrift[3] = fRunEnd;
273
dde68d36 274 fArrayDz=new TObjArray();
817766d5 275 fAlignITSTPC = new TObjArray; //alignemnt array ITS TPC match
276 fAlignTRDTPC = new TObjArray; //alignemnt array ITS TPC match
277 fAlignTOFTPC = new TObjArray; //alignemnt array ITS TPC match
3e55050f 278 fAlignITSTPC->SetOwner(kTRUE);
279 fAlignTRDTPC->SetOwner(kTRUE);
280 fAlignTOFTPC->SetOwner(kTRUE);
281
c74c5f6c 282
d3ce44cb 283 fCosmiMatchingHisto[0]=new TH1F("Cosmics matching","p0-all" ,100,-10*0.5356 ,10*0.5356 );
284 fCosmiMatchingHisto[1]=new TH1F("Cosmics matching","p1-all" ,100,-10*4.541 ,10*4.541 );
285 fCosmiMatchingHisto[2]=new TH1F("Cosmics matching","p2-all" ,100,-10*0.01134 ,10*0.01134 );
286 fCosmiMatchingHisto[3]=new TH1F("Cosmics matching","p3-all" ,100,-10*0.004644,10*0.004644);
287 fCosmiMatchingHisto[4]=new TH1F("Cosmics matching","p4-all" ,100,-10*0.03773 ,10*0.03773 );
288 fCosmiMatchingHisto[5]=new TH1F("Cosmics matching","p0-isPair",100,-10*0.5356 ,10*0.5356 );
289 fCosmiMatchingHisto[6]=new TH1F("Cosmics matching","p1-isPair",100,-10*4.541 ,10*4.541 );
290 fCosmiMatchingHisto[7]=new TH1F("Cosmics matching","p2-isPair",100,-10*0.01134 ,10*0.01134 );
291 fCosmiMatchingHisto[8]=new TH1F("Cosmics matching","p3-isPair",100,-10*0.004644,10*0.004644);
292 fCosmiMatchingHisto[9]=new TH1F("Cosmics matching","p4-isPair",100,-10*0.03773 ,10*0.03773 );
70036ef3 293 for (Int_t i=0;i<12;i++) {
294 fTPCVertex[i]=0;
295 }
b1fea14f 296 for (Int_t i=0;i<5;i++) {
297 fTPCVertexCorrelation[i]=0;
298 }
391ffdb2 299 BookDistortionMaps();
70036ef3 300
74235403 301}
c74c5f6c 302
303AliTPCcalibTime::~AliTPCcalibTime(){
304 //
6ba68367 305 // Virtual Destructor
c74c5f6c 306 //
660f0543 307 static Int_t counter=0;
308 if (1) {
309 TTimeStamp s;
310 Int_t time=s;
311 AliInfo(Form("Counter Destructor\t%s\t%d\t%d",GetName(),counter,time));
312 counter++;
313 }
d3ce44cb 314 for(Int_t i=0;i<3;i++){
315 if(fHistVdriftLaserA[i]){
316 delete fHistVdriftLaserA[i];
317 fHistVdriftLaserA[i]=NULL;
318 }
319 if(fHistVdriftLaserC[i]){
320 delete fHistVdriftLaserC[i];
321 fHistVdriftLaserC[i]=NULL;
322 }
323 }
dde68d36 324 if(fArrayDz){
325 fArrayDz->SetOwner();
326 fArrayDz->Delete();
327 delete fArrayDz;
328 fArrayDz=NULL;
2be25cc9 329 }
d3ce44cb 330 for(Int_t i=0;i<5;i++){
331 if(fCosmiMatchingHisto[i]){
332 delete fCosmiMatchingHisto[i];
333 fCosmiMatchingHisto[i]=NULL;
334 }
74235403 335 }
391ffdb2 336
337 for (Int_t i=0;i<5;i++) {
b9908d0b 338 delete fResHistoTPCCE[i];
391ffdb2 339 delete fResHistoTPCITS[i];
340 delete fResHistoTPCTRD[i];
0b736a46 341 delete fResHistoTPCTOF[i];
391ffdb2 342 delete fResHistoTPCvertex[i];
b9908d0b 343 fResHistoTPCCE[i]=0;
391ffdb2 344 fResHistoTPCITS[i]=0;
345 fResHistoTPCTRD[i]=0;
0b736a46 346 fResHistoTPCTOF[i]=0;
391ffdb2 347 fResHistoTPCvertex[i]=0;
348 }
349
472f0066 350 if (fTPCVertex[0]) {
70036ef3 351 for (Int_t i=0;i<12;i++) delete fTPCVertex[i];
352 }
472f0066 353 if (fTPCVertexCorrelation[0]) {
b1fea14f 354 for (Int_t i=0;i<5;i++) delete fTPCVertexCorrelation[i];
355 }
70036ef3 356
3e55050f 357 fAlignITSTPC->SetOwner(kTRUE);
358 fAlignTRDTPC->SetOwner(kTRUE);
359 fAlignTOFTPC->SetOwner(kTRUE);
360
817766d5 361 fAlignITSTPC->Delete();
362 fAlignTRDTPC->Delete();
363 fAlignTOFTPC->Delete();
3e55050f 364 delete fAlignITSTPC;
365 delete fAlignTRDTPC;
366 delete fAlignTOFTPC;
c74c5f6c 367}
2be25cc9 368
9170f0f0 369// Bool_t AliTPCcalibTime::IsLaser(const AliESDEvent *const /*event*/) const{
370// //
371// // Indicator is laser event not yet implemented - to be done using trigger info or event specie
372// //
373// return kTRUE; //More accurate creteria to be added
374// }
375// Bool_t AliTPCcalibTime::IsCosmics(const AliESDEvent *const /*event*/){
376// //
377// // Indicator is cosmic event not yet implemented - to be done using trigger info or event specie
378// //
379
380// return kTRUE; //More accurate creteria to be added
381// }
382// Bool_t AliTPCcalibTime::IsBeam(const AliESDEvent *const /*event*/) const{
383// //
384// // Indicator is physic event not yet implemented - to be done using trigger info or event specie
385// //
386
387// return kTRUE; //More accurate creteria to be added
388// }
389void AliTPCcalibTime::ResetCurrent(){
6ba68367 390 //
9170f0f0 391 //ResetCurrent
6ba68367 392 //
dde68d36 393 fDz=0; //Reset current dz
2be25cc9 394}
6ba68367 395
396
397
2be25cc9 398void AliTPCcalibTime::Process(AliESDEvent *event){
6ba68367 399 //
400 // main function to make calibration
401 //
2be25cc9 402 if(!event) return;
403 if (event->GetNumberOfTracks()<2) return;
404 ResetCurrent();
9170f0f0 405 //if(IsLaser (event))
406 ProcessLaser (event);
407 //if(IsCosmics(event))
408 ProcessCosmic(event);
409 //if(IsBeam (event))
410 ProcessBeam (event);
2be25cc9 411}
412
413void AliTPCcalibTime::ProcessLaser(AliESDEvent *event){
c74c5f6c 414 //
2be25cc9 415 // Fit drift velocity using laser
416 //
417 // 0. cuts
dde68d36 418 const Int_t kMinTracks = 40; // minimal number of laser tracks
419 const Int_t kMinTracksSide = 20; // minimal number of tracks per side
420 const Float_t kMaxDeltaZ = 30.; // maximal trigger delay
421 const Float_t kMaxDeltaV = 0.05; // maximal deltaV
2be25cc9 422 const Float_t kMaxRMS = 0.1; // maximal RMS of tracks
dde68d36 423 //
2be25cc9 424 /*
425 TCut cutRMS("sqrt(laserA.fElements[4])<0.1&&sqrt(laserC.fElements[4])<0.1");
426 TCut cutZ("abs(laserA.fElements[0]-laserC.fElements[0])<3");
427 TCut cutV("abs(laserA.fElements[1]-laserC.fElements[1])<0.01");
428 TCut cutY("abs(laserA.fElements[2]-laserC.fElements[2])<2");
429 TCut cutAll = cutRMS+cutZ+cutV+cutY;
430 */
431 if (event->GetNumberOfTracks()<kMinTracks) return;
c74c5f6c 432 //
2be25cc9 433 if(!fLaser) fLaser = new AliTPCcalibLaser("laserTPC","laserTPC",kFALSE);
434 fLaser->Process(event);
435 if (fLaser->GetNtracks()<kMinTracks) return; // small amount of tracks cut
dde68d36 436 if (fLaser->fFitAside->GetNrows()==0 && fLaser->fFitCside->GetNrows()==0) return; // no fit neither a or C side
c74c5f6c 437 //
2be25cc9 438 // debug streamer - activate stream level
439 // Use it for tuning of the cuts
da6c0bc9 440 //
dde68d36 441 // cuts to be applied
442 //
443 Int_t isReject[2]={0,0};
444 //
445 // not enough tracks
446 if (TMath::Abs((*fLaser->fFitAside)[3]) < kMinTracksSide) isReject[0]|=1;
447 if (TMath::Abs((*fLaser->fFitCside)[3]) < kMinTracksSide) isReject[1]|=1;
448 // unreasonable z offset
449 if (TMath::Abs((*fLaser->fFitAside)[0])>kMaxDeltaZ) isReject[0]|=2;
450 if (TMath::Abs((*fLaser->fFitCside)[0])>kMaxDeltaZ) isReject[1]|=2;
451 // unreasonable drift velocity
452 if (TMath::Abs((*fLaser->fFitAside)[1]-1)>kMaxDeltaV) isReject[0]|=4;
453 if (TMath::Abs((*fLaser->fFitCside)[1]-1)>kMaxDeltaV) isReject[1]|=4;
454 // big chi2
455 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitAside)[4]))>kMaxRMS ) isReject[0]|=8;
456 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitCside)[4]))>kMaxRMS ) isReject[1]|=8;
457
458
459
460
2be25cc9 461 if (fStreamLevel>0){
462 printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
463
da6c0bc9 464 TTreeSRedirector *cstream = GetDebugStreamer();
465 if (cstream){
466 TTimeStamp tstamp(fTime);
2be25cc9 467 (*cstream)<<"laserInfo"<<
da6c0bc9 468 "run="<<fRun<< // run number
469 "event="<<fEvent<< // event number
470 "time="<<fTime<< // time stamp of event
471 "trigger="<<fTrigger<< // trigger
472 "mag="<<fMagF<< // magnetic field
da6c0bc9 473 //laser
dde68d36 474 "rejectA="<<isReject[0]<<
475 "rejectC="<<isReject[1]<<
2be25cc9 476 "laserA.="<<fLaser->fFitAside<<
477 "laserC.="<<fLaser->fFitCside<<
478 "laserAC.="<<fLaser->fFitACside<<
479 "trigger="<<event->GetFiredTriggerClasses()<<
da6c0bc9 480 "\n";
481 }
482 }
2be25cc9 483 //
2be25cc9 484 // fill histos
485 //
012c4694 486 TVectorD vdriftA(5), vdriftC(5),vdriftAC(6);
2be25cc9 487 vdriftA=*(fLaser->fFitAside);
488 vdriftC=*(fLaser->fFitCside);
489 vdriftAC=*(fLaser->fFitACside);
490 Int_t npointsA=0, npointsC=0;
491 Float_t chi2A=0, chi2C=0;
492 npointsA= TMath::Nint(vdriftA[3]);
493 chi2A= vdriftA[4];
494 npointsC= TMath::Nint(vdriftC[3]);
495 chi2C= vdriftC[4];
2be25cc9 496
97d3fea0 497 if (npointsA>kMinTracksSide || npointsC>kMinTracksSide){
498 TVectorD *fitA = new TVectorD(6);
499 TVectorD *fitC = new TVectorD(6);
500 for (Int_t ipar=0; ipar<5; ipar++){
501 (*fitA)[ipar]=vdriftA[ipar];
502 (*fitC)[ipar]=vdriftC[ipar];
503 }
504 (*fitA)[5]=fTime;
505 (*fitC)[5]=fTime;
506 fArrayLaserA->AddLast(fitA);
507 fArrayLaserC->AddLast(fitC);
508 }
509 //
510
2be25cc9 511 TTimeStamp tstamp(fTime);
512 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
513 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
dde68d36 514 Double_t driftA=0, driftC=0;
515 if (vdriftA[1]>1.-kMaxDeltaV) driftA = 1./vdriftA[1]-1.;
516 if (vdriftC[1]>1.-kMaxDeltaV) driftC = 1./vdriftC[1]-1.;
517 //
518 Double_t vecDriftLaserA[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftA,event->GetRunNumber()};
519 Double_t vecDriftLaserC[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftC,event->GetRunNumber()};
520 // Double_t vecDrift[4] ={fTime,(ptrelative0+ptrelative1)/2.0,1./((*(fLaser->fFitACside))[1])-1,event->GetRunNumber()};
521
522 for (Int_t icalib=0;icalib<3;icalib++){
523 if (icalib==0){ //z0 shift
524 vecDriftLaserA[2]=vdriftA[0]/250.;
525 vecDriftLaserC[2]=vdriftC[0]/250.;
2be25cc9 526 }
dde68d36 527 if (icalib==1){ //vdrel shift
528 vecDriftLaserA[2]=driftA;
529 vecDriftLaserC[2]=driftC;
2be25cc9 530 }
dde68d36 531 if (icalib==2){ //gy shift - full gy - full drift
532 vecDriftLaserA[2]=vdriftA[2]/250.;
533 vecDriftLaserC[2]=vdriftC[2]/250.;
2be25cc9 534 }
391ffdb2 535 //if (isReject[0]==0) fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
536 //if (isReject[1]==0) fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
537 fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
538 fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
2be25cc9 539 }
2be25cc9 540}
c74c5f6c 541
2811495d 542void AliTPCcalibTime::ProcessCosmic(const AliESDEvent *const event){
6ba68367 543 //
544 // process Cosmic event - track matching A side C side
545 //
c74c5f6c 546 if (!event) {
547 Printf("ERROR: ESD not available");
548 return;
549 }
550 if (event->GetTimeStamp() == 0 ) {
551 Printf("no time stamp!");
552 return;
553 }
74235403 554
2be25cc9 555 //fd
c74c5f6c 556 // Find cosmic pairs
557 //
558 // Track0 is choosen in upper TPC part
559 // Track1 is choosen in lower TPC part
560 //
cc65e4f5 561 const Int_t kMinClustersCross =30;
562 const Int_t kMinClusters =80;
c74c5f6c 563 Int_t ntracks=event->GetNumberOfTracks();
564 if (ntracks==0) return;
74235403 565 if (ntracks > fCutTracks) return;
566
a390f11f 567 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
3e0d9569 568 AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
c74c5f6c 569
570 TObjArray tpcSeeds(ntracks);
571 Double_t vtxx[3]={0,0,0};
572 Double_t svtxx[3]={0.000001,0.000001,100.};
573 AliESDVertex vtx(vtxx,svtxx);
574 //
575 // track loop
576 //
cc65e4f5 577 TArrayI clusterSideA(ntracks);
578 TArrayI clusterSideC(ntracks);
c74c5f6c 579 for (Int_t i=0;i<ntracks;++i) {
cc65e4f5 580 clusterSideA[i]=0;
581 clusterSideC[i]=0;
74235403 582 AliESDtrack *track = event->GetTrack(i);
583
584 const AliExternalTrackParam * trackIn = track->GetInnerParam();
585 const AliExternalTrackParam * trackOut = track->GetOuterParam();
586 if (!trackIn) continue;
587 if (!trackOut) continue;
588
6ba68367 589 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
b1c27e4f 590 if (!friendTrack) continue;
3e55050f 591 if (friendTrack) ProcessSame(track,friendTrack,event);
6ba68367 592 if (friendTrack) ProcessAlignITS(track,friendTrack,event,esdFriend);
0dac7b3a 593 if (friendTrack) ProcessAlignTRD(track,friendTrack);
817766d5 594 if (friendTrack) ProcessAlignTOF(track,friendTrack);
74235403 595 TObject *calibObject;
596 AliTPCseed *seed = 0;
597 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
cc65e4f5 598 if (seed) {
599 tpcSeeds.AddAt(seed,i);
600 Int_t nA=0, nC=0;
601 for (Int_t irow=159;irow>0;irow--) {
602 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
603 if (!cl) continue;
604 if ((cl->GetDetector()%36)<18) nA++;
605 if ((cl->GetDetector()%36)>=18) nC++;
606 }
607 clusterSideA[i]=nA;
608 clusterSideC[i]=nC;
609 }
c74c5f6c 610 }
c74c5f6c 611 if (ntracks<2) return;
612 //
613 // Find pairs
614 //
d3ce44cb 615
c74c5f6c 616 for (Int_t i=0;i<ntracks;++i) {
74235403 617 AliESDtrack *track0 = event->GetTrack(i);
c74c5f6c 618 // track0 - choosen upper part
619 if (!track0) continue;
620 if (!track0->GetOuterParam()) continue;
621 if (track0->GetOuterParam()->GetAlpha()<0) continue;
622 Double_t d1[3];
623 track0->GetDirection(d1);
624 for (Int_t j=0;j<ntracks;++j) {
625 if (i==j) continue;
74235403 626 AliESDtrack *track1 = event->GetTrack(j);
627 //track 1 lower part
628 if (!track1) continue;
629 if (!track1->GetOuterParam()) continue;
cc65e4f5 630 if (track0->GetTPCNcls()+ track1->GetTPCNcls()< kMinClusters) continue;
631 Int_t nAC = TMath::Max( TMath::Min(clusterSideA[i], clusterSideC[j]),
632 TMath::Min(clusterSideC[i], clusterSideA[j]));
633 if (nAC<kMinClustersCross) continue;
634 Int_t nA0=clusterSideA[i];
635 Int_t nC0=clusterSideC[i];
636 Int_t nA1=clusterSideA[j];
637 Int_t nC1=clusterSideC[j];
3e55050f 638 // if (track1->GetOuterParam()->GetAlpha()>0) continue;
74235403 639 //
640 Double_t d2[3];
641 track1->GetDirection(d2);
642
643 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
644 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
645 if (! seed0) continue;
646 if (! seed1) continue;
647 Float_t dir = (d1[0]*d2[0] + d1[1]*d2[1] + d1[2]*d2[2]);
648 Float_t dist0 = track0->GetLinearD(0,0);
649 Float_t dist1 = track1->GetLinearD(0,0);
650 //
651 // conservative cuts - convergence to be guarantied
652 // applying before track propagation
3e55050f 653 if (TMath::Abs(TMath::Abs(dist0)-TMath::Abs(dist1))>fCutMaxD) continue; // distance to the 0,0
654 if (TMath::Abs(dir)<TMath::Abs(fCutMinDir)) continue; // direction vector product
74235403 655 Float_t bz = AliTracker::GetBz();
656 Float_t dvertex0[2]; //distance to 0,0
657 Float_t dvertex1[2]; //distance to 0,0
658 track0->GetDZ(0,0,0,bz,dvertex0);
659 track1->GetDZ(0,0,0,bz,dvertex1);
660 if (TMath::Abs(dvertex0[1])>250) continue;
661 if (TMath::Abs(dvertex1[1])>250) continue;
662 //
663 //
664 //
665 Float_t dmax = TMath::Max(TMath::Abs(dist0),TMath::Abs(dist1));
666 AliExternalTrackParam param0(*track0);
667 AliExternalTrackParam param1(*track1);
668 //
669 // Propagate using Magnetic field and correct fo material budget
670 //
9963b5e2 671 AliTracker::PropagateTrackTo(&param0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
672 AliTracker::PropagateTrackTo(&param1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
74235403 673 //
674 // Propagate rest to the 0,0 DCA - z should be ignored
675 //
676 //Bool_t b0 = ;
677 param0.PropagateToDCA(&vtx,bz,1000);
678 //Bool_t b1 =
679 param1.PropagateToDCA(&vtx,bz,1000);
74235403 680 param0.GetDZ(0,0,0,bz,dvertex0);
681 param1.GetDZ(0,0,0,bz,dvertex1);
dde68d36 682 Double_t xyz0[3];
683 Double_t xyz1[3];
74235403 684 param0.GetXYZ(xyz0);
685 param1.GetXYZ(xyz1);
686 Bool_t isPair = IsPair(&param0,&param1);
d3ce44cb 687 Bool_t isCross = IsCross(track0, track1);
dde68d36 688 Bool_t isSame = IsSame(track0, track1);
d3ce44cb 689
dde68d36 690 THnSparse* hist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
691 TString shortName=hist->ClassName();
692 shortName+="_MEAN_VDRIFT_COSMICS_";
693 delete hist;
694 hist=NULL;
695
3e55050f 696 if((isSame) || (isCross && isPair)){
697 if (track0->GetTPCNcls()+ track1->GetTPCNcls()> 80) {
74235403 698 fDz = param0.GetZ() - param1.GetZ();
cc65e4f5 699 Double_t sign=(nA0>nA1)? 1:-1;
700 fDz*=sign;
74235403 701 TTimeStamp tstamp(fTime);
702 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
703 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
dde68d36 704 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
705 THnSparse* curHist=NULL;
706 TString name="";
707
708 name=shortName;
709 name+=event->GetFiredTriggerClasses();
710 name.ToUpper();
711 curHist=(THnSparseF*)fArrayDz->FindObject(name);
74235403 712 if(!curHist){
dde68d36 713 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
714 fArrayDz->AddLast(curHist);
74235403 715 }
dde68d36 716// curHist=(THnSparseF*)(fMapDz->GetValue(event->GetFiredTriggerClasses()));
717// if(!curHist){
718// curHist=new THnSparseF(event->GetFiredTriggerClasses(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
719// fMapDz->Add(new TObjString(event->GetFiredTriggerClasses()),curHist);
720// }
721 curHist->Fill(vecDrift);
74235403 722
dde68d36 723 name=shortName;
724 name+="ALL";
725 name.ToUpper();
726 curHist=(THnSparseF*)fArrayDz->FindObject(name);
74235403 727 if(!curHist){
dde68d36 728 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
729 fArrayDz->AddLast(curHist);
74235403 730 }
dde68d36 731// curHist=(THnSparseF*)(fMapDz->GetValue("all"));
732// if(!curHist){
733// curHist=new THnSparseF("all","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
734// fMapDz->Add(new TObjString("all"),curHist);
735// }
736 curHist->Fill(vecDrift);
74235403 737 }
738 }
d3ce44cb 739 TTreeSRedirector *cstream = GetDebugStreamer();
740 if (fStreamLevel>0){
741 if (cstream){
742 (*cstream)<<"trackInfo"<<
743 "tr0.="<<track0<<
744 "tr1.="<<track1<<
745 "p0.="<<&param0<<
746 "p1.="<<&param1<<
cc65e4f5 747 "nAC="<<nAC<<
748 "nA0="<<nA0<<
749 "nA1="<<nA1<<
750 "nC0="<<nC0<<
751 "nC1="<<nC1<<
d3ce44cb 752 "isPair="<<isPair<<
753 "isCross="<<isCross<<
dde68d36 754 "isSame="<<isSame<<
d3ce44cb 755 "fDz="<<fDz<<
756 "fRun="<<fRun<<
757 "fTime="<<fTime<<
758 "\n";
759 }
760 }
c74c5f6c 761 } // end 2nd order loop
762 } // end 1st order loop
74235403 763
2be25cc9 764 if (fStreamLevel>0){
765 TTreeSRedirector *cstream = GetDebugStreamer();
766 if (cstream){
2be25cc9 767 (*cstream)<<"timeInfo"<<
768 "run="<<fRun<< // run number
769 "event="<<fEvent<< // event number
770 "time="<<fTime<< // time stamp of event
771 "trigger="<<fTrigger<< // trigger
772 "mag="<<fMagF<< // magnetic field
773 // Environment values
2be25cc9 774 //
775 // accumulated values
776 //
777 "fDz="<<fDz<< //! current delta z
2be25cc9 778 "trigger="<<event->GetFiredTriggerClasses()<<
779 "\n";
780 }
781 }
a390f11f 782 if (GetDebugLevel()>20) printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
c74c5f6c 783}
784
70036ef3 785void AliTPCcalibTime::ProcessBeam(const AliESDEvent *const event){
786 //
787 //
788 //
789 const Int_t kMinClusters =80;
b1fea14f 790 const Int_t kMinTracks =2; // minimal number of tracks to define the vertex
791 const Int_t kMinTracksVertex=30; // minimal number of tracks to define the cumulative vertex
70036ef3 792 const Double_t kMaxTgl =1.2; // maximal Tgl (z angle)
b1fea14f 793 const Double_t kMinPt =0.2; // minimal pt
794 const Double_t kMaxD0 =5.; // cut on distance to the primary vertex first guess
70036ef3 795 const Double_t kMaxZ0 =20;
b1fea14f 796 const Double_t kMaxD =2.5; // cut on distance to the primary vertex
797 const Double_t kMaxZ =4; // maximal z distance between tracks form the same side
798 const Double_t kMaxChi2 =15; // maximal chi2 of the TPCvertex
799 const Double_t kCumulCovarXY=0.003; //increase the error of cumul vertex 30 microns profile
800 const Double_t kCumulCovarZ=250.; //increase the error of cumul vertex
801 const Double_t kMaxDvertex = 1.0; // cut to accept the vertex;
802 //
803 Int_t flags=0;
804 const Int_t kBuffSize=100;
805 static Double_t deltaZ[kBuffSize]={0};
806 static Int_t counterZ=0;
807 static AliKFVertex cumulVertexA, cumulVertexC, cumulVertexAC; // cumulative vertex
808 AliKFVertex vertexA, vertexC;
70036ef3 809
810 Float_t dca0[2]={0,0};
811 Double_t dcaVertex[2]={0,0};
70036ef3 812 Int_t ntracks=event->GetNumberOfTracks();
813 if (ntracks==0) return;
012c4694 814 if (ntracks > fCutTracks) return;
70036ef3 815 //
816 AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
817 //
818 // Divide tracks to A and C side tracks - using the cluster indexes
819 TObjArray tracksA(ntracks);
820 TObjArray tracksC(ntracks);
6ba68367 821 //
70036ef3 822 AliESDVertex *vertexSPD = (AliESDVertex *)event->GetPrimaryVertexSPD();
823 AliESDVertex *vertex = (AliESDVertex *)event->GetPrimaryVertex();
824 AliESDVertex *vertexTracks = (AliESDVertex *)event->GetPrimaryVertexTracks();
825 Double_t vertexZA[10000], vertexZC[10000];
6ba68367 826 //
70036ef3 827 Int_t ntracksA= 0;
828 Int_t ntracksC= 0;
829 //
830 for (Int_t itrack=0;itrack<ntracks;itrack++) {
831 AliESDtrack *track = event->GetTrack(itrack);
832 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(itrack);
833 if (!friendTrack) continue;
834 if (TMath::Abs(track->GetTgl())>kMaxTgl) continue;
835 if (TMath::Abs(track->Pt())<kMinPt) continue;
836 const AliExternalTrackParam * trackIn = track->GetInnerParam();
837 TObject *calibObject=0;
838 AliTPCseed *seed = 0;
839 Int_t nA=0, nC=0;
840 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
841 if (seed) {
842 for (Int_t irow=159;irow>0;irow--) {
843 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
844 if (!cl) continue;
845 if ((cl->GetDetector()%36)<18) nA++;
846 if ((cl->GetDetector()%36)>=18) nC++;
847 }
848 if ((nA>kMinClusters || nC>kMinClusters) && (nA*nC==0) ){
849 track->GetImpactParameters(dca0[0],dca0[1]);
850 if (TMath::Abs(dca0[0])>kMaxD0) continue;
851 if (TMath::Abs(dca0[1])>kMaxZ0) continue;
852 AliExternalTrackParam pTPCvertex(*trackIn);
853 if (!AliTracker::PropagateTrackToBxByBz(&pTPCvertex,4.+4.*TMath::Abs(dca0[0]),0.1,2,kTRUE)) continue;
854 pTPCvertex.PropagateToDCA(vertex,AliTracker::GetBz(), kMaxD, dcaVertex,0);
855 if (TMath::Abs(dcaVertex[0])>kMaxD) continue;
856 if (nA>kMinClusters &&nC==0) { tracksA.AddLast(pTPCvertex.Clone()); vertexZA[ntracksA++] = pTPCvertex.GetZ();}
857 if (nC>kMinClusters &&nA==0) {tracksC.AddLast(pTPCvertex.Clone()); vertexZC[ntracksC++] = pTPCvertex.GetZ();}
858 }
859 }
860 }
b1fea14f 861 Double_t medianZA=TMath::Median(ntracksA, vertexZA); // tracks median
862 Double_t medianZC=TMath::Median(ntracksC, vertexZC); // tracks median
70036ef3 863 //
864 ntracksA= tracksA.GetEntriesFast();
865 ntracksC= tracksC.GetEntriesFast();
866 if (ntracksA>kMinTracks && ntracksC>kMinTracks){
b1fea14f 867 deltaZ[counterZ%kBuffSize]=medianZA-medianZC;
70036ef3 868 counterZ+=1;
b1fea14f 869 Double_t medianDelta=(counterZ>=kBuffSize)? TMath::Median(kBuffSize, deltaZ): TMath::Median(counterZ, deltaZ);
70036ef3 870 if (TMath::Abs(medianDelta-(medianZA-medianZC))>kMaxZ) flags+=16;
871 // increse the error of cumulative vertex at the beginning of event
872 cumulVertexA.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
873 cumulVertexA.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
874 cumulVertexA.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
875 cumulVertexC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
876 cumulVertexC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
877 cumulVertexC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
878 cumulVertexAC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
879 cumulVertexAC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
880 cumulVertexAC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
881 //
882 for (Int_t iA=0; iA<ntracksA; iA++){
883 if (flags!=0) continue;
884 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksA.At(iA);
885 if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
886 AliKFParticle part(*aliTrack,211);
887 vertexA+=part;
70036ef3 888 }
889 for (Int_t iC=0; iC<ntracksC; iC++){
890 if (flags!=0) continue;
891 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksC.At(iC);
892 if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
893 AliKFParticle part(*aliTrack,211);
894 vertexC+=part;
70036ef3 895 }
b1fea14f 896 //
70036ef3 897 if (vertexA.GetNDF()<kMinTracks) flags+=32;
898 if (vertexC.GetNDF()<kMinTracks) flags+=32;
899 if (TMath::Abs(vertexA.Z()-medianZA)>kMaxZ) flags+=1; //apply cuts
900 if (TMath::Abs(vertexC.Z()-medianZC)>kMaxZ) flags+=2;
901 if (TMath::Abs(vertexA.GetChi2()/vertexA.GetNDF()+vertexC.GetChi2()/vertexC.GetNDF())> kMaxChi2) flags+=4;
b1fea14f 902 //
903 if (flags==0){
904 for (Int_t iA=0; iA<ntracksA; iA++){
905 if (flags!=0) continue;
906 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksA.At(iA);
907 if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
908 AliKFParticle part(*aliTrack,211);
909 cumulVertexA+=part;
910 cumulVertexAC+=part;
911 }
912 for (Int_t iC=0; iC<ntracksC; iC++){
913 if (flags!=0) continue;
914 AliExternalTrackParam *aliTrack = (AliExternalTrackParam *)tracksC.At(iC);
915 if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
916 AliKFParticle part(*aliTrack,211);
917 cumulVertexC+=part;
918 cumulVertexAC+=part;
919 }
70036ef3 920 //
b1fea14f 921 if (TMath::Abs(cumulVertexA.X()-vertexA.X())>kMaxDvertex) flags+=64;
922 if (TMath::Abs(cumulVertexA.Y()-vertexA.Y())>kMaxDvertex) flags+=64;
923 if (TMath::Abs(cumulVertexA.Z()-vertexA.Z())>kMaxDvertex) flags+=64;
70036ef3 924 //
b1fea14f 925 if (TMath::Abs(cumulVertexC.X()-vertexC.X())>kMaxDvertex) flags+=64;
926 if (TMath::Abs(cumulVertexC.Y()-vertexC.Y())>kMaxDvertex) flags+=64;
927 if (TMath::Abs(cumulVertexC.Z()-vertexC.Z())>kMaxDvertex) flags+=64;
928
929
930 if ( flags==0 && cumulVertexC.GetNDF()>kMinTracksVertex&&cumulVertexA.GetNDF()>kMinTracksVertex){
931 Double_t cont[2]={0,fTime};
932 //
933 cont[0]= cumulVertexA.X();
934 fTPCVertex[0]->Fill(cont);
935 cont[0]= cumulVertexC.X();
936 fTPCVertex[1]->Fill(cont);
937 cont[0]= 0.5*(cumulVertexA.X()-cumulVertexC.X());
938 fTPCVertex[2]->Fill(cont);
939 cont[0]= 0.5*(cumulVertexA.X()+cumulVertexC.X())-vertexSPD->GetX();
940 fTPCVertex[3]->Fill(cont);
941 //
942 cont[0]= cumulVertexA.Y();
943 fTPCVertex[4]->Fill(cont);
944 cont[0]= cumulVertexC.Y();
945 fTPCVertex[5]->Fill(cont);
946 cont[0]= 0.5*(cumulVertexA.Y()-cumulVertexC.Y());
947 fTPCVertex[6]->Fill(cont);
948 cont[0]= 0.5*(cumulVertexA.Y()+cumulVertexC.Y())-vertexSPD->GetY();
949 fTPCVertex[7]->Fill(cont);
950 //
951 //
952 cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z());
953 fTPCVertex[8]->Fill(cont);
954 cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
955 fTPCVertex[9]->Fill(cont);
956 cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
957 fTPCVertex[10]->Fill(cont);
958 cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z())-vertexSPD->GetZ();
959 fTPCVertex[11]->Fill(cont);
960 //
961 Double_t correl[2]={0,0};
962 //
963 correl[0]=cumulVertexC.Z();
964 correl[1]=cumulVertexA.Z();
965 fTPCVertexCorrelation[0]->Fill(correl); // fill A side :TPC
966 correl[0]=cumulVertexA.Z();
967 correl[1]=cumulVertexC.Z();
968 fTPCVertexCorrelation[1]->Fill(correl); // fill C side :TPC
969 //
970 correl[0]=vertexSPD->GetZ();
971 correl[1]=cumulVertexA.Z()-correl[0];
972 fTPCVertexCorrelation[2]->Fill(correl); // fill A side :ITS
973 correl[1]=cumulVertexC.Z()-correl[0];
974 fTPCVertexCorrelation[3]->Fill(correl); // fill C side :ITS
975 correl[1]=0.5*(cumulVertexA.Z()+cumulVertexC.Z())-correl[0];
976 fTPCVertexCorrelation[4]->Fill(correl); // fill C side :ITS
977 }
978 }
70036ef3 979 TTreeSRedirector *cstream = GetDebugStreamer();
980 if (cstream){
981 /*
982 TCut cutChi2= "sqrt(vA.fChi2/vA.fNDF+vC.fChi2/vC.fNDF)<10"; // chi2 Cut e.g 10
983 TCut cutXY= "sqrt((vA.fP[0]-vC.fP[0])^2+(vA.fP[0]-vC.fP[1])^2)<5"; // vertex Cut
984 TCut cutZ= "abs(vA.fP[2]-mZA)<3&&abs(vC.fP[2]-mZC)<5"; // vertex Cut
985 tree->Draw("sqrt(vA.fChi2/vA.fNDF)","sqrt(vA.fChi2/vA.fNDF)<100","")
986
987 */
660f0543 988 //vertexA.Print();
989 //vertexC.Print();
70036ef3 990 (*cstream)<<"vertexTPC"<<
991 "flags="<<flags<< // rejection flags
992 "vSPD.="<<vertexSPD<< // SPD vertex
993 "vT.="<<vertexTracks<< // track vertex
994 "v.="<<vertex<< // esd vertex
995 "mZA="<<medianZA<< // median Z position at vertex A side
996 "mZC="<<medianZC<< // median Z position at vertex C side
997 "mDelta="<<medianDelta<< // median delta A side -C side
998 "counter="<<counterZ<< // counter Z
999 //
1000 "vA.="<<&vertexA<< // vertex A side
1001 "vC.="<<&vertexC<< // vertex C side
1002 "cvA.="<<&cumulVertexA<< // cumulative vertex A side
1003 "cvC.="<<&cumulVertexC<< // cumulative vertex C side
1004 "cvAC.="<<&cumulVertexAC<< // cumulative vertex A+C side
1005 "nA="<<ntracksA<< // contributors
1006 "nC="<<ntracksC<< // contributors
1007 "\n";
1008 }
1009 }
1010 tracksA.Delete();
1011 tracksC.Delete();
dde68d36 1012}
1013
6ba68367 1014void AliTPCcalibTime::Analyze(){
1015 //
1016 // Special macro to analyze result of calibration and extract calibration entries
1017 // Not yet ported to the Analyze function yet
1018 //
1019}
74235403 1020
2811495d 1021THnSparse* AliTPCcalibTime::GetHistoDrift(const char* name) const
1022{
6ba68367 1023 //
1024 // Get histogram for given trigger mask
1025 //
dde68d36 1026 TIterator* iterator = fArrayDz->MakeIterator();
1027 iterator->Reset();
1028 TString newName=name;
1029 newName.ToUpper();
1030 THnSparse* newHist=new THnSparseF(newName,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1031 THnSparse* addHist=NULL;
1032 while((addHist=(THnSparseF*)iterator->Next())){
472f0066 1033 // if(!addHist) continue;
dde68d36 1034 TString histName=addHist->GetName();
1035 if(!histName.Contains(newName)) continue;
1036 addHist->Print();
1037 newHist->Add(addHist);
74235403 1038 }
dde68d36 1039 return newHist;
74235403 1040}
1041
2811495d 1042TObjArray* AliTPCcalibTime::GetHistoDrift() const
1043{
6ba68367 1044 //
1045 // return array of histograms
1046 //
dde68d36 1047 return fArrayDz;
74235403 1048}
1049
1050TGraphErrors* AliTPCcalibTime::GetGraphDrift(const char* name){
6ba68367 1051 //
1052 // Make a drift velocity (delta Z) graph
1053 //
dde68d36 1054 THnSparse* histoDrift=GetHistoDrift(name);
1055 TGraphErrors* graphDrift=NULL;
1056 if(histoDrift){
1057 graphDrift=FitSlices(histoDrift,2,0,400,100,0.05,0.95, kTRUE);
1058 TString end=histoDrift->GetName();
1059 Int_t pos=end.Index("_");
1060 end=end(pos,end.Capacity()-pos);
1061 TString graphName=graphDrift->ClassName();
1062 graphName+=end;
1063 graphName.ToUpper();
1064 graphDrift->SetName(graphName);
74235403 1065 }
1066 return graphDrift;
1067}
1068
dde68d36 1069TObjArray* AliTPCcalibTime::GetGraphDrift(){
6ba68367 1070 //
1071 // make a array of drift graphs
1072 //
dde68d36 1073 TObjArray* arrayGraphDrift=new TObjArray();
1074 TIterator* iterator=fArrayDz->MakeIterator();
74235403 1075 iterator->Reset();
dde68d36 1076 THnSparse* addHist=NULL;
1077 while((addHist=(THnSparseF*)iterator->Next())) arrayGraphDrift->AddLast(GetGraphDrift(addHist->GetName()));
1078 return arrayGraphDrift;
74235403 1079}
1080
dde68d36 1081AliSplineFit* AliTPCcalibTime::GetFitDrift(const char* name){
6ba68367 1082 //
1083 // Make a fit AliSplinefit of drift velocity
1084 //
d3ce44cb 1085 TGraph* graphDrift=GetGraphDrift(name);
dde68d36 1086 AliSplineFit* fitDrift=NULL;
74235403 1087 if(graphDrift && graphDrift->GetN()){
dde68d36 1088 fitDrift=new AliSplineFit();
1089 fitDrift->SetGraph(graphDrift);
1090 fitDrift->SetMinPoints(graphDrift->GetN()+1);
1091 fitDrift->InitKnots(graphDrift,2,0,0.001);
1092 fitDrift->SplineFit(0);
1093 TString end=graphDrift->GetName();
1094 Int_t pos=end.Index("_");
1095 end=end(pos,end.Capacity()-pos);
1096 TString fitName=fitDrift->ClassName();
1097 fitName+=end;
1098 fitName.ToUpper();
1099 //fitDrift->SetName(fitName);
74235403 1100 delete graphDrift;
dde68d36 1101 graphDrift=NULL;
74235403 1102 }
1103 return fitDrift;
1104}
1105
74235403 1106
2811495d 1107Long64_t AliTPCcalibTime::Merge(TCollection *const li) {
6ba68367 1108 //
1109 // Object specific merging procedure
1110 //
c74c5f6c 1111 TIterator* iter = li->MakeIterator();
1112 AliTPCcalibTime* cal = 0;
660f0543 1113 //
c74c5f6c 1114 while ((cal = (AliTPCcalibTime*)iter->Next())) {
1115 if (!cal->InheritsFrom(AliTPCcalibTime::Class())) {
1116 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
1117 return -1;
1118 }
2be25cc9 1119 for (Int_t imeas=0; imeas<3; imeas++){
1120 if (cal->GetHistVdriftLaserA(imeas) && cal->GetHistVdriftLaserA(imeas)){
1121 fHistVdriftLaserA[imeas]->Add(cal->GetHistVdriftLaserA(imeas));
1122 fHistVdriftLaserC[imeas]->Add(cal->GetHistVdriftLaserC(imeas));
1123 }
1124 }
9112627b 1125 //
b1fea14f 1126 if (fTPCVertexCorrelation && cal->fTPCVertexCorrelation){
1127 for (Int_t imeas=0; imeas<5; imeas++){
1128 if (fTPCVertexCorrelation[imeas] && cal->fTPCVertexCorrelation[imeas]) fTPCVertexCorrelation[imeas]->Add(cal->fTPCVertexCorrelation[imeas]);
1129 }
1130 }
1131
472f0066 1132 if (fTPCVertex[0] && cal->fTPCVertex[0])
70036ef3 1133 for (Int_t imeas=0; imeas<12; imeas++){
660f0543 1134 if (fTPCVertex[imeas] && cal->fTPCVertex[imeas]) fTPCVertex[imeas]->Add(cal->fTPCVertex[imeas]);
70036ef3 1135 }
1136
660f0543 1137 if (fMemoryMode>0) for (Int_t imeas=0; imeas<5; imeas++){
1138 if (fMemoryMode>1){
1139 if ( cal->GetResHistoTPCCE(imeas) && cal->GetResHistoTPCCE(imeas)){
1140 fResHistoTPCCE[imeas]->Add(cal->fResHistoTPCCE[imeas]);
1141 }else{
1142 fResHistoTPCCE[imeas]=(THnSparse*)cal->fResHistoTPCCE[imeas]->Clone();
1143 }
b9908d0b 1144 }
660f0543 1145 //
1146 if ((fMemoryMode>0) &&cal->GetResHistoTPCITS(imeas) && cal->GetResHistoTPCITS(imeas)){
1147 if (fMemoryMode>1 || (imeas%2)==1) fResHistoTPCITS[imeas]->Add(cal->fResHistoTPCITS[imeas]);
1148 if (fMemoryMode>1) fResHistoTPCvertex[imeas]->Add(cal->fResHistoTPCvertex[imeas]);
33fc6b9e 1149 }
660f0543 1150 //
1151 if ((fMemoryMode>1) && cal->fResHistoTPCTRD[imeas]){
33fc6b9e 1152 if (fResHistoTPCTRD[imeas])
1153 fResHistoTPCTRD[imeas]->Add(cal->fResHistoTPCTRD[imeas]);
1154 else
1155 fResHistoTPCTRD[imeas]=(THnSparse*)cal->fResHistoTPCTRD[imeas]->Clone();
1156 }
660f0543 1157 //
1158 if ((fMemoryMode>1) && cal->fResHistoTPCTOF[imeas]){
33fc6b9e 1159 if (fResHistoTPCTOF[imeas])
1160 fResHistoTPCTOF[imeas]->Add(cal->fResHistoTPCTOF[imeas]);
1161 else
1162 fResHistoTPCTOF[imeas]=(THnSparse*)cal->fResHistoTPCTOF[imeas]->Clone();
9112627b 1163 }
660f0543 1164 //
97d3fea0 1165 if (cal->fArrayLaserA){
1166 fArrayLaserA->Expand(fArrayLaserA->GetEntriesFast()+cal->fArrayLaserA->GetEntriesFast());
1167 fArrayLaserC->Expand(fArrayLaserC->GetEntriesFast()+cal->fArrayLaserC->GetEntriesFast());
1168 for (Int_t ical=0; ical<cal->fArrayLaserA->GetEntriesFast(); ical++){
1169 if (cal->fArrayLaserA->UncheckedAt(ical)) fArrayLaserA->AddLast(cal->fArrayLaserA->UncheckedAt(ical)->Clone());
1170 if (cal->fArrayLaserC->UncheckedAt(ical)) fArrayLaserC->AddLast(cal->fArrayLaserC->UncheckedAt(ical)->Clone());
1171 }
1172 }
1173
9112627b 1174 }
dde68d36 1175 TObjArray* addArray=cal->GetHistoDrift();
1176 if(!addArray) return 0;
1177 TIterator* iterator = addArray->MakeIterator();
2be25cc9 1178 iterator->Reset();
dde68d36 1179 THnSparse* addHist=NULL;
660f0543 1180 if ((fMemoryMode>1)) while((addHist=(THnSparseF*)iterator->Next())){
472f0066 1181 // if(!addHist) continue;
2be25cc9 1182 addHist->Print();
dde68d36 1183 THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
2be25cc9 1184 if(!localHist){
1185 localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
dde68d36 1186 fArrayDz->AddLast(localHist);
2be25cc9 1187 }
1188 localHist->Add(addHist);
1189 }
391ffdb2 1190
d3ce44cb 1191 for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
3e55050f 1192 //
1193 // Merge alignment
1194 //
1195 for (Int_t itype=0; itype<3; itype++){
1196 //
1197 //
1198 TObjArray *arr0= 0;
1199 TObjArray *arr1= 0;
1200 if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
1201 if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
1202 if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
1203 if (!arr1) continue;
1204 if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
1205 if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
1206 arr0->Expand(arr1->GetEntriesFast());
1207 }
1208 for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
1209 AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
1210 AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
1211 if (!kalman1) continue;
1212 if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
1213 kalman0->SetRejectOutliers(kFALSE);
1214 kalman0->Merge(kalman1);
1215 }
1216 }
1217
c74c5f6c 1218 }
c74c5f6c 1219 return 0;
c74c5f6c 1220}
1221
c74c5f6c 1222Bool_t AliTPCcalibTime::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
c74c5f6c 1223 /*
1224 // 0. Same direction - OPOSITE - cutDir +cutT
1225 TCut cutDir("cutDir","dir<-0.99")
1226 // 1.
1227 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
1228 //
1229 // 2. The same rphi
1230 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
1231 //
c74c5f6c 1232 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
1233 // 1/Pt diff cut
1234 */
1235 const Double_t *p0 = tr0->GetParameter();
1236 const Double_t *p1 = tr1->GetParameter();
74235403 1237 fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
1238 fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
d3ce44cb 1239 fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
74235403 1240 fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
1241 fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
1242
c74c5f6c 1243 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
1244 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
da6c0bc9 1245 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
c74c5f6c 1246 Double_t d0[3], d1[3];
1247 tr0->GetDirection(d0);
1248 tr1->GetDirection(d1);
1249 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
74235403 1250
d3ce44cb 1251 fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
1252 fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
1253 fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
1254 fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
1255 fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
1256
c74c5f6c 1257 return kTRUE;
1258}
2811495d 1259Bool_t AliTPCcalibTime::IsCross(AliESDtrack *const tr0, AliESDtrack *const tr1){
cc65e4f5 1260 //
1261 // check if the cosmic pair of tracks crossed A/C side
1262 //
1263 Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
1264 if (result==kFALSE) return result;
1265 result=kTRUE;
1266 return result;
d3ce44cb 1267}
1268
2811495d 1269Bool_t AliTPCcalibTime::IsSame(AliESDtrack *const tr0, AliESDtrack *const tr1){
3e55050f 1270 //
1271 // track crossing the CE
1272 // 0. minimal number of clusters
1273 // 1. Same sector +-1
1274 // 2. Inner and outer track param on opposite side
1275 // 3. Outer and inner track parameter close each to other
1276 // 3.
1277 Bool_t result=kTRUE;
1278 //
1279 // inner and outer on opposite sides in z
1280 //
1281 const Int_t knclCut0 = 30;
1282 const Double_t kalphaCut = 0.4;
1283 //
1284 // 0. minimal number of clusters
1285 //
1286 if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
1287 if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
1288 //
1289 // 1. alpha cut - sector+-1
1290 //
1291 if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
1292 //
1293 // 2. Z crossing
1294 //
1295 if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
1296 if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
1297 if (result==kFALSE){
1298 return result;
1299 }
1300 //
1301 //
1302 const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
1303 const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
1304 const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
1305 const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
1306 //
1307 if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD) result&=kFALSE;
1308 if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
1309 if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1310 if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1311 if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD) result&=kFALSE;
1312 if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1313 if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1314 if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1315 if (result==kTRUE){
1316 result=kTRUE; // just to put break point here
1317 }
1318 return result;
dde68d36 1319}
1320
31aa7c5c 1321
2811495d 1322void AliTPCcalibTime::ProcessSame(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event){
3e55050f 1323 //
1324 // Process TPC tracks crossing CE
1325 //
1326 // 0. Select only track crossing the CE
1327 // 1. Cut on the track length
b9908d0b 1328 // 2. Refit the the track on A and C side separatelly
3e55050f 1329 // 3. Fill time histograms
1330 const Int_t kMinNcl=100;
1331 const Int_t kMinNclS=25; // minimul number of clusters on the sides
b9908d0b 1332 const Double_t pimass=TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1333 const Double_t kMaxDy=1; // maximal distance in y
1334 const Double_t kMaxDsnp=0.05; // maximal distance in snp
1335 const Double_t kMaxDtheta=0.05; // maximal distance in theta
1336
3e55050f 1337 if (!friendTrack->GetTPCOut()) return;
1338 //
1339 // 0. Select only track crossing the CE
1340 //
1341 if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1342 //
1343 // 1. cut on track length
1344 //
1345 if (track->GetTPCNcls()<kMinNcl) return;
1346 //
1347 // 2. Refit track sepparatel on A and C side
1348 //
1349 TObject *calibObject;
1350 AliTPCseed *seed = 0;
1351 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1352 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1353 }
1354 if (!seed) return;
1355 //
1356 AliExternalTrackParam trackIn(*track->GetInnerParam());
1357 AliExternalTrackParam trackOut(*track->GetOuterParam());
1358 Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1359 Double_t xyz[3]={0,0.,0.0};
1360 Double_t bz =0;
1361 Int_t nclIn=0,nclOut=0;
70036ef3 1362 trackIn.ResetCovariance(1000.);
1363 trackOut.ResetCovariance(1000.);
3e55050f 1364 //
1365 //2.a Refit inner
1366 //
b9908d0b 1367 Int_t sideIn=0;
3e55050f 1368 for (Int_t irow=0;irow<159;irow++) {
1369 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1370 if (!cl) continue;
1371 if (cl->GetX()<80) continue;
b9908d0b 1372 if (sideIn==0){
1373 if (cl->GetDetector()%36<18) sideIn=1;
1374 if (cl->GetDetector()%36>=18) sideIn=-1;
1375 }
1376 if (sideIn== -1 && (cl->GetDetector()%36)<18) break;
1377 if (sideIn== 1 &&(cl->GetDetector()%36)>=18) break;
3e55050f 1378 Int_t sector = cl->GetDetector();
1379 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1380 if (TMath::Abs(dalpha)>0.01){
1381 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1382 }
1383 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1384 trackIn.GetXYZ(xyz);
1385 bz = AliTracker::GetBz(xyz);
b9908d0b 1386 AliTracker::PropagateTrackToBxByBz(&trackIn,r[0],1.,pimass,kFALSE);
3e55050f 1387 if (!trackIn.PropagateTo(r[0],bz)) break;
1388 nclIn++;
1389 trackIn.Update(&r[1],cov);
1390 }
1391 //
1392 //2.b Refit outer
b9908d0b 1393 //
1394 Int_t sideOut=0;
3e55050f 1395 for (Int_t irow=159;irow>0;irow--) {
1396 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1397 if (!cl) continue;
1398 if (cl->GetX()<80) continue;
b9908d0b 1399 if (sideOut==0){
1400 if (cl->GetDetector()%36<18) sideOut=1;
1401 if (cl->GetDetector()%36>=18) sideOut=-1;
1402 if (sideIn==sideOut) break;
1403 }
1404 if (sideOut== -1 && (cl->GetDetector()%36)<18) break;
1405 if (sideOut== 1 &&(cl->GetDetector()%36)>=18) break;
1406 //
3e55050f 1407 Int_t sector = cl->GetDetector();
1408 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1409 if (TMath::Abs(dalpha)>0.01){
1410 if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1411 }
1412 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1413 trackOut.GetXYZ(xyz);
1414 bz = AliTracker::GetBz(xyz);
b9908d0b 1415 AliTracker::PropagateTrackToBxByBz(&trackOut,r[0],1.,pimass,kFALSE);
3e55050f 1416 if (!trackOut.PropagateTo(r[0],bz)) break;
1417 nclOut++;
1418 trackOut.Update(&r[1],cov);
1419 }
1420 trackOut.Rotate(trackIn.GetAlpha());
1421 Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1422 trackIn.PropagateTo(meanX,bz);
1423 trackOut.PropagateTo(meanX,bz);
b9908d0b 1424 if (TMath::Abs(trackIn.GetY()-trackOut.GetY())>kMaxDy) return;
1425 if (TMath::Abs(trackIn.GetSnp()-trackOut.GetSnp())>kMaxDsnp) return;
1426 if (TMath::Abs(trackIn.GetTgl()-trackOut.GetTgl())>kMaxDtheta) return;
1427 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1428 FillResHistoTPCCE(&trackIn,&trackOut);
1429 }
3e55050f 1430 TTreeSRedirector *cstream = GetDebugStreamer();
1431 if (cstream){
1432 TVectorD gxyz(3);
1433 trackIn.GetXYZ(gxyz.GetMatrixArray());
1434 TTimeStamp tstamp(fTime);
3e55050f 1435 (*cstream)<<"tpctpc"<<
1436 "run="<<fRun<< // run number
1437 "event="<<fEvent<< // event number
1438 "time="<<fTime<< // time stamp of event
1439 "trigger="<<fTrigger<< // trigger
1440 "mag="<<fMagF<< // magnetic field
3e55050f 1441 //
b9908d0b 1442 "sideIn="<<sideIn<< // side at inner part
1443 "sideOut="<<sideOut<< // side at puter part
3e55050f 1444 "xyz.="<<&gxyz<< // global position
1445 "tIn.="<<&trackIn<< // refitterd track in
1446 "tOut.="<<&trackOut<< // refitter track out
1447 "nclIn="<<nclIn<< //
1448 "nclOut="<<nclOut<< //
1449 "\n";
1450 }
1451 //
1452 // 3. Fill time histograms
1453 // Debug stremaer expression
1454 // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1455 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1456 fDz = trackOut.GetZ()-trackIn.GetZ();
1457 if (trackOut.GetTgl()<0) fDz*=-1.;
1458 TTimeStamp tstamp(fTime);
1459 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1460 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1461 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1462 //
1463 // fill histograms per trigger class and itegrated
1464 //
1465 THnSparse* curHist=NULL;
1466 for (Int_t itype=0; itype<2; itype++){
1467 TString name="MEAN_VDRIFT_CROSS_";
1468 if (itype==0){
1469 name+=event->GetFiredTriggerClasses();
1470 name.ToUpper();
1471 }else{
1472 name+="ALL";
1473 }
1474 curHist=(THnSparseF*)fArrayDz->FindObject(name);
1475 if(!curHist){
1476 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1477 fArrayDz->AddLast(curHist);
1478 }
1479 curHist->Fill(vecDrift);
1480 }
1481 }
1482
1483}
817766d5 1484
2811495d 1485void AliTPCcalibTime::ProcessAlignITS(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event, AliESDfriend *const esdFriend){
817766d5 1486 //
3e55050f 1487 // Process track - Update TPC-ITS alignment
1488 // Updates:
1489 // 0. Apply standartd cuts
1490 // 1. Recalucluate the current statistic median/RMS
1491 // 2. Apply median+-rms cut
1492 // 3. Update kalman filter
1493 //
1494 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1495 const Int_t kMinITS = 3; // minimal number of ITS cluster
1496 const Double_t kMinZ = 10; // maximal dz distance
b96c3aef 1497 const Double_t kMaxDy = 2.; // maximal dy distance
b1fea14f 1498 const Double_t kMaxAngle= 0.07; // maximal angular distance
3e55050f 1499 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1500 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
97d3fea0 1501 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
3e55050f 1502 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1503 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
b96c3aef 1504 const Double_t kMinPt = 0.3; // minimal pt
b1fea14f 1505 const Double_t kMax1Pt=0.5; //maximal 1/pt distance
be67055b 1506 const Int_t kN=50; // deepnes of history
3e55050f 1507 static Int_t kglast=0;
1508 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
817766d5 1509 //
3e55050f 1510 // 0. Apply standard cuts
b1fea14f 1511 //
817766d5 1512 Int_t dummycl[1000];
817766d5 1513 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
b96c3aef 1514 if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
817766d5 1515 if (!track->GetInnerParam()) return;
1516 if (!track->GetOuterParam()) return;
b842d904 1517 if (track->GetInnerParam()->Pt()<kMinPt) return;
817766d5 1518 // exclude crossing track
1519 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
70036ef3 1520 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ/3.) return;
b96c3aef 1521 if (track->GetInnerParam()->GetX()>90) return;
817766d5 1522 //
1523 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
b1fea14f 1524 //
1525 AliExternalTrackParam pITS; // ITS standalone if possible
1526 AliExternalTrackParam pITS2; //TPC-ITS track
1527 if (friendTrack->GetITSOut()){
1528 pITS2=(*(friendTrack->GetITSOut())); //TPC-ITS track - snapshot ITS out
1529 pITS2.Rotate(pTPC.GetAlpha());
1530 AliTracker::PropagateTrackToBxByBz(&pITS2,pTPC.GetX(),0.1,0.1,kFALSE);
1531 }
be67055b 1532
b842d904 1533 AliESDfriendTrack *itsfriendTrack=0;
1534 //
1535 // try to find standalone ITS track corresponing to the TPC if possible
1536 //
1537 Bool_t hasAlone=kFALSE;
1538 Int_t ntracks=event->GetNumberOfTracks();
1539 for (Int_t i=0; i<ntracks; i++){
b1fea14f 1540 AliESDtrack * trackITS = event->GetTrack(i);
1541 if (!trackITS) continue;
1542 if (trackITS->GetITSclusters(dummycl)<kMinITS) continue; // minimal amount of clusters
6ba68367 1543 itsfriendTrack = esdFriend->GetTrack(i);
b842d904 1544 if (!itsfriendTrack) continue;
1545 if (!itsfriendTrack->GetITSOut()) continue;
b1fea14f 1546
1547 if (TMath::Abs(pTPC.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1548 if (TMath::Abs(pTPC.GetSigned1Pt()-itsfriendTrack->GetITSOut()->GetSigned1Pt())> kMax1Pt) continue;
b842d904 1549 pITS=(*(itsfriendTrack->GetITSOut()));
1550 //
1551 pITS.Rotate(pTPC.GetAlpha());
be67055b 1552 AliTracker::PropagateTrackToBxByBz(&pITS,pTPC.GetX(),0.1,0.1,kFALSE);
b1fea14f 1553 if (TMath::Abs(pTPC.GetY()-pITS.GetY())> kMaxDy) continue;
1554 if (TMath::Abs(pTPC.GetSnp()-pITS.GetSnp())> kMaxAngle) continue;
b842d904 1555 hasAlone=kTRUE;
1556 }
b1fea14f 1557 if (!hasAlone) {
1558 if (track->GetITSclusters(dummycl)<kMinITS) return;
1559 pITS=pITS2; // use combined track if it has ITS
1560 }
b842d904 1561 //
3e55050f 1562 if (TMath::Abs(pITS.GetY()-pTPC.GetY()) >kMaxDy) return;
1563 if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1564 if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1565 //
1566 // 1. Update median and RMS info
1567 //
b322e06a 1568 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
3e55050f 1569 TVectorD vecDeltaN(5);
1570 Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1571 vecDelta[4]=0;
1572 for (Int_t i=0;i<4;i++){
1573 vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1574 kgdP[i][kglast%kN]=vecDelta[i];
1575 }
1576 kglast=(kglast+1);
1577 Int_t entries=(kglast<kN)?kglast:kN;
1578 for (Int_t i=0;i<4;i++){
1579 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1580 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1581 vecDeltaN[i] = 0;
1582 if (vecRMS[i]>0.){
1583 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1584 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1585 }
1586 }
817766d5 1587 //
3e55050f 1588 // 2. Apply median+-rms cut
817766d5 1589 //
3e55050f 1590 if (kglast<3) return; //median and RMS to be defined
1591 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1592 //
1593 // 3. Update alignment
817766d5 1594 //
97d3fea0 1595 Int_t htime = (fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time bins number
4486a91f 1596 if (fAlignITSTPC->GetEntriesFast()<htime){
817766d5 1597 fAlignITSTPC->Expand(htime*2+20);
1598 }
1599 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1600 if (!align){
3e55050f 1601 // make Alignment object if doesn't exist
817766d5 1602 align=new AliRelAlignerKalman();
3e55050f 1603 align->SetRunNumber(fRun);
1604 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
97d3fea0 1605 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
3e55050f 1606 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1607 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1608 align->SetRejectOutliers(kFALSE);
1609
1610 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
817766d5 1611 align->SetMagField(fMagF);
1612 fAlignITSTPC->AddAt(align,htime);
1613 }
817766d5 1614 align->AddTrackParams(&pITS,&pTPC);
5647625c 1615 Double_t averageTime = fTime;
1616 if (align->GetTimeStamp()>0&&align->GetNUpdates()>0){
97d3fea0 1617 averageTime=((Double_t(align->GetTimeStamp())*Double_t(align->GetNUpdates())+Double_t(fTime)))/(Double_t(align->GetNUpdates())+1.);
5647625c 1618 }
1619 align->SetTimeStamp(Int_t(averageTime));
1620
3e55050f 1621 align->SetRunNumber(fRun );
391ffdb2 1622 Float_t dca[2],cov[3];
1623 track->GetImpactParameters(dca,cov);
882f5c06 1624 if (TMath::Abs(dca[0])<kMaxDy){
391ffdb2 1625 FillResHistoTPCITS(&pTPC,&pITS);
1626 FillResHistoTPC(track);
1627 }
3e55050f 1628 //
1629 Int_t nupdates=align->GetNUpdates();
1630 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1631 align->SetRejectOutliers(kFALSE);
817766d5 1632 TTreeSRedirector *cstream = GetDebugStreamer();
1633 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
817766d5 1634 TVectorD gpTPC(3), gdTPC(3);
1635 TVectorD gpITS(3), gdITS(3);
1636 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1637 pTPC.GetDirection(gdTPC.GetMatrixArray());
1638 pITS.GetXYZ(gpITS.GetMatrixArray());
1639 pITS.GetDirection(gdITS.GetMatrixArray());
1640 (*cstream)<<"itstpc"<<
1641 "run="<<fRun<< // run number
1642 "event="<<fEvent<< // event number
1643 "time="<<fTime<< // time stamp of event
1644 "trigger="<<fTrigger<< // trigger
1645 "mag="<<fMagF<< // magnetic field
817766d5 1646 //
b842d904 1647 "hasAlone="<<hasAlone<< // has ITS standalone ?
1648 "track.="<<track<< // track info
3e55050f 1649 "nmed="<<kglast<< // number of entries to define median and RMS
1650 "vMed.="<<&vecMedian<< // median of deltas
1651 "vRMS.="<<&vecRMS<< // rms of deltas
1652 "vDelta.="<<&vecDelta<< // delta in respect to median
1653 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1654 "t.="<<track<< // ful track - find proper cuts
1655 "a.="<<align<< // current alignment
b842d904 1656 "pITS.="<<&pITS<< // track param ITS
1657 "pITS2.="<<&pITS2<< // track param ITS+TPC
3e55050f 1658 "pTPC.="<<&pTPC<< // track param TPC
1659 "gpTPC.="<<&gpTPC<< // global position TPC
1660 "gdTPC.="<<&gdTPC<< // global direction TPC
1661 "gpITS.="<<&gpITS<< // global position ITS
1662 "gdITS.="<<&gdITS<< // global position ITS
817766d5 1663 "\n";
1664 }
817766d5 1665}
3e55050f 1666
1667
1668
1669
2811495d 1670void AliTPCcalibTime::ProcessAlignTRD(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
0dac7b3a 1671 //
3e55050f 1672 // Process track - Update TPC-TRD alignment
1673 // Updates:
1674 // 0. Apply standartd cuts
1675 // 1. Recalucluate the current statistic median/RMS
1676 // 2. Apply median+-rms cut
1677 // 3. Update kalman filter
1678 //
1679 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1680 const Int_t kMinTRD = 50; // minimal number of TRD cluster
1681 const Double_t kMinZ = 20; // maximal dz distance
9112627b 1682 const Double_t kMaxDy = 5.; // maximal dy distance
1683 const Double_t kMaxAngle= 0.1; // maximal angular distance
1684 const Double_t kSigmaCut= 10; // maximal sigma distance to median
3e55050f 1685 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
97d3fea0 1686 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
3e55050f 1687 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1688 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
be67055b 1689 const Double_t kRefX = 275; // reference X
1690 const Int_t kN=50; // deepnes of history
3e55050f 1691 static Int_t kglast=0;
1692 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
0dac7b3a 1693 //
3e55050f 1694 // 0. Apply standard cuts
0dac7b3a 1695 //
1696 Int_t dummycl[1000];
1697 if (track->GetTRDclusters(dummycl)<kMinTRD) return; // minimal amount of clusters
1698 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
3e55050f 1699 if (!friendTrack->GetTRDIn()) return;
be67055b 1700 if (!track->IsOn(AliESDtrack::kTRDrefit)) return;
1701 if (!track->IsOn(AliESDtrack::kTRDout)) return;
0dac7b3a 1702 if (!track->GetInnerParam()) return;
be67055b 1703 if (!friendTrack->GetTPCOut()) return;
0dac7b3a 1704 // exclude crossing track
be67055b 1705 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
3e55050f 1706 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
0dac7b3a 1707 //
be67055b 1708 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(friendTrack->GetTPCOut()));
1709 AliTracker::PropagateTrackToBxByBz(&pTPC,kRefX,0.1,0.1,kFALSE);
0dac7b3a 1710 AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
3e55050f 1711 pTRD.Rotate(pTPC.GetAlpha());
be67055b 1712 // pTRD.PropagateTo(pTPC.GetX(),fMagF);
1713 AliTracker::PropagateTrackToBxByBz(&pTRD,pTPC.GetX(),0.1,0.1,kFALSE);
1714
3e55050f 1715 ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.; // increas sys errors
1716 ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1717
1718 if (TMath::Abs(pTRD.GetY()-pTPC.GetY()) >kMaxDy) return;
1719 if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
9112627b 1720 // if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
0dac7b3a 1721 //
3e55050f 1722 // 1. Update median and RMS info
0dac7b3a 1723 //
b322e06a 1724 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
3e55050f 1725 TVectorD vecDeltaN(5);
1726 Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1727 vecDelta[4]=0;
1728 for (Int_t i=0;i<4;i++){
1729 vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1730 kgdP[i][kglast%kN]=vecDelta[i];
1731 }
1732 kglast=(kglast+1);
1733 Int_t entries=(kglast<kN)?kglast:kN;
1734 for (Int_t i=0;i<4;i++){
1735 vecMedian[i] = TMath::Median(entries,kgdP[i]);
be67055b 1736
3e55050f 1737 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1738 vecDeltaN[i] = 0;
1739 if (vecRMS[i]>0.){
1740 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1741 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1742 }
1743 }
1744 //
1745 // 2. Apply median+-rms cut
1746 //
1747 if (kglast<3) return; //median and RMS to be defined
1748 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1749 //
1750 // 3. Update alignment
0dac7b3a 1751 //
5647625c 1752 //Int_t htime = fTime/3600; //time in hours
97d3fea0 1753 Int_t htime = (Int_t)(fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time in half hour
4486a91f 1754 if (fAlignTRDTPC->GetEntriesFast()<htime){
0dac7b3a 1755 fAlignTRDTPC->Expand(htime*2+20);
1756 }
1757 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1758 if (!align){
3e55050f 1759 // make Alignment object if doesn't exist
0dac7b3a 1760 align=new AliRelAlignerKalman();
3e55050f 1761 align->SetRunNumber(fRun);
1762 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
97d3fea0 1763 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
3e55050f 1764 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1765 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1766 align->SetRejectOutliers(kFALSE);
1767 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
0dac7b3a 1768 align->SetMagField(fMagF);
1769 fAlignTRDTPC->AddAt(align,htime);
1770 }
0dac7b3a 1771 align->AddTrackParams(&pTRD,&pTPC);
5647625c 1772 //align->SetTimeStamp(fTime);
1773 Double_t averageTime = fTime;
1774 if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1775 averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1776 //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1777 }
1778 align->SetTimeStamp((Int_t)averageTime);
1779
1780 //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1781
3e55050f 1782 align->SetRunNumber(fRun );
0b736a46 1783 Float_t dca[2],cov[3];
1784 track->GetImpactParameters(dca,cov);
1785 if (TMath::Abs(dca[0])<kMaxDy){
1786 FillResHistoTPCTRD(&pTPC,&pTRD); //only primaries
1787 }
3e55050f 1788 //
1789 Int_t nupdates=align->GetNUpdates();
1790 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1791 align->SetRejectOutliers(kFALSE);
0dac7b3a 1792 TTreeSRedirector *cstream = GetDebugStreamer();
1793 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
0dac7b3a 1794 TVectorD gpTPC(3), gdTPC(3);
1795 TVectorD gpTRD(3), gdTRD(3);
1796 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1797 pTPC.GetDirection(gdTPC.GetMatrixArray());
1798 pTRD.GetXYZ(gpTRD.GetMatrixArray());
1799 pTRD.GetDirection(gdTRD.GetMatrixArray());
3e55050f 1800 (*cstream)<<"trdtpc"<<
0dac7b3a 1801 "run="<<fRun<< // run number
1802 "event="<<fEvent<< // event number
1803 "time="<<fTime<< // time stamp of event
1804 "trigger="<<fTrigger<< // trigger
1805 "mag="<<fMagF<< // magnetic field
0dac7b3a 1806 //
3e55050f 1807 "nmed="<<kglast<< // number of entries to define median and RMS
1808 "vMed.="<<&vecMedian<< // median of deltas
1809 "vRMS.="<<&vecRMS<< // rms of deltas
1810 "vDelta.="<<&vecDelta<< // delta in respect to median
1811 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1812 "t.="<<track<< // ful track - find proper cuts
1813 "a.="<<align<< // current alignment
1814 "pTRD.="<<&pTRD<< // track param TRD
1815 "pTPC.="<<&pTPC<< // track param TPC
1816 "gpTPC.="<<&gpTPC<< // global position TPC
1817 "gdTPC.="<<&gdTPC<< // global direction TPC
1818 "gpTRD.="<<&gpTRD<< // global position TRD
1819 "gdTRD.="<<&gdTRD<< // global position TRD
0dac7b3a 1820 "\n";
1821 }
0dac7b3a 1822}
817766d5 1823
1824
2811495d 1825void AliTPCcalibTime::ProcessAlignTOF(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
817766d5 1826 //
817766d5 1827 //
3e55050f 1828 // Process track - Update TPC-TOF alignment
1829 // Updates:
1830 // -1. Make a TOF "track"
1831 // 0. Apply standartd cuts
1832 // 1. Recalucluate the current statistic median/RMS
1833 // 2. Apply median+-rms cut
1834 // 3. Update kalman filter
1835 //
1836 const Int_t kMinTPC = 80; // minimal number of TPC cluster
b842d904 1837 // const Double_t kMinZ = 10; // maximal dz distance
3e55050f 1838 const Double_t kMaxDy = 5.; // maximal dy distance
70036ef3 1839 const Double_t kMaxAngle= 0.05; // maximal angular distance
3e55050f 1840 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1841 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
97d3fea0 1842 const Double_t kT0Err = 3.; // initial uncertainty of the T0 time
3e55050f 1843 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1844
1845 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
be67055b 1846 const Int_t kN=50; // deepnes of history
3e55050f 1847 static Int_t kglast=0;
1848 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
817766d5 1849 //
3e55050f 1850 // -1. Make a TOF track-
1851 // Clusters are not in friends - use alingment points
1852 //
1853 if (track->GetTOFsignal()<=0) return;
1854 if (!friendTrack->GetTPCOut()) return;
1855 if (!track->GetInnerParam()) return;
be67055b 1856 if (!friendTrack->GetTPCOut()) return;
817766d5 1857 const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1858 if (!points) return;
be67055b 1859 AliExternalTrackParam pTPC(*(friendTrack->GetTPCOut()));
3e55050f 1860 AliExternalTrackParam pTOF(pTPC);
1861 Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
817766d5 1862 Int_t npoints = points->GetNPoints();
1863 AliTrackPoint point;
3e55050f 1864 Int_t naccept=0;
817766d5 1865 //
1866 for (Int_t ipoint=0;ipoint<npoints;ipoint++){
817766d5 1867 points->GetPoint(point,ipoint);
1868 Float_t xyz[3];
1869 point.GetXYZ(xyz);
1870 Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
3e55050f 1871 if (r<350) continue;
817766d5 1872 if (r>400) continue;
3e55050f 1873 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1874 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);
1875 AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1876 pTPC.PropagateTo(lpoint.GetX(),fMagF);
1877 pTOF=pTPC;
817766d5 1878 ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1879 ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
3e55050f 1880 ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1881 ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
0dac7b3a 1882 ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1883 ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
3e55050f 1884 naccept++;
1885 }
1886 if (naccept==0) return; // no tof match clusters
1887 //
1888 // 0. Apply standard cuts
1889 //
1890 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1891 // exclude crossing track
be67055b 1892 if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
3e55050f 1893 //
1894 if (TMath::Abs(pTOF.GetY()-pTPC.GetY()) >kMaxDy) return;
1895 if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1896 if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1897 //
1898 // 1. Update median and RMS info
1899 //
b322e06a 1900 TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
3e55050f 1901 TVectorD vecDeltaN(5);
1902 Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1903 vecDelta[4]=0;
1904 for (Int_t i=0;i<4;i++){
1905 vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1906 kgdP[i][kglast%kN]=vecDelta[i];
1907 }
1908 kglast=(kglast+1);
1909 Int_t entries=(kglast<kN)?kglast:kN;
1910 Bool_t isOK=kTRUE;
1911 for (Int_t i=0;i<4;i++){
1912 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1913 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1914 vecDeltaN[i] = 0;
1915 if (vecRMS[i]>0.){
1916 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1917 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1918 if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
817766d5 1919 }
3e55050f 1920 }
1921 //
1922 // 2. Apply median+-rms cut
1923 //
1924 if (kglast<10) return; //median and RMS to be defined
1925 if (!isOK) return;
1926 //
1927 // 3. Update alignment
1928 //
5647625c 1929 //Int_t htime = fTime/3600; //time in hours
97d3fea0 1930 Int_t htime = (Int_t)(fTime-fTimeKalmanBin)/fTimeKalmanBin; //time bin
4486a91f 1931 if (fAlignTOFTPC->GetEntriesFast()<htime){
3e55050f 1932 fAlignTOFTPC->Expand(htime*2+20);
1933 }
1934 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1935 if (!align){
1936 // make Alignment object if doesn't exist
1937 align=new AliRelAlignerKalman();
1938 align->SetRunNumber(fRun);
1939 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
97d3fea0 1940 (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
3e55050f 1941 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1942 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1943 align->SetRejectOutliers(kFALSE);
1944 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1945 align->SetMagField(fMagF);
1946 fAlignTOFTPC->AddAt(align,htime);
1947 }
1948 align->AddTrackParams(&pTOF,&pTPC);
0b736a46 1949 Float_t dca[2],cov[3];
1950 track->GetImpactParameters(dca,cov);
1951 if (TMath::Abs(dca[0])<kMaxDy){
1952 FillResHistoTPCTOF(&pTPC,&pTOF);
1953 }
5647625c 1954 //align->SetTimeStamp(fTime);
1955 Double_t averageTime = fTime;
1956 if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1957 averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1958 //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1959 }
1960 align->SetTimeStamp((Int_t)averageTime);
1961
1962 //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1963
3e55050f 1964 align->SetRunNumber(fRun );
1965 //
1966 Int_t nupdates=align->GetNUpdates();
1967 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1968 align->SetRejectOutliers(kFALSE);
1969 TTreeSRedirector *cstream = GetDebugStreamer();
1970 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
3e55050f 1971 TVectorD gpTPC(3), gdTPC(3);
1972 TVectorD gpTOF(3), gdTOF(3);
1973 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1974 pTPC.GetDirection(gdTPC.GetMatrixArray());
1975 pTOF.GetXYZ(gpTOF.GetMatrixArray());
1976 pTOF.GetDirection(gdTOF.GetMatrixArray());
1977 (*cstream)<<"toftpc"<<
1978 "run="<<fRun<< // run number
1979 "event="<<fEvent<< // event number
1980 "time="<<fTime<< // time stamp of event
1981 "trigger="<<fTrigger<< // trigger
1982 "mag="<<fMagF<< // magnetic field
3e55050f 1983 //
1984 "nmed="<<kglast<< // number of entries to define median and RMS
1985 "vMed.="<<&vecMedian<< // median of deltas
1986 "vRMS.="<<&vecRMS<< // rms of deltas
1987 "vDelta.="<<&vecDelta<< // delta in respect to median
1988 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1989 "t.="<<track<< // ful track - find proper cuts
1990 "a.="<<align<< // current alignment
1991 "pTOF.="<<&pTOF<< // track param TOF
1992 "pTPC.="<<&pTPC<< // track param TPC
1993 "gpTPC.="<<&gpTPC<< // global position TPC
1994 "gdTPC.="<<&gdTPC<< // global direction TPC
1995 "gpTOF.="<<&gpTOF<< // global position TOF
1996 "gdTOF.="<<&gdTOF<< // global position TOF
1997 "\n";
817766d5 1998 }
817766d5 1999}
3e55050f 2000
2001
391ffdb2 2002void AliTPCcalibTime::BookDistortionMaps(){
2003 //
2004 // Book ndimensional histograms of distortions/residuals
2005 // Only primary tracks are selected for analysis
2006 //
2007
b9908d0b 2008 Double_t xminTrack[5], xmaxTrack[5];
2009 Int_t binsTrack[5];
2010 TString axisName[5];
2011 TString axisTitle[5];
391ffdb2 2012 //
882f5c06 2013 binsTrack[0] =50;
391ffdb2 2014 axisName[0] ="#Delta";
2015 axisTitle[0] ="#Delta";
2016 //
70036ef3 2017 binsTrack[1] =44;
2018 xminTrack[1] =-1.1; xmaxTrack[1]=1.1;
391ffdb2 2019 axisName[1] ="tanTheta";
2020 axisTitle[1] ="tan(#Theta)";
2021 //
70036ef3 2022 binsTrack[2] =180;
391ffdb2 2023 xminTrack[2] =-TMath::Pi(); xmaxTrack[2]=TMath::Pi();
2024 axisName[2] ="phi";
2025 axisTitle[2] ="#phi";
2026 //
2027 binsTrack[3] =20;
2028 xminTrack[3] =-1.; xmaxTrack[3]=1.; // 0.33 GeV cut
2029 axisName[3] ="snp";
b9908d0b 2030 axisTitle[3] ="snp";
2031 //
2032 binsTrack[4] =10;
2033 xminTrack[4] =120.; xmaxTrack[4]=215.; // crossing radius for CE only
2034 axisName[4] ="r";
2035 axisTitle[4] ="r(cm)";
391ffdb2 2036 //
2037 // delta y
be67055b 2038 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
b9908d0b 2039 fResHistoTPCCE[0] = new THnSparseS("TPCCE#Delta_{Y} (cm)","#Delta_{Y} (cm)", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2040 fResHistoTPCITS[0] = new THnSparseS("TPCITS#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2041 fResHistoTPCvertex[0] = new THnSparseS("TPCVertex#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
be67055b 2042 xminTrack[0] =-1.5; xmaxTrack[0]=1.5; //
391ffdb2 2043 fResHistoTPCTRD[0] = new THnSparseS("TPCTRD#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2044 xminTrack[0] =-5; xmaxTrack[0]=5; //
2045 fResHistoTPCTOF[0] = new THnSparseS("TPCTOF#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2046 //
2047 // delta z
9112627b 2048 xminTrack[0] =-3.; xmaxTrack[0]=3.; //
b9908d0b 2049 fResHistoTPCCE[1] = new THnSparseS("TPCCE#Delta_{Z} (cm)","#Delta_{Z} (cm)", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2050 fResHistoTPCITS[1] = new THnSparseS("TPCITS#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2051 fResHistoTPCvertex[1] = new THnSparseS("TPCVertex#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2052 fResHistoTPCTRD[1] = new THnSparseS("TPCTRD#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2053 xminTrack[0] =-5.; xmaxTrack[0]=5.; //
2054 fResHistoTPCTOF[1] = new THnSparseS("TPCTOF#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2055 //
2056 // delta snp-P2
882f5c06 2057 xminTrack[0] =-0.015; xmaxTrack[0]=0.015; //
b9908d0b 2058 fResHistoTPCCE[2] = new THnSparseS("TPCCE#Delta_{#phi}","#Delta_{#phi}", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2059 fResHistoTPCITS[2] = new THnSparseS("TPCITS#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
9112627b 2060 fResHistoTPCvertex[2] = new THnSparseS("TPCITSv#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2061 fResHistoTPCTRD[2] = new THnSparseS("TPCTRD#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2062 fResHistoTPCTOF[2] = new THnSparseS("TPCTOF#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2063 //
2064 // delta theta-P3
9112627b 2065 xminTrack[0] =-0.025; xmaxTrack[0]=0.025; //
b9908d0b 2066 fResHistoTPCCE[3] = new THnSparseS("TPCCE#Delta_{#theta}","#Delta_{#theta}", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2067 fResHistoTPCITS[3] = new THnSparseS("TPCITS#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
9112627b 2068 fResHistoTPCvertex[3] = new THnSparseS("TPCITSv#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2069 fResHistoTPCTRD[3] = new THnSparseS("TPCTRD#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2070 fResHistoTPCTOF[3] = new THnSparseS("TPCTOF#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2071 //
2072 // delta theta-P4
9112627b 2073 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
b9908d0b 2074 fResHistoTPCCE[4] = new THnSparseS("TPCCE#Delta_{1/pt}","#Delta_{1/pt}", 5, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2075 fResHistoTPCITS[4] = new THnSparseS("TPCITS#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
9112627b 2076 fResHistoTPCvertex[4] = new THnSparseS("TPCITSv#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2077 fResHistoTPCTRD[4] = new THnSparseS("TPCTRD#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
0b736a46 2078 fResHistoTPCTOF[4] = new THnSparseS("TPCTOF#Delta_{1/pt}","#Delta_{1/pt}", 4, binsTrack,xminTrack, xmaxTrack);
391ffdb2 2079 //
2080 for (Int_t ivar=0;ivar<4;ivar++){
b9908d0b 2081 for (Int_t ivar2=0;ivar2<5;ivar2++){
2082 fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2083 fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2084 if (ivar2<4){
2085 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2086 fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2087 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2088 fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2089 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2090 fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2091 }
391ffdb2 2092 }
2093 }
70036ef3 2094 //
2095 // Book vertex: time histograms
2096 //
2097 Int_t binsVertex[2]={500, fTimeBins};
2098 Double_t aminVertex[2]={-5,fTimeStart};
2099 Double_t amaxVertex[2]={5, fTimeEnd};
2100 const char* hnames[12]={"TPCXAside", "TPCXCside","TPCXACdiff","TPCXAPCdiff",
2101 "TPCYAside", "TPCYCside","TPCYACdiff","TPCYAPCdiff",
2102 "TPCZAPCside", "TPCZAMCside","TPCZACdiff","TPCZAPCdiff"};
2103 const char* anames[12]={"x (cm) - A side ", "x (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{x} (cm) - TPC-Common",
2104 "y (cm) - A side ", "y (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{y} (cm) - TPC-Common",
2105 "z (cm)", "#Delta_{Z} (cm) A-C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{Z} (cm) TPC-common"};
2106 for (Int_t ihis=0; ihis<12; ihis++) {
2107 if (ihis>=8) aminVertex[0]=-20.;
2108 if (ihis>=8) amaxVertex[0]=20.;
2109 fTPCVertex[ihis]=new THnSparseF(hnames[ihis],hnames[ihis],2,binsVertex,aminVertex,amaxVertex);
2110 fTPCVertex[ihis]->GetAxis(1)->SetTitle("Time");
2111 fTPCVertex[ihis]->GetAxis(0)->SetTitle(anames[ihis]);
2112 }
b1fea14f 2113
2114 Int_t binsVertexC[2]={40, 300};
2115 Double_t aminVertexC[2]={-20,-30};
2116 Double_t amaxVertexC[2]={20,30};
2117 const char* hnamesC[5]={"TPCA_TPC","TPCC_TPC","TPCA_ITS","TPCC_ITS","TPC_ITS"};
2118 for (Int_t ihis=0; ihis<5; ihis++) {
2119 fTPCVertexCorrelation[ihis]=new THnSparseF(hnamesC[ihis],hnamesC[ihis],2,binsVertexC,aminVertexC,amaxVertexC);
2120 fTPCVertexCorrelation[ihis]->GetAxis(1)->SetTitle("z (cm)");
2121 fTPCVertexCorrelation[ihis]->GetAxis(0)->SetTitle("z (cm)");
2122 }
391ffdb2 2123}
2124
2125
b9908d0b 2126void AliTPCcalibTime::FillResHistoTPCCE(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pTPCOut ){
2127 //
2128 // fill residual histograms pTPCOut-pTPCin - trac crossing CE
2129 // Histogram
2130 //
660f0543 2131 if (fMemoryMode<2) return;
b9908d0b 2132 Double_t histoX[5];
2133 Double_t xyz[3];
2134 pTPCIn->GetXYZ(xyz);
2135 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2136 histoX[1]= pTPCIn->GetTgl();
2137 histoX[2]= phi;
2138 histoX[3]= pTPCIn->GetSnp();
2139 histoX[4]= pTPCIn->GetX();
2140 AliExternalTrackParam lout(*pTPCOut);
2141 lout.Rotate(pTPCIn->GetAlpha());
2142 lout.PropagateTo(pTPCIn->GetX(),fMagF);
2143 //
2144 for (Int_t ihisto=0; ihisto<5; ihisto++){
2145 histoX[0]=lout.GetParameter()[ihisto]-pTPCIn->GetParameter()[ihisto];
2146 fResHistoTPCCE[ihisto]->Fill(histoX);
2147 }
2148}
391ffdb2 2149void AliTPCcalibTime::FillResHistoTPCITS(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pITSOut ){
2150 //
2151 // fill residual histograms pTPCIn-pITSOut
2152 // Histogram is filled only for primary tracks
2153 //
2154 Double_t histoX[4];
2155 Double_t xyz[3];
2156 pTPCIn->GetXYZ(xyz);
2157 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2158 histoX[1]= pTPCIn->GetTgl();
2159 histoX[2]= phi;
2160 histoX[3]= pTPCIn->GetSnp();
9112627b 2161 AliExternalTrackParam lits(*pITSOut);
2162 lits.Rotate(pTPCIn->GetAlpha());
2163 lits.PropagateTo(pTPCIn->GetX(),fMagF);
391ffdb2 2164 //
2165 for (Int_t ihisto=0; ihisto<5; ihisto++){
9112627b 2166 histoX[0]=pTPCIn->GetParameter()[ihisto]-lits.GetParameter()[ihisto];
391ffdb2 2167 fResHistoTPCITS[ihisto]->Fill(histoX);
2168 }
2169}
2170
2171
2172void AliTPCcalibTime::FillResHistoTPC(const AliESDtrack * pTrack){
2173 //
2174 // fill residual histograms pTPC - vertex
2175 // Histogram is filled only for primary tracks
2176 //
660f0543 2177 if (fMemoryMode<2) return;
391ffdb2 2178 Double_t histoX[4];
2179 const AliExternalTrackParam * pTPCIn = pTrack->GetInnerParam();
be67055b 2180 AliExternalTrackParam pTPCvertex(*(pTrack->GetInnerParam()));
9112627b 2181 //
2182 AliExternalTrackParam lits(*pTrack);
be67055b 2183 if (TMath::Abs(pTrack->GetY())>3) return; // beam pipe
2184 pTPCvertex.Rotate(lits.GetAlpha());
2185 //pTPCvertex.PropagateTo(pTPCvertex->GetX(),fMagF);
2186 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,2,kFALSE);
2187 AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,0.1,kFALSE);
391ffdb2 2188 Double_t xyz[3];
2189 pTPCIn->GetXYZ(xyz);
2190 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2191 histoX[1]= pTPCIn->GetTgl();
2192 histoX[2]= phi;
2193 histoX[3]= pTPCIn->GetSnp();
2194 //
2195 Float_t dca[2], cov[3];
2196 pTrack->GetImpactParametersTPC(dca,cov);
882f5c06 2197 for (Int_t ihisto=0; ihisto<5; ihisto++){
be67055b 2198 histoX[0]=pTPCvertex.GetParameter()[ihisto]-lits.GetParameter()[ihisto];
2199 // if (ihisto<2) histoX[0]=dca[ihisto];
9112627b 2200 fResHistoTPCvertex[ihisto]->Fill(histoX);
391ffdb2 2201 }
2202}
2203
2204
2205void AliTPCcalibTime::FillResHistoTPCTRD(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTRDIn ){
2206 //
2207 // fill resuidual histogram TPCout-TRDin
2208 //
660f0543 2209 if (fMemoryMode<2) return;
391ffdb2 2210 Double_t histoX[4];
2211 Double_t xyz[3];
2212 pTPCOut->GetXYZ(xyz);
2213 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2214 histoX[1]= pTPCOut->GetTgl();
2215 histoX[2]= phi;
2216 histoX[3]= pTPCOut->GetSnp();
2217 //
882f5c06 2218 AliExternalTrackParam ltrd(*pTRDIn);
2219 ltrd.Rotate(pTPCOut->GetAlpha());
be67055b 2220 // ltrd.PropagateTo(pTPCOut->GetX(),fMagF);
2221 AliTracker::PropagateTrackToBxByBz(&ltrd,pTPCOut->GetX(),0.1,0.1,kFALSE);
882f5c06 2222
391ffdb2 2223 for (Int_t ihisto=0; ihisto<5; ihisto++){
882f5c06 2224 histoX[0]=pTPCOut->GetParameter()[ihisto]-ltrd.GetParameter()[ihisto];
391ffdb2 2225 fResHistoTPCTRD[ihisto]->Fill(histoX);
2226 }
2227
2228}
0b736a46 2229
2230void AliTPCcalibTime::FillResHistoTPCTOF(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTOFIn ){
2231 //
2232 // fill resuidual histogram TPCout-TOFin
2233 // track propagated to the TOF position
660f0543 2234 if (fMemoryMode<2) return;
0b736a46 2235 Double_t histoX[4];
2236 Double_t xyz[3];
2237
2238 AliExternalTrackParam ltpc(*pTPCOut);
2239 ltpc.Rotate(pTOFIn->GetAlpha());
2240 AliTracker::PropagateTrackToBxByBz(&ltpc,pTOFIn->GetX(),0.1,0.1,kFALSE);
2241 //
2242 ltpc.GetXYZ(xyz);
2243 Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2244 histoX[1]= ltpc.GetTgl();
2245 histoX[2]= phi;
2246 histoX[3]= ltpc.GetSnp();
2247 //
2248 for (Int_t ihisto=0; ihisto<2; ihisto++){
2249 histoX[0]=ltpc.GetParameter()[ihisto]-pTOFIn->GetParameter()[ihisto];
2250 fResHistoTPCTOF[ihisto]->Fill(histoX);
2251 }
2252
2253}