]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/macros/mini/ConfigLStar.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / macros / mini / ConfigLStar.C
1 /***************************************************************************
2   rama.chandra.baral@cern.ch & sarita.sahoo@cern.ch - last modified on 20/08/2014
3
4 // *** Configuration script for L*, anti-L*, syst. analysis for pp and p-Pb runs ***
5 // 
6 // A configuration script for RSN package needs to define the followings:
7 //
8 // (1) decay tree of each resonance to be studied, which is needed to select
9 //     true pairs and to assign the right mass to all candidate daughters
10 // (2) cuts at all levels: single daughters, tracks, events
11 // (3) output objects: histograms or trees
12 ****************************************************************************/
13
14 Bool_t ConfigLStar
15 (  
16     AliRsnMiniAnalysisTask *task, 
17     Bool_t                 isMC, 
18     Bool_t                 isPP,
19     const char             *suffix,
20     AliRsnCutSet           *cutsPair,
21     Int_t                  aodFilterBit = 5,
22     AliRsnCutSetDaughterParticle::ERsnDaughterCutSet cutPrCandidate = AliRsnCutSetDaughterParticle::kFastTOFpidNsigma,
23     AliRsnCutSetDaughterParticle::ERsnDaughterCutSet cutKaCandidate = AliRsnCutSetDaughterParticle::kFastTOFpidNsigma,
24     Float_t                nsigmaPr = 2.0,
25     Float_t                nsigmaKa = 2.0,
26     Bool_t                 enableTrkSyst = kFALSE,
27     Char_t                 DCAxyFormula[100] = "0.0182+0.035/pt^1.01",
28     Double_t               dcazmax = 3.2,
29     Double_t               minNcls = 70,
30     Double_t               maxX2cls = 4.0,
31     Double_t               minCrossedRows = 50.0,
32     Double_t               maxClsCrossedRows = 0.8,
33
34     Bool_t                 enableMonitor = kTRUE,
35     Bool_t                 IsMcTrueOnly = kFALSE,
36     Int_t                  signedPdg = 3124,
37
38     TString                monitorOpt = "",  //Flag for AddMonitorOutput.C e.g."NoSIGN"
39     Bool_t                 useCrossedRows = kFALSE,
40     AliRsnMiniValue::EType yaxisVar = AliRsnMiniValue::kPt,
41     Bool_t                 useMixLS = 0
42 )
43 {
44   // manage suffix
45   if (strlen(suffix) > 0) suffix = Form("_%s", suffix);
46   
47   // set daughter cuts
48   AliRsnCutSetDaughterParticle * cutSetQ;
49   AliRsnCutSetDaughterParticle * cutSetP;
50   AliRsnCutSetDaughterParticle * cutSetK;
51
52
53   //vary track quality cuts for systematic checks
54   if(enableTrkSyst) {
55     AliRsnCutTrackQuality * trkQualityCut =  new AliRsnCutTrackQuality("QualityCut");
56     
57     //trkQualityCut->DisableAll();//disable all cuts, filter bit, pT, eta, and DCAxy cuts will be reset later
58     trkQualityCut->SetAODTestFilterBit(aodFilterBit);//reset the filter bit cut 
59     trkQualityCut->SetCheckOnlyFilterBit(kFALSE);//tells the cut object to check all other cuts individually,
60     trkQualityCut->SetDCARPtFormula(DCAxyFormula);
61     trkQualityCut->SetDCAZmax(dcazmax);
62
63     if(useCrossedRows) {
64       trkQualityCut->SetMinNCrossedRowsTPC(minCrossedRows, kTRUE);
65       trkQualityCut->SetMinNCrossedRowsOverFindableClsTPC(maxClsCrossedRows, kTRUE); }
66     else trkQualityCut->SetTPCminNClusters(minNcls);
67
68     trkQualityCut->SetTPCmaxChi2(maxX2cls);
69     trkQualityCut->SetRejectKinkDaughters(kTRUE);
70     //trkQualityCut->SetSPDminNClusters(AliESDtrackCuts::kAny);
71     trkQualityCut->SetSPDminNClusters(1);
72     trkQualityCut->SetITSminNClusters(0);
73     trkQualityCut->SetITSmaxChi2(36);
74     trkQualityCut->SetMaxChi2TPCConstrainedGlobal(36);
75     trkQualityCut->AddStatusFlag(AliESDtrack::kTPCin   , kTRUE);//already in defaults 2011
76     trkQualityCut->AddStatusFlag(AliESDtrack::kTPCrefit, kTRUE);//already in defaults 2011
77     trkQualityCut->AddStatusFlag(AliESDtrack::kITSrefit, kTRUE);//already in defaults 2011
78     trkQualityCut->SetPtRange(0.15, 20.0);
79     trkQualityCut->SetEtaRange(-0.8, 0.8);
80     
81     trkQualityCut->Print();
82
83     if(isPP) cutSetQ  = new AliRsnCutSetDaughterParticle(Form("cutQ_bit%i",aodFilterBit), trkQualityCut, AliRsnCutSetDaughterParticle::kQualityStd2010, AliPID::kPion, -1.0);
84     else     cutSetQ  = new AliRsnCutSetDaughterParticle(Form("cutQ_bit%i",aodFilterBit), trkQualityCut, AliRsnCutSetDaughterParticle::kQualityStd2011, AliPID::kPion, -1.0);
85     cutSetP           = new AliRsnCutSetDaughterParticle(Form("cutP%i_%2.1fsigma",cutPrCandidate, nsigmaPr), trkQualityCut, cutPrCandidate, AliPID::kProton, nsigmaPr);
86     cutSetK           = new AliRsnCutSetDaughterParticle(Form("cutK%i_%2.1fsigma",cutKaCandidate, nsigmaKa), trkQualityCut, cutKaCandidate, AliPID::kKaon, nsigmaKa);
87
88   }
89
90   else
91     {
92       //default cuts 2010 for pp and 2011 for p-Pb
93       if(isPP) {
94         cutSetQ  = new AliRsnCutSetDaughterParticle("cutQuality", AliRsnCutSetDaughterParticle::kQualityStd2010, AliPID::kPion, -1.0, aodFilterBit);
95         cutSetP  = new AliRsnCutSetDaughterParticle(Form("cutProton_%2.1fsigma",nsigmaPr), cutPrCandidate, AliPID::kProton, nsigmaPr, aodFilterBit);
96         cutSetK  = new AliRsnCutSetDaughterParticle(Form("cutKaon_%2.1f2sigma",nsigmaKa), cutKaCandidate, AliPID::kKaon, nsigmaKa, aodFilterBit);
97       }
98       else {
99         cutSetQ  = new AliRsnCutSetDaughterParticle("cutQuality", AliRsnCutSetDaughterParticle::kQualityStd2011, AliPID::kPion, -1.0, aodFilterBit);
100         cutSetQ->SetUse2011StdQualityCuts(kTRUE);
101         cutSetP = new AliRsnCutSetDaughterParticle(Form("cutProton2011_%2.1fsigma",nsigmaPr), cutPrCandidate, AliPID::kProton, nsigmaPr, aodFilterBit);
102         cutSetP->SetUse2011StdQualityCuts(kTRUE);
103         cutSetK  = new AliRsnCutSetDaughterParticle(Form("cutKaon2011_%2.1f2sigma",nsigmaKa), cutKaCandidate, AliPID::kKaon, nsigmaKa, aodFilterBit);
104         cutSetK->SetUse2011StdQualityCuts(kTRUE);
105       }
106     }
107   
108   Int_t iCutQ = task->AddTrackCuts(cutSetQ);
109   Int_t iCutP = task->AddTrackCuts(cutSetP);
110   Int_t iCutK = task->AddTrackCuts(cutSetK);
111   
112   if (enableMonitor){
113     Printf("======== Cut monitoring enabled");
114     gROOT->LoadMacro("$ALICE_ROOT/PWGLF/RESONANCES/macros/mini/AddMonitorOutput.C");
115     AddMonitorOutput(isMC, cutSetQ->GetMonitorOutput(), monitorOpt.Data());
116     AddMonitorOutput(isMC, cutSetP->GetMonitorOutput(), monitorOpt.Data());
117     AddMonitorOutput(isMC, cutSetK->GetMonitorOutput(), monitorOpt.Data());
118   }  
119   
120   // -- Values ------------------------------------------------------------------------------------
121   /* invariant mass   */ Int_t imID   = task->CreateValue(AliRsnMiniValue::kInvMass, kFALSE);
122   /* IM resolution    */ Int_t resID  = task->CreateValue(AliRsnMiniValue::kInvMassRes, kTRUE);
123   /* transv. momentum */ Int_t ptID   = task->CreateValue(AliRsnMiniValue::kPt, kFALSE);
124   /* centrality       */ Int_t centID = task->CreateValue(AliRsnMiniValue::kMult, kFALSE);
125   /* pseudorapidity   */ Int_t etaID  = task->CreateValue(AliRsnMiniValue::kEta, kFALSE);
126   /* rapidity         */ Int_t yID    = task->CreateValue(AliRsnMiniValue::kY, kFALSE);
127   /* 1st daughter pt  */ Int_t fdpt   = task->CreateValue(AliRsnMiniValue::kFirstDaughterPt, kFALSE);
128   /* 2nd daughter pt  */ Int_t sdpt   = task->CreateValue(AliRsnMiniValue::kSecondDaughterPt, kFALSE);
129   /* 1st daughter p   */ Int_t fdp    = task->CreateValue(AliRsnMiniValue::kFirstDaughterP, kFALSE);
130   /* 2nd daughter p   */ Int_t sdp    = task->CreateValue(AliRsnMiniValue::kSecondDaughterP, kFALSE);
131
132   // -- Create all needed outputs -----------------------------------------------------------------
133   // use an array for more compact writing, which are different on mixing and charges
134   // [0] = unlike
135   // [1] = mixing
136   // [2] = like ++
137   // [3] = like --
138
139   Bool_t  use     [12] = { !IsMcTrueOnly, !IsMcTrueOnly, !IsMcTrueOnly,!IsMcTrueOnly,!IsMcTrueOnly,!IsMcTrueOnly, isMC    ,  isMC   , isMC   ,  isMC  , useMixLS , useMixLS};
140   Bool_t  useIM   [12] = { 1            ,  1           ,  1           ,  1          ,  1          , 1           ,  1      ,   1     , 0      ,   0    , 1        , 1       };
141   TString name    [12] = {"UnlikePM"    , "UnlikeMP"   , "MixingPM"   , "MixingMP"  , "LikePP"    , "LikeMM"    ,"TruesPM","TruesMP","ResPM" ,"ResMP" ,"MixingPP","MixingMM"};
142   TString comp    [12] = {"PAIR"        , "PAIR"       , "MIX"        , "MIX"       , "PAIR"      , "PAIR"      , "TRUE"  , "TRUE"  , "TRUE" , "TRUE" , "MIX"    ,"MIX"   };
143   TString output  [12] = {"SPARSE"      , "SPARSE"     , "SPARSE"     , "SPARSE"    , "SPARSE"    , "SPARSE"    , "SPARSE","SPARSE" ,"SPARSE","SPARSE","SPARSE"  ,"SPARSE"};
144   Char_t  charge1 [12] = {'+'           , '-'          , '+'          , '-'         , '+'         , '-'         , '+'     ,  '-'    , '+'    ,  '-'   , '+'      , '-'    };
145   Char_t  charge2 [12] = {'-'           , '+'          , '-'          , '+'         , '+'         , '-'         , '-'     ,  '+'    , '-'    ,  '+'   ,'+'       , '-'    };
146   Int_t   cutID1  [12] = { iCutP        ,  iCutP       ,  iCutP       ,  iCutP      ,  iCutP      ,  iCutP      ,  iCutP  ,  iCutP  ,  iCutP ,  iCutP , iCutP    , iCutP };
147   Int_t   cutID2  [12] = { iCutK        ,  iCutK       ,  iCutK       ,  iCutK      ,  iCutK      ,  iCutK      ,  iCutK  ,  iCutK  ,  iCutK ,  iCutK , iCutK    , iCutK };
148   
149   //TString output  [10] = {"HIST"   , "HIST"   , "HIST"   , "HIST"   , "HIST"  , "HIST"  , "HIST"  ,  "HIST"  , "HIST"  ,  "HIST"  };
150
151   for (Int_t i = 0; i < 12; i++) {
152     if (!use[i]) continue;
153     AliRsnMiniOutput *out = task->CreateOutput(Form("Lstar_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
154     out->SetCutID(0, cutID1[i]);
155     out->SetCutID(1, cutID2[i]);
156     out->SetDaughter(0, AliRsnDaughter::kProton);
157     out->SetDaughter(1, AliRsnDaughter::kKaon);
158     out->SetCharge(0, charge1[i]);
159     out->SetCharge(1, charge2[i]);
160     out->SetMotherPDG(signedPdg);
161     out->SetMotherMass(1.51953);
162     out->SetPairCuts(cutsPair);
163
164     // axis X: invmass (or resolution)
165     if (useIM[i]) 
166       out->AddAxis(imID, 800, 1.4, 2.2);
167     else
168       out->AddAxis(resID, 200, -0.02, 0.02);
169     
170     // axis Y: transverse momentum of pair as default - else chosen value
171     if (yaxisVar==AliRsnMiniValue::kFirstDaughterPt) {
172       out->AddAxis(fdpt, 100, 0.0, 10.0);
173       out->AddAxis(sdpt, 100, 0.0, 10.0);  }
174     else
175       if (yaxisVar==AliRsnMiniValue::kFirstDaughterP) {
176         out->AddAxis(fdp, 100, 0.0, 10.0);
177         out->AddAxis(sdp, 100, 0.0, 10.0);  }
178       else 
179         out->AddAxis(ptID, 100, 0.0, 10.0); //default use mother pt
180
181     // axis Z: centrality-multiplicity
182     if (!isPP)
183       out->AddAxis(centID, 100, 0.0, 100.0);
184     else 
185       out->AddAxis(centID, 400, 0.0, 400.0);
186     
187     // axis W: pseudorapidity
188     // out->AddAxis(etaID, 20, -1.0, 1.0);
189     // axis J: rapidity
190     // out->AddAxis(yID, 10, -0.5, 0.5);
191     
192   }   
193   
194   if (isMC){   
195     // create output
196     AliRsnMiniOutput *outm = task->CreateOutput(Form("Lstar_Mother%s", suffix), "SPARSE", "MOTHER");
197     outm->SetDaughter(0, AliRsnDaughter::kProton);
198     outm->SetDaughter(1, AliRsnDaughter::kKaon);
199     outm->SetMotherPDG(signedPdg);
200     outm->SetMotherMass(1.51953);
201     // pair cuts
202     outm->SetPairCuts(cutsPair);
203     // binnings
204     outm->AddAxis(imID, 800, 1.4, 2.2);
205
206     // axis Y: transverse momentum of pair as default - else chosen value
207     if (yaxisVar==AliRsnMiniValue::kFirstDaughterPt) {
208       outm->AddAxis(fdpt, 100, 0.0, 10.0);
209       outm->AddAxis(sdpt, 100, 0.0, 10.0);  }
210     else
211       if (yaxisVar==AliRsnMiniValue::kFirstDaughterP) {
212         outm->AddAxis(fdp, 100, 0.0, 10.0);
213         outm->AddAxis(sdp, 100, 0.0, 10.0);  }
214       else 
215         outm->AddAxis(ptID, 100, 0.0, 10.0); //default use mother pt
216
217
218     if (!isPP){
219       outm->AddAxis(centID, 100, 0.0, 100.0);
220     }   else    { 
221       outm->AddAxis(centID, 400, 0.0, 400.0);
222     }
223   }
224   return kTRUE;
225 }