1 /**************************************************************************
2 * Copyright(c) 1998-2010, 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 **************************************************************************/
18 /////////////////////////////////////////////////////////////
20 // Class for cuts on AOD reconstructed Lc->V0+X
22 // Modified by A.De Caro - decaro@sa.infn.it
24 /////////////////////////////////////////////////////////////
26 #include <Riostream.h>
28 #include <TDatabasePDG.h>
31 #include "AliAnalysisManager.h"
32 #include "AliInputEventHandler.h"
33 #include "AliPIDResponse.h"
34 #include "AliRDHFCutsLctoV0.h"
35 #include "AliAODRecoCascadeHF.h"
36 #include "AliAODTrack.h"
37 #include "AliESDtrack.h"
38 #include "AliESDVertex.h"
39 #include "AliAODVertex.h"
46 ClassImp(AliRDHFCutsLctoV0)
48 //--------------------------------------------------------------------------
49 AliRDHFCutsLctoV0::AliRDHFCutsLctoV0(const char* name, Short_t /*v0channel*/) :
58 // Default Constructor
63 TString varNames[nvars]={"Lc inv. mass if K0S [GeV/c2]", // 0
64 "Lc inv. mass if Lambda [GeV/c2]", // 1
65 "K0S inv. mass [GeV/c2]", // 2
66 "Lambda/LambdaBar inv. mass[GeV/c2]", // 3
67 "pT min bachelor track [GeV/c]", // 4
68 "pT min V0-positive track [GeV/c]", // 5
69 "pT min V0-negative track [GeV/c]", // 6
70 "dca cascade (prong-to-prong) cut [cm]", // 7
71 "dca V0 (prong-to-prong) cut [number of sigmas]", // 8
72 "V0 cosPA min wrt PV", // 9
73 "d0 max bachelor wrt PV [cm]", // 10
74 "d0 max V0 wrt PV [cm]", // 11
75 "mass K0S veto [GeV/c2]", // 12
76 "mass Lambda/LambdaBar veto [GeV/c2]", // 13
77 "mass Gamma veto [GeV/c2]", // 14
78 "pT min V0 track [GeV/c]", // 15
82 Bool_t isUpperCut[nvars]={kTRUE, // 0
100 SetVarNames(nvars,varNames,isUpperCut);
101 Bool_t forOpt[nvars]={kFALSE, // 0
119 SetVarsForOpt(nvars,forOpt);
121 Float_t limits[2]={0,999999999.};
139 //--------------------------------------------------------------------------
140 AliRDHFCutsLctoV0::AliRDHFCutsLctoV0(const AliRDHFCutsLctoV0 &source) :
142 fPidSelectionFlag(source.fPidSelectionFlag),
146 fV0Type(source.fV0Type)
147 /*fV0channel(source.fV0channel)*/
153 if (source.fPidHFV0pos) fPidHFV0pos = new AliAODPidHF(*(source.fPidHFV0pos));
154 else fPidHFV0pos = new AliAODPidHF();
155 if (source.fPidHFV0neg) fPidHFV0neg = new AliAODPidHF(*(source.fPidHFV0neg));
156 else fPidHFV0neg = new AliAODPidHF();
158 if (source.fV0daughtersCuts) AddTrackCutsV0daughters(source.fV0daughtersCuts);
159 else fV0daughtersCuts = new AliESDtrackCuts();
162 //--------------------------------------------------------------------------
163 AliRDHFCutsLctoV0 &AliRDHFCutsLctoV0::operator=(const AliRDHFCutsLctoV0 &source)
166 // assignment operator
169 if (this != &source) {
171 AliRDHFCuts::operator=(source);
172 fPidSelectionFlag = source.fPidSelectionFlag;
174 fPidHFV0pos = new AliAODPidHF(*(source.fPidHFV0pos));
176 fPidHFV0neg = new AliAODPidHF(*(source.fPidHFV0neg));
178 delete fV0daughtersCuts;
179 fV0daughtersCuts = new AliESDtrackCuts(*(source.fV0daughtersCuts));
181 fV0Type = source.fV0Type;
189 //---------------------------------------------------------------------------
190 AliRDHFCutsLctoV0::~AliRDHFCutsLctoV0() {
192 // Default Destructor
204 if (fV0daughtersCuts) {
205 delete fV0daughtersCuts;
211 //---------------------------------------------------------------------------
212 void AliRDHFCutsLctoV0::GetCutVarsForOpt(AliAODRecoDecayHF *d,Float_t *vars,Int_t nvars,Int_t *pdgdaughters) {
214 // Fills in vars the values of the variables
217 if (pdgdaughters[0]==-9999) return; // dummy
219 if (nvars!=fnVarsForOpt) {
220 printf("AliRDHFCutsLctoV0::GetCutVarsForOpt: wrong number of variables\n");
224 Double_t mLcPDG = TDatabasePDG::Instance()->GetParticle(4122)->Mass();
225 Double_t mk0sPDG = TDatabasePDG::Instance()->GetParticle(310)->Mass();
226 Double_t mLPDG = TDatabasePDG::Instance()->GetParticle(3122)->Mass();
228 AliAODRecoCascadeHF *dd = (AliAODRecoCascadeHF*)d;
230 // Get the v0 and all daughter tracks
231 AliAODTrack *bachelorTrack = (AliAODTrack*)dd->GetBachelor();
232 AliAODv0 *v0 = (AliAODv0*)dd->Getv0();
233 AliAODTrack *v0positiveTrack = (AliAODTrack*)dd->Getv0PositiveTrack();
234 AliAODTrack *v0negativeTrack = (AliAODTrack*)dd->Getv0NegativeTrack();
237 // cut on cascade mass, if K0S + p
238 if (fVarsForOpt[0]) {
240 vars[iter]=TMath::Abs(dd->InvMassLctoK0sP()-mLcPDG);
242 // cut on cascade mass, if Lambda/LambdaBar + pi
243 if (fVarsForOpt[1]) {
245 vars[iter]=TMath::Abs(dd->InvMassLctoLambdaPi()-mLcPDG);
248 // cut on V0 mass if K0S
249 if (fVarsForOpt[2]) {
251 vars[iter]=TMath::Abs(v0->MassK0Short()-mk0sPDG);
254 // cut on V0 mass if Lambda/LambdaBar
255 if (fVarsForOpt[3]) {
257 if (bachelorTrack->Charge()==1) {
259 vars[iter]=TMath::Abs(v0->MassLambda()-mLPDG);
260 } else if (bachelorTrack->Charge()==-1) {
262 vars[iter]=TMath::Abs(v0->MassAntiLambda()-mLPDG);
267 // cut bachelor min pt
268 if (fVarsForOpt[4]) {
270 vars[iter]=bachelorTrack->Pt();
273 // cut on V0-positive min pt
274 if (fVarsForOpt[5]) {
276 vars[iter]=v0positiveTrack->Pt();
279 // cut on V0-negative min pt
280 if (fVarsForOpt[6]) {
282 vars[iter]=v0negativeTrack->Pt();
285 // cut on cascade dca (prong-to-prong)
286 if (fVarsForOpt[7]) {
288 vars[iter]=dd->GetDCA(); // prong-to-prong cascade DCA
291 // cut on V0 dca (prong-to-prong)
292 if (fVarsForOpt[8]) {
294 vars[iter]=v0->GetDCA(); // prong-to-prong V0 DCA
297 // cut on V0 cosPA wrt PV
298 if (fVarsForOpt[9]) {
300 vars[iter]=dd->CosV0PointingAngle(); // cosine of V0 pointing angle wrt primary vertex
303 // cut on bachelor transverse impact parameter wrt PV
304 if (fVarsForOpt[10]) {
306 vars[iter]=dd->Getd0Prong(0); // bachelor transverse impact parameter wrt primary vertex
309 // cut on V0 transverse impact parameter wrt PV
310 if (fVarsForOpt[11]) {
312 vars[iter]=dd->Getd0Prong(1); // V0 transverse impact parameter wrt primary vertex
315 // cut on K0S invariant mass veto
316 if (fVarsForOpt[12]) {
318 vars[iter]=TMath::Abs(v0->MassK0Short()-mk0sPDG); // K0S invariant mass veto
321 // cut on Lambda/LambdaBar invariant mass veto
322 if (fVarsForOpt[13]) {
324 if (bachelorTrack->Charge()==1) {
326 vars[iter]=TMath::Abs(v0->MassLambda()-mLPDG);
327 } else if (bachelorTrack->Charge()==-1) {
329 vars[iter]=TMath::Abs(v0->MassAntiLambda()-mLPDG);
334 // cut on gamma invariant mass veto
335 if (fVarsForOpt[14]) {
337 vars[iter]= v0->InvMass2Prongs(0,1,11,11); // gamma invariant mass veto
342 if (fVarsForOpt[15]) {
344 vars[iter]= v0->Pt(); // V0 pT min
350 //---------------------------------------------------------------------------
351 Int_t AliRDHFCutsLctoV0::IsSelected(TObject* obj,Int_t selectionLevel) {
357 AliFatal("Cut matrice not inizialized. Exit...");
361 AliAODRecoCascadeHF* d = (AliAODRecoCascadeHF*)obj;
363 AliDebug(2,"AliAODRecoCascadeHF null");
367 if (!d->GetSecondaryVtx()) {
368 AliDebug(2,"No secondary vertex for cascade");
372 if (d->GetNDaughters()!=2) {
373 AliDebug(2,Form("No 2 daughters for current cascade (nDaughters=%d)",d->GetNDaughters()));
377 AliAODv0 * v0 = dynamic_cast<AliAODv0*>(d->Getv0());
378 AliAODTrack * bachelorTrack = dynamic_cast<AliAODTrack*>(d->GetBachelor());
379 if (!v0 || !bachelorTrack) {
380 AliDebug(2,"No V0 or no bachelor for current cascade");
384 if (bachelorTrack->GetID()<0) {
385 AliDebug(2,Form("Bachelor has negative ID %d",bachelorTrack->GetID()));
389 if (!v0->GetSecondaryVtx()) {
390 AliDebug(2,"No secondary vertex for V0 by cascade");
394 if (v0->GetNDaughters()!=2) {
395 AliDebug(2,Form("No 2 daughters for V0 of current cascade (onTheFly=%d, nDaughters=%d)",v0->GetOnFlyStatus(),v0->GetNDaughters()));
400 // Get the V0 daughter tracks
401 AliAODTrack *v0positiveTrack = dynamic_cast<AliAODTrack*>(d->Getv0PositiveTrack());
402 AliAODTrack *v0negativeTrack = dynamic_cast<AliAODTrack*>(d->Getv0NegativeTrack());
403 if (!v0positiveTrack || !v0negativeTrack ) {
404 AliDebug(2,"No V0 daughters' objects");
408 if (v0positiveTrack->GetID()<0 || v0negativeTrack->GetID()<0) {
409 AliDebug(2,Form("At least one of V0 daughters has negative ID %d %d",v0positiveTrack->GetID(),v0negativeTrack->GetID()));
413 //if(fUseTrackSelectionWithFilterBits && d->HasBadDaughters()) return 0;
414 if ( fUseTrackSelectionWithFilterBits && !(bachelorTrack->TestFilterMask(BIT(4))) ) {
415 AliDebug(2,"Check on the bachelor FilterBit: no BIT(4). Candidate rejected.");
419 Int_t returnvalueTrack = 7;
421 // selection on daughter tracks
422 if (selectionLevel==AliRDHFCuts::kAll ||
423 selectionLevel==AliRDHFCuts::kTracks) {
425 if (!AreLctoV0DaughtersSelected(d)) return 0;
429 Bool_t okLck0sp=kTRUE, okLcLpi=kTRUE, okLcLBarpi=kTRUE;
431 // selection on candidate
432 if (selectionLevel==AliRDHFCuts::kAll ||
433 selectionLevel==AliRDHFCuts::kCandidate) {
435 Double_t pt = d->Pt();
436 Int_t ptbin = PtBin(pt);
438 Double_t mLcPDG = TDatabasePDG::Instance()->GetParticle(4122)->Mass();
439 Double_t mk0sPDG = TDatabasePDG::Instance()->GetParticle(310)->Mass();
440 Double_t mLPDG = TDatabasePDG::Instance()->GetParticle(3122)->Mass();
443 Double_t mk0s = v0->MassK0Short();
444 Double_t mLck0sp = d->InvMassLctoK0sP();
447 Double_t mlambda = v0->MassLambda();
448 Double_t malambda = v0->MassAntiLambda();
449 Double_t mLcLpi = d->InvMassLctoLambdaPi();
451 Bool_t okK0spipi=kTRUE, okLppi=kTRUE, okLBarpip=kTRUE;
452 Bool_t isNotK0S = kTRUE, isNotLambda = kTRUE, isNotLambdaBar = kTRUE, isNotGamma = kTRUE;
454 // cut on Lc mass with K0S+p hypothesis
455 if (TMath::Abs(mLck0sp-mLcPDG) > fCutsRD[GetGlobalIndex(0,ptbin)]) {
457 AliDebug(4,Form(" cascade mass is %2.2e and does not correspond to Lambda_c into K0S+p cut",mLck0sp));
460 // cuts on the V0 mass: K0S case
461 if (TMath::Abs(mk0s-mk0sPDG) > fCutsRD[GetGlobalIndex(2,ptbin)]) {
463 AliDebug(4,Form(" V0 mass is %2.2e and does not correspond to K0S cut",mk0s));
466 // cut on Lc mass with Lambda+pi hypothesis
467 if (TMath::Abs(mLcLpi-mLcPDG) > fCutsRD[GetGlobalIndex(1,ptbin)]) {
470 AliDebug(4,Form(" cascade mass is %2.2e and does not correspond to Lambda_c into Lambda+pi cut",mLcLpi));
473 // cuts on the V0 mass: Lambda/LambdaBar case
474 //if ( !(bachelorTrack->Charge()==+1 && TMath::Abs(mlambda-mLPDG) <= fCutsRD[GetGlobalIndex(3,ptbin)] ) ) {
475 if ( TMath::Abs(mlambda-mLPDG) > fCutsRD[GetGlobalIndex(3,ptbin)] ) {
477 AliDebug(4,Form(" V0 mass is %2.2e and does not correspond to LambdaBar cut",mlambda));
480 //if ( !(bachelorTrack->Charge()==-1 && TMath::Abs(malambda-mLPDG) <= fCutsRD[GetGlobalIndex(3,ptbin)] ) ) {
481 if ( TMath::Abs(malambda-mLPDG) > fCutsRD[GetGlobalIndex(3,ptbin)] ) {
483 AliDebug(4,Form(" V0 mass is %2.2e and does not correspond to LambdaBar cut",malambda));
486 // cut on K0S invariant mass veto
487 if (TMath::Abs(v0->MassK0Short()-mk0sPDG) < fCutsRD[GetGlobalIndex(12,ptbin)]) { // K0S invariant mass veto
488 AliDebug(4,Form(" veto on K0S invariant mass doesn't pass the cut"));
493 // cut on Lambda/LambdaBar invariant mass veto
494 if (TMath::Abs(v0->MassLambda()-mLPDG) < fCutsRD[GetGlobalIndex(13,ptbin)]) { // Lambda invariant mass veto
495 AliDebug(4,Form(" veto on Lambda invariant mass doesn't pass the cut"));
499 if (TMath::Abs(v0->MassAntiLambda()-mLPDG) < fCutsRD[GetGlobalIndex(13,ptbin)] ) { // LambdaBar invariant mass veto
500 AliDebug(4,Form(" veto on LambdaBar invariant mass doesn't pass the cut"));
501 isNotLambdaBar=kFALSE;
505 // cut on gamma invariant mass veto
506 if (v0->InvMass2Prongs(0,1,11,11) < fCutsRD[GetGlobalIndex(14,ptbin)]) { // K0S invariant mass veto
507 AliDebug(4,Form(" veto on gamma invariant mass doesn't pass the cut"));
512 okLck0sp = okLck0sp && okK0spipi && isNotLambda && isNotLambdaBar && isNotGamma;
513 okLcLpi = okLcLpi && okLppi && isNotK0S && isNotLambdaBar && isNotGamma;
514 okLcLBarpi = okLcLBarpi && okLBarpip && isNotK0S && isNotLambda && isNotGamma;
516 if (!okLck0sp && !okLcLpi && !okLcLBarpi) return 0;
518 // cuts on the minimum pt of the tracks
519 if (TMath::Abs(bachelorTrack->Pt()) < fCutsRD[GetGlobalIndex(4,ptbin)]) {
520 AliDebug(4,Form(" bachelor track Pt=%2.2e > %2.2e",bachelorTrack->Pt(),fCutsRD[GetGlobalIndex(4,ptbin)]));
523 if (TMath::Abs(v0positiveTrack->Pt()) < fCutsRD[GetGlobalIndex(5,ptbin)]) {
524 AliDebug(4,Form(" V0-positive track Pt=%2.2e > %2.2e",v0positiveTrack->Pt(),fCutsRD[GetGlobalIndex(5,ptbin)]));
527 if (TMath::Abs(v0negativeTrack->Pt()) < fCutsRD[GetGlobalIndex(6,ptbin)]) {
528 AliDebug(4,Form(" V0-negative track Pt=%2.2e > %2.2e",v0negativeTrack->Pt(),fCutsRD[GetGlobalIndex(6,ptbin)]));
532 // cut on cascade dca (prong-to-prong)
533 if ( TMath::Abs(d->GetDCA()) > fCutsRD[GetGlobalIndex(7,ptbin)] ) { // prong-to-prong cascade DCA
534 AliDebug(4,Form(" cascade tracks DCA don't pass the cut"));
538 // cut on V0 dca (prong-to-prong)
539 if ( TMath::Abs(v0->GetDCA()) > fCutsRD[GetGlobalIndex(8,ptbin)] ) { // prong-to-prong V0 DCA
540 AliDebug(4,Form(" V0 DCA don't pass the cut"));
544 // cut on V0 cosine of pointing angle wrt PV
545 if (d->CosV0PointingAngle() < fCutsRD[GetGlobalIndex(9,ptbin)]) { // cosine of V0 pointing angle wrt primary vertex
546 AliDebug(4,Form(" V0 cosine of pointing angle doesn't pass the cut"));
550 // cut on bachelor transverse impact parameter wrt PV
551 if (TMath::Abs(d->Getd0Prong(0)) > fCutsRD[GetGlobalIndex(10,ptbin)]) { // bachelor transverse impact parameter wrt PV
552 AliDebug(4,Form(" bachelor transverse impact parameter doesn't pass the cut"));
556 // cut on V0 transverse impact parameter wrt PV
557 if (TMath::Abs(d->Getd0Prong(1)) > fCutsRD[GetGlobalIndex(11,ptbin)]) { // V0 transverse impact parameter wrt PV
558 AliDebug(4,Form(" V0 transverse impact parameter doesn't pass the cut"));
563 if (v0->Pt() < fCutsRD[GetGlobalIndex(15,ptbin)]) { // V0 pT min
564 AliDebug(4,Form(" V0 track Pt=%2.2e > %2.2e",v0->Pt(),fCutsRD[GetGlobalIndex(15,ptbin)]));
570 Int_t returnvalue = okLck0sp+2*okLcLBarpi+4*okLcLpi;
575 3 Lc->K0S + p AND Lc->LambdaBar + pi
577 5 Lc->K0S + p AND Lc->Lambda + pi
578 6 Lc->LambdaBar + pi AND Lc->Lambda + pi
579 7 Lc->K0S + p AND Lc->LambdaBar + pi AND Lc->Lambda + pi
582 Int_t returnvaluePID = 7;
584 // selection on candidate
585 if (selectionLevel==AliRDHFCuts::kAll ||
586 selectionLevel==AliRDHFCuts::kCandidate ||
587 selectionLevel==AliRDHFCuts::kPID )
588 returnvaluePID = IsSelectedPID(d);
590 //if (fUsePID && returnvaluePID==0) return 0;
592 Int_t returnvalueTot = 0;
594 returnvalueTot = CombineCuts(returnvalueTrack,returnvalue,returnvaluePID);
596 returnvalueTot = CombineCuts(returnvalueTrack,returnvalue,7);
598 return returnvalueTot;
601 //---------------------------------------------------------------------------
602 Int_t AliRDHFCutsLctoV0::IsSelectedPID(AliAODRecoDecayHF* obj) {
604 // fPidHF -> PID object for bachelor
605 // fPidHFV0pos -> PID object for positive V0 daughter
606 // fPidHFV0neg -> PID object for negative V0 daughter
608 if (!fUsePID || !obj) {
609 AliDebug(2,"PID selection inactive. Candidate accepted.");
610 return 7; // all hypothesis are valid
613 if (fPidHF->GetPidResponse()==0x0 ||
614 fPidHFV0pos->GetPidResponse()==0x0 ||
615 fPidHFV0neg->GetPidResponse()==0x0) {
616 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
617 AliInputEventHandler *inputHandler=(AliInputEventHandler*)mgr->GetInputEventHandler();
618 AliPIDResponse *pidResp=inputHandler->GetPIDResponse();
619 fPidHF->SetPidResponse(pidResp);
620 fPidHFV0pos->SetPidResponse(pidResp);
621 fPidHFV0neg->SetPidResponse(pidResp);
622 fPidHF->SetOldPid(kFALSE);
623 fPidHFV0pos->SetOldPid(kFALSE);
624 fPidHFV0neg->SetOldPid(kFALSE);
627 AliAODRecoCascadeHF *objD = (AliAODRecoCascadeHF*)obj;
629 Bool_t isPeriodd = fPidHF->GetOnePad();
630 Bool_t isMC = fPidHF->GetMC();
633 fPidHFV0pos->SetOnePad(kTRUE);
634 fPidHFV0neg->SetOnePad(kTRUE);
637 fPidHFV0neg->SetMC(kTRUE);
638 fPidHFV0pos->SetMC(kTRUE);
641 AliAODTrack *bachelor = (AliAODTrack*)objD->GetBachelor();
642 AliAODTrack *v0Pos = (AliAODTrack*)objD->Getv0PositiveTrack();
643 AliAODTrack *v0Neg = (AliAODTrack*)objD->Getv0NegativeTrack();
645 if (!bachelor || !v0Pos || !v0Neg) return 0;
647 Bool_t okLcK0Sp = kTRUE; // K0S case
648 Bool_t okLcLambdaBarPi = kTRUE; // LambdaBar case
649 Bool_t okLcLambdaPi = kTRUE; // Lambda case
651 CheckPID(bachelor,v0Neg,v0Pos,okLcK0Sp,okLcLambdaBarPi,okLcLambdaPi);
653 Int_t returnvalue = okLcK0Sp+2*okLcLambdaBarPi+4*okLcLambdaPi;
657 //-----------------------
658 void AliRDHFCutsLctoV0::CheckPID(AliAODTrack *bachelor, AliAODTrack *v0Neg, AliAODTrack *v0Pos,
659 Bool_t &isBachelorID1, Bool_t &isV0NegID2, Bool_t &isV0PosID4) {
660 // identification strategy
662 Int_t trackIDtof = -1;
663 Int_t trackIDtpc = -1;
665 Bool_t dummy1 = kFALSE;
666 Bool_t dummy2 = kFALSE;
667 Bool_t dummy4 = kFALSE;
671 Double_t nTPCsigmasPr = -999;
672 Double_t nTOFsigmasPr = -999;
674 Bool_t trackIDtofB = -1;
675 Bool_t trackIDtpcB = -1;
677 switch (fPidSelectionFlag) {
682 trackIDtof = fPidHF->ApplyPidTOFRaw(bachelor,4);
683 trackIDtpc = fPidHF->ApplyPidTPCRaw(bachelor,4);
684 AliDebug(1,Form(" fPidHF->ApplyPidTOFRaw(bachelor,4)=%d fPidHF->ApplyPidTPCRaw(bachelor,4)=%d",trackIDtof,trackIDtpc));
685 isBachelorID1 = (trackIDtof==4) && (trackIDtpc==4); // K0S case
686 //isBachelorID2 = (fPidHF->ApplyPidTPCRaw(bachelor,2)==2) && (fPidHF->ApplyPidTOFRaw(bachelor,2)==2); // LambdaBar case
687 //isBachelorID4 = isBachelorID2; // Lambda case
690 trackIDtof = fPidHFV0neg->ApplyPidTOFRaw(v0Neg,4);
691 trackIDtpc = fPidHFV0neg->ApplyPidTPCRaw(v0Neg,4);
692 AliDebug(1,Form(" fPidHFV0neg->ApplyPidTOFRaw(v0Neg,4)=%d fPidHFV0neg->ApplyPidTPCRaw(v0Neg,4)=%d",trackIDtof,trackIDtpc));
693 //isV0NegID1 = (fPidHFV0neg->ApplyPidTPCRaw(v0Neg,2)==2) && (fPidHFV0neg->ApplyPidTOFRaw(v0Neg,2)==2); // K0S case
694 isV0NegID2 = (trackIDtof==4) && (trackIDtpc==4); // LambdaBar case
695 //isV0NegID4 = isV0NegID1; // Lambda case
698 trackIDtof = fPidHFV0pos->ApplyPidTOFRaw(v0Pos,4);
699 trackIDtpc = fPidHFV0pos->ApplyPidTPCRaw(v0Pos,4);
700 AliDebug(1,Form(" fPidHFV0pos->ApplyPidTOFRaw(v0Pos,4)=%d fPidHFV0pos->ApplyPidTPCRaw(v0POS,4)=%d",trackIDtof,trackIDtpc));
701 //isV0PosID1 = (fPidHFV0pos->ApplyPidTPCRaw(v0Pos,2)==2) && (fPidHFV0pos->ApplyPidTOFRaw(v0Pos,2)==2); // K0S case
702 //isV0PosID2 = isV0PosID1; // LambdaBar case
703 isV0PosID4 = (trackIDtof==4) && (trackIDtpc==4); // Lambda case
709 trackIDtof = fPidHF->ApplyPidTOFRaw(bachelor,4);
710 trackIDtpc = fPidHF->ApplyPidTPCRaw(bachelor,4);
711 AliDebug(1,Form(" fPidHF->ApplyPidTOFRaw(bachelor,4)=%d fPidHFV0->ApplyPidTPCRaw(bachelor,4)=%d",trackIDtof,trackIDtpc));
712 isBachelorID1 = ( trackIDtof==4 );
713 dummy1 = ( !(fPidHF->CheckTOFPIDStatus(bachelor)) && (trackIDtpc==4) &&
714 fPidHF->IsExcluded(bachelor,2,2.,"TPC") && fPidHF->IsExcluded(bachelor,3,2.,"TPC") ); // K0S case
715 isBachelorID1 = isBachelorID1 || dummy1;
719 trackIDtof = fPidHFV0neg->ApplyPidTOFRaw(v0Neg,4);
720 trackIDtpc = fPidHFV0neg->ApplyPidTPCRaw(v0Neg,4);
721 AliDebug(1,Form(" fPidHFV0neg->ApplyPidTOFRaw(v0Neg,4)=%d fPidHFV0neg->ApplyPidTPCRaw(v0Neg,4)=%d",trackIDtof,trackIDtpc));
722 isV0NegID2 = ( trackIDtof==4 );
723 dummy2 = ( !(fPidHFV0neg->CheckTOFPIDStatus(v0Neg)) && (trackIDtpc==4) &&
724 fPidHFV0neg->IsExcluded(v0Neg,2,2.,"TPC") && fPidHFV0neg->IsExcluded(v0Neg,3,2.,"TPC") ); // LambdaBar case
725 isV0NegID2 = isV0NegID2 || dummy2;
729 trackIDtof = fPidHFV0pos->ApplyPidTOFRaw(v0Pos,4);
730 trackIDtpc = fPidHFV0pos->ApplyPidTPCRaw(v0Pos,4);
731 AliDebug(1,Form(" fPidHFV0pos->ApplyPidTOFRaw(v0Pos,4)=%d fPidHFV0pos->ApplyPidTPCRaw(v0Pos,4)=%d",trackIDtof,trackIDtpc));
732 isV0PosID4 = ( trackIDtof==4 );
733 dummy4 = ( !(fPidHFV0pos->CheckTOFPIDStatus(v0Pos)) && (trackIDtpc==4) &&
734 fPidHFV0pos->IsExcluded(v0Pos,2,2.,"TPC") && fPidHFV0pos->IsExcluded(v0Pos,3,2.,"TPC") ); // Lambda case
735 isV0PosID4 = isV0PosID4 || dummy4;
743 tofID = fPidHF->GetnSigmaTOF(bachelor,4,nTOFsigmasPr);
745 tpcID = fPidHF->GetnSigmaTPC(bachelor,4,nTPCsigmasPr);
746 trackIDtofB = (tofID==1) && ( (bachelor->P()>=1.0 && bachelor->P()<2.5 && TMath::Abs(nTOFsigmasPr)<3) ||
747 (bachelor->P()>=2.5 && nTOFsigmasPr>-2 && nTOFsigmasPr<3) );
748 trackIDtpcB = (tpcID==1) && ( (bachelor->P()<1.0 && TMath::Abs(nTPCsigmasPr)<2) ||
749 (bachelor->P()>=1.0 && TMath::Abs(nTPCsigmasPr)<3) );
750 AliDebug(1,Form(" trackIDtofB=%d trackIDtpcB=%d",trackIDtofB,trackIDtpcB));
751 isBachelorID1 = (bachelor->P()<1 && trackIDtpcB) || (bachelor->P()>=1 && trackIDtpcB && trackIDtofB); // K0S case
755 tofID = fPidHFV0neg->GetnSigmaTOF(v0Neg,4,nTOFsigmasPr);
757 tpcID = fPidHFV0neg->GetnSigmaTPC(v0Neg,4,nTPCsigmasPr);
758 trackIDtofB = (tofID==1) && ( (v0Neg->P()>=1.0 && v0Neg->P()<2.5 && TMath::Abs(nTOFsigmasPr)<3) ||
759 (v0Neg->P()>=2.5 && nTOFsigmasPr>-2 && nTOFsigmasPr<3) );
760 trackIDtpcB = (tpcID==1) && ( (v0Neg->P()<1.0 && TMath::Abs(nTPCsigmasPr)<2) ||
761 (v0Neg->P()>=1.0 && TMath::Abs(nTPCsigmasPr)<3) );
762 AliDebug(1,Form(" trackIDtofB=%d trackIDtpcB=%d",trackIDtofB,trackIDtpcB));
763 isV0NegID2 = (v0Neg->P()<1 && trackIDtpcB) || (v0Neg->P()>=1 && trackIDtpcB && trackIDtofB); // LambdaBar case
767 tofID = fPidHFV0pos->GetnSigmaTOF(v0Pos,4,nTOFsigmasPr);
769 tpcID = fPidHFV0pos->GetnSigmaTPC(v0Pos,4,nTPCsigmasPr);
770 trackIDtofB = (tofID==1) && ( (v0Pos->P()>=1.0 && v0Pos->P()<2.5 && TMath::Abs(nTOFsigmasPr)<3) ||
771 (v0Pos->P()>=2.5 && nTOFsigmasPr>-2 && nTOFsigmasPr<3) );
772 trackIDtpcB = (tpcID==1) && ( (v0Pos->P()<1.0 && TMath::Abs(nTPCsigmasPr)<2) ||
773 (v0Pos->P()>=1.0 && TMath::Abs(nTPCsigmasPr)<3) );
774 AliDebug(1,Form(" trackIDtofB=%d trackIDtpcB=%d",trackIDtofB,trackIDtpcB));
775 isV0PosID4 = (v0Pos->P()<1 && trackIDtpcB) || (v0Pos->P()>=1 && trackIDtpcB && trackIDtofB); // Lambda case
782 tofID = fPidHF->GetnSigmaTOF(bachelor,4,nTOFsigmasPr);
784 tpcID = fPidHF->GetnSigmaTPC(bachelor,4,nTPCsigmasPr);
785 trackIDtofB = (tofID==1) && ( (bachelor->P()>=1.0 && bachelor->P()<2.5 && TMath::Abs(nTOFsigmasPr)<3) ||
786 (bachelor->P()>=2.5 && nTOFsigmasPr>-2 && nTOFsigmasPr<3) );
787 trackIDtpcB = (tpcID==1) && ( (bachelor->P()<1.0 && TMath::Abs(nTPCsigmasPr)<2) ||
788 (bachelor->P()>=1.0 && bachelor->P()<3.0 && TMath::Abs(nTPCsigmasPr)<3) ||
789 (bachelor->P()>=3.0 && nTPCsigmasPr>-3 && nTPCsigmasPr<2) );
790 AliDebug(1,Form(" trackIDtofB=%d trackIDtpcB=%d",trackIDtofB,trackIDtpcB));
791 isBachelorID1 = (bachelor->P()<1 && trackIDtpcB) || (bachelor->P()>=1 && trackIDtpcB && trackIDtofB); // K0S case
795 tofID = fPidHFV0neg->GetnSigmaTOF(v0Neg,4,nTOFsigmasPr);
797 tpcID = fPidHFV0neg->GetnSigmaTPC(v0Neg,4,nTPCsigmasPr);
798 trackIDtofB = (tofID==1) && ( (v0Neg->P()>=1.0 && v0Neg->P()<2.5 && TMath::Abs(nTOFsigmasPr)<3) ||
799 (v0Neg->P()>=2.5 && nTOFsigmasPr>-2 && nTOFsigmasPr<3) );
800 trackIDtpcB = (tpcID==1) && ( (v0Neg->P()<1.0 && TMath::Abs(nTPCsigmasPr)<2) ||
801 (v0Neg->P()>=1.0 && v0Neg->P()<3.0 && TMath::Abs(nTPCsigmasPr)<3) ||
802 (v0Neg->P()>=3.0 && nTPCsigmasPr>-3 && nTPCsigmasPr<2) );
803 AliDebug(1,Form(" trackIDtofB=%d trackIDtpcB=%d",trackIDtofB,trackIDtpcB));
804 isV0NegID2 = (v0Neg->P()<1 && trackIDtpcB) || (v0Neg->P()>=1 && trackIDtpcB && trackIDtofB); // LambdaBar case
808 tofID = fPidHFV0pos->GetnSigmaTOF(v0Pos,4,nTOFsigmasPr);
810 tpcID = fPidHFV0pos->GetnSigmaTPC(v0Pos,4,nTPCsigmasPr);
811 trackIDtofB = (tofID==1) && ( (v0Pos->P()>=1.0 && v0Pos->P()<2.5 && TMath::Abs(nTOFsigmasPr)<3) ||
812 (v0Pos->P()>=2.5 && nTOFsigmasPr>-2 && nTOFsigmasPr<3) );
813 trackIDtpcB = (tpcID==1) && ( (v0Pos->P()<1.0 && TMath::Abs(nTPCsigmasPr)<2) ||
814 (v0Pos->P()>=1.0 && v0Pos->P()<3.0 && TMath::Abs(nTPCsigmasPr)<3) ||
815 (v0Pos->P()>=3.0 && nTPCsigmasPr>-3 && nTPCsigmasPr<2) );
816 AliDebug(1,Form(" trackIDtofB=%d trackIDtpcB=%d",trackIDtofB,trackIDtpcB));
817 isV0PosID4 = (v0Pos->P()<1 && trackIDtpcB) || (v0Pos->P()>=1 && trackIDtpcB && trackIDtofB); // Lambda case
825 Int_t AliRDHFCutsLctoV0::CombineCuts(Int_t returnvalueTrack, Int_t returnvalue, Int_t returnvaluePID) const {
827 // combine track selection, topological cuts and PID
830 Int_t returnvalueTot=returnvalueTrack&returnvalue;
831 returnvalueTot=returnvalueTot&returnvaluePID;
833 return returnvalueTot;
836 //----------------------------------
837 Int_t AliRDHFCutsLctoV0::IsSelectedSingleCut(TObject* obj, Int_t selectionLevel, Int_t cutIndex) {
839 // Apply selection on single cut
843 AliDebug(2,"Cut matrice not inizialized. Exit...");
847 AliAODRecoCascadeHF* d = (AliAODRecoCascadeHF*)obj;
849 AliDebug(2,"AliAODRecoCascadeHF null");
853 if (!d->GetSecondaryVtx()) {
854 AliDebug(2,"No secondary vertex for cascade");
858 if (d->GetNDaughters()!=2) {
859 AliDebug(2,Form("No 2 daughters for current cascade (nDaughters=%d)",d->GetNDaughters()));
863 AliAODv0 * v0 = dynamic_cast<AliAODv0*>(d->Getv0());
864 AliAODTrack * bachelorTrack = dynamic_cast<AliAODTrack*>(d->GetBachelor());
865 if (!v0 || !bachelorTrack) {
866 AliDebug(2,"No V0 or no bachelor for current cascade");
870 if (bachelorTrack->GetID()<0) {
871 AliDebug(2,Form("Bachelor has negative ID %d",bachelorTrack->GetID()));
875 if (!v0->GetSecondaryVtx()) {
876 AliDebug(2,"No secondary vertex for V0 by cascade");
880 if (v0->GetNDaughters()!=2) {
881 AliDebug(2,Form("No 2 daughters for V0 of current cascade (onTheFly=%d, nDaughters=%d)",v0->GetOnFlyStatus(),v0->GetNDaughters()));
886 // Get the V0 daughter tracks
887 AliAODTrack *v0positiveTrack = dynamic_cast<AliAODTrack*>(d->Getv0PositiveTrack());
888 AliAODTrack *v0negativeTrack = dynamic_cast<AliAODTrack*>(d->Getv0NegativeTrack());
889 if (!v0positiveTrack || !v0negativeTrack ) {
890 AliDebug(2,"No V0 daughters' objects");
894 if (v0positiveTrack->GetID()<0 || v0negativeTrack->GetID()<0) {
895 AliDebug(2,Form("At least one of V0 daughters has negative ID %d %d",v0positiveTrack->GetID(),v0negativeTrack->GetID()));
899 //if(fUseTrackSelectionWithFilterBits && d->HasBadDaughters()) return 0;
900 if ( fUseTrackSelectionWithFilterBits && !(bachelorTrack->TestFilterMask(BIT(4))) ) {
901 AliDebug(2,"Check on the bachelor FilterBit: no BIT(4). Candidate rejected.");
906 // selection on daughter tracks
907 if (selectionLevel==AliRDHFCuts::kAll ||
908 selectionLevel==AliRDHFCuts::kTracks) {
910 if (!AreLctoV0DaughtersSelected(d)) return 0;
914 Bool_t okLck0sp=kFALSE, okLcLpi=kFALSE, okLcLBarpi=kFALSE;
916 // selection on candidate
917 if (selectionLevel==AliRDHFCuts::kAll ||
918 selectionLevel==AliRDHFCuts::kCandidate) {
920 Double_t pt = d->Pt();
921 Int_t ptbin = PtBin(pt);
923 Double_t mLcPDG = TDatabasePDG::Instance()->GetParticle(4122)->Mass();
924 Double_t mk0sPDG = TDatabasePDG::Instance()->GetParticle(310)->Mass();
925 Double_t mLPDG = TDatabasePDG::Instance()->GetParticle(3122)->Mass();
928 Double_t mk0s = v0->MassK0Short();
929 Double_t mLck0sp = d->InvMassLctoK0sP();
932 Double_t mlambda = v0->MassLambda();
933 Double_t malambda = v0->MassAntiLambda();
934 Double_t mLcLpi = d->InvMassLctoLambdaPi();
938 // cut on Lc mass with K0S+p hypothesis
939 okLck0sp = TMath::Abs(mLck0sp-mLcPDG)<=fCutsRD[GetGlobalIndex(0,ptbin)];
944 // cut on Lc mass with Lambda+pi hypothesis
946 okLcLpi = TMath::Abs(mLcLpi-mLcPDG)<=fCutsRD[GetGlobalIndex(1,ptbin)];
947 okLcLBarpi = okLcLpi;
950 // cuts on the V0 mass: K0S case
951 okLck0sp = TMath::Abs(mk0s-mk0sPDG)<=fCutsRD[GetGlobalIndex(2,ptbin)];
956 // cuts on the V0 mass: Lambda/LambdaBar case
958 okLcLpi = TMath::Abs(mlambda-mLPDG)<=fCutsRD[GetGlobalIndex(3,ptbin)];
959 //okLcLpi = okLcLpi && (bachelorTrack->Charge()==+1);
960 okLcLBarpi = TMath::Abs(malambda-mLPDG)<=fCutsRD[GetGlobalIndex(3,ptbin)];
961 //okLcLBarpi = okLcLBarpi && (bachelorTrack->Charge()==-1);
964 // cuts on the minimum pt of bachelor
965 okLck0sp = TMath::Abs(bachelorTrack->Pt())>=fCutsRD[GetGlobalIndex(4,ptbin)];
967 okLcLBarpi = okLck0sp;
970 // cuts on the minimum pt of positive V0daughter
971 okLck0sp = TMath::Abs(v0positiveTrack->Pt())>=fCutsRD[GetGlobalIndex(5,ptbin)];
973 okLcLBarpi = okLck0sp;
976 // cuts on the minimum pt of negative V0daughter
977 okLck0sp = TMath::Abs(v0negativeTrack->Pt())>=fCutsRD[GetGlobalIndex(6,ptbin)];
979 okLcLBarpi = okLck0sp;
982 // cut on cascade dca
983 okLck0sp = TMath::Abs(d->GetDCA())<=fCutsRD[GetGlobalIndex(7,ptbin)];
985 okLcLBarpi = okLck0sp;
989 okLck0sp = TMath::Abs(v0->GetDCA())<=fCutsRD[GetGlobalIndex(8,ptbin)];
991 okLcLBarpi = okLck0sp;
994 // cut on V0 cosine of pointing angle wrt PV
995 okLck0sp = d->CosV0PointingAngle()>=fCutsRD[GetGlobalIndex(9,ptbin)];
997 okLcLBarpi = okLck0sp;
1000 // cut on bachelor transverse impact parameter wrt PV
1001 okLck0sp = TMath::Abs(d->Getd0Prong(0))<=fCutsRD[GetGlobalIndex(10,ptbin)];
1003 okLcLBarpi = okLck0sp;
1006 // cut on V0 transverse impact parameter wrt PV
1007 okLck0sp = TMath::Abs(d->Getd0Prong(1))<=fCutsRD[GetGlobalIndex(11,ptbin)];
1009 okLcLBarpi = okLck0sp;
1012 // cut on K0S invariant mass veto
1013 okLcLpi = TMath::Abs(mk0s-mk0sPDG)>=fCutsRD[GetGlobalIndex(12,ptbin)];
1014 okLcLBarpi = TMath::Abs(mk0s-mk0sPDG)>=fCutsRD[GetGlobalIndex(12,ptbin)];
1017 // cut on Lambda/LambdaBar invariant mass veto
1018 okLck0sp = (TMath::Abs(mlambda-mLPDG)>=fCutsRD[GetGlobalIndex(13,ptbin)] &&
1019 TMath::Abs(malambda-mLPDG)>=fCutsRD[GetGlobalIndex(13,ptbin)]);
1022 // cut on gamma invariant mass veto
1023 okLck0sp = v0->InvMass2Prongs(0,1,11,11)>=fCutsRD[GetGlobalIndex(14,ptbin)];
1025 okLcLBarpi = okLck0sp;
1029 okLck0sp = v0->Pt()>=fCutsRD[GetGlobalIndex(15,ptbin)];
1031 okLcLBarpi = okLck0sp;
1036 Int_t returnvalue = okLck0sp+2*okLcLBarpi+4*okLcLpi;
1040 2 Lc->LambdaBar + pi
1041 3 Lc->K0S + p AND Lc->LambdaBar + pi
1043 5 Lc->K0S + p AND Lc->Lambda + pi
1044 6 Lc->LambdaBar + pi AND Lc->Lambda + pi
1045 7 Lc->K0S + p AND Lc->LambdaBar + pi AND Lc->Lambda + pi
1050 Int_t returnvaluePID = 7;
1053 if (selectionLevel==AliRDHFCuts::kAll ||
1054 selectionLevel==AliRDHFCuts::kCandidate ||
1055 selectionLevel==AliRDHFCuts::kPID )
1056 returnvaluePID = IsSelectedPID(d);
1059 Int_t returnvalueTot = 0;
1061 //returnvalueTot = CombineCuts(returnvalue,returnvaluePID);
1063 returnvalueTot = returnvalue;
1065 return returnvalueTot;
1068 //----------------------------------
1069 void AliRDHFCutsLctoV0::SetStandardCutsPP2010() {
1071 SetName("LctoV0ProductionCuts");
1072 SetTitle("Production cuts for Lc->V0+bachelor analysis");
1074 AliESDtrackCuts* esdTrackCuts=new AliESDtrackCuts();
1075 esdTrackCuts->SetRequireSigmaToVertex(kFALSE);
1077 esdTrackCuts->SetRequireTPCRefit(kTRUE);
1078 esdTrackCuts->SetRequireITSRefit(kTRUE);
1079 esdTrackCuts->SetMinNClustersITS(0);//(4); // default is 5
1080 esdTrackCuts->SetMinNClustersTPC(70);
1081 //esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
1082 // AliESDtrackCuts::kAny);
1083 // default is kBoth, otherwise kAny
1084 esdTrackCuts->SetMinDCAToVertexXY(0.);
1085 esdTrackCuts->SetPtRange(0.3,1.e10);
1086 //esdTrackCuts->SetEtaRange(-0.8,+0.8);
1087 esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
1088 AddTrackCuts(esdTrackCuts);
1091 AliESDtrackCuts* esdTrackCutsV0daughters=new AliESDtrackCuts();
1092 esdTrackCutsV0daughters->SetRequireSigmaToVertex(kFALSE);
1094 esdTrackCutsV0daughters->SetRequireTPCRefit(kTRUE);
1095 esdTrackCutsV0daughters->SetRequireITSRefit(kFALSE);//(kTRUE);
1096 esdTrackCutsV0daughters->SetMinNClustersITS(0);//(4); // default is 5
1097 esdTrackCutsV0daughters->SetMinNClustersTPC(70);
1098 //esdTrackCutsV0daughters->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
1099 // AliESDtrackCuts::kAny);
1100 // default is kBoth, otherwise kAny
1101 esdTrackCutsV0daughters->SetMinDCAToVertexXY(0.);
1102 esdTrackCutsV0daughters->SetPtRange(0.,1.e10);
1103 esdTrackCutsV0daughters->SetAcceptKinkDaughters(kFALSE);
1104 AddTrackCutsV0daughters(esdTrackCutsV0daughters);
1106 const Int_t nptbins=1;
1108 ptbins=new Float_t[nptbins+1];
1110 ptbins[1]=99999999.;
1112 SetPtBins(nptbins+1,ptbins);
1113 SetPtBins(nptbins+1,ptbins);
1115 const Int_t nvars=17;
1117 Float_t** prodcutsval;
1118 prodcutsval=new Float_t*[nvars];
1119 for(Int_t ic=0;ic<nvars;ic++){prodcutsval[ic]=new Float_t[nptbins];}
1120 for(Int_t ipt2=0;ipt2<nptbins;ipt2++){
1121 prodcutsval[0][ipt2]=1.; // inv. mass if K0S [GeV/c2]
1122 prodcutsval[1][ipt2]=1.; // inv. mass if Lambda [GeV/c2]
1123 prodcutsval[2][ipt2]=0.05; // inv. mass V0 if K0S [GeV/c2]
1124 prodcutsval[3][ipt2]=0.05; // inv. mass V0 if Lambda [GeV/c2]
1125 prodcutsval[4][ipt2]=0.3; // pT min bachelor track [GeV/c] // AOD by construction
1126 prodcutsval[5][ipt2]=0.; // pT min V0-positive track [GeV/c]
1127 prodcutsval[6][ipt2]=0.; // pT min V0-negative track [GeV/c]
1128 prodcutsval[7][ipt2]=1000.; // dca cascade cut [cm]
1129 prodcutsval[8][ipt2]=1.5; // dca V0 cut [nSigma] // it's 1.5 x offline V0s
1130 prodcutsval[9][ipt2]=-1.; // cosPA V0 cut // it's 0.90 x offline V0s at reconstruction level, 0.99 at filtering level
1131 prodcutsval[10][ipt2]=3.; // d0 max bachelor wrt PV [cm]
1132 prodcutsval[11][ipt2]=1000.;// d0 max V0 wrt PV [cm]
1133 prodcutsval[12][ipt2]=0.; // mass K0S veto [GeV/c2]
1134 prodcutsval[13][ipt2]=0.; // mass Lambda/LambdaBar veto [GeV/c2]
1135 prodcutsval[14][ipt2]=0.; // mass Gamma veto [GeV/c2]
1136 prodcutsval[15][ipt2]=0.; // pT min V0 track [GeV/c]
1137 prodcutsval[16][ipt2]=0.; // V0 type cut
1139 SetCuts(nvars,nptbins,prodcutsval);
1141 SetGlobalIndex(nvars,nptbins);
1142 SetPtBins(nptbins+1,ptbins);
1146 //1. bachelor: default one
1147 AliAODPidHF* pidObjBachelor = new AliAODPidHF();
1148 Double_t sigmasBac[5]={3.,1.,1.,3.,3.}; // 0, 1(A), 2(A) -> TPC; 3 -> TOF; 4 -> ITS
1149 pidObjBachelor->SetSigma(sigmasBac);
1150 pidObjBachelor->SetAsym(kFALSE);
1151 pidObjBachelor->SetMatch(1);
1152 pidObjBachelor->SetTPC(kTRUE);
1153 pidObjBachelor->SetTOF(kTRUE);
1154 pidObjBachelor->SetTOFdecide(kFALSE);
1155 SetPidHF(pidObjBachelor);
1158 AliAODPidHF* pidObjV0pos = new AliAODPidHF();
1159 Double_t sigmasV0pos[5]={3.,1.,1.,3.,3.}; // 0, 1(A), 2(A) -> TPC; 3 -> TOF; 4 -> ITS
1160 pidObjV0pos->SetSigma(sigmasV0pos);
1161 pidObjV0pos->SetAsym(kFALSE);
1162 pidObjV0pos->SetMatch(1);
1163 pidObjV0pos->SetTPC(kTRUE);
1164 pidObjV0pos->SetTOF(kTRUE);
1165 pidObjV0pos->SetTOFdecide(kFALSE);
1166 SetPidV0pos(pidObjV0pos);
1169 AliAODPidHF* pidObjV0neg = new AliAODPidHF();
1170 Double_t sigmasV0neg[5]={3.,1.,1.,3.,3.}; // 0, 1(A), 2(A) -> TPC; 3 -> TOF; 4 -> ITS
1171 pidObjV0neg->SetSigma(sigmasV0neg);
1172 pidObjV0neg->SetAsym(kFALSE);
1173 pidObjV0neg->SetMatch(1);
1174 pidObjV0neg->SetTPC(kTRUE);
1175 pidObjV0neg->SetTOF(kTRUE);
1176 pidObjV0neg->SetTOFdecide(kFALSE);
1177 SetPidV0neg(pidObjV0neg);
1179 SetUsePID(kFALSE);//(kTRUE);
1183 for(Int_t iiv=0;iiv<nvars;iiv++){
1184 delete [] prodcutsval[iiv];
1186 delete [] prodcutsval;
1192 delete pidObjBachelor;
1193 pidObjBachelor=NULL;
1203 //------------------
1204 void AliRDHFCutsLctoV0::SetStandardCutsPbPb2010() {
1206 SetName("LctoV0ProductionCuts");
1207 SetTitle("Production cuts for Lc->V0+bachelor analysis");
1209 SetStandardCutsPP2010();
1213 //------------------
1214 void AliRDHFCutsLctoV0::SetStandardCutsPbPb2011() {
1216 // Default 2010 PbPb cut object
1217 SetStandardCutsPbPb2010();
1220 // Enable all 2011 PbPb run triggers
1222 SetTriggerClass("");
1223 ResetMaskAndEnableMBTrigger();
1224 EnableCentralTrigger();
1225 EnableSemiCentralTrigger();
1228 //-----------------------
1229 Int_t AliRDHFCutsLctoV0::GetV0Type(){
1231 const Int_t nvars = this->GetNVars() ;
1232 //Float_t *vArray =GetCuts();
1233 //fV0Type = vArray[nvars-1];
1234 fV0Type = (this->GetCuts())[nvars-1];
1235 //this->GetCuts(vArray);
1236 TString *sVarNames=GetVarNames();
1238 if(sVarNames[nvars-1].Contains("V0 type")) return (Int_t)fV0Type;
1239 else {AliInfo("AliRDHFCutsLctoV0 Last variable is not the V0 type!!!"); return -999;}
1242 //---------------------------------------------------------------------------
1243 void AliRDHFCutsLctoV0::PrintAll() const {
1245 // print all cuts values
1248 printf("Minimum vtx type %d\n",fMinVtxType);
1249 printf("Minimum vtx contr %d\n",fMinVtxContr);
1250 printf("Max vtx red chi2 %f\n",fMaxVtxRedChi2);
1251 printf("Min SPD mult %d\n",fMinSPDMultiplicity);
1252 printf("Use PID %d (PID selection flag = %d) OldPid=%d\n",(Int_t)fUsePID,(Int_t)fPidSelectionFlag,fPidHF ? fPidHF->GetOldPid() : -1);
1253 printf("Remove daughters from vtx %d\n",(Int_t)fRemoveDaughtersFromPrimary);
1254 printf("Recompute primary vertex %d\n",(Int_t)fRecomputePrimVertex);
1255 printf("Physics selection: %s\n",fUsePhysicsSelection ? "Yes" : "No");
1256 printf("Pileup rejection: %s\n",(fOptPileup > 0) ? "Yes" : "No");
1257 printf("UseTrackSelectionWithFilterBits: %s\n",fUseTrackSelectionWithFilterBits ? "Yes" : "No");
1258 printf("Reject kink: %s\n",fKinkReject ? "Yes" : "No");
1259 if(fOptPileup==1) printf(" -- Reject pileup event");
1260 if(fOptPileup==2) printf(" -- Reject tracks from pileup vtx");
1261 if(fUseCentrality>0) {
1262 TString estimator="";
1263 if(fUseCentrality==1) estimator = "V0";
1264 if(fUseCentrality==2) estimator = "Tracks";
1265 if(fUseCentrality==3) estimator = "Tracklets";
1266 if(fUseCentrality==4) estimator = "SPD clusters outer";
1267 printf("Centrality class considered: %.1f-%.1f, estimated with %s",fMinCentrality,fMaxCentrality,estimator.Data());
1269 if(fIsCandTrackSPDFirst) printf("Check for candidates with pt < %2.2f, that daughters fullfill kFirst criteria\n",fMaxPtCandTrackSPDFirst);
1272 cout<<"Array of variables"<<endl;
1273 for(Int_t iv=0;iv<fnVars;iv++){
1274 cout<<fVarNames[iv]<<"\t";
1279 cout<<"Array of optimization"<<endl;
1280 for(Int_t iv=0;iv<fnVars;iv++){
1281 cout<<fVarsForOpt[iv]<<"\t";
1286 cout<<"Array of upper/lower cut"<<endl;
1287 for(Int_t iv=0;iv<fnVars;iv++){
1288 cout<<fIsUpperCut[iv]<<"\t";
1293 cout<<"Array of ptbin limits"<<endl;
1294 for(Int_t ib=0;ib<fnPtBinLimits;ib++){
1295 cout<<fPtBinLimits[ib]<<"\t";
1300 cout<<"Matrix of cuts"<<endl;
1301 for(Int_t iv=0;iv<fnVars;iv++){
1302 for(Int_t ib=0;ib<fnPtBins;ib++){
1303 cout<<"fCutsRD["<<iv<<"]["<<ib<<"] = "<<fCutsRD[GetGlobalIndex(iv,ib)]<<"\t";
1311 Float_t eta1=0, eta2=0; fTrackCuts->GetEtaRange(eta1,eta2);
1312 cout << " etaRange for Bachelor: [" << eta1 << "," << eta2 << "]\n";
1314 if (fV0daughtersCuts) {
1315 Float_t eta3=0, eta4=0; fV0daughtersCuts->GetEtaRange(eta3,eta4);
1316 cout << " etaRange for V0daughters: [" << eta3 << "," << eta4 << "]\n";
1322 //-------------------------
1323 Bool_t AliRDHFCutsLctoV0::IsInFiducialAcceptance(Double_t pt, Double_t y) const
1327 // Checking if Lc is in fiducial acceptance region
1331 if(fMaxRapidityCand>-998.){
1332 if(TMath::Abs(y) > fMaxRapidityCand) return kFALSE;
1337 // applying cut for pt > 5 GeV
1338 AliDebug(2,Form("pt of Lambda_c = %f (> 5), cutting at |y| < 0.8",pt));
1339 if (TMath::Abs(y) > 0.8) return kFALSE;
1342 // appliying smooth cut for pt < 5 GeV
1343 Double_t maxFiducialY = -0.2/15*pt*pt+1.9/15*pt+0.5;
1344 Double_t minFiducialY = 0.2/15*pt*pt-1.9/15*pt-0.5;
1345 AliDebug(2,Form("pt of Lambda_c = %f (< 5), cutting according to the fiducial zone [%f, %f]\n",pt,minFiducialY,maxFiducialY));
1346 if (y < minFiducialY || y > maxFiducialY) return kFALSE;
1351 //---------------------------------------------------------------------------
1352 Bool_t AliRDHFCutsLctoV0::AreLctoV0DaughtersSelected(AliAODRecoDecayHF *dd) const{
1354 // Daughter tracks selection
1357 AliAODRecoCascadeHF* d = (AliAODRecoCascadeHF*)dd;
1359 AliDebug(2,"AliAODRecoCascadeHF null");
1364 AliFatal("Cut object is not defined for bachelor. Candidate accepted.");
1368 AliAODTrack * bachelorTrack = dynamic_cast<AliAODTrack*>(d->GetBachelor());
1369 if (!bachelorTrack) return kFALSE;
1371 if (fIsCandTrackSPDFirst && d->Pt()<fMaxPtCandTrackSPDFirst) {
1372 if(!bachelorTrack->HasPointOnITSLayer(0)) return kFALSE;
1375 if (fKinkReject != (!(fTrackCuts->GetAcceptKinkDaughters())) ) {
1376 AliError(Form("Not compatible setting: fKinkReject=%1d - fTrackCuts->GetAcceptKinkDaughters()=%1d",fKinkReject, fTrackCuts->GetAcceptKinkDaughters()));
1380 AliAODVertex *vAOD = d->GetPrimaryVtx();
1381 Double_t pos[3]; vAOD->GetXYZ(pos);
1382 Double_t cov[6]; vAOD->GetCovarianceMatrix(cov);
1383 const AliESDVertex vESD(pos,cov,100.,100);
1385 if (!IsDaughterSelected(bachelorTrack,&vESD,fTrackCuts)) return kFALSE;
1387 if (!fV0daughtersCuts) {
1388 AliFatal("Cut object is not defined for V0daughters. Candidate accepted.");
1392 AliAODv0 * v0 = dynamic_cast<AliAODv0*>(d->Getv0());
1393 if (!v0) return kFALSE;
1394 AliAODTrack *v0positiveTrack = dynamic_cast<AliAODTrack*>(d->Getv0PositiveTrack());
1395 if (!v0positiveTrack) return kFALSE;
1396 AliAODTrack *v0negativeTrack = dynamic_cast<AliAODTrack*>(d->Getv0NegativeTrack());
1397 if (!v0negativeTrack) return kFALSE;
1400 Float_t etaMin=0, etaMax=0; fV0daughtersCuts->GetEtaRange(etaMin,etaMax);
1401 if ( (v0positiveTrack->Eta()<=etaMin || v0positiveTrack->Eta()>=etaMax) ||
1402 (v0negativeTrack->Eta()<=etaMin || v0negativeTrack->Eta()>=etaMax) ) return kFALSE;
1403 Float_t ptMin=0, ptMax=0; fV0daughtersCuts->GetPtRange(ptMin,ptMax);
1404 if ( (v0positiveTrack->Pt()<=ptMin || v0positiveTrack->Pt()>=ptMax) ||
1405 (v0negativeTrack->Pt()<=ptMin || v0negativeTrack->Pt()>=ptMax) ) return kFALSE;
1407 // Condition on nTPCclusters
1408 if (fV0daughtersCuts->GetMinNClusterTPC()>0) {
1409 if ( ( ( v0positiveTrack->GetTPCClusterInfo(2,1) ) < fV0daughtersCuts->GetMinNClusterTPC() ) ||
1410 ( ( v0negativeTrack->GetTPCClusterInfo(2,1) ) < fV0daughtersCuts->GetMinNClusterTPC() ) ) return kFALSE;
1414 if (v0->GetOnFlyStatus()==kFALSE) { // only for offline V0s
1415 if (fV0daughtersCuts->GetRequireTPCRefit()) {
1416 if( !(v0positiveTrack->GetStatus() & AliESDtrack::kTPCrefit)) return kFALSE;
1417 if( !(v0negativeTrack->GetStatus() & AliESDtrack::kTPCrefit)) return kFALSE;
1421 if (!fV0daughtersCuts->GetAcceptKinkDaughters()) {
1422 AliAODVertex *maybeKinkPos = (AliAODVertex*)v0positiveTrack->GetProdVertex();
1423 AliAODVertex *maybeKinkNeg = (AliAODVertex*)v0negativeTrack->GetProdVertex();
1424 if (maybeKinkPos->GetType()==AliAODVertex::kKink ||
1425 maybeKinkNeg->GetType()==AliAODVertex::kKink) return kFALSE;
1427 // Findable clusters > 0 condition - from V0 analysis
1428 //if( v0positiveTrack->GetTPCNclsF()<=0 || v0negativeTrack->GetTPCNclsF()<=0 ) return kFALSE;
1430 Float_t lPosTrackCrossedRows = v0positiveTrack->GetTPCClusterInfo(2,1);
1431 Float_t lNegTrackCrossedRows = v0positiveTrack->GetTPCClusterInfo(2,1);
1432 fTreeVariableLeastNbrCrossedRows = (Int_t) lPosTrackCrossedRows;
1433 if( lNegTrackCrossedRows < fTreeVariableLeastNbrCrossedRows )
1434 fTreeVariableLeastNbrCrossedRows = (Int_t) lNegTrackCrossedRows;
1435 //Compute ratio Crossed Rows / Findable clusters
1436 //Note: above test avoids division by zero!
1437 Float_t lPosTrackCrossedRowsOverFindable = lPosTrackCrossedRows / ((double)(pTrack->GetTPCNclsF()));
1438 Float_t lNegTrackCrossedRowsOverFindable = lNegTrackCrossedRows / ((double)(nTrack->GetTPCNclsF()));
1439 fTreeVariableLeastRatioCrossedRowsOverFindable = lPosTrackCrossedRowsOverFindable;
1440 if( lNegTrackCrossedRowsOverFindable < fTreeVariableLeastRatioCrossedRowsOverFindable )
1441 fTreeVariableLeastRatioCrossedRowsOverFindable = lNegTrackCrossedRowsOverFindable;
1442 //Lowest Cut Level for Ratio Crossed Rows / Findable = 0.8, set here
1443 if ( fTreeVariableLeastRatioCrossedRowsOverFindable < 0.8 ) return kFALSE;