1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 // Secondary vertexing construction Class
17 // Construct secondary vertex from Beauty hadron with electron and
18 // hadrons, then apply selection criteria
21 // MinJung Kweon <minjung@physi.uni-heidelberg.de>
25 #include <TIterator.h>
26 #include <TParticle.h>
28 #include <AliESDVertex.h>
29 #include <AliESDEvent.h>
30 #include <AliAODEvent.h>
31 #include <AliVTrack.h>
32 #include <AliESDtrack.h>
33 #include <AliAODTrack.h>
34 #include "AliHFEsecVtx.h"
35 #include <AliKFParticle.h>
36 #include <AliKFVertex.h>
38 #include <AliMCEvent.h>
39 #include <AliAODMCParticle.h>
40 #include "AliHFEpairs.h"
41 #include "AliHFEsecVtxs.h"
42 #include "AliHFEtrackFilter.h"
43 #include "AliHFEmcQA.h"
44 #include "AliHFEtools.h"
46 ClassImp(AliHFEsecVtx)
48 //_______________________________________________________________________________________________
49 AliHFEsecVtx::AliHFEsecVtx():
87 // Default constructor
93 //_______________________________________________________________________________________________
94 AliHFEsecVtx::AliHFEsecVtx(const AliHFEsecVtx &p):
101 ,fUseMCPID(p.fUseMCPID)
103 ,fNparents(p.fNparents)
107 ,fNoOfHFEpairs(p.fNoOfHFEpairs)
108 ,fNoOfHFEsecvtxs(p.fNoOfHFEsecvtxs)
109 ,fArethereSecVtx(p.fArethereSecVtx)
118 ,fSignedLxy(p.fSignedLxy)
119 ,fSignedLxy2(p.fSignedLxy2)
121 ,fInvmass(p.fInvmass)
122 ,fInvmassSigma(p.fInvmassSigma)
125 ,fNsectrk2prim(p.fNsectrk2prim)
126 ,fVtxchi2Tightcut(p.fVtxchi2Tightcut)
127 ,fVtxchi2Loosecut(p.fVtxchi2Loosecut)
135 fFilter = new AliHFEtrackFilter(*p.fFilter);
138 //_______________________________________________________________________________________________
140 AliHFEsecVtx::operator=(const AliHFEsecVtx &)
143 // Assignment operator
146 AliInfo("Not yet implemented.");
150 //_______________________________________________________________________________________________
151 AliHFEsecVtx::~AliHFEsecVtx()
157 //cout << "Analysis Done." << endl;
161 //__________________________________________
162 void AliHFEsecVtx::Init()
165 // set pdg code and index
170 fParentSelect[0][0] = 411;
171 fParentSelect[0][1] = 421;
172 fParentSelect[0][2] = 431;
173 fParentSelect[0][3] = 4122;
174 fParentSelect[0][4] = 4132;
175 fParentSelect[0][5] = 4232;
176 fParentSelect[0][6] = 4332;
178 fParentSelect[1][0] = 511;
179 fParentSelect[1][1] = 521;
180 fParentSelect[1][2] = 531;
181 fParentSelect[1][3] = 5122;
182 fParentSelect[1][4] = 5132;
183 fParentSelect[1][5] = 5232;
184 fParentSelect[1][6] = 5332;
186 // momentum ranges to apply pt dependent cuts
195 // momentum dependent DCA cut values (preliminary)
203 fFilter = new AliHFEtrackFilter("SecVtxFilter");
204 fFilter->MakeCutStepRecKineITSTPC();
205 fFilter->MakeCutStepPrimary();
208 Bool_t AliHFEsecVtx::Process(AliVTrack *signalTrack){
212 //if(signalTrack->Pt() < 1.0) return;
215 AliESDtrack *track = dynamic_cast<AliESDtrack *>(signalTrack);
218 AliDebug(1, "no esd track pointer, return\n");
222 Bool_t bTagged = kFALSE;
224 FillHistos(0,track); // wo any cuts
228 AliESDtrack *htrack = 0x0;
230 fFilter->FilterTracks(fESD1);
231 TObjArray *candidates = fFilter->GetFilteredTracks();
232 TIterator *trackIter = candidates->MakeIterator();
233 while((htrack = dynamic_cast<AliESDtrack *>(trackIter->Next()))){
234 if(track->GetID() == htrack->GetID()) continue; // since it is for tagging single electron, don't need additional condition
235 if (htrack->Pt()<1.0) continue;
236 PairAnalysis(track, htrack, htrack->GetID()); // e-h pairing
239 for(int ip=0; ip<HFEpairs()->GetEntriesFast(); ip++){
241 AliHFEpairs *pair = (AliHFEpairs*) (HFEpairs()->UncheckedAt(ip));
242 //if(!(pair->GetPairCode()>1. && pair->GetPairCode()<4.)) // apply various cuts
243 if(pair->GetKFChi2()>5.) // only apply vertex chi2 cut for the moment
244 HFEpairs()->RemoveAt(ip);
247 HFEpairs()->Compress();
248 if(HFEpairs()->GetEntriesFast()) FillHistos(1,track); // after paired
249 if(HFEpairs()->GetEntriesFast()) RunSECVTX(track); // secondary vertexing with e,h1,h2,.. tracks
250 if(HFEsecvtxs()->GetEntriesFast()) FillHistos(2,track); // after secvtxing
251 for(int ip=0; ip<HFEsecvtxs()->GetEntriesFast(); ip++){
252 AliHFEsecVtxs *secvtx=0x0;
253 secvtx = (AliHFEsecVtxs*) (HFEsecvtxs()->UncheckedAt(ip));
254 // here you apply cuts, then if it doesn't pass the cut, remove it from the HFEsecvtxs()
255 if(!(secvtx->GetInvmass()>2.0 && secvtx->GetInvmass()<5.2) || !(secvtx->GetSignedLxy2()>0.08 && secvtx->GetSignedLxy2()<1.5) || !(secvtx->GetKFIP2()>-0.1 && secvtx->GetKFIP2()<0.1))
256 HFEsecvtxs()->RemoveAt(ip);
259 // fill histos for raw spectra
260 if(HFEsecvtxs()->GetEntriesFast()) {
261 FillHistos(3,track); //after secvtx cut
271 //_______________________________________________________________________________________________
272 void AliHFEsecVtx::CreateHistograms(TList * const qaList)
280 fSecVtxList = qaList;
281 fSecVtxList->SetName("SecVtx");
290 fkSourceLabel[kAll]="all";
291 fkSourceLabel[kDirectCharm]="directCharm";
292 fkSourceLabel[kDirectBeauty]="directBeauty";
293 fkSourceLabel[kBeautyCharm]="beauty2charm";
294 fkSourceLabel[kGamma]="gamma";
295 fkSourceLabel[kPi0]="pi0";
296 fkSourceLabel[kElse]="others";
297 fkSourceLabel[kBeautyGamma]="beauty22gamma";
298 fkSourceLabel[kBeautyPi0]="beauty22pi0";
299 fkSourceLabel[kBeautyElse]="beauty22others";
302 TString hnopt="secvtx_";
303 for (Int_t isource = 0; isource < 10; isource++ ){
307 //qaList->AddLast(fSecVtxList);
310 //_______________________________________________________________________________________________
311 void AliHFEsecVtx::GetPrimaryCondition()
314 // get primary characteristics and set
318 AliKFVertex primVtxCopy(*(fESD1->GetPrimaryVertex()));
319 if( primVtxCopy.GetNDF() <1 ) return;
320 fPVx = primVtxCopy.GetX();
321 fPVy = primVtxCopy.GetY();
324 AliKFVertex primVtxCopy(*(fAOD1->GetPrimaryVertex()));
325 if( primVtxCopy.GetNDF() <1 ) return;
326 fPVx = primVtxCopy.GetX();
327 fPVy = primVtxCopy.GetY();
331 //_______________________________________________________________________________________________
332 void AliHFEsecVtx::PairAnalysis(AliVTrack* track1, AliVTrack* track2, Int_t index2)
335 // calculate e-h pair characteristics and tag pair
338 //later consider the below
339 Float_t dca1[2]={-999.,-999.}, dca2[2]={-999.,-999.};
340 Float_t cov1[3]={-999.,-999.,-999.}, cov2[3]={-999.,-999.,-999.};
342 Double_t dca1aod[2]={-999.,-999.}, dca2aod[2]={-999.,-999.};
343 Double_t cov1aod[3]={-999.,-999.,-999.}, cov2aod[3]={-999.,-999.,-999.};
345 if (IsAODanalysis()){
346 const AliAODVertex *primVtx = fAOD1->GetPrimaryVertex();
347 AliESDtrack esdTrk1(track1);
348 AliESDtrack esdTrk2(track2);
349 esdTrk1.PropagateToDCA(primVtx,0.,10000.,dca1aod,cov1aod);
350 esdTrk2.PropagateToDCA(primVtx,0.,10000.,dca2aod,cov2aod);
353 ((AliESDtrack*)track1)->GetImpactParameters(dca1,cov1);
354 ((AliESDtrack*)track2)->GetImpactParameters(dca2,cov2);
357 // apply pt dependent dca cut on hadrons
358 /*for(int ibin=0; ibin<6; ibin++){
359 if((track2->Pt()>fPtRng[ibin] && track2->Pt()<fPtRng[ibin+1]) && TMath::Abs(dca2[0])<fDcaCut[ibin]) return;
362 // get KF particle input pid
363 Int_t pdg1 = GetPDG(track1);
364 Int_t pdg2 = GetPDG(track2);
366 if(pdg1==-1 || pdg2==-1) {
367 //printf("out if considered pid range \n");
371 // create KF particle of pair
372 if(IsAODanalysis()) AliKFParticle::SetField(fAOD1->GetMagneticField());
373 else AliKFParticle::SetField(fESD1->GetMagneticField());
375 AliKFParticle kfTrack[2];
376 kfTrack[0] = AliKFParticle(*track1, pdg1);
377 kfTrack[1] = AliKFParticle(*track2, pdg2);
379 AliKFParticle kfSecondary(kfTrack[0],kfTrack[1]);
381 //secondary vertex point from kf particle
382 Double_t kfx = kfSecondary.GetX();
383 Double_t kfy = kfSecondary.GetY();
384 //Double_t kfz = kfSecondary.GetZ();
386 //momentum at the decay point from kf particle
387 Double_t kfpx = kfSecondary.GetPx();
388 Double_t kfpy = kfSecondary.GetPy();
389 //Double_t kfpz = kfSecondary.GetPz();
392 /* //directly use of ESD vertex
393 const AliESDVertex *pvertex = fESD1->GetPrimaryVertex();
395 pvertex->GetXYZ(xyzVtx);
397 Double_t dx = kfx-xyzVtx[0];
398 Double_t dy = kfy-xyzVtx[1];*/
400 AliKFVertex primVtxCopy(*(fESD1->GetPrimaryVertex()));
401 if( primVtxCopy.GetNDF() <1 ) return;
402 fPVx = primVtxCopy.GetX();
403 fPVy = primVtxCopy.GetY();
405 // printf("esdx= %lf kfx= %lf esdy= %lf kfy= %lf\n",xyzVtx[0],fPVx,xyzVtx[1],fPVy);
407 Double_t dx = kfx-fPVx;
408 Double_t dy = kfy-fPVy;
410 // discriminating variables
411 // invariant mass of the KF particle
412 Double_t invmass = -1;
413 Double_t invmassSigma = -1;
414 kfSecondary.GetMass(invmass,invmassSigma);
416 // chi2 of the KF particle
417 Double_t kfchi2 = -1;
418 if(kfSecondary.GetNDF()>0) kfchi2 = TMath::Sqrt(TMath::Abs(kfSecondary.GetChi2()/kfSecondary.GetNDF()));
420 // opening angle between two particles in XY plane
421 Double_t phi = kfTrack[0].GetAngleXY(kfTrack[1]);
422 Double_t cosphi = TMath::Cos(phi);
424 // DCA from primary to e-h KF particle (impact parameter of KF particle)
425 Double_t vtx[2]={fPVx, fPVy};
426 Double_t kfip = kfSecondary.GetDistanceFromVertexXY(vtx);
428 // projection of kf vertex vector to the kf momentum direction
429 Double_t signedLxy=-999.;
430 if((dx*kfpx+dy*kfpy)>0) signedLxy = TMath::Sqrt(dx*dx+dy*dy);
431 if((dx*kfpx+dy*kfpy)<0) signedLxy = -1*TMath::Sqrt(dx*dx+dy*dy);
432 //[the other way to think about]
433 //Double_t psqr = kfpx*kfpx+kfpy*kfpy;
434 //if(psqr>0) signedLxy=(dx*kfpx+dy*kfpy)/TMath::Sqrt(psqr);
436 //recalculating primary vertex after removing secvtx tracks --------------------------
438 trkid[0] = track1->GetID();
439 trkid[1] = track2->GetID();
441 RecalcPrimvtx(2, trkid, kfTrack);
442 Double_t dx2 = kfx-fPVx2;
443 Double_t dy2 = kfy-fPVy2;
445 // IP of sec particle recalculated based on recalculated primary vertex
446 Double_t vtx2[2]={fPVx2, fPVy2};
447 Double_t kfip2 = kfSecondary.GetDistanceFromVertexXY(vtx2);
448 // signed Lxy recalculated based on recalculated primary vertex
449 Double_t signedLxy2=-999.;
450 if((dx2*kfpx+dy2*kfpy)>0) signedLxy2 = TMath::Sqrt(dx2*dx2+dy2*dy2);
451 if((dx2*kfpx+dy2*kfpy)<0) signedLxy2 = -1*TMath::Sqrt(dx2*dx2+dy2*dy2);
452 //------------------------------------------------------------------------------------
456 if (HasMCData()) paircode = GetPairCode(track1,track2);
459 hfepair.SetTrkLabel(index2);
460 hfepair.SetInvmass(invmass);
461 hfepair.SetKFChi2(kfchi2);
462 hfepair.SetOpenangle(phi);
463 hfepair.SetCosOpenangle(cosphi);
464 hfepair.SetSignedLxy(signedLxy);
465 hfepair.SetSignedLxy2(signedLxy2);
466 hfepair.SetKFIP(kfip);
467 hfepair.SetKFIP2(kfip2);
468 hfepair.SetPairCode(paircode);
469 AddHFEpairToArray(&hfepair);
472 // fill into container for later QA
480 //if(cov1[0]>0) dataE[6]=Double_t(dca1[0]/cov1[0]);
481 //if(cov2[0]>0) dataE[7]=Double_t(dca2[0]/cov2[0]);
482 //dataE[6]=track1->Pt();
483 //dataE[7]=track2->Pt();
484 //dataE[6]=dca1[0]; //mjtmp
485 //dataE[7]=dca2[0]; //mjtmp
486 //dataE[8]=TMath::Abs(dca1[0]);
487 //dataE[9]=TMath::Abs(dca2[0]);
489 fPairQA->Fill(dataE);
493 //_______________________________________________________________________________________________
494 void AliHFEsecVtx::RunSECVTX(AliVTrack *track)
497 // run secondary vertexing algorithm and do tagging
500 //printf("number of considered pairs= %d\n",HFEpairs()->GetEntriesFast());
501 FindSECVTXCandid(track);
504 //_______________________________________________________________________________________________
505 void AliHFEsecVtx::FindSECVTXCandid(AliVTrack *track)
508 // find secondary vertex candidate and store them into container
511 AliVTrack *htrack[20];
512 //Int_t htracklabel[20];
513 //Int_t paircode[20];
514 //Double_t vtxchi2[20];
515 //Double_t dataE[7]={-999.,-999.,-999.,-999.,-1.,0,0};
517 fVtxchi2Tightcut=3.; // tight cut for pair
518 fVtxchi2Loosecut=5.; // loose cut for secvtx
520 if (HFEpairs()->GetEntriesFast()>20){
521 AliDebug(3, "number of paired hadron is over maximum(20)");
525 // get paired track objects
526 AliHFEpairs *pair=0x0;
527 if (IsAODanalysis()){
528 for (int ip=0; ip<HFEpairs()->GetEntriesFast(); ip++){
529 pair = (AliHFEpairs*) (HFEpairs()->UncheckedAt(ip));
530 //htracklabel[ip] = pair->GetTrkLabel();
531 htrack[ip] = fAOD1->GetTrack(pair->GetTrkLabel());
532 //if(pair->GetPairCode()==2 || pair->GetPairCode()==3) paircode[ip]=1;
533 //else paircode[ip]=0;
534 //vtxchi2[ip] = pair->GetKFChi2();
538 for (int ip=0; ip<HFEpairs()->GetEntriesFast(); ip++){
539 pair = (AliHFEpairs*) (HFEpairs()->UncheckedAt(ip));
540 //htracklabel[ip] = pair->GetTrkLabel();
541 htrack[ip] = fESD1->GetTrack(pair->GetTrkLabel());
542 //if(pair->GetPairCode()==2 || pair->GetPairCode()==3) paircode[ip]=1;
543 //else paircode[ip]=0;
544 //vtxchi2[ip] = pair->GetKFChi2();
548 Int_t nPairs = HFEpairs()->GetEntriesFast();
551 // 1 electron candidate + 1 track
553 if (pair->GetKFChi2() < fVtxchi2Tightcut) { // you can also put single track cut -> here you apply very tight cut for the pair
554 Fill2TrkSECVTX(track, pair);
558 //--------------------------------------------------------------
560 // 1 electron candidate + 2 tracks
562 CalcSECVTXProperty(track, htrack[0], htrack[1]); // calculate secondary vertex property
564 if (fKFchi2 < fVtxchi2Loosecut) { // -> here you apply rather loose cut
565 Fill3TrkSECVTX(track, 0, 1);
567 else{ // if doesn't pass the sec vtx chi2 cut
568 for(int jp=0; jp<2; jp++){
569 pair = (AliHFEpairs*) (HFEpairs()->UncheckedAt(jp));
570 if (pair->GetKFChi2() < fVtxchi2Tightcut){
571 Fill2TrkSECVTX(track, pair);
577 //--------------------------------------------------------------
579 // 1 electron candidate + 3 tracks
581 CalcSECVTXProperty(track, htrack[0], htrack[1], htrack[2]); // calculate secondary vertex property
583 if (fKFchi2 < fVtxchi2Loosecut) {
584 Fill4TrkSECVTX(track, 0, 1, 2);
588 for (int i=0; i<nPairs-1; i++){
589 for (int j=i+1; j<nPairs; j++){
590 CalcSECVTXProperty(track, htrack[i], htrack[j]);
591 if (fKFchi2 < fVtxchi2Loosecut) {
593 Fill3TrkSECVTX(track, i, j);
597 if(!fArethereSecVtx){
598 for(int jp=0; jp<nPairs; jp++){
599 pair = (AliHFEpairs*) (HFEpairs()->UncheckedAt(jp));
600 if (pair->GetKFChi2() < fVtxchi2Tightcut){
601 Fill2TrkSECVTX(track, pair);
608 //--------------------------------------------------------------
610 // 1 electron candidate + more than 3 tracks
613 for (int ih1=0; ih1<nPairs-2; ih1++){
614 for (int ih2=ih1+1; ih2<nPairs-1; ih2++){
615 for (int ih3=ih2+1; ih3<nPairs; ih3++){
616 CalcSECVTXProperty(track, htrack[ih1], htrack[ih2], htrack[ih3]); // calculate secondary vertex property
617 if (fKFchi2 < fVtxchi2Loosecut) {
619 Fill4TrkSECVTX(track, ih1, ih2, ih3);
624 if (!fArethereSecVtx){
626 for (int i=0; i<nPairs-1; i++){
627 for (int j=i+1; j<nPairs; j++){
628 CalcSECVTXProperty(track, htrack[i], htrack[j]);
629 if (fKFchi2 < fVtxchi2Loosecut) {
631 Fill3TrkSECVTX(track, i, j);
636 if (!fArethereSecVtx){
637 for(int jp=0; jp<nPairs; jp++){
638 pair = (AliHFEpairs*) (HFEpairs()->UncheckedAt(jp));
639 if (pair->GetKFChi2() < fVtxchi2Tightcut){
640 Fill2TrkSECVTX(track, pair);
646 //--------------------------------------------------------------
650 //_______________________________________________________________________________________________
651 void AliHFEsecVtx::Fill4TrkSECVTX(AliVTrack* track, Int_t ipair, Int_t jpair, Int_t kpair)
654 // fill 3 tracks' secondary vertex properties
657 Double_t dataE[9]={-999.,-999.,-999.,-999.,-1.,0,0,-999.,-999.};
659 Int_t paircode1 = 0, paircode2 = 0, paircode3 = 0;
660 Int_t htracklabel1 = 0, htracklabel2= 0;
663 AliHFEpairs *pair1=0x0;
664 AliHFEpairs *pair2=0x0;
665 AliHFEpairs *pair3=0x0;
666 pair1 = (AliHFEpairs*) (HFEpairs()->UncheckedAt(ipair));
667 pair2 = (AliHFEpairs*) (HFEpairs()->UncheckedAt(jpair));
668 pair3 = (AliHFEpairs*) (HFEpairs()->UncheckedAt(kpair));
670 htracklabel1 = pair1->GetTrkLabel();
671 htracklabel2 = pair2->GetTrkLabel();
673 if (pair1->GetPairCode()==2 || pair1->GetPairCode()==3) paircode1=1;
675 if (pair2->GetPairCode()==2 || pair2->GetPairCode()==3) paircode2=1;
677 if (pair3->GetPairCode()==2 || pair3->GetPairCode()==3) paircode3=1;
681 AliHFEsecVtxs hfesecvtx;
682 hfesecvtx.SetTrkLabel1(htracklabel1); // mj: not much meaningful for the moment
683 hfesecvtx.SetTrkLabel2(htracklabel2); // mj: not much meaningful for the moment
684 if(HasMCData()) hfesecvtx.SetMCCode(GetElectronSource(TMath::Abs(track->GetLabel())));
685 hfesecvtx.SetKFChi2(fKFchi2);
686 hfesecvtx.SetInvmass(fInvmass);
687 hfesecvtx.SetSignedLxy(fSignedLxy);
688 hfesecvtx.SetSignedLxy2(fSignedLxy2);
689 hfesecvtx.SetKFIP(fKFip);
690 hfesecvtx.SetKFIP2(fKFip2);
691 AddHFEsecvtxToArray(&hfesecvtx);
698 if(HasMCData()) dataE[4]=GetElectronSource(TMath::Abs(track->GetLabel()));
699 dataE[5]=4; //# of associated tracks
700 if(paircode1 & paircode2 & paircode3) dataE[6]=1;
701 else if(!paircode1 & !paircode2 & !paircode3) dataE[6]=0;
703 dataE[7]=fSignedLxy2;
704 dataE[8]=track->Pt();
705 fSecvtxQA->Fill(dataE);
708 //_______________________________________________________________________________________________
709 void AliHFEsecVtx::Fill3TrkSECVTX(AliVTrack* track, Int_t ipair, Int_t jpair)
712 // fill 3 tracks' secondary vertex properties
715 Double_t dataE[9]={-999.,-999.,-999.,-999.,-1.,0,0,-999.,-999.};
717 Int_t paircode1 = 0, paircode2 = 0;
718 Int_t htracklabel1 = 0, htracklabel2 = 0;
721 AliHFEpairs *pair1=0x0;
722 AliHFEpairs *pair2=0x0;
723 pair1 = (AliHFEpairs*) (HFEpairs()->UncheckedAt(ipair));
724 pair2 = (AliHFEpairs*) (HFEpairs()->UncheckedAt(jpair));
726 htracklabel1 = pair1->GetTrkLabel();
727 htracklabel2 = pair2->GetTrkLabel();
729 if (pair1->GetPairCode()==2 || pair1->GetPairCode()==3) paircode1=1;
731 if (pair2->GetPairCode()==2 || pair2->GetPairCode()==3) paircode2=1;
735 // fill secondary vertex container
736 AliHFEsecVtxs hfesecvtx;
737 hfesecvtx.SetTrkLabel1(htracklabel1);
738 hfesecvtx.SetTrkLabel2(htracklabel2);
739 if(HasMCData()) hfesecvtx.SetMCCode(GetElectronSource(TMath::Abs(track->GetLabel())));
740 hfesecvtx.SetKFChi2(fKFchi2);
741 hfesecvtx.SetInvmass(fInvmass);
742 hfesecvtx.SetSignedLxy(fSignedLxy);
743 hfesecvtx.SetSignedLxy2(fSignedLxy2);
744 hfesecvtx.SetKFIP(fKFip);
745 hfesecvtx.SetKFIP2(fKFip2);
746 AddHFEsecvtxToArray(&hfesecvtx);
749 // fill debugging THnSparse
754 if(HasMCData()) dataE[4]=GetElectronSource(TMath::Abs(track->GetLabel()));
756 if(paircode1 & paircode2) dataE[6]=1;
757 else if(!paircode1 & !paircode2) dataE[6]=0;
759 dataE[7]=fSignedLxy2;
760 dataE[8]=track->Pt();
761 fSecvtxQA->Fill(dataE);
765 //_______________________________________________________________________________________________
766 void AliHFEsecVtx::Fill2TrkSECVTX(AliVTrack* track, const AliHFEpairs *pair)
769 // fill 2 tracks' secondary vertex properties
772 Double_t dataE[9]={-999.,-999.,-999.,-999.,-1.,0,0,-999.,-999.};
775 if (pair->GetPairCode()==2 || pair->GetPairCode()==3) paircode=1;
778 // fill secondary vertex container
779 AliHFEsecVtxs hfesecvtx;
780 hfesecvtx.SetTrkLabel1(pair->GetTrkLabel());
781 hfesecvtx.SetTrkLabel2(-999);
782 if(HasMCData()) hfesecvtx.SetMCCode(GetElectronSource(TMath::Abs(track->GetLabel())));
783 hfesecvtx.SetInvmass(pair->GetInvmass());
784 hfesecvtx.SetKFChi2(pair->GetKFChi2());
785 hfesecvtx.SetSignedLxy(pair->GetSignedLxy());
786 hfesecvtx.SetSignedLxy2(pair->GetSignedLxy2());
787 hfesecvtx.SetKFIP(pair->GetKFIP());
788 hfesecvtx.SetKFIP2(pair->GetKFIP2());
789 AddHFEsecvtxToArray(&hfesecvtx);
792 // fill debugging THnSparse
793 dataE[0]=pair->GetInvmass();
794 dataE[1]=pair->GetKFChi2();
795 dataE[2]=pair->GetSignedLxy();
796 dataE[3]=pair->GetKFIP();
797 if (HasMCData()) dataE[4]=GetElectronSource(TMath::Abs(track->GetLabel()));
798 dataE[5]=2; //# of associated tracks
800 dataE[7]=pair->GetSignedLxy2();
801 dataE[8]=track->Pt();
802 fSecvtxQA->Fill(dataE);
806 //_______________________________________________________________________________________________
807 void AliHFEsecVtx::CalcSECVTXProperty(AliVTrack* track1, AliVTrack* track2, AliVTrack* track3)
810 // calculate secondary vertex properties
813 // get KF particle input pid
814 Int_t pdg1 = GetPDG(track1);
815 Int_t pdg2 = GetPDG(track2);
816 Int_t pdg3 = GetPDG(track3);
818 if(pdg1==-1 || pdg2==-1 || pdg3==-1) {
819 //printf("out if considered pid range \n");
823 // create KF particle of pair
824 if(IsAODanalysis()) AliKFParticle::SetField(fAOD1->GetMagneticField());
825 else AliKFParticle::SetField(fESD1->GetMagneticField());
826 AliKFParticle kfTrack[3];
827 kfTrack[0] = AliKFParticle(*track1, pdg1);
828 kfTrack[1] = AliKFParticle(*track2, pdg2);
829 kfTrack[2] = AliKFParticle(*track3, pdg3);
831 AliKFParticle kfSecondary(kfTrack[0],kfTrack[1],kfTrack[2]);
832 //AliKFParticle kfSecondary(kfTrack1,kfTrack2,kfTrack3);
834 //secondary vertex point from kf particle
835 Double_t kfx = kfSecondary.GetX();
836 Double_t kfy = kfSecondary.GetY();
837 //Double_t kfz = kfSecondary.GetZ();
839 //momentum at the decay point from kf particle
840 Double_t kfpx = kfSecondary.GetPx();
841 Double_t kfpy = kfSecondary.GetPy();
842 //Double_t kfpz = kfSecondary.GetPz();
844 Double_t dx = kfx-fPVx;
845 Double_t dy = kfy-fPVy;
847 // discriminating variables ----------------------------------------------------------
849 if(kfSecondary.GetNDF()>0) fKFchi2 = TMath::Sqrt(TMath::Abs(kfSecondary.GetChi2()/kfSecondary.GetNDF()));
851 // invariant mass of the KF particle
852 kfSecondary.GetMass(fInvmass,fInvmassSigma);
854 // DCA from primary to e-h KF particle (impact parameter of KF particle)
855 Double_t vtx[2]={fPVx, fPVy};
856 fKFip = kfSecondary.GetDistanceFromVertexXY(vtx);
858 if((dx*kfpx+dy*kfpy)>0) fSignedLxy= TMath::Sqrt(dx*dx+dy*dy);
859 if((dx*kfpx+dy*kfpy)<0) fSignedLxy= -1*TMath::Sqrt(dx*dx+dy*dy);
860 //[the other way to think about] - projection of kf vertex vector to the kf momentum direction
861 //Double_t psqr = kfpx*kfpx+kfpy*kfpy;
862 //if(psqr>0) fSignedLxy=(dx*kfpx+dy*kfpy)/TMath::Sqrt(psqr);
865 //recalculating primary vertex after removing secvtx tracks --------------------------
867 trkid[0] = track1->GetID();
868 trkid[1] = track2->GetID();
869 trkid[2] = track3->GetID();
871 RecalcPrimvtx(3, trkid, kfTrack);
872 Double_t dx2 = kfx-fPVx2;
873 Double_t dy2 = kfy-fPVy2;
875 // IP of sec particle recalculated based on recalculated primary vertex
876 Double_t vtx2[2]={fPVx2, fPVy2};
877 fKFip2 = kfSecondary.GetDistanceFromVertexXY(vtx2);
878 // signed Lxy recalculated based on recalculated primary vertex
879 if((dx2*kfpx+dy2*kfpy)>0) fSignedLxy2= TMath::Sqrt(dx2*dx2+dy2*dy2);
880 if((dx2*kfpx+dy2*kfpy)<0) fSignedLxy2= -1*TMath::Sqrt(dx2*dx2+dy2*dy2);
881 //------------------------------------------------------------------------------------
885 //_______________________________________________________________________________________________
886 void AliHFEsecVtx::CalcSECVTXProperty(AliVTrack* track1, AliVTrack* track2, AliVTrack* track3, AliVTrack* track4)
889 // calculate secondary vertex properties
892 // get KF particle input pid
893 Int_t pdg1 = GetPDG(track1);
894 Int_t pdg2 = GetPDG(track2);
895 Int_t pdg3 = GetPDG(track3);
896 Int_t pdg4 = GetPDG(track4);
898 if(pdg1==-1 || pdg2==-1 || pdg3==-1 || pdg4==-1) {
899 //printf("out if considered pid range \n");
903 // create KF particle of pair
904 if(IsAODanalysis()) AliKFParticle::SetField(fAOD1->GetMagneticField());
905 else AliKFParticle::SetField(fESD1->GetMagneticField());
907 AliKFParticle kfTrack[4];
908 kfTrack[0] = AliKFParticle(*track1, pdg1);
909 kfTrack[1] = AliKFParticle(*track2, pdg2);
910 kfTrack[2] = AliKFParticle(*track3, pdg3);
911 kfTrack[3] = AliKFParticle(*track4, pdg4);
913 AliKFParticle kfSecondary(kfTrack[0],kfTrack[1],kfTrack[2],kfTrack[3]);
915 //secondary vertex point from kf particle
916 Double_t kfx = kfSecondary.GetX();
917 Double_t kfy = kfSecondary.GetY();
918 //Double_t kfz = kfSecondary.GetZ();
920 //momentum at the decay point from kf particle
921 Double_t kfpx = kfSecondary.GetPx();
922 Double_t kfpy = kfSecondary.GetPy();
923 //Double_t kfpz = kfSecondary.GetPz();
925 Double_t dx = kfx-fPVx;
926 Double_t dy = kfy-fPVy;
928 // discriminating variables ----------------------------------------------------------
930 if(kfSecondary.GetNDF()>0) fKFchi2 = TMath::Sqrt(TMath::Abs(kfSecondary.GetChi2()/kfSecondary.GetNDF()));
932 // invariant mass of the KF particle
933 kfSecondary.GetMass(fInvmass,fInvmassSigma);
935 // DCA from primary to e-h KF particle (impact parameter of KF particle)
936 Double_t vtx[2]={fPVx, fPVy};
937 fKFip = kfSecondary.GetDistanceFromVertexXY(vtx);
939 if((dx*kfpx+dy*kfpy)>0) fSignedLxy= TMath::Sqrt(dx*dx+dy*dy);
940 if((dx*kfpx+dy*kfpy)<0) fSignedLxy= -1*TMath::Sqrt(dx*dx+dy*dy);
941 //[the other way to think about] - projection of kf vertex vector to the kf momentum direction
942 //Double_t psqr = kfpx*kfpx+kfpy*kfpy;
943 //if(psqr>0) fSignedLxy=(dx*kfpx+dy*kfpy)/TMath::Sqrt(psqr);
945 //recalculating primary vertex after removing secvtx tracks --------------------------
947 trkid[0] = track1->GetID();
948 trkid[1] = track2->GetID();
949 trkid[2] = track3->GetID();
950 trkid[3] = track4->GetID();
952 RecalcPrimvtx(4, trkid, kfTrack);
953 Double_t dx2 = kfx-fPVx2;
954 Double_t dy2 = kfy-fPVy2;
956 // IP of sec particle recalculated based on recalculated primary vertex
957 Double_t vtx2[2]={fPVx2, fPVy2};
958 fKFip2 = kfSecondary.GetDistanceFromVertexXY(vtx2);
959 // signed Lxy recalculated based on recalculated primary vertex
960 if((dx2*kfpx+dy2*kfpy)>0) fSignedLxy2= TMath::Sqrt(dx2*dx2+dy2*dy2);
961 if((dx2*kfpx+dy2*kfpy)<0) fSignedLxy2= -1*TMath::Sqrt(dx2*dx2+dy2*dy2);
962 //------------------------------------------------------------------------------------
966 //_______________________________________________________________________________________________
967 void AliHFEsecVtx::RecalcPrimvtx(Int_t nkftrk, const Int_t * const trkid, const AliKFParticle * const kftrk){
969 // reccalculate primary vertex after removing considering track in the calculation
972 const AliESDVertex *primvtx = fESD1->GetPrimaryVertex();
974 AliKFVertex kfESDprimary;
975 Int_t n = primvtx->GetNIndices();
980 if (n>0 && primvtx->GetStatus()){
981 kfESDprimary = AliKFVertex(*primvtx);
982 UShort_t *priIndex = primvtx->GetIndices();
983 for(Int_t j=0; j<nkftrk; j++){
984 for (Int_t i=0;i<n;i++){
985 Int_t idx = Int_t(priIndex[i]);
986 if (idx == trkid[j]){
987 kfESDprimary -= kftrk[j];
994 fPVx2 = kfESDprimary.GetX();
995 fPVy2 = kfESDprimary.GetY();
999 //_______________________________________________________________________________________________
1000 Int_t AliHFEsecVtx::GetMCPID(const AliESDtrack *track)
1006 AliMCParticle *mctrack = NULL;
1007 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(track->GetLabel()))))) return 0;
1008 TParticle *mcpart = mctrack->Particle();
1010 if ( !mcpart ) return 0;
1011 Int_t pdgCode = mcpart->GetPdgCode();
1016 //_______________________________________________________________________________________________
1017 Int_t AliHFEsecVtx::GetPairOriginESD(AliESDtrack* trk1, AliESDtrack* trk2)
1020 // return pdg code of the origin(source) of the pair
1023 // ---*---*---*-----ancester A----- track1
1026 // => if they originated from same ancester,
1027 // then return "the absolute value of pdg code of ancester A"
1029 // ---*---*---B hadron-----ancester A----- track1
1032 // => if they originated from same ancester, and this ancester originally comes from B hadrons
1033 // then return -1*"the absolute value of pdg code of ancester A"
1035 // caution : it can also return parton pdg code if it originated from same string or gluon spliting.
1038 if (trk1->GetLabel()<0 || trk2->GetLabel()<0) return 0;
1040 AliMCParticle *mctrack = NULL;
1041 AliMCParticle *mctrack1 = NULL;
1042 AliMCParticle *mctrack2 = NULL;
1043 if(!(mctrack1 = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(trk1->GetLabel()))))) return 0;
1044 if(!(mctrack2 = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(trk2->GetLabel()))))) return 0;
1045 TParticle *part1 = mctrack1->Particle();
1046 TParticle *part2 = mctrack2->Particle();
1048 TParticle* part2cp = part2;
1049 if (!(part1) || !(part2)) return 0;
1053 //if the two tracks' mother's label is same, get the mother info
1054 //in case of charm, check if it originated from beauty
1055 for (Int_t i=0; i<10; i++){ //check up to 10 ancesters
1056 Int_t label1 = part1->GetFirstMother();
1057 if (label1 < 0) return 0;
1059 for (Int_t j=0; j<10; j++){ //check up to 10 ancesters
1060 Int_t label2 = part2->GetFirstMother();
1061 if (label2 < 0) break;
1063 if (label1 == label2){ //check if two tracks are originated from same mother
1064 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(label2))))) return 0;
1065 TParticle* commonmom = mctrack->Particle();
1067 srcpdg = abs(commonmom->GetPdgCode());
1069 //check ancester to see if it is originally from beauty
1070 for (Int_t k=0; k<10; k++){ //check up to 10 ancesters
1071 Int_t ancesterlabel = commonmom->GetFirstMother();
1072 if (ancesterlabel < 0) return srcpdg; // if there is no more commonancester, return commonmom's pdg
1074 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(ancesterlabel))))) return 0;
1075 TParticle* commonancester = mctrack->Particle();
1077 Int_t ancesterpdg = abs(commonancester->GetPdgCode());
1079 for (Int_t l=0; l<fNparents; l++){
1080 if (abs(ancesterpdg)==fParentSelect[1][l]){
1081 srcpdg = -1*srcpdg; //multiply -1 for hadron from bottom
1085 commonmom = commonancester;
1088 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(label2))))) return 0;
1089 part2 = mctrack->Particle(); //if their mother is different, go to earlier generation of 2nd particle
1091 if (!(part2)) break;
1093 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(label1))))) return 0;
1094 part1 = mctrack->Particle(); //if their mother is different, go to earlier generation of 1st particle
1096 if (!(part1)) return 0;
1102 //_______________________________________________________________________________________________
1103 Int_t AliHFEsecVtx::GetPairOriginAOD(AliAODTrack* trk1, AliAODTrack* trk2)
1107 // return pdg code of the origin(source) of the pair
1110 // ---*---*---*-----ancester A----- track1
1113 // => if they originated from same ancester,
1114 // then return "the absolute value of pdg code of ancester A"
1116 // ---*---*---B hadron-----ancester A----- track1
1119 // => if they originated from same ancester, and this ancester originally comes from B hadrons
1120 // then return -1*"the absolute value of pdg code of ancester A"
1122 // caution : it can also return parton pdg code if it originated from same string or gluon spliting.
1125 if (trk1->GetLabel()<0 || trk2->GetLabel()<0) return 0;
1126 AliAODMCParticle *part1 = (AliAODMCParticle*)fMCArray->At(trk1->GetLabel());
1127 AliAODMCParticle *part2 = (AliAODMCParticle*)fMCArray->At(trk2->GetLabel());
1128 AliAODMCParticle *part2cp = part2;
1129 if (!(part1) || !(part2)) return 0;
1133 //if the two tracks' mother's label is same, get the mother info
1134 //in case of charm, check if it originated from beauty
1135 for (Int_t i=0; i<10; i++){ //check up to 10 ancesters
1136 Int_t label1 = part1->GetMother();
1137 if (label1 < 0) return 0;
1139 for (Int_t j=0; j<10; j++){ //check up to 10 ancesters
1140 Int_t label2 = part2->GetMother();
1141 if (label2 < 0) return 0;
1143 if (label1 == label2){ //check if two tracks are originated from same mother
1144 AliAODMCParticle *commonmom = (AliAODMCParticle*)fMCArray->At(label1);
1145 srcpdg = abs(commonmom->GetPdgCode());
1147 //check ancester to see if it is originally from beauty
1148 for (Int_t k=0; k<10; k++){ //check up to 10 ancesters
1149 Int_t ancesterlabel = commonmom->GetMother();
1150 if (ancesterlabel < 0) return srcpdg; // if there is no more commonancester, return commonmom's pdg
1152 AliAODMCParticle *commonancester = (AliAODMCParticle*)fMCArray->At(ancesterlabel);
1153 Int_t ancesterpdg = abs(commonancester->GetPdgCode());
1155 for (Int_t l=0; l<fNparents; l++){
1156 if (abs(ancesterpdg)==fParentSelect[1][l]){
1157 srcpdg = -1*srcpdg; //multiply -1 for charm from bottom
1161 commonmom = commonancester;
1164 part2 = (AliAODMCParticle*)fMCArray->At(label2); //if their mother is different, go to earlier generation of 2nd particle
1165 if (!(part2)) break;
1167 part1 = (AliAODMCParticle*)fMCArray->At(label1); //if their mother is different, go to earlier generation of 2nd particle
1169 if (!(part1)) return 0;
1175 //_______________________________________________________________________________________________
1176 Int_t AliHFEsecVtx::GetPairCode(const AliVTrack* const trk1, const AliVTrack* const trk2)
1179 // return pair code which is predefinded as:
1180 // kDirectCharm, kDirectBeauty, kBeautyCharm, kGamma, kPi0,
1181 // kElse, kBeautyGamma, kBeautyPi0, kBeautyElse
1185 Int_t srccode = kElse;
1187 if(IsAODanalysis()) srcpdg = GetPairOriginAOD((AliAODTrack*)trk1,(AliAODTrack*)trk2);
1188 else srcpdg = GetPairOriginESD((AliESDtrack*)trk1,(AliESDtrack*)trk2);
1190 if (srcpdg < 0) srccode = kBeautyElse;
1191 for (Int_t i=0; i<fNparents; i++){
1192 if (abs(srcpdg)==fParentSelect[0][i]) {
1193 if (srcpdg>0) srccode = kDirectCharm;
1194 if (srcpdg<0) srccode = kBeautyCharm;
1196 if (abs(srcpdg)==fParentSelect[1][i]) {
1197 if (srcpdg>0) srccode = kDirectBeauty;
1198 if (srcpdg<0) return kElse;
1201 if (srcpdg == 22) srccode = kGamma;
1202 if (srcpdg == -22) srccode = kBeautyGamma;
1203 if (srcpdg == 111) srccode = kPi0;
1204 if (srcpdg == -111) srccode = kBeautyPi0;
1209 //_______________________________________________________________________________________________
1210 Int_t AliHFEsecVtx::GetElectronSource(Int_t iTrack)
1213 // return decay electron's origin
1217 AliDebug(1, "Stack label is negative, return\n");
1221 AliMCParticle *mctrack = NULL;
1222 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(iTrack))))) return -1;
1223 TParticle *mcpart = mctrack->Particle();
1226 AliDebug(1, "no mc particle, return\n");
1230 if ( abs(mcpart->GetPdgCode()) != 11 ) return kMisID;
1232 // if ( abs(mcpart->GetPdgCode()) != 11 ) return -1; // check if it is electron !
1234 Int_t iLabel = mcpart->GetFirstMother();
1236 AliDebug(1, "Stack label is negative, return\n");
1241 Bool_t isFinalOpenCharm = kFALSE;
1243 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(iLabel))))) return -1;
1244 TParticle *partMother = mctrack->Particle();
1246 Int_t maPdgcode = partMother->GetPdgCode();
1248 // if the mother is charmed hadron
1249 if ( int(abs(maPdgcode)/100.) == kCharm || int(abs(maPdgcode)/1000.) == kCharm ) {
1251 for (Int_t i=0; i<fNparents; i++){
1252 if (abs(maPdgcode)==fParentSelect[0][i]){
1253 isFinalOpenCharm = kTRUE;
1256 if (!isFinalOpenCharm) return -1;
1258 // iterate until you find B hadron as a mother or become top ancester
1259 for (Int_t i=1; i<100; i++){ // check back to the 100 generation older
1261 Int_t jLabel = partMother->GetFirstMother();
1263 origin = kDirectCharm;
1266 if (jLabel < 0){ // safety protection
1267 AliDebug(1, "Stack label is negative, return\n");
1271 // if there is an ancester
1272 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(jLabel))))) return -1;
1273 TParticle *grandMa = mctrack->Particle();
1275 Int_t grandMaPDG = grandMa->GetPdgCode();
1277 for (Int_t j=0; j<fNparents; j++){
1278 if (abs(grandMaPDG)==fParentSelect[1][j]){
1279 origin = kBeautyCharm;
1284 partMother = grandMa;
1285 } // end of iteration
1287 else if ( int(abs(maPdgcode)/100.) == kBeauty || int(abs(maPdgcode)/1000.) == kBeauty ) {
1288 for (Int_t i=0; i<fNparents; i++){
1289 if (abs(maPdgcode)==fParentSelect[1][i]){
1290 origin = kDirectBeauty;
1296 //============ gamma ================
1297 else if ( abs(maPdgcode) == 22 ) {
1300 // iterate until you find B hadron as a mother or become top ancester
1301 for (Int_t i=1; i<100; i++){ // check back to the 100 generation older
1303 Int_t jLabel = partMother->GetFirstMother();
1308 if (jLabel < 0){ // safety protection
1309 AliDebug(1, "Stack label is negative, return\n");
1313 // if there is an ancester
1314 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(jLabel))))) return -1;
1315 TParticle *grandMa = mctrack->Particle();
1317 Int_t grandMaPDG = grandMa->GetPdgCode();
1319 for (Int_t j=0; j<fNparents; j++){
1320 if (abs(grandMaPDG)==fParentSelect[1][j]){
1321 origin = kBeautyGamma;
1326 partMother = grandMa;
1327 } // end of iteration
1332 //============ pi0 ================
1333 else if ( abs(maPdgcode) == 111 ) {
1336 // iterate until you find B hadron as a mother or become top ancester
1337 for (Int_t i=1; i<100; i++){ // check back to the 100 generation older
1339 Int_t jLabel = partMother->GetFirstMother();
1344 if (jLabel < 0){ // safety protection
1345 AliDebug(1, "Stack label is negative, return\n");
1349 // if there is an ancester
1350 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(jLabel))))) return -1;
1351 TParticle *grandMa = mctrack->Particle();
1353 Int_t grandMaPDG = grandMa->GetPdgCode();
1355 for (Int_t j=0; j<fNparents; j++){
1356 if (abs(grandMaPDG)==fParentSelect[1][j]){
1357 origin = kBeautyPi0;
1362 partMother = grandMa;
1363 } // end of iteration
1370 // iterate until you find B hadron as a mother or become top ancester
1371 for (Int_t i=1; i<100; i++){ // check back to the 100 generation older
1373 Int_t jLabel = partMother->GetFirstMother();
1378 if (jLabel < 0){ // safety protection
1379 AliDebug(1, "Stack label is negative, return\n");
1383 // if there is an ancester
1384 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(jLabel))))) return -1;
1385 TParticle *grandMa = mctrack->Particle();
1387 Int_t grandMaPDG = grandMa->GetPdgCode();
1389 for (Int_t j=0; j<fNparents; j++){
1390 if (abs(grandMaPDG)==fParentSelect[1][j]){
1391 origin = kBeautyElse;
1396 partMother = grandMa;
1397 } // end of iteration
1403 //_______________________________________________________________________________________________
1404 Int_t AliHFEsecVtx::GetPDG(const AliVTrack *track)
1407 // get KF particle input pdg for mass hypothesis
1412 if (fUseMCPID && HasMCData()){
1413 pdgCode = GetMCPDG(track);
1418 GetESDPID((AliESDtrack*)track, pid, prob);
1420 case 0: pdgCode = 11; break;
1421 case 1: pdgCode = 13; break;
1422 case 2: pdgCode = 211; break;
1423 case 3: pdgCode = 321; break;
1424 case 4: pdgCode = 2212; break;
1425 default: pdgCode = -1;
1429 Int_t pid = ((AliAODTrack*)track)->GetMostProbablePID();
1431 case 0: pdgCode = 11; break;
1432 case 1: pdgCode = 13; break;
1433 case 2: pdgCode = 211; break;
1434 case 3: pdgCode = 321; break;
1435 case 4: pdgCode = 2212; break;
1436 default: pdgCode = -1;
1443 //_______________________________________________________________________________________________
1444 Int_t AliHFEsecVtx::GetMCPDG(const AliVTrack *track)
1447 // return mc pdg code
1450 Int_t label = TMath::Abs(track->GetLabel());
1452 AliMCParticle *mctrack = NULL;
1454 if (IsAODanalysis()) {
1455 AliAODMCParticle *mcpart = (AliAODMCParticle*)fMCArray->At(label);
1456 if ( !mcpart ) return 0;
1457 pdgCode = mcpart->GetPdgCode();
1460 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(label))))) return 0;
1461 TParticle *mcpart = mctrack->Particle();
1463 if ( !mcpart ) return 0;
1464 pdgCode = mcpart->GetPdgCode();
1470 //______________________________________________________________________________
1471 void AliHFEsecVtx::GetESDPID(const AliESDtrack *track, Int_t &recpid, Double_t &recprob)
1474 // calculate likehood for esd pid
1483 Double_t probability[5];
1485 // get probability of the diffenrent particle types
1486 track->GetESDpid(probability);
1488 // find most probable particle in ESD pid
1489 // 0:Electron - 1:Muon - 2:Pion - 3:Kaon - 4:Proton
1490 ipid = TMath::LocMax(5,probability);
1491 max = TMath::MaxElement(5,probability);
1497 //_____________________________________________________________________________
1498 void AliHFEsecVtx::AddHFEpairToArray(const AliHFEpairs* const pair)
1501 // Add a HFE pair to the array
1504 Int_t n = HFEpairs()->GetEntriesFast();
1505 if(n!=fNoOfHFEpairs)AliError(Form("fNoOfHFEpairs != HFEpairs()->GetEntriesFast %i != %i \n", fNoOfHFEpairs, n));
1506 new((*HFEpairs())[n]) AliHFEpairs(*pair);
1509 //_____________________________________________________________________________
1510 TClonesArray *AliHFEsecVtx::HFEpairs()
1513 // Returns the list of HFE pairs
1517 fHFEpairs = new TClonesArray("AliHFEpairs", 200);
1522 //_____________________________________________________________________________
1523 void AliHFEsecVtx::DeleteHFEpairs()
1526 // Delete the list of HFE pairs
1530 fHFEpairs->Delete();
1535 //_____________________________________________________________________________
1536 void AliHFEsecVtx::InitHFEpairs()
1539 // Initialization should be done before make all possible pairs for a given electron candidate
1545 //_____________________________________________________________________________
1546 void AliHFEsecVtx::AddHFEsecvtxToArray(const AliHFEsecVtxs* const secvtx)
1549 // Add a HFE secondary vertex to the array
1552 Int_t n = HFEsecvtxs()->GetEntriesFast();
1553 if(n!=fNoOfHFEsecvtxs)AliError(Form("fNoOfHFEsecvtxs != HFEsecvtxs()->GetEntriesFast %i != %i \n", fNoOfHFEsecvtxs, n));
1554 new((*HFEsecvtxs())[n]) AliHFEsecVtxs(*secvtx);
1557 //_____________________________________________________________________________
1558 TClonesArray *AliHFEsecVtx::HFEsecvtxs()
1561 // Returns the list of HFE secvtx
1565 fHFEsecvtxs = new TClonesArray("AliHFEsecVtxs", 200);
1570 //_____________________________________________________________________________
1571 void AliHFEsecVtx::DeleteHFEsecvtxs()
1574 // Delete the list of HFE pairs
1578 fHFEsecvtxs->Delete();
1579 //delete fHFEsecvtx;
1583 //_____________________________________________________________________________
1584 void AliHFEsecVtx::InitHFEsecvtxs()
1587 // Initialization should be done
1590 fNoOfHFEsecvtxs = 0;
1593 //____________________________________________________________
1594 void AliHFEsecVtx::MakeContainer(){
1600 const Int_t nDimPair=6;
1601 Int_t nBinPair[nDimPair] = {200, 500, 314, 2000, 2000, 13};
1602 //Int_t nBinPair[nDimPair] = {200, 500, 314, 2000, 2000, 13, 60, 60, 2000, 2000};
1603 const Double_t kInvmassmin = 0., kInvmassmax = 20.;
1604 const Double_t kKFChi2min = 0, kKFChi2max= 50;
1605 const Double_t kOpenanglemin = 0, kOpenanglemax = 3.14;
1606 const Double_t kSignedLxymin = -10, kSignedLxymax= 10;
1607 const Double_t kKFIPmin = -10, kKFIPmax= 10;
1608 const Double_t kPairCodemin = -1, kPairCodemax= 12;
1609 //const Double_t kPtmin = 0, kPtmax= 30;
1610 //const Double_t kDCAsigmin = -5, kDCAsigmax= 5;
1612 Double_t* binEdgesPair[nDimPair];
1613 for(Int_t ivar = 0; ivar < nDimPair; ivar++)
1614 binEdgesPair[ivar] = new Double_t[nBinPair[ivar] + 1];
1616 for(Int_t i=0; i<=nBinPair[0]; i++) binEdgesPair[0][i]=(Double_t)kInvmassmin + (kInvmassmax - kInvmassmin)/nBinPair[0]*(Double_t)i;
1617 for(Int_t i=0; i<=nBinPair[1]; i++) binEdgesPair[1][i]=(Double_t)kKFChi2min + (kKFChi2max - kKFChi2min)/nBinPair[1]*(Double_t)i;
1618 for(Int_t i=0; i<=nBinPair[2]; i++) binEdgesPair[2][i]=(Double_t)kOpenanglemin + (kOpenanglemax - kOpenanglemin)/nBinPair[2]*(Double_t)i;
1619 for(Int_t i=0; i<=nBinPair[3]; i++) binEdgesPair[3][i]=(Double_t)kSignedLxymin + (kSignedLxymax - kSignedLxymin)/nBinPair[3]*(Double_t)i;
1620 for(Int_t i=0; i<=nBinPair[4]; i++) binEdgesPair[4][i]=(Double_t)kKFIPmin + (kKFIPmax - kKFIPmin)/nBinPair[4]*(Double_t)i;
1621 for(Int_t i=0; i<=nBinPair[5]; i++) binEdgesPair[5][i]=(Double_t)kPairCodemin + (kPairCodemax - kPairCodemin)/nBinPair[5]*(Double_t)i;
1622 //for(Int_t i=0; i<=nBinPair[6]; i++) binEdgesPair[6][i]=(Double_t)kPtmin + (kPtmax - kPtmin)/nBinPair[6]*(Double_t)i;
1623 //for(Int_t i=0; i<=nBinPair[7]; i++) binEdgesPair[7][i]=binEdgesPair[6][i];
1624 //for(Int_t i=0; i<=nBinPair[6]; i++) binEdgesPair[6][i]=(Double_t)kDCAsigmin + (kDCAsigmax - kDCAsigmin)/nBinPair[6]*(Double_t)i;
1625 //for(Int_t i=0; i<=nBinPair[7]; i++) binEdgesPair[7][i]=binEdgesPair[6][i];
1627 fPairQA = new THnSparseF("pairQA", "QA for Pair; invmass[GeV/c^2]; KF chi2; opening angle; signed Lxy; KF ip; pair code; dca1; dca2", nDimPair, nBinPair);
1628 //fPairQA = new THnSparseF("pairQA", "QA for Pair; invmass[GeV/c^2]; KF chi2; opening angle; signed Lxy; KF ip; pair code; pt1; pt2; dca1; dca2", nDimPair, nBinPair);
1629 for(Int_t idim = 0; idim < nDimPair; idim++){
1630 fPairQA->SetBinEdges(idim, binEdgesPair[idim]);
1633 fSecVtxList->AddAt(fPairQA,0);
1635 const Int_t nDimSecvtx=9;
1636 Double_t* binEdgesSecvtx[nDimSecvtx];
1637 Int_t nBinSecvtx[nDimSecvtx] = {200, 500, 2000, 2000, 13, 10, 4, 2000, 500};
1638 const Double_t kNtrksmin = 0, kNtrksmax= 10;
1639 const Double_t kTrueBmin = 0, kTrueBmax= 4;
1640 const Double_t kPtmin = 0, kPtmax= 50;
1641 for(Int_t ivar = 0; ivar < nDimSecvtx; ivar++)
1642 binEdgesSecvtx[ivar] = new Double_t[nBinSecvtx[ivar] + 1];
1644 for(Int_t i=0; i<=nBinSecvtx[0]; i++) binEdgesSecvtx[0][i]=binEdgesPair[0][i];
1645 for(Int_t i=0; i<=nBinSecvtx[1]; i++) binEdgesSecvtx[1][i]=binEdgesPair[1][i];
1646 for(Int_t i=0; i<=nBinSecvtx[2]; i++) binEdgesSecvtx[2][i]=binEdgesPair[3][i];
1647 for(Int_t i=0; i<=nBinSecvtx[3]; i++) binEdgesSecvtx[3][i]=binEdgesPair[4][i];
1648 for(Int_t i=0; i<=nBinSecvtx[4]; i++) binEdgesSecvtx[4][i]=binEdgesPair[5][i];
1649 for(Int_t i=0; i<=nBinSecvtx[5]; i++) binEdgesSecvtx[5][i]=(Double_t)kNtrksmin + (kNtrksmax - kNtrksmin)/nBinSecvtx[5]*(Double_t)i;
1650 for(Int_t i=0; i<=nBinSecvtx[6]; i++) binEdgesSecvtx[6][i]=(Double_t)kTrueBmin + (kTrueBmax - kTrueBmin)/nBinSecvtx[6]*(Double_t)i;
1651 for(Int_t i=0; i<=nBinSecvtx[7]; i++) binEdgesSecvtx[7][i]=binEdgesPair[3][i];
1652 for(Int_t i=0; i<=nBinSecvtx[8]; i++) binEdgesSecvtx[8][i]=(Double_t)kPtmin + (kPtmax - kPtmin)/nBinSecvtx[8]*(Double_t)i;
1654 fSecvtxQA = new THnSparseF("secvtxQA", "QA for Secvtx; invmass[GeV/c^2]; KF chi2; signed Lxy; KF ip; pair code; n tracks ", nDimSecvtx, nBinSecvtx);
1655 for(Int_t idim = 0; idim < nDimSecvtx; idim++){
1656 fSecvtxQA->SetBinEdges(idim, binEdgesSecvtx[idim]);
1659 fSecVtxList->AddAt(fSecvtxQA,1);
1661 for(Int_t ivar = 0; ivar < nDimPair; ivar++)
1662 delete[] binEdgesPair[ivar];
1663 for(Int_t ivar = 0; ivar < nDimSecvtx; ivar++)
1664 delete[] binEdgesSecvtx[ivar];
1667 //____________________________________________________________
1668 void AliHFEsecVtx::MakeHistos(Int_t step){
1674 TString hname=Form("step%d",step);
1677 const Double_t kPtbound[2] = {0.1, 20.};
1679 iBin[0] = 44; // bins in pt
1680 Double_t* binEdges[1];
1681 binEdges[0] = AliHFEtools::MakeLogarithmicBinning(iBin[0], kPtbound[0], kPtbound[1]);
1683 fSecVtxList->AddAt(new TH1F(hname+"taggedElec", "pT of e", iBin[0],binEdges[0]), step);
1684 fSecVtxList->AddAt(new TH1F(hname+"charmElec", "pT of e", iBin[0],binEdges[0]), step+1);
1685 fSecVtxList->AddAt(new TH1F(hname+"beautyElec", "pT of e", iBin[0],binEdges[0]), step+2);
1686 fSecVtxList->AddAt(new TH1F(hname+"conversionElec", "pT of e", iBin[0],binEdges[0]), step+3);
1687 fSecVtxList->AddAt(new TH1F(hname+"ebgElec", "pT of e", iBin[0],binEdges[0]), step+4);
1688 fSecVtxList->AddAt(new TH1F(hname+"hcontaminElec", "pT of e", iBin[0],binEdges[0]), step+5);
1689 fSecVtxList->AddAt(new TH1F(hname+"elseElec", "pT of e", iBin[0],binEdges[0]), step+6);
1690 delete[] binEdges[0];
1693 //____________________________________________________________
1694 void AliHFEsecVtx::FillHistos(Int_t step, const AliESDtrack *track){
1702 AliMCParticle *mctrack = NULL;
1703 TParticle* mcpart = NULL;
1705 (static_cast<TH1F *>(fSecVtxList->At(step)))->Fill(track->Pt()); // electrons tagged
1707 if(HasMCData() && fMCQA){
1708 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(track->GetLabel()))))) return;
1709 mcpart = mctrack->Particle();
1711 Int_t esource=fMCQA->GetElecSource(mcpart);
1713 //if(!(dynamic_cast<TH1F *>(fSecVtxList->At(step+1)))) return;
1714 (static_cast<TH1F *>(fSecVtxList->At(step+1)))->Fill(mcpart->Pt()); //charm
1716 else if(esource==2 || esource==3) {
1717 //if(!(dynamic_cast<TH1F *>(fSecVtxList->At(step+2)))) return;
1718 (static_cast<TH1F *>(fSecVtxList->At(step+2)))->Fill(mcpart->Pt()); //beauty
1720 else if(esource>12 && esource<19) {
1721 //else if(esource==4) {
1722 //if(!(dynamic_cast<TH1F *>(fSecVtxList->At(step+3)))) return;
1723 (static_cast<TH1F *>(fSecVtxList->At(step+3)))->Fill(mcpart->Pt()); //conversion
1725 else if(esource==7) {
1726 //if(!(dynamic_cast<TH1F *>(fSecVtxList->At(step+5)))) return;
1727 (static_cast<TH1F *>(fSecVtxList->At(step+5)))->Fill(mcpart->Pt()); //contamination
1729 else if(!(esource<0)) {
1730 //if(!(dynamic_cast<TH1F *>(fSecVtxList->At(step+4)))) return;
1731 (static_cast<TH1F *>(fSecVtxList->At(step+4)))->Fill(mcpart->Pt()); //e backgrounds
1734 //if(!(dynamic_cast<TH1F *>(fSecVtxList->At(step+6)))) return;
1735 (static_cast<TH1F *>(fSecVtxList->At(step+6)))->Fill(mcpart->Pt()); //something else?