3 /**************************************************************************
4 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
6 * Author: The ALICE Off-line Project. *
7 * Contributors are mentioned in the code where appropriate. *
9 * Permission to use, copy, modify and distribute this software and its *
10 * documentation strictly for non-commercial purposes is hereby granted *
11 * without fee, provided that the above copyright notice appears in all *
12 * copies and that both the copyright notice and this permission notice *
13 * appear in the supporting documentation. The authors make no claims *
14 * about the suitability of this software for any purpose. It is *
15 * provided "as is" without express or implied warranty. *
16 **************************************************************************/
19 Comments to be written here:
20 1. What do we calibrate.
21 2. How to interpret results
23 4. Analysis using debug streamers.
28 // To make cosmic scan the user interaction neccessary
35 #include "Riostream.h"
45 #include "THnSparse.h"
46 #include "TDatabasePDG.h"
48 #include "AliTPCclusterMI.h"
49 #include "AliTPCseed.h"
50 #include "AliESDVertex.h"
51 #include "AliESDEvent.h"
52 #include "AliESDfriend.h"
53 #include "AliESDInputHandler.h"
54 #include "AliAnalysisManager.h"
56 #include "AliTracker.h"
58 #include "AliTPCCalROC.h"
59 #include "AliTPCParam.h"
62 #include "AliTPCcalibCosmic.h"
63 #include "TTreeStream.h"
64 #include "AliTPCTracklet.h"
65 //#include "AliESDcosmic.h"
66 #include "AliRecoParam.h"
67 #include "AliMultiplicity.h"
68 #include "AliTPCTransform.h"
69 #include "AliTPCcalibDB.h"
70 #include "AliTPCseed.h"
71 #include "AliGRPObject.h"
72 #include "AliTPCCorrection.h"
73 ClassImp(AliTPCcalibCosmic)
76 AliTPCcalibCosmic::AliTPCcalibCosmic()
85 fCutMaxD(5), // maximal distance in rfi ditection
86 fCutMaxDz(40), // maximal distance in z ditection
87 fCutTheta(0.03), // maximal distan theta
88 fCutMinDir(-0.99), // direction vector products
89 fCosmicTree(0) // tree with cosmic data
92 // CONSTRUCTOR - SEE COMMENTS ABOVE
94 AliInfo("Default Constructor");
95 for (Int_t ihis=0; ihis<6;ihis++){
99 for (Int_t ihis=0; ihis<4;ihis++){
100 fHistodEdxMax[ihis] =0;
101 fHistodEdxTot[ihis] =0;
106 AliTPCcalibCosmic::AliTPCcalibCosmic(const Text_t *name, const Text_t *title)
115 fCutMaxD(5), // maximal distance in rfi ditection
116 fCutMaxDz(40), // maximal distance in z ditection
117 fCutTheta(0.03), // maximal distan theta
118 fCutMinDir(-0.99), // direction vector products
119 fCosmicTree(0) // tree with cosmic data
122 // cONSTRUCTOR - SEE COMENTS ABOVE
127 fHistNTracks = new TH1F("ntracks","Number of Tracks per Event; number of tracks per event; number of tracks",501,-0.5,500.5);
128 fClusters = new TH1F("signal","Number of Clusters per track; number of clusters per track n_{cl}; counts",160,0,160);
129 fModules = new TH2F("sector","Acorde hits; z (cm); x(cm)",1200,-650,650,600,-700,700);
130 fHistPt = new TH1F("Pt","Pt distribution; p_{T} (GeV); counts",2000,0,50);
131 fDeDx = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",500,0.01,100.,500,2.,1000);
133 fDeDxMIP = new TH1F("DeDxMIP","MIP region; TPC signal (a.u.);counts ",500,2.,1000);
135 AliInfo("Non Default Constructor");
139 AliTPCcalibCosmic::~AliTPCcalibCosmic(){
143 for (Int_t ihis=0; ihis<6;ihis++){
144 delete fHistoDelta[ihis];
145 delete fHistoPull[ihis];
147 for (Int_t ihis=0; ihis<4;ihis++){
148 delete fHistodEdxTot[ihis];
149 delete fHistodEdxMax[ihis];
152 delete fHistNTracks; // histogram showing number of ESD tracks per event
153 delete fClusters; // histogram showing the number of clusters per track
154 delete fModules; // 2d histogram of tracks which are propagated to the ACORDE scintillator array
155 delete fHistPt; // Pt histogram of reconstructed tracks
156 delete fDeDx; // dEdx spectrum showing the different particle types
157 delete fDeDxMIP; // TPC signal close to the MIP region of muons 0.4 < p < 0.45 GeV
161 void AliTPCcalibCosmic::Init(){
164 // Make performance histograms
167 // tracking performance bins
168 // 0 - delta of interest
169 // 1 - min (track0, track1) number of clusters
170 // 2 - R - vertex radius
172 // 4 - P2 - snp(phi) at inner wall of TPC
173 // 5 - P3 - tan(theta) at inner wall of TPC
174 // 6 - P4 - 1/pt mean
177 // 9 - is corss indicator
179 Double_t xminTrack[10], xmaxTrack[10];
181 TString axisName[10];
184 axisName[0] ="#Delta";
187 xminTrack[1] =80; xmaxTrack[1]=160;
188 axisName[1] ="N_{cl}";
191 xminTrack[2] =0; xmaxTrack[2]=90; //
192 axisName[2] ="dca_{r} (cm)";
195 xminTrack[3] =-250; xmaxTrack[3]=250; //
196 axisName[3] ="z (cm)";
199 xminTrack[4] =-0.8; xmaxTrack[4]=0.8; //
200 axisName[4] ="sin(#phi)";
203 xminTrack[5] =-1; xmaxTrack[5]=1; //
204 axisName[5] ="tan(#theta)";
207 xminTrack[6] =-2; xmaxTrack[6]=2; //
208 axisName[6] ="1/pt (1/GeV)";
211 xminTrack[7] =1; xmaxTrack[7]=1000; //
212 axisName[7] ="pt (GeV)";
215 xminTrack[8] =0; xmaxTrack[8]=TMath::Pi(); //
216 axisName[8] ="alpha";
219 xminTrack[9] =-0.1; xmaxTrack[9]=2.1; //
220 axisName[9] ="cross";
223 xminTrack[0] =-1; xmaxTrack[0]=1; //
224 fHistoDelta[0] = new THnSparseS("#Delta_{Y} (cm)","#Delta_{Y} (cm)", ndim, binsTrack,xminTrack, xmaxTrack);
225 xminTrack[0] =-5; xmaxTrack[0]=5; //
226 fHistoPull[0] = new THnSparseS("#Delta_{Y} (unit)","#Delta_{Y} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
229 xminTrack[0] =-1; xmaxTrack[0]=1; //
230 fHistoDelta[1] = new THnSparseS("#Delta_{Z} (cm)","#Delta_{Z} (cm)", ndim, binsTrack,xminTrack, xmaxTrack);
231 xminTrack[0] =-5; xmaxTrack[0]=5; //
232 fHistoPull[1] = new THnSparseS("#Delta_{Z} (unit)","#Delta_{Z} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
235 xminTrack[0] =-10; xmaxTrack[0]=10; //
236 fHistoDelta[2] = new THnSparseS("#Delta_{#phi} (mrad)","#Delta_{#phi} (mrad)", ndim, binsTrack,xminTrack, xmaxTrack);
237 xminTrack[0] =-5; xmaxTrack[0]=5; //
238 fHistoPull[2] = new THnSparseS("#Delta_{#phi} (unit)","#Delta_{#phi} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
241 xminTrack[0] =-10; xmaxTrack[0]=10; //
242 fHistoDelta[3] = new THnSparseS("#Delta_{#theta} (mrad)","#Delta_{#theta} (mrad)", ndim, binsTrack,xminTrack, xmaxTrack);
243 xminTrack[0] =-5; xmaxTrack[0]=5; //
244 fHistoPull[3] = new THnSparseS("#Delta_{#theta} (unit)","#Delta_{#theta} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
247 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
248 fHistoDelta[4] = new THnSparseS("#Delta_{1/pt} (1/GeV)","#Delta_{1/pt} (1/GeV)", ndim, binsTrack,xminTrack, xmaxTrack);
249 xminTrack[0] =-5; xmaxTrack[0]=5; //
250 fHistoPull[4] = new THnSparseS("#Delta_{1/pt} (unit)","#Delta_{1/pt} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
254 xminTrack[0] =-0.5; xmaxTrack[0]=0.5; //
255 fHistoDelta[5] = new THnSparseS("#Delta_{pt}/p_{t}","#Delta_{pt}/p_{t}", ndim, binsTrack,xminTrack, xmaxTrack);
256 xminTrack[0] =-5; xmaxTrack[0]=5; //
257 fHistoPull[5] = new THnSparseS("#Delta_{pt}/p_{t} (unit)","#Delta_{pt}/p_{t} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
260 for (Int_t idedx=0;idedx<4;idedx++){
261 xminTrack[0] =0.5; xmaxTrack[0]=1.5; //
263 xminTrack[1] =10; xmaxTrack[1]=160;
265 fHistodEdxMax[idedx] = new THnSparseS(Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
266 Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
267 ndim, binsTrack,xminTrack, xmaxTrack);
268 fHistodEdxTot[idedx] = new THnSparseS(Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
269 Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
270 ndim, binsTrack,xminTrack, xmaxTrack);
275 for (Int_t ivar=0;ivar<6;ivar++){
276 for (Int_t ivar2=0;ivar2<ndim;ivar2++){
277 fHistoDelta[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
278 fHistoDelta[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
279 fHistoPull[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
280 fHistoPull[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
281 BinLogX(fHistoDelta[ivar],7);
282 BinLogX(fHistoPull[ivar],7);
284 fHistodEdxMax[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
285 fHistodEdxMax[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
286 fHistodEdxTot[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
287 fHistodEdxTot[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
288 BinLogX(fHistodEdxMax[ivar],7);
289 BinLogX(fHistodEdxTot[ivar],7);
296 void AliTPCcalibCosmic::Add(const AliTPCcalibCosmic* cosmic){
298 // merge the content of the cosmic componentnts
300 for (Int_t ivar=0; ivar<6;ivar++){
301 if (fHistoDelta[ivar] && cosmic->fHistoDelta[ivar]){
302 fHistoDelta[ivar]->Add(cosmic->fHistoDelta[ivar]);
304 if (fHistoPull[ivar] && cosmic->fHistoPull[ivar]){
305 fHistoPull[ivar]->Add(cosmic->fHistoPull[ivar]);
308 for (Int_t ivar=0; ivar<4;ivar++){
309 if (fHistodEdxMax[ivar] && cosmic->fHistodEdxMax[ivar]){
310 fHistodEdxMax[ivar]->Add(cosmic->fHistodEdxMax[ivar]);
312 if (fHistodEdxTot[ivar] && cosmic->fHistodEdxTot[ivar]){
313 fHistodEdxTot[ivar]->Add(cosmic->fHistodEdxTot[ivar]);
316 if (cosmic->fCosmicTree){
318 fCosmicTree = new TTree("pairs","pairs");
319 fCosmicTree->SetDirectory(0);
321 AliTPCcalibCosmic::AddTree(fCosmicTree,cosmic->fCosmicTree);
328 void AliTPCcalibCosmic::Process(AliESDEvent *event) {
330 // Process of the ESD event - fill calibration components
333 Printf("ERROR: ESD not available");
336 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
338 Printf("ERROR: esdFriend not available");
344 // COSMIC not signed properly
345 // UInt_t specie = event->GetEventSpecie(); // select only cosmic events
346 //if (specie==AliRecoParam::kCosmic || specie==AliRecoParam::kCalib) {
351 FindCosmicPairs(event);
352 const AliMultiplicity *multiplicity = event->GetMultiplicity();
353 Int_t ntracklets = multiplicity->GetNumberOfTracklets();
354 if (ntracklets>6) return; // filter out "normal" event with high multiplicity
355 const TString &trigger = event->GetFiredTriggerClasses();
356 if (trigger.Contains("C0OB0")==0) return;
359 FindPairs(event); // nearly everything takes place in find pairs...
361 if (GetDebugLevel()>20) printf("Hallo world: Im here and processing an event\n");
362 Int_t ntracks=event->GetNumberOfTracks();
363 fHistNTracks->Fill(ntracks);
368 void AliTPCcalibCosmic::FillHistoPerformance(const AliExternalTrackParam *par0, const AliExternalTrackParam *par1, const AliExternalTrackParam *inner0, const AliExternalTrackParam */*inner1*/, AliTPCseed *seed0, AliTPCseed *seed1, const AliExternalTrackParam *param0Combined , Int_t cross){
370 // par0,par1 - parameter of tracks at DCA 0
371 // inner0,inner1 - parameter of tracks at the TPC entrance
372 // seed0, seed1 - detailed track information
373 // param0Combined - Use combined track parameters for binning
375 Int_t kMinCldEdx =20;
376 Int_t ncl0 = seed0->GetNumberOfClusters();
377 Int_t ncl1 = seed1->GetNumberOfClusters();
378 const Double_t kpullCut = 10;
384 Double_t radius0 = TMath::Sqrt(xyz0[0]*xyz0[0]+xyz0[1]*xyz0[1]);
385 Double_t radius1 = TMath::Sqrt(xyz1[0]*xyz1[0]+xyz1[1]*xyz1[1]);
386 inner0->GetXYZ(xyz0);
387 Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
389 x[1] = TMath::Min(ncl0,ncl1);
390 x[2] = (radius0+radius1)*0.5;
391 x[3] = param0Combined->GetZ();
392 x[4] = inner0->GetSnp();
393 x[5] = param0Combined->GetTgl();
394 x[6] = param0Combined->GetSigned1Pt();
395 x[7] = param0Combined->Pt();
401 delta[0] = (par0->GetY()+par1->GetY());
402 delta[1] = (par0->GetZ()-par1->GetZ());
403 delta[2] = (par0->GetAlpha()-par1->GetAlpha()-TMath::Pi());
404 delta[3] = (par0->GetTgl()+par1->GetTgl());
405 delta[4] = (par0->GetParameter()[4]+par1->GetParameter()[4]);
406 delta[5] = (par0->Pt()-par1->Pt())/((par0->Pt()+par1->Pt())*0.5);
408 sigma[0] = TMath::Sqrt(par0->GetSigmaY2()+par1->GetSigmaY2());
409 sigma[1] = TMath::Sqrt(par0->GetSigmaZ2()+par1->GetSigmaZ2());
410 sigma[2] = TMath::Sqrt(par0->GetSigmaSnp2()+par1->GetSigmaSnp2());
411 sigma[3] = TMath::Sqrt(par0->GetSigmaTgl2()+par1->GetSigmaTgl2());
412 sigma[4] = TMath::Sqrt(par0->GetSigma1Pt2()+par1->GetSigma1Pt2());
413 sigma[5] = sigma[4]*((par0->Pt()+par1->Pt())*0.5);
416 for (Int_t ivar=0;ivar<6;ivar++){
417 if (sigma[ivar]==0) isOK=kFALSE;
418 x[0]= delta[ivar]/sigma[ivar];
419 if (TMath::Abs(x[0])>kpullCut) isOK = kFALSE;
423 if (isOK) for (Int_t ivar=0;ivar<6;ivar++){
424 x[0]= delta[ivar]; // Modifiation 10.10 use not normalized deltas
425 if (ivar==2 || ivar ==3) x[0]*=1000; // angles in radian
426 fHistoDelta[ivar]->Fill(x);
428 x[0]= delta[ivar]/sigma[ivar];
429 fHistoPull[ivar]->Fill(x);
434 // Fill dedx performance
436 for (Int_t ipad=0; ipad<4;ipad++){
442 if (ipad==0) row1=63;
443 if (ipad==1) {row0=63; row1=63+64;}
444 if (ipad==2) {row0=128;}
445 Int_t nclUp = TMath::Nint(seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
446 Int_t nclDown = TMath::Nint(seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
447 Int_t minCl = TMath::Min(nclUp,nclDown);
448 if (minCl<kMinCldEdx) continue;
451 Float_t dEdxTotUp = seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
452 Float_t dEdxTotDown = seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
453 Float_t dEdxMaxUp = seed0->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
454 Float_t dEdxMaxDown = seed1->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
456 if (dEdxTotDown<=0) continue;
457 if (dEdxMaxDown<=0) continue;
458 x[0]=dEdxTotUp/dEdxTotDown;
459 fHistodEdxTot[ipad]->Fill(x);
460 x[0]=dEdxMaxUp/dEdxMaxDown;
461 fHistodEdxMax[ipad]->Fill(x);
468 void AliTPCcalibCosmic::FindPairs(const AliESDEvent *event){
472 // Track0 is choosen in upper TPC part
473 // Track1 is choosen in lower TPC part
475 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
476 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
477 Int_t ntracks=event->GetNumberOfTracks();
478 TObjArray tpcSeeds(ntracks);
479 if (ntracks==0) return;
480 Double_t vtxx[3]={0,0,0};
481 Double_t svtxx[3]={0.000001,0.000001,100.};
482 AliESDVertex vtx(vtxx,svtxx);
486 for (Int_t i=0;i<ntracks;++i) {
487 AliESDtrack *track = event->GetTrack(i);
488 fClusters->Fill(track->GetTPCNcls());
490 const AliExternalTrackParam * trackIn = track->GetInnerParam();
491 const AliExternalTrackParam * trackOut = track->GetOuterParam();
492 if (!trackIn) continue;
493 if (!trackOut) continue;
494 if (ntracks>4 && TMath::Abs(trackIn->GetTgl())<0.0015) continue; // filter laser
497 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
498 if (!friendTrack) continue;
499 TObject *calibObject;
500 AliTPCseed *seed = 0;
501 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
502 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
504 if (seed) tpcSeeds.AddAt(seed,i);
506 Double_t meanP = 0.5*(trackIn->GetP() + trackOut->GetP());
507 if (seed && track->GetTPCNcls() > 80 + 60/(1+TMath::Exp(-meanP+5))) {
508 fDeDx->Fill(meanP, seed->CookdEdxNorm(0.0,0.45,0,0,159));
510 if (meanP > 0.4 && meanP < 0.45) fDeDxMIP->Fill(seed->CookdEdxNorm(0.0,0.45,0,0,159));
512 // if (GetDebugLevel()>0&&meanP>0.2&&seed->CookdEdxNorm(0.0,0.45,0,0,159)>300) {
513 // //TFile *curfile = AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile();
514 // //if (curfile) printf(">>> p+ in file: %s \t event: %i \t Number of ESD tracks: %i \n", curfile->GetName(), (int)event->GetEventNumberInFile(), (int)ntracks);
515 // // if (track->GetOuterParam()->GetAlpha()<0) cout << " Polartiy: " << track->GetSign() << endl;
522 if (ntracks<2) return;
526 for (Int_t i=0;i<ntracks;++i) {
527 AliESDtrack *track0 = event->GetTrack(i);
528 // track0 - choosen upper part
529 if (!track0) continue;
530 if (!track0->GetOuterParam()) continue;
531 if (track0->GetOuterParam()->GetAlpha()<0) continue;
533 track0->GetDirection(dir0);
534 for (Int_t j=0;j<ntracks;++j) {
536 AliESDtrack *track1 = event->GetTrack(j);
538 if (!track1) continue;
539 if (!track1->GetOuterParam()) continue;
540 if (track1->GetOuterParam()->GetAlpha()>0) continue;
543 track1->GetDirection(dir1);
545 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
546 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
547 if (! seed0) continue;
548 if (! seed1) continue;
549 Float_t dedx0 = seed0->CookdEdxNorm(0.05,0.55,0,0,159);
550 Float_t dedx1 = seed1->CookdEdxNorm(0.05,0.55,0,0,159);
552 Float_t dedx0I = seed0->CookdEdxNorm(0.05,0.55,0,0,63);
553 Float_t dedx1I = seed1->CookdEdxNorm(0.05,0.55,0,0,63);
555 Float_t dedx0O = seed0->CookdEdxNorm(0.05,0.55,0,64,159);
556 Float_t dedx1O = seed1->CookdEdxNorm(0.05,0.55,0,64,159);
558 Float_t dir = (dir0[0]*dir1[0] + dir0[1]*dir1[1] + dir0[2]*dir1[2]);
559 Float_t d0 = track0->GetLinearD(0,0);
560 Float_t d1 = track1->GetLinearD(0,0);
562 // conservative cuts - convergence to be guarantied
563 // applying before track propagation
564 if (TMath::Abs(d0+d1)>fCutMaxD) continue; // distance to the 0,0
565 if (dir>fCutMinDir) continue; // direction vector product
566 Float_t bz = AliTracker::GetBz();
567 Float_t dvertex0[2]; //distance to 0,0
568 Float_t dvertex1[2]; //distance to 0,0
569 track0->GetDZ(0,0,0,bz,dvertex0);
570 track1->GetDZ(0,0,0,bz,dvertex1);
571 if (TMath::Abs(dvertex0[1])>250) continue;
572 if (TMath::Abs(dvertex1[1])>250) continue;
576 Float_t dmax = TMath::Max(TMath::Abs(d0),TMath::Abs(d1));
577 AliExternalTrackParam param0(*track0);
578 AliExternalTrackParam param1(*track1);
580 // Propagate using Magnetic field and correct fo material budget
584 Double_t maxsnp=0.90;
585 AliTracker::PropagateTrackToBxByBz(¶m0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE,maxsnp,sign0);
586 AliTracker::PropagateTrackToBxByBz(¶m1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE,maxsnp,sign1);
588 // Propagate rest to the 0,0 DCA - z should be ignored
590 Bool_t b0 = param0.PropagateToDCA(&vtx,bz,1000);
591 Bool_t b1 = param1.PropagateToDCA(&vtx,bz,1000);
593 param0.GetDZ(0,0,0,bz,dvertex0);
594 param1.GetDZ(0,0,0,bz,dvertex1);
595 if (TMath::Abs(param0.GetZ()-param1.GetZ())>fCutMaxDz) continue;
597 Double_t xyz0[3];//,pxyz0[3];
598 Double_t xyz1[3];//,pxyz1[3];
601 Bool_t isPair = IsPair(¶m0,¶m1);
603 if (isPair) FillAcordeHist(track0);
604 if (isPair &¶m0.Pt()>1) {
605 const TString &trigger = event->GetFiredTriggerClasses();
606 UInt_t specie = event->GetEventSpecie();
607 printf("COSMIC ?\t%s\t%d\t%f\t%f\n", trigger.Data(),specie, param0.GetZ(), param1.GetZ());
610 // combined track params
612 AliExternalTrackParam *par0U=MakeCombinedTrack(¶m0,¶m1);
613 AliExternalTrackParam *par1U=MakeCombinedTrack(¶m1,¶m0);
617 TTreeSRedirector * cstream = GetDebugStreamer();
618 //printf("My stream=%p\n",(void*)cstream);
619 AliExternalTrackParam *ip0 = (AliExternalTrackParam *)track0->GetInnerParam();
620 AliExternalTrackParam *ip1 = (AliExternalTrackParam *)track1->GetInnerParam();
621 AliExternalTrackParam *op0 = (AliExternalTrackParam *)track0->GetOuterParam();
622 AliExternalTrackParam *op1 = (AliExternalTrackParam *)track1->GetOuterParam();
623 Bool_t isCrossI = ip0->GetZ()*ip1->GetZ()<0;
624 Bool_t isCrossO = op0->GetZ()*op1->GetZ()<0;
625 Double_t alpha0 = TMath::ATan2(dir0[1],dir0[0]);
626 Double_t alpha1 = TMath::ATan2(dir1[1],dir1[0]);
630 Int_t cross =0; // 0 no cross, 2 cross on both sides
631 if (isCrossI) cross+=1;
632 if (isCrossO) cross+=1;
633 FillHistoPerformance(¶m0, ¶m1, ip0, ip1, seed0, seed1,par0U, cross);
635 (*cstream) << "Track0" <<
636 "run="<<fRun<< // run number
637 "event="<<fEvent<< // event number
638 "time="<<fTime<< // time stamp of event
639 "trigger="<<fTrigger<< // trigger
640 "triggerClass="<<&fTriggerClass<< // trigger
641 "mag="<<fMagF<< // magnetic field
642 "dir="<<dir<< // direction
643 "OK="<<isPair<< // will be accepted
644 "b0="<<b0<< // propagate status
645 "b1="<<b1<< // propagate status
646 "crossI="<<isCrossI<< // cross inner
647 "crossO="<<isCrossO<< // cross outer
649 "Orig0.=" << track0 << // original track 0
650 "Orig1.=" << track1 << // original track 1
651 "Tr0.="<<¶m0<< // track propagated to the DCA 0,0
652 "Tr1.="<<¶m1<< // track propagated to the DCA 0,0
653 "Ip0.="<<ip0<< // inner param - upper
654 "Ip1.="<<ip1<< // inner param - lower
655 "Op0.="<<op0<< // outer param - upper
656 "Op1.="<<op1<< // outer param - lower
657 "Up0.="<<par0U<< // combined track 0
658 "Up1.="<<par1U<< // combined track 1
660 "v00="<<dvertex0[0]<< // distance using kalman
661 "v01="<<dvertex0[1]<< //
662 "v10="<<dvertex1[0]<< //
663 "v11="<<dvertex1[1]<< //
664 "d0="<<d0<< // linear distance to 0,0
665 "d1="<<d1<< // linear distance to 0,0
669 "x00="<<xyz0[0]<< // global position close to vertex
673 "x10="<<xyz1[0]<< // global position close to vertex
679 "dir00="<<dir0[0]<< // direction upper
683 "dir10="<<dir1[0]<< // direction lower
688 "Seed0.=" << seed0 << // original seed 0
689 "Seed1.=" << seed1 << // original seed 1
691 "dedx0="<<dedx0<< // dedx0 - all
692 "dedx1="<<dedx1<< // dedx1 - all
694 "dedx0I="<<dedx0I<< // dedx0 - inner ROC
695 "dedx1I="<<dedx1I<< // dedx1 - inner ROC
697 "dedx0O="<<dedx0O<< // dedx0 - outer ROC
698 "dedx1O="<<dedx1O<< // dedx1 - outer ROC
711 void AliTPCcalibCosmic::FillAcordeHist(AliESDtrack *upperTrack) {
713 // Pt cut to select straight tracks which can be easily propagated to ACORDE which is outside the magnetic field
714 if (upperTrack->Pt() < 10 || upperTrack->GetTPCNcls() < 80) return;
716 const Double_t acordePlane = 850.; // distance of the central Acorde detectors to the beam line at y =0
717 const Double_t roof = 210.5; // distance from x =0 to end of magnet roof
720 upperTrack->GetXYZ(r);
722 upperTrack->GetDirection(d);
724 z = r[2] + (d[2]/d[1])*(acordePlane - r[1]);
725 x = r[0] + (d[0]/d[1])*(acordePlane - r[1]);
728 x = r[0] + (d[0]/(d[0]+d[1]))*(acordePlane+roof-r[0]-r[1]);
729 z = r[2] + (d[2]/(d[0]+d[1]))*(acordePlane+roof-r[0]-r[1]);
732 x = r[0] + (d[0]/(d[1]-d[0]))*(acordePlane+roof+r[0]-r[1]);
733 z = r[2] + (d[2]/(d[1]-d[0]))*(acordePlane+roof+r[0]-r[1]);
736 fModules->Fill(z, x);
742 Long64_t AliTPCcalibCosmic::Merge(TCollection *const li) {
747 TIterator* iter = li->MakeIterator();
748 AliTPCcalibCosmic* cal = 0;
750 while ((cal = (AliTPCcalibCosmic*)iter->Next())) {
751 if (!cal->InheritsFrom(AliTPCcalibCosmic::Class())) {
752 //Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
756 fHistNTracks->Add(cal->GetHistNTracks());
757 fClusters->Add(cal-> GetHistClusters());
758 fModules->Add(cal->GetHistAcorde());
759 fHistPt->Add(cal->GetHistPt());
760 fDeDx->Add(cal->GetHistDeDx());
761 fDeDxMIP->Add(cal->GetHistMIP());
769 Bool_t AliTPCcalibCosmic::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1) const{
773 // 0. Same direction - OPOSITE - cutDir +cutT
774 TCut cutDir("cutDir","dir<-0.99")
776 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
779 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
783 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
786 const Double_t *p0 = tr0->GetParameter();
787 const Double_t *p1 = tr1->GetParameter();
788 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
789 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
790 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
792 Double_t d0[3], d1[3];
793 tr0->GetDirection(d0);
794 tr1->GetDirection(d1);
795 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
802 Double_t AliTPCcalibCosmic::CalculateMIPvalue(TH1F * hist) {
804 // Calculate the MIP value - gaussian fit used
807 TF1 * funcDoubleGaus = new TF1("funcDoubleGaus", "gaus(0)+gaus(3)",0,1000);
808 funcDoubleGaus->SetParameters(hist->GetEntries()*0.75,hist->GetMean()/1.3,hist->GetMean()*0.10,
809 hist->GetEntries()*0.25,hist->GetMean()*1.3,hist->GetMean()*0.10);
810 hist->Fit(funcDoubleGaus);
811 Double_t aMIPvalue = TMath::Min(funcDoubleGaus->GetParameter(1),funcDoubleGaus->GetParameter(4));
813 delete funcDoubleGaus;
822 void AliTPCcalibCosmic::CalculateBetheParams(TH2F */*hist*/, Double_t * /*initialParam*/) {
824 // Not implemented yet
831 void AliTPCcalibCosmic::BinLogX(THnSparse *const h, Int_t axisDim) {
833 // Method for the correct logarithmic binning of histograms
835 TAxis *axis = h->GetAxis(axisDim);
836 int bins = axis->GetNbins();
838 Double_t from = axis->GetXmin();
839 Double_t to = axis->GetXmax();
840 Double_t *newBins = new Double_t[bins + 1];
843 Double_t factor = pow(to/from, 1./bins);
845 for (int i = 1; i <= bins; i++) {
846 newBins[i] = factor * newBins[i-1];
848 axis->Set(bins, newBins);
854 void AliTPCcalibCosmic::BinLogX(TH1 *const h) {
856 // Method for the correct logarithmic binning of histograms
858 TAxis *axis = h->GetXaxis();
859 int bins = axis->GetNbins();
861 Double_t from = axis->GetXmin();
862 Double_t to = axis->GetXmax();
863 Double_t *newBins = new Double_t[bins + 1];
866 Double_t factor = pow(to/from, 1./bins);
868 for (int i = 1; i <= bins; i++) {
869 newBins[i] = factor * newBins[i-1];
871 axis->Set(bins, newBins);
877 AliExternalTrackParam *AliTPCcalibCosmic::MakeTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
879 // Make a atrack using the kalman update of track0 and track1
881 AliExternalTrackParam *par1R= new AliExternalTrackParam(*track1);
882 par1R->Rotate(track0->GetAlpha());
883 par1R->PropagateTo(track0->GetX(),AliTracker::GetBz());
886 Double_t * param = (Double_t*)par1R->GetParameter();
887 Double_t * covar = (Double_t*)par1R->GetCovariance();
895 covar[6] *=-1.; covar[7] *=-1.; covar[8] *=-1.;
896 //covar[10]*=-1.; covar[11]*=-1.; covar[12]*=-1.;
901 AliExternalTrackParam *AliTPCcalibCosmic::MakeCombinedTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
903 // Make combined track
906 AliExternalTrackParam * par1T = MakeTrack(track0,track1);
907 AliExternalTrackParam * par0U = new AliExternalTrackParam(*track0);
909 UpdateTrack(*par0U,*par1T);
915 void AliTPCcalibCosmic::UpdateTrack(AliExternalTrackParam &track1, const AliExternalTrackParam &track2){
917 // Update track 1 with track 2
921 TMatrixD vecXk(5,1); // X vector
922 TMatrixD covXk(5,5); // X covariance
923 TMatrixD matHk(5,5); // vector to mesurement
924 TMatrixD measR(5,5); // measurement error
925 TMatrixD vecZk(5,1); // measurement
927 TMatrixD vecYk(5,1); // Innovation or measurement residual
928 TMatrixD matHkT(5,5);
929 TMatrixD matSk(5,5); // Innovation (or residual) covariance
930 TMatrixD matKk(5,5); // Optimal Kalman gain
931 TMatrixD mat1(5,5); // update covariance matrix
932 TMatrixD covXk2(5,5); //
933 TMatrixD covOut(5,5);
935 Double_t *param1=(Double_t*) track1.GetParameter();
936 Double_t *covar1=(Double_t*) track1.GetCovariance();
937 Double_t *param2=(Double_t*) track2.GetParameter();
938 Double_t *covar2=(Double_t*) track2.GetCovariance();
940 // copy data to the matrix
941 for (Int_t ipar=0; ipar<5; ipar++){
942 for (Int_t jpar=0; jpar<5; jpar++){
943 covXk(ipar,jpar) = covar1[track1.GetIndex(ipar, jpar)];
944 measR(ipar,jpar) = covar2[track2.GetIndex(ipar, jpar)];
948 vecXk(ipar,0) = param1[ipar];
949 vecZk(ipar,0) = param2[ipar];
958 vecYk = vecZk-matHk*vecXk; // Innovation or measurement residual
959 matHkT=matHk.T(); matHk.T();
960 matSk = (matHk*(covXk*matHkT))+measR; // Innovation (or residual) covariance
962 matKk = (covXk*matHkT)*matSk; // Optimal Kalman gain
963 vecXk += matKk*vecYk; // updated vector
964 covXk2 = (mat1-(matKk*matHk));
965 covOut = covXk2*covXk;
969 // copy from matrix to parameters
982 for (Int_t ipar=0; ipar<5; ipar++){
983 param1[ipar]= vecXk(ipar,0) ;
984 for (Int_t jpar=0; jpar<5; jpar++){
985 covar1[track1.GetIndex(ipar, jpar)]=covOut(ipar,jpar);
992 void AliTPCcalibCosmic::FindCosmicPairs(const AliESDEvent * event) {
994 // find cosmic pairs trigger by random trigger
997 AliESDVertex *vertexSPD = (AliESDVertex *)event->GetPrimaryVertexSPD();
998 AliESDVertex *vertexTPC = (AliESDVertex *)event->GetPrimaryVertexTPC();
999 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
1000 const Double_t kMinPt=1;
1001 const Double_t kMinPtMax=0.8;
1002 const Double_t kMinNcl=50;
1003 const Double_t kMaxDelta[5]={2,600,0.02,0.02,0.1};
1004 Int_t ntracks=event->GetNumberOfTracks();
1005 // Float_t dcaTPC[2]={0,0};
1006 // Float_t covTPC[3]={0,0,0};
1008 UInt_t specie = event->GetEventSpecie(); // skip laser events
1009 if (specie==AliRecoParam::kCalib) return;
1013 for (Int_t itrack0=0;itrack0<ntracks;itrack0++) {
1014 AliESDtrack *track0 = event->GetTrack(itrack0);
1015 if (!track0) continue;
1016 if (!track0->IsOn(AliESDtrack::kTPCrefit)) continue;
1018 if (TMath::Abs(AliTracker::GetBz())>1&&track0->Pt()<kMinPt) continue;
1019 if (track0->GetTPCncls()<kMinNcl) continue;
1020 if (TMath::Abs(track0->GetY())<kMaxDelta[0]) continue;
1021 if (track0->GetKinkIndex(0)>0) continue;
1022 const Double_t * par0=track0->GetParameter(); //track param at rhe DCA
1025 //track0->GetImpactParametersTPC(dcaTPC,covTPC);
1026 //if (TMath::Abs(dcaTPC[0])<kMaxDelta[0]) continue;
1027 //if (TMath::Abs(dcaTPC[1])<kMaxDelta[0]*2) continue;
1028 // const AliExternalTrackParam * trackIn0 = track0->GetInnerParam();
1029 for (Int_t itrack1=itrack0+1;itrack1<ntracks;itrack1++) {
1030 AliESDtrack *track1 = event->GetTrack(itrack1);
1031 if (!track1) continue;
1032 if (!track1->IsOn(AliESDtrack::kTPCrefit)) continue;
1033 if (track1->GetKinkIndex(0)>0) continue;
1034 if (TMath::Abs(AliTracker::GetBz())>1&&track1->Pt()<kMinPt) continue;
1035 if (track1->GetTPCncls()<kMinNcl) continue;
1036 if (TMath::Abs(AliTracker::GetBz())>1&&TMath::Max(track1->Pt(), track0->Pt())<kMinPtMax) continue;
1037 if (TMath::Abs(track1->GetY())<kMaxDelta[0]) continue;
1038 //track1->GetImpactParametersTPC(dcaTPC,covTPC);
1039 // if (TMath::Abs(dcaTPC[0])<kMaxDelta[0]) continue;
1040 //if (TMath::Abs(dcaTPC[1])<kMaxDelta[0]*2) continue;
1042 const Double_t* par1=track1->GetParameter(); //track param at rhe DCA
1044 Bool_t isPair=kTRUE;
1045 for (Int_t ipar=0; ipar<5; ipar++){
1046 if (ipar==4&&TMath::Abs(AliTracker::GetBz())<1) continue; // 1/pt not defined for B field off
1047 if (TMath::Abs(TMath::Abs(par0[ipar])-TMath::Abs(par1[ipar]))>kMaxDelta[ipar]) isPair=kFALSE;
1049 if (!isPair) continue;
1050 if (TMath::Abs(TMath::Abs(track0->GetAlpha()-track1->GetAlpha())-TMath::Pi())>kMaxDelta[2]) isPair=kFALSE;
1051 //delta with correct sign
1053 TCut cut0="abs(t1.fP[0]+t0.fP[0])<2"
1054 TCut cut3="abs(t1.fP[3]+t0.fP[3])<0.02"
1055 TCut cut4="abs(t1.fP[4]+t0.fP[4])<0.2"
1057 if (TMath::Abs(par0[0]+par1[0])>kMaxDelta[0]) isPair=kFALSE; //delta y opposite sign
1058 if (TMath::Abs(par0[3]+par1[3])>kMaxDelta[3]) isPair=kFALSE; //delta tgl opposite sign
1059 if (TMath::Abs(AliTracker::GetBz())>1 && TMath::Abs(par0[4]+par1[4])>kMaxDelta[4]) isPair=kFALSE; //delta 1/pt opposite sign
1060 if (!isPair) continue;
1061 // const AliExternalTrackParam * trackIn1 = track1->GetInnerParam();
1064 TTreeSRedirector * pcstream = GetDebugStreamer();
1065 Int_t ntracksSPD = vertexSPD->GetNContributors();
1066 Int_t ntracksTPC = vertexTPC->GetNContributors();
1068 AliESDfriendTrack *friendTrack0 = esdFriend->GetTrack(itrack0);
1069 if (!friendTrack0) continue;
1070 AliESDfriendTrack *friendTrack1 = esdFriend->GetTrack(itrack1);
1071 if (!friendTrack1) continue;
1072 TObject *calibObject;
1073 AliTPCseed *seed0 = 0;
1074 AliTPCseed *seed1 = 0;
1076 for (Int_t l=0;(calibObject=friendTrack0->GetCalibObject(l));++l) {
1077 if ((seed0=dynamic_cast<AliTPCseed*>(calibObject))) break;
1079 for (Int_t l=0;(calibObject=friendTrack1->GetCalibObject(l));++l) {
1080 if ((seed1=dynamic_cast<AliTPCseed*>(calibObject))) break;
1084 (*pcstream)<<"pairs"<<
1085 "run="<<fRun<< // run number
1086 "event="<<fEvent<< // event number
1087 "time="<<fTime<< // time stamp of event
1088 "trigger="<<fTrigger<< // trigger
1089 "triggerClass="<<&fTriggerClass<< // trigger
1090 "mag="<<fMagF<< // magnetic field
1092 "nSPD="<<ntracksSPD<<
1093 "nTPC="<<ntracksTPC<<
1094 "vSPD.="<<vertexSPD<< //primary vertex -SPD
1095 "vTPC.="<<vertexTPC<< //primary vertex -TPC
1096 "t0.="<<track0<< //track0
1097 "t1.="<<track1<< //track1
1098 "ft0.="<<friendTrack0<< //track0
1099 "ft1.="<<friendTrack1<< //track1
1100 "s0.="<<seed0<< //track0
1101 "s1.="<<seed1<< //track1
1105 fCosmicTree = new TTree("pairs","pairs");
1106 fCosmicTree->SetDirectory(0);
1108 if (fCosmicTree->GetEntries()==0){
1110 fCosmicTree->SetDirectory(0);
1111 fCosmicTree->Branch("t0.",&track0);
1112 fCosmicTree->Branch("t1.",&track1);
1113 fCosmicTree->Branch("ft0.",&friendTrack0);
1114 fCosmicTree->Branch("ft1.",&friendTrack1);
1116 fCosmicTree->SetBranchAddress("t0.",&track0);
1117 fCosmicTree->SetBranchAddress("t1.",&track1);
1118 fCosmicTree->SetBranchAddress("ft0.",&friendTrack0);
1119 fCosmicTree->SetBranchAddress("ft1.",&friendTrack1);
1121 fCosmicTree->Fill();
1127 void AliTPCcalibCosmic::Terminate(){
1129 // copy the cosmic tree to memory resident tree
1131 static Int_t counter=0;
1132 printf("AliTPCcalibCosmic::Terminate\t%d\n",counter);
1134 AliTPCcalibBase::Terminate();
1138 void AliTPCcalibCosmic::AddTree(TTree* treeOutput, TTree * treeInput){
1140 // Add the content of tree:
1141 // Notice automatic copy of tree in ROOT does not work for such complicated tree
1143 AliESDtrack *track0=new AliESDtrack;
1144 AliESDtrack *track1=new AliESDtrack;
1145 AliESDfriendTrack *ftrack0=new AliESDfriendTrack;
1146 AliESDfriendTrack *ftrack1=new AliESDfriendTrack;
1147 treeInput->SetBranchAddress("t0.",&track0);
1148 treeInput->SetBranchAddress("t1.",&track1);
1149 treeInput->SetBranchAddress("ft0.",&ftrack0);
1150 treeInput->SetBranchAddress("ft1.",&ftrack1);
1151 if (treeOutput->GetEntries()==0){
1153 treeOutput->SetDirectory(0);
1154 treeOutput->Branch("t0.",&track0);
1155 treeOutput->Branch("t1.",&track1);
1156 treeOutput->Branch("ft0.",&ftrack0);
1157 treeOutput->Branch("ft1.",&ftrack1);
1159 treeOutput->SetBranchAddress("t0.",&track0);
1160 treeOutput->SetBranchAddress("t1.",&track1);
1161 treeOutput->SetBranchAddress("ft0.",&ftrack0);
1162 treeOutput->SetBranchAddress("ft1.",&ftrack1);
1164 Int_t entries= treeInput->GetEntries();
1165 for (Int_t i=0; i<entries; i++){
1166 treeInput->GetEntry(i);
1173 void AliTPCcalibCosmic::MakeFitTree(TTree * treeInput, TTreeSRedirector *pcstream, const TObjArray * corrArray, Int_t step, Int_t run){
1176 // refit the tracks with original points + corrected points for each correction
1178 // treeInput - tree with cosmic tracks
1179 // pcstream - debug output
1182 // Loop over pair of cosmic tracks:
1183 // 1. Find trigger offset between cosmic event and "physic" trigger
1184 // 2. Refit tracks with current transformation
1185 // 3. Refit tracks using additional "primitive" distortion on top
1186 // Best correction estimated as linear combination of corrections
1187 // minimizing the observed distortions
1188 // Observed distortions - matching in the y,z, snp, theta and 1/pt
1190 const Double_t kResetCov=20.;
1191 const Double_t kMaxDelta[5]={1,40,0.03,0.01,0.2};
1192 const Double_t kSigma=2.;
1193 const Double_t kMaxTime=1050;
1194 const Double_t kMaxSnp=0.8;
1195 Int_t ncorr=corrArray->GetEntries();
1196 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
1197 AliTPCParam *param = AliTPCcalibDB::Instance()->GetParameters();
1198 AliGRPObject* grp = AliTPCcalibDB::Instance()->GetGRP(run);
1199 Double_t time=0.5*(grp->GetTimeStart() +grp->GetTimeEnd());
1200 transform->SetCurrentRun(run);
1201 transform->SetCurrentTimeStamp(TMath::Nint(time));
1203 for (Int_t i=0;i<15;i++) covar[i]=0;
1204 covar[0]=kSigma*kSigma;
1205 covar[2]=kSigma*kSigma;
1206 covar[5]=kSigma*kSigma/Float_t(150*150);
1207 covar[9]=kSigma*kSigma/Float_t(150*150);
1209 Double_t *distortions = new Double_t[ncorr+1];
1211 AliESDtrack *track0=new AliESDtrack;
1212 AliESDtrack *track1=new AliESDtrack;
1213 AliESDfriendTrack *ftrack0=new AliESDfriendTrack;
1214 AliESDfriendTrack *ftrack1=new AliESDfriendTrack;
1215 treeInput->SetBranchAddress("t0.",&track0);
1216 treeInput->SetBranchAddress("t1.",&track1);
1217 treeInput->SetBranchAddress("ft0.",&ftrack0);
1218 treeInput->SetBranchAddress("ft1.",&ftrack1);
1219 Int_t entries= treeInput->GetEntries();
1220 for (Int_t i=0; i<entries; i+=step){
1221 treeInput->GetEntry(i);
1222 if (i%20==0) printf("%d\n",i);
1223 if (!ftrack0->GetTPCOut()) continue;
1224 if (!ftrack1->GetTPCOut()) continue;
1225 AliTPCseed *seed0=0;
1226 AliTPCseed *seed1=0;
1227 TObject *calibObject;
1228 for (Int_t l=0;(calibObject=ftrack0->GetCalibObject(l));++l) {
1229 if ((seed0=dynamic_cast<AliTPCseed*>(calibObject))) break;
1231 for (Int_t l=0;(calibObject=ftrack1->GetCalibObject(l));++l) {
1232 if ((seed1=dynamic_cast<AliTPCseed*>(calibObject))) break;
1234 if (!seed0) continue;
1235 if (!seed1) continue;
1236 if (TMath::Abs(seed0->GetSnp())>kMaxSnp) continue;
1237 if (TMath::Abs(seed1->GetSnp())>kMaxSnp) continue;
1240 Int_t nclA0=0, nclC0=0; // number of clusters
1241 Int_t nclA1=0, nclC1=0; // number of clusters
1242 Int_t ncl0=0,ncl1=0;
1243 Double_t rmin0=300, rmax0=-300; // variables to estimate the time0 - trigger offset
1244 Double_t rmin1=300, rmax1=-300;
1245 Double_t tmin0=2000, tmax0=-2000;
1246 Double_t tmin1=2000, tmax1=-2000;
1249 // calculate trigger offset usig "missing clusters"
1250 for (Int_t irow=0; irow<159; irow++){
1251 AliTPCclusterMI *cluster0=seed0->GetClusterPointer(irow);
1252 if (cluster0 &&cluster0->GetX()>10){
1253 if (cluster0->GetX()<rmin0) { rmin0=cluster0->GetX(); tmin0=cluster0->GetTimeBin();}
1254 if (cluster0->GetX()>rmax0) { rmax0=cluster0->GetX(); tmax0=cluster0->GetTimeBin();}
1256 if (cluster0->GetDetector()%36<18) nclA0++;
1257 if (cluster0->GetDetector()%36>=18) nclC0++;
1259 AliTPCclusterMI *cluster1=seed1->GetClusterPointer(irow);
1260 if (cluster1&&cluster1->GetX()>10){
1261 if (cluster1->GetX()<rmin1) { rmin1=cluster1->GetX(); tmin1=cluster1->GetTimeBin();}
1262 if (cluster1->GetX()>rmax1) { rmax1=cluster1->GetX(); tmax1=cluster1->GetTimeBin();}
1264 if (cluster1->GetDetector()%36<18) nclA1++;
1265 if (cluster1->GetDetector()%36>=18) nclC1++;
1268 Int_t cosmicType=0; // types of cosmic topology
1269 if ((nclA0>nclC0) && (nclA1>nclC1)) cosmicType=0; // AA side
1270 if ((nclA0<nclC0) && (nclA1<nclC1)) cosmicType=1; // CC side
1271 if ((nclA0>nclC0) && (nclA1<nclC1)) cosmicType=2; // AC side
1272 if ((nclA0<nclC0) && (nclA1>nclC1)) cosmicType=3; // CA side
1273 //if ((nclA0>nclC0) && (nclA1<nclC1)) cosmicType=6; // AC side out of time
1274 //if ((nclA0>nclC0) && (nclA1<nclC1)) cosmicType=7; // CA side out of time
1276 Double_t deltaTime = 0; // correction for trigger not in time - equalizing the time arival
1277 if ((cosmicType>1)&&TMath::Abs(track1->GetZ()-track0->GetZ())>4){
1279 deltaTime=0.5*(track1->GetZ()-track0->GetZ())/param->GetZWidth();
1280 if (nclA0>nclC0) deltaTime*=-1; // if A side track
1283 TVectorD vectorDT(8);
1284 Int_t crossCounter=0;
1285 Double_t deltaTimeCross = AliTPCcalibCosmic::GetDeltaTime(rmin0, rmax0, rmin1, rmax1, tmin0, tmax0, tmin1, tmax1, TMath::Abs(track0->GetY()),vectorDT);
1286 Bool_t isOKTrigger=kTRUE;
1287 for (Int_t ic=0; ic<6;ic++) {
1288 if (TMath::Abs(vectorDT[ic])>0) {
1289 if (vectorDT[ic]+vectorDT[6]<0) isOKTrigger=kFALSE;
1290 if (vectorDT[ic]+vectorDT[7]>kMaxTime) isOKTrigger=kFALSE;
1296 Double_t deltaTimeCluster=deltaTime;
1297 if ((cosmicType==0 || cosmicType==1) && crossCounter>0){
1298 deltaTimeCluster=deltaTimeCross;
1301 if (nclA0*nclC0>0 || nclA1*nclC1>0) cosmicType+=16; // mixed A side C side - bad for visualization
1303 // Apply current transformation
1306 for (Int_t irow=0; irow<159; irow++){
1307 AliTPCclusterMI *cluster0=seed0->GetClusterPointer(irow);
1308 if (cluster0 &&cluster0->GetX()>10){
1309 Double_t x0[3]={cluster0->GetRow(),cluster0->GetPad(),cluster0->GetTimeBin()+deltaTimeCluster};
1310 Int_t index0[1]={cluster0->GetDetector()};
1311 transform->Transform(x0,index0,0,1);
1312 cluster0->SetX(x0[0]);
1313 cluster0->SetY(x0[1]);
1314 cluster0->SetZ(x0[2]);
1317 AliTPCclusterMI *cluster1=seed1->GetClusterPointer(irow);
1318 if (cluster1&&cluster1->GetX()>10){
1319 Double_t x1[3]={cluster1->GetRow(),cluster1->GetPad(),cluster1->GetTimeBin()+deltaTimeCluster};
1320 Int_t index1[1]={cluster1->GetDetector()};
1321 transform->Transform(x1,index1,0,1);
1322 cluster1->SetX(x1[0]);
1323 cluster1->SetY(x1[1]);
1324 cluster1->SetZ(x1[2]);
1329 Double_t alpha=track0->GetAlpha(); // rotation frame
1330 Double_t cos = TMath::Cos(alpha);
1331 Double_t sin = TMath::Sin(alpha);
1332 Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1333 AliExternalTrackParam btrack0=*(ftrack0->GetTPCOut());
1334 AliExternalTrackParam btrack1=*(ftrack1->GetTPCOut());
1335 btrack0.Rotate(alpha);
1336 btrack1.Rotate(alpha);
1337 // change the sign for track 1
1338 Double_t* par1=(Double_t*)btrack0.GetParameter();
1341 btrack0.AddCovariance(covar);
1342 btrack1.AddCovariance(covar);
1343 btrack0.ResetCovariance(kResetCov);
1344 btrack1.ResetCovariance(kResetCov);
1347 TObjArray tracks0(ncorr+1);
1348 TObjArray tracks1(ncorr+1);
1350 Double_t dEdx0Tot=seed0->CookdEdxAnalytical(0.02,0.6,kTRUE);
1351 Double_t dEdx1Tot=seed1->CookdEdxAnalytical(0.02,0.6,kTRUE);
1352 Double_t dEdx0Max=seed0->CookdEdxAnalytical(0.02,0.6,kFALSE);
1353 Double_t dEdx1Max=seed1->CookdEdxAnalytical(0.02,0.6,kFALSE);
1354 //if (TMath::Abs((dEdx0Max+1)/(dEdx0Tot+1)-1.)>0.1) isOK=kFALSE;
1355 //if (TMath::Abs((dEdx1Max+1)/(dEdx1Tot+1)-1.)>0.1) isOK=kFALSE;
1357 for (Int_t icorr=-1; icorr<ncorr; icorr++){
1358 AliExternalTrackParam rtrack0=btrack0;
1359 AliExternalTrackParam rtrack1=btrack1;
1360 AliTPCCorrection *corr = 0;
1361 if (icorr>=0) corr = (AliTPCCorrection*)corrArray->At(icorr);
1363 for (Int_t irow=159; irow>0; irow--){
1364 AliTPCclusterMI *cluster=seed0->GetClusterPointer(irow);
1365 if (!cluster) continue;
1367 Double_t rD[3]={cluster->GetX(),cluster->GetY(),cluster->GetZ()};
1368 transform->RotatedGlobal2Global(cluster->GetDetector()%36,rD); // transform to global
1369 Float_t r[3]={rD[0],rD[1],rD[2]};
1371 corr->DistortPoint(r, cluster->GetDetector());
1374 Double_t cov[3]={0.01,0.,0.01};
1375 Double_t lx =cos*r[0]+sin*r[1];
1376 Double_t ly =-sin*r[0]+cos*r[1];
1377 rD[1]=ly; rD[0]=lx; rD[2]=r[2]; //transform to track local
1378 if (!AliTracker::PropagateTrackToBxByBz(&rtrack0, lx,mass,1.,kFALSE)) isOKT=kFALSE;;
1379 if (!rtrack0.Update(&rD[1],cov)) isOKT =kFALSE;
1380 if (icorr<0) ncl0++;
1383 for (Int_t irow=159; irow>0; irow--){
1384 AliTPCclusterMI *cluster=seed1->GetClusterPointer(irow);
1385 if (!cluster) continue;
1387 Double_t rD[3]={cluster->GetX(),cluster->GetY(),cluster->GetZ()};
1388 transform->RotatedGlobal2Global(cluster->GetDetector()%36,rD);
1389 Float_t r[3]={rD[0],rD[1],rD[2]};
1391 corr->DistortPoint(r, cluster->GetDetector());
1394 Double_t cov[3]={0.01,0.,0.01};
1395 Double_t lx =cos*r[0]+sin*r[1];
1396 Double_t ly =-sin*r[0]+cos*r[1];
1397 rD[1]=ly; rD[0]=lx; rD[2]=r[2];
1398 if (!AliTracker::PropagateTrackToBxByBz(&rtrack1, lx,mass,1.,kFALSE)) isOKT=kFALSE;
1399 if (!rtrack1.Update(&rD[1],cov)) isOKT=kFALSE;
1400 if (icorr<0) ncl1++;
1402 if (!AliTracker::PropagateTrackToBxByBz(&rtrack0, 0,mass,10.,kFALSE)) isOKT=kFALSE;
1403 if (!AliTracker::PropagateTrackToBxByBz(&rtrack1, 0,mass,10.,kFALSE)) isOKT=kFALSE;
1404 if (!AliTracker::PropagateTrackToBxByBz(&rtrack0, 0,mass,1.,kFALSE)) isOKT=kFALSE;
1405 if (!AliTracker::PropagateTrackToBxByBz(&rtrack1, 0,mass,1.,kFALSE)) isOKT=kFALSE;
1406 const Double_t *param0=rtrack0.GetParameter();
1407 const Double_t *param1=rtrack1.GetParameter();
1408 for (Int_t ipar=0; ipar<4; ipar++){
1409 if (TMath::Abs(param1[ipar]-param0[ipar])>kMaxDelta[ipar]) isOK=kFALSE;
1411 tracks0.AddAt(rtrack0.Clone(), icorr+1);
1412 tracks1.AddAt(rtrack1.Clone(), icorr+1);
1413 AliExternalTrackParam out0=*(ftrack0->GetTPCOut());
1414 AliExternalTrackParam out1=*(ftrack1->GetTPCOut());
1415 Int_t nentries=TMath::Min(ncl0,ncl1);
1418 (*pcstream)<<"cosmic"<<
1419 "isOK="<<isOK<< // correct all propagation update and also residuals
1420 "isOKT="<<isOKT<< // correct all propagation update
1421 "isOKTrigger="<<isOKTrigger<< // correct? estimate of trigger offset
1425 "cross="<<crossCounter<<
1426 "vDT.="<<&vectorDT<<
1428 "dTime="<<deltaTime<< // delta time using the A-c side cross
1429 "dTimeCross="<<deltaTimeCross<< // delta time using missing clusters
1431 "dEdx0Max="<<dEdx0Max<<
1432 "dEdx0Tot="<<dEdx0Tot<<
1433 "dEdx1Max="<<dEdx1Max<<
1434 "dEdx1Tot="<<dEdx1Tot<<
1436 "track0.="<<track0<< // original track 0
1437 "track1.="<<track1<< // original track 1
1438 "out0.="<<&out0<< // outer track 0
1439 "out1.="<<&out1<< // outer track 1
1440 "rtrack0.="<<&rtrack0<< // refitted track with current transform
1441 "rtrack1.="<<&rtrack1<< //
1444 "entries="<<nentries<< // number of clusters
1451 Int_t nentries=TMath::Min(ncl0,ncl1);
1452 for (Int_t ipar=0; ipar<5; ipar++){
1453 for (Int_t icorr=-1; icorr<ncorr; icorr++){
1454 AliTPCCorrection *corr = 0;
1455 if (icorr>=0) corr = (AliTPCCorrection*)corrArray->At(icorr);
1457 AliExternalTrackParam *param0=(AliExternalTrackParam *) tracks0.At(icorr+1);
1458 AliExternalTrackParam *param1=(AliExternalTrackParam *) tracks1.At(icorr+1);
1459 distortions[icorr+1]=param1->GetParameter()[ipar]-param0->GetParameter()[ipar];
1461 distortions[icorr+1]-=distortions[0];
1465 Double_t bz=AliTrackerBase::GetBz();
1467 param0->GetXYZ(gxyz);
1469 Double_t theta=param0->GetParameter()[3];
1470 Double_t phi = alpha;
1471 Double_t snp = track0->GetInnerParam()->GetSnp();
1472 Double_t mean= distortions[0];
1473 Int_t index = param0->GetIndex(ipar,ipar);
1474 Double_t rms=TMath::Sqrt(param1->GetCovariance()[index]+param1->GetCovariance()[index]);
1475 if (crossCounter<1) rms*=1;
1476 Double_t sector=9*phi/TMath::Pi();
1477 Double_t dsec = sector-TMath::Nint(sector+0.5);
1478 Double_t gx=gxyz[0],gy=gxyz[1],gz=gxyz[2];
1479 Double_t refX=TMath::Sqrt(gx*gx+gy*gy);
1481 // Double_t pt=(param0->GetSigned1Pt()+param1->GetSigned1Pt())*0.5;
1482 Double_t pt=(param0->GetSigned1Pt()+param1->GetSigned1Pt())*0.5;
1484 (*pcstream)<<"fit"<< // dump valus for fit
1485 "run="<<run<< //run number
1486 "bz="<<bz<< // magnetic filed used
1487 "dtype="<<dtype<< // detector match type 20
1488 "ptype="<<ipar<< // parameter type
1489 "theta="<<theta<< // theta
1490 "phi="<<phi<< // phi
1491 "snp="<<snp<< // snp
1492 "mean="<<mean<< // mean dist value
1493 "rms="<<rms<< // rms
1497 "refX="<<refX<< // reference radius
1498 "gx="<<gx<< // global position
1499 "gy="<<gy<< // global position
1500 "gz="<<gz<< // global position
1501 "dRrec="<<dRrec<< // delta Radius in reconstruction
1503 "id="<<cosmicType<< //type of the cosmic used
1504 "entries="<<nentries;// number of entries in bin
1505 (*pcstream)<<"cosmicDebug"<<
1506 "p0.="<<param0<< // dump distorted track 0
1507 "p1.="<<param1; // dump distorted track 1
1510 (*pcstream)<<"fit"<<
1511 Form("%s=",corr->GetName())<<distortions[icorr+1]; // dump correction value
1512 (*pcstream)<<"cosmicDebug"<<
1513 Form("%s=",corr->GetName())<<distortions[icorr+1]<< // dump correction value
1514 Form("%sp0.=",corr->GetName())<<param0<< // dump distorted track 0
1515 Form("%sp1.=",corr->GetName())<<param1; // dump distorted track 1
1517 } //loop corrections
1518 (*pcstream)<<"fit"<<"isOK="<<isOK<<"\n";
1519 (*pcstream)<<"cosmicDebug"<<"isOK="<<isOK<<"\n";
1520 } //loop over parameters
1523 delete [] distortions;
1528 Double_t AliTPCcalibCosmic::GetDeltaTime(Double_t rmin0, Double_t rmax0, Double_t rmin1, Double_t rmax1, Double_t tmin0, Double_t tmax0, Double_t tmin1, Double_t tmax1, Double_t dcaR, TVectorD &vectorDT)
1531 // Estimate trigger offset between random cosmic event and "physics" trigger
1532 // Efficiency about 50 % of cases:
1534 // 0. Tracks crossing A side and C side - no match in z - 30 % of cases
1535 // 1. Track only on one side and dissapear at small or at high radiuses - 50 % of cases
1536 // 1.a) rmax<Rc && tmax<Tcmax && tmax>tmin => deltaT=Tcmax-tmax
1537 // 1.b) rmin>Rcmin && tmin<Tcmax && tmin>tmax => deltaT=Tcmax-tmin
1538 // // case the z matching gives proper time
1539 // 1.c) rmax<Rc && tmax>Tcmin && tmax<tmin => deltaT=-tmax
1542 // TCut cutStop = "(min(rmax0,rmax1)<235||abs(rmin0-rmin1)>10)"; // tracks not registered for full time
1547 // 0-3 - occur - wrong correlation
1548 // 1-2 - occur - wrong correlation
1550 // 2-3 - occur - small number of outlyers -20%
1557 const Double_t kMaxRCut=235; // max radius
1558 const Double_t kMinRCut=TMath::Max(dcaR,90.); // min radius
1559 const Double_t kMaxDCut=30; // max distance for minimal radius
1560 const Double_t kMinTime=110;
1561 const Double_t kMaxTime=950;
1564 vectorDT[6]=TMath::Min(TMath::Min(tmin0,tmin1),TMath::Min(tmax0,tmax1));
1565 vectorDT[7]=TMath::Max(TMath::Max(tmin0,tmin1),TMath::Max(tmax0,tmax1));
1566 if (TMath::Min(rmax0,rmax1)<kMaxRCut){
1567 // max cross - deltaT>0
1568 if (rmax0<kMaxRCut && tmax0 <kMaxTime && tmax0>tmin0) vectorDT[0]=kMaxTime-tmax0; // disapear at CE
1569 if (rmax1<kMaxRCut && tmax1 <kMaxTime && tmax1>tmin1) vectorDT[1]=kMaxTime-tmax1; // disapear at CE
1570 // min cross - deltaT<0 - OK they are correlated
1571 if (rmax0<kMaxRCut && tmax0 >kMinTime && tmax0<tmin0) vectorDT[2]=-tmax0; // disapear at ROC
1572 if (rmax1<kMaxRCut && tmax1 >kMinTime && tmax1<tmin1) vectorDT[3]=-tmax1; // disapear at ROC
1574 if (rmin0> kMinRCut+kMaxDCut && tmin0 <kMaxTime && tmin0>tmax0) vectorDT[4]=kMaxTime-tmin0;
1575 if (rmin1> kMinRCut+kMaxDCut && tmin1 <kMaxTime && tmin1>tmax1) vectorDT[5]=kMaxTime-tmin1;
1577 for (Int_t i=0; i<6;i++) {
1578 if (TMath::Abs(vectorDT[i])>0) {
1580 if (vectorDT[i]+vectorDT[6]<0) isOK=kFALSE;
1581 if (vectorDT[i]+vectorDT[7]>kMaxTime) isOK=kFALSE;
1582 if (isOK) deltaT=vectorDT[i];