]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/macros/mini/ConfigKStarPPb.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / macros / mini / ConfigKStarPPb.C
1 /***************************************************************************
2               fbellini@cern.ch - last modified on 17/02/2014
3
4 // *** Configuration script for K*, anti-K* analysis with 2013 pPb 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 Bool_t ConfigKStarPPb
14 (  
15     AliRsnMiniAnalysisTask *task, 
16     Bool_t                 isMC, 
17     Bool_t                 isPP,
18     const char             *suffix,
19     AliRsnCutSet           *cutsPair,
20     Int_t                  aodFilterBit = 5,
21     Int_t                  customQualityCutsID = AliRsnCutSetDaughterParticle::kDisableCustom,
22     AliRsnCutSetDaughterParticle::ERsnDaughterCutSet cutPiCandidate = AliRsnCutSetDaughterParticle::kTOFpidKstarPbPb2010,
23     AliRsnCutSetDaughterParticle::ERsnDaughterCutSet cutKaCandidate = AliRsnCutSetDaughterParticle::kTOFpidKstarPbPb2010,
24     Float_t                nsigmaPi = 2.0,
25     Float_t                nsigmaKa = 2.0,
26     Bool_t                 enableMonitor = kTRUE,
27     Bool_t                 IsMcTrueOnly = kFALSE,
28     TString                monitorOpt = "",
29     Bool_t                 useMixLS = 0,
30     Bool_t                 checkReflex = 0,
31     AliRsnMiniValue::EType yaxisVar = AliRsnMiniValue::kPt
32 )
33 {
34   // manage suffix
35   if (strlen(suffix) > 0) suffix = Form("_%s", suffix);
36   
37   // set daughter cuts
38   AliRsnCutSetDaughterParticle * cutSetQ;
39   AliRsnCutSetDaughterParticle * cutSetPi;
40   AliRsnCutSetDaughterParticle * cutSetK;
41   
42   AliRsnCutTrackQuality * trkQualityCut =  new AliRsnCutTrackQuality("myQualityCut");
43   if (SetCustomQualityCut(trkQualityCut, customQualityCutsID, aodFilterBit)) {
44     //Set custom quality cuts for systematic checks
45     cutSetQ  = new AliRsnCutSetDaughterParticle(Form("cutQ_bit%i",aodFilterBit), trkQualityCut, AliRsnCutSetDaughterParticle::kQualityStd2011, AliPID::kPion, -1.0);
46     cutSetPi = new AliRsnCutSetDaughterParticle(Form("cutPi%i_%2.1fsigma",cutPiCandidate, nsigmaPi), trkQualityCut, cutPiCandidate, AliPID::kPion, nsigmaPi);
47     cutSetK  = new AliRsnCutSetDaughterParticle(Form("cutK%i_%2.1fsigma",cutPiCandidate, nsigmaKa), trkQualityCut, cutKaCandidate, AliPID::kKaon, nsigmaKa);
48   } else {
49     //use defult quality cuts (std 2010 or 2011)
50     cutSetQ  = new AliRsnCutSetDaughterParticle(Form("cutQ_bit%i",aodFilterBit), AliRsnCutSetDaughterParticle::kQualityStd2011, AliPID::kPion, -1.0, aodFilterBit, kTRUE);
51     cutSetQ->SetUse2011StdQualityCuts(kTRUE);
52     cutSetPi = new AliRsnCutSetDaughterParticle(Form("cutPi%i_%2.1fsigma",cutPiCandidate, nsigmaPi), cutPiCandidate, AliPID::kPion, nsigmaPi, aodFilterBit,kTRUE);
53     cutSetPi->SetUse2011StdQualityCuts(kTRUE);
54     cutSetK  = new AliRsnCutSetDaughterParticle(Form("cutK%i_%2.1fsigma",cutPiCandidate, nsigmaKa), cutKaCandidate, AliPID::kKaon, nsigmaKa, aodFilterBit,kTRUE);
55     cutSetK->SetUse2011StdQualityCuts(kTRUE);
56    }
57   
58   Int_t iCutQ = task->AddTrackCuts(cutSetQ);
59   Int_t iCutPi = task->AddTrackCuts(cutSetPi);
60   Int_t iCutK = task->AddTrackCuts(cutSetK);
61   
62   if (enableMonitor){
63     Printf("======== Cut monitoring enabled");
64     gROOT->LoadMacro("$ALICE_ROOT/PWGLF/RESONANCES/macros/mini/AddMonitorOutput.C");
65     AddMonitorOutput(isMC, cutSetQ->GetMonitorOutput(), monitorOpt.Data());
66     AddMonitorOutput(isMC, cutSetPi->GetMonitorOutput(), monitorOpt.Data());
67     AddMonitorOutput(isMC, cutSetK->GetMonitorOutput()), monitorOpt.Data();
68   }  
69   
70   // -- Values ------------------------------------------------------------------------------------
71   /* invariant mass   */ Int_t imID   = task->CreateValue(AliRsnMiniValue::kInvMass, kFALSE);
72   /* IM resolution    */ Int_t resID  = task->CreateValue(AliRsnMiniValue::kInvMassRes, kTRUE);
73   /* transv. momentum */ Int_t ptID   = task->CreateValue(AliRsnMiniValue::kPt, kFALSE);
74   /* centrality       */ Int_t centID = task->CreateValue(AliRsnMiniValue::kMult, kFALSE);
75   /* pseudorapidity   */ Int_t etaID  = task->CreateValue(AliRsnMiniValue::kEta, kFALSE);
76   /* rapidity         */ Int_t yID    = task->CreateValue(AliRsnMiniValue::kY, kFALSE);
77   /* 1st daughter pt  */ Int_t fdpt   = task->CreateValue(AliRsnMiniValue::kFirstDaughterPt, kFALSE);
78   /* 2nd daughter pt  */ Int_t sdpt   = task->CreateValue(AliRsnMiniValue::kSecondDaughterPt, kFALSE);
79   /* 1st daughter p   */ Int_t fdp    = task->CreateValue(AliRsnMiniValue::kFirstDaughterP, kFALSE);
80   /* 2nd daughter p   */ Int_t sdp    = task->CreateValue(AliRsnMiniValue::kSecondDaughterP, kFALSE);
81   
82   // -- Create all needed outputs -----------------------------------------------------------------
83   // use an array for more compact writing, which are different on mixing and charges
84   // [0] = unlike
85   // [1] = mixing
86   // [2] = like ++
87   // [3] = like --
88
89   Bool_t  use     [12] = {!IsMcTrueOnly,!IsMcTrueOnly,!IsMcTrueOnly,!IsMcTrueOnly,!IsMcTrueOnly,!IsMcTrueOnly ,isMC,isMC,isMC,isMC, useMixLS , useMixLS    };
90   Bool_t  useIM   [12] = { 1        , 1         ,  1       ,  1       ,  1      ,  1      ,  1      ,  1      ,  0      , 0      , 1         , 1           };
91   TString name    [12] = {"UnlikePM", "UnlikeMP","MixingPM","MixingMP", "LikePP", "LikeMM","TruesPM","TruesMP", "ResPM" ,"ResMP" , "MixingPP", "MixingMM"  };
92   TString comp    [12] = {"PAIR"    , "PAIR"    , "MIX"    , "MIX"    , "PAIR"  , "PAIR"  , "TRUE"  , "TRUE"  , "TRUE"  ,"TRUE"  , "MIX"     , "MIX"       };
93   TString output  [12] = {"SPARSE"  , "SPARSE"  , "SPARSE" , "SPARSE" , "SPARSE", "SPARSE", "SPARSE","SPARSE" , "SPARSE","SPARSE", "SPARSE"  , "SPARSE"    };
94   Int_t   pdgCode [12] = {313       , 313       , 313      , 313      , 313     , 313     , 313     , -313    ,  313    , -313   ,  313      , 313         };
95   Char_t  charge1 [12] = {'+'       , '-'       , '+'      , '-'      , '+'     , '-'     , '+'     ,  '-'    , '+'     , '-'    , '+'       , '-'         };
96   Char_t  charge2 [12] = {'-'       , '+'       , '-'      , '+'      , '+'     , '-'     , '-'     ,  '+'    , '-'     , '+'    , '+'       , '-'         };
97   Int_t   cutID1  [12] = { iCutK    ,  iCutK    ,  iCutK   ,  iCutK   ,  iCutK  ,  iCutK  ,  iCutK  ,   iCutK ,  iCutK  , iCutK  , iCutK     , iCutK       };
98   Int_t   cutID2  [12] = { iCutPi   ,  iCutPi   ,  iCutPi  ,  iCutPi  ,  iCutPi ,  iCutPi ,  iCutPi ,   iCutPi,  iCutPi , iCutPi , iCutPi    , iCutPi      };
99   
100   for (Int_t i = 0; i < 12; i++) {
101     if (!use[i]) continue;
102     AliRsnMiniOutput *out = task->CreateOutput(Form("kstar_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
103     out->SetCutID(0, cutID1[i]);
104     out->SetCutID(1, cutID2[i]);
105     out->SetDaughter(0, AliRsnDaughter::kKaon);
106     out->SetDaughter(1, AliRsnDaughter::kPion);
107     out->SetCharge(0, charge1[i]);
108     out->SetCharge(1, charge2[i]);
109     out->SetMotherPDG(pdgCode[i]);
110     out->SetMotherMass(0.89594);
111     out->SetPairCuts(cutsPair);
112
113     // axis X: invmass (or resolution)
114     if (useIM[i]) 
115       out->AddAxis(imID, 90, 0.6, 1.5);
116     else
117       out->AddAxis(resID, 200, -0.02, 0.02);
118     
119     // axis Y: transverse momentum of pair as default - else chosen value
120     if (yaxisVar==AliRsnMiniValue::kFirstDaughterPt)
121       out->AddAxis(fdpt, 100, 0.0, 10.0);
122     else
123       if (yaxisVar==AliRsnMiniValue::kSecondDaughterPt)
124         out->AddAxis(sdpt, 100, 0.0, 10.0);
125       else
126         if (yaxisVar==AliRsnMiniValue::kFirstDaughterP)
127           out->AddAxis(fdp, 100, 0.0, 10.0);
128         else
129           if (yaxisVar==AliRsnMiniValue::kSecondDaughterP)
130             out->AddAxis(sdp, 100, 0.0, 10.0);
131           else 
132             out->AddAxis(ptID, 200, 0.0, 20.0); //default use mother pt
133     
134     // axis Z: centrality-multiplicity
135     if (!isPP)
136       out->AddAxis(centID, 100, 0.0, 100.0);
137     else 
138       out->AddAxis(centID, 400, 0.0, 400.0);
139     // axis W: pseudorapidity
140     // out->AddAxis(etaID, 20, -1.0, 1.0);
141     // axis J: rapidity
142     // out->AddAxis(yID, 10, -0.5, 0.5);
143   }   
144   
145   if (isMC){   
146     //get mothers for K* PDG = 313
147     AliRsnMiniOutput *outm = task->CreateOutput(Form("Ks_Mother%s", suffix), "SPARSE", "MOTHER");
148     outm->SetDaughter(0, AliRsnDaughter::kKaon);
149     outm->SetDaughter(1, AliRsnDaughter::kPion);
150     outm->SetMotherPDG(313);
151     outm->SetMotherMass(0.89594);
152     outm->SetPairCuts(cutsPair);
153     outm->AddAxis(imID, 90, 0.6, 1.5);
154     outm->AddAxis(ptID, 200, 0.0, 20.0);
155     if (!isPP){
156       outm->AddAxis(centID, 100, 0.0, 100.0);
157     }   else    { 
158       outm->AddAxis(centID, 400, 0.0, 400.0);
159     }
160     
161     //get mothers for antiK* PDG = -313
162     AliRsnMiniOutput *outam = task->CreateOutput(Form("antiKs_Mother%s", suffix), "SPARSE", "MOTHER");
163     outam->SetDaughter(0, AliRsnDaughter::kKaon);
164     outam->SetDaughter(1, AliRsnDaughter::kPion);
165     outam->SetMotherPDG(-313);
166     outam->SetMotherMass(0.89594);
167     outam->SetPairCuts(cutsPair);
168     outam->AddAxis(imID, 90, 0.6, 1.5);
169     outam->AddAxis(ptID, 200, 0.0, 20.0);
170     if (!isPP){
171       outam->AddAxis(centID, 100, 0.0, 100.0);
172     }   else    { 
173       outam->AddAxis(centID, 400, 0.0, 400.0);
174     }
175     
176     //get phase space of the decay from mothers
177     AliRsnMiniOutput *outps = task->CreateOutput(Form("Ks_phaseSpace%s", suffix), "HIST", "TRUE");
178     outps->SetDaughter(0, AliRsnDaughter::kKaon);
179     outps->SetDaughter(1, AliRsnDaughter::kPion);
180     outps->SetCutID(0, iCutK);
181     outps->SetCutID(1, iCutPi);
182     outps->SetMotherPDG(313);
183     outps->SetMotherMass(0.89594);
184     outps->SetPairCuts(cutsPair);
185     outps->AddAxis(fdpt, 50, 0.0, 5.0);
186     outps->AddAxis(sdpt, 50, 0.0, 5.0);
187     outps->AddAxis(ptID, 100, 0.0, 10.0);
188     
189     AliRsnMiniOutput *outaps = task->CreateOutput(Form("antiKs_phaseSpace%s", suffix), "HIST", "TRUE");
190     outaps->SetDaughter(0, AliRsnDaughter::kKaon);
191     outaps->SetDaughter(1, AliRsnDaughter::kPion);
192     outaps->SetCutID(0, iCutK);
193     outaps->SetCutID(1, iCutPi);
194     outaps->SetMotherPDG(-313);
195     outaps->SetMotherMass(0.89594);
196     outaps->SetPairCuts(cutsPair);
197     outaps->AddAxis(fdpt, 50, 0.0, 5.0);
198     outaps->AddAxis(sdpt, 50, 0.0, 5.0);
199     outaps->AddAxis(ptID, 100, 0.0, 10.0);
200    
201     //get reflections
202     if (checkReflex) { 
203
204       AliRsnMiniOutput *outreflex = task->CreateOutput(Form("Ks_reflex%s", suffix), "SPARSE", "TRUE");
205       outreflex->SetDaughter(0, AliRsnDaughter::kKaon);
206       outreflex->SetDaughter(1, AliRsnDaughter::kPion);
207       outreflex->SetCutID(0, iCutPi);
208       outreflex->SetCutID(1, iCutK);
209       outreflex->SetMotherPDG(313);
210       outreflex->SetMotherMass(0.89594);
211       outreflex->SetPairCuts(cutsPair);
212       outreflex->AddAxis(imID, 90, 0.6, 1.5);
213       outreflex->AddAxis(ptID, 200, 0.0, 20.0);
214       if (!isPP){
215         outreflex->AddAxis(centID, 100, 0.0, 100.0);
216       }   else    { 
217         outreflex->AddAxis(centID, 400, 0.0, 400.0);
218       }
219       
220       AliRsnMiniOutput *outareflex = task->CreateOutput(Form("antiKs_reflex%s", suffix), "SPARSE", "TRUE");
221       outareflex->SetDaughter(0, AliRsnDaughter::kKaon);
222       outareflex->SetDaughter(1, AliRsnDaughter::kPion);
223       outareflex->SetCutID(0, iCutPi);
224       outareflex->SetCutID(1, iCutK);
225       outareflex->SetMotherPDG(-313);
226       outareflex->SetMotherMass(0.89594);
227       outareflex->SetPairCuts(cutsPair);
228       outareflex->AddAxis(imID, 90, 0.6, 1.5);
229       outareflex->AddAxis(ptID, 100, 0.0, 10.0);
230       if (!isPP){
231         outareflex->AddAxis(centID, 100, 0.0, 100.0);
232       }   else    { 
233         outareflex->AddAxis(centID, 400, 0.0, 400.0);
234       }
235
236     }//end reflections
237   }//end MC
238   
239    return kTRUE;
240 }
241
242 //-------------------------------------------------------  
243 Bool_t SetCustomQualityCut(AliRsnCutTrackQuality * trkQualityCut, Int_t customQualityCutsID = 0, Int_t customFilterBit = 0)
244 {
245   //Sets configuration for track quality object different from std quality cuts.
246   //Returns kTRUE if track quality cut object is successfully defined,
247   //returns kFALSE if an invalid set of cuts (customQualityCutsID) is chosen or if the
248   //object to be configured does not exist.
249   
250   /* NOTES FROM PRODUCTION LHC13b pass3 - AOD filtered with v5-03-Rev-20
251   //(http://svnweb.cern.ch/world/wsvn/AliRoot/tags/v5-03-Rev-20/ANALYSIS/macros/AddTaskESDFilter.C)
252
253   //filter bit 0: Cuts on primary tracks
254   // AliESDtrackCuts* esdTrackCutsL = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
255
256   //filter bit 4: std but looser dca cut
257   // AliESDtrackCuts* esdTrackCutsH = AliESDtrackCuts::GetStandardITSTPCTrackCuts2011(kFALSE);
258   // esdTrackCutsH->SetMaxDCAToVertexXY(2.4);
259   // esdTrackCutsH->SetMaxDCAToVertexZ(3.2);
260   // esdTrackCutsH->SetDCAToVertex2D(kTRUE);
261
262   //filter bit 5:  AliESDtrackCuts::GetStandardITSTPCTrackCuts2011();
263
264    //filter bit 10: standard cuts with tight DCA cut, using cluster cut instead of crossed rows (a la 2010 default)
265    //AliESDtrackCuts* esdTrackCutsH2Cluster = AliESDtrackCuts::GetStandardITSTPCTrackCuts2011(kTRUE, 0);
266    */
267
268   if ((!trkQualityCut) || (customQualityCutsID<=0) || (customQualityCutsID>=AliRsnCutSetDaughterParticle::kNcustomQualityCuts)){
269     Printf("::::: SetCustomQualityCut:: use default quality cuts specified in task configuration.");
270     return kFALSE;
271   }
272   //for pA 2013
273   //trkQualityCut->SetDefaults2011();//with filter bit=10
274   //reset filter bit to very loose cuts 
275   trkQualityCut->SetAODTestFilterBit(customFilterBit); 
276   //apply all other cuts "by hand"
277   trkQualityCut->SetCheckOnlyFilterBit(kFALSE);
278   trkQualityCut->SetMinNCrossedRowsTPC(70, kTRUE);
279   trkQualityCut->SetMinNCrossedRowsOverFindableClsTPC(0.8, kTRUE);
280   trkQualityCut->SetMaxChi2TPCConstrainedGlobal(36);//used for ESD only - for AOD does not correspond to any cut
281   trkQualityCut->SetTPCmaxChi2(4.0); //already in filter bit 0
282   trkQualityCut->SetRejectKinkDaughters(kTRUE); //already in filter bit 0
283   trkQualityCut->SetSPDminNClusters(AliESDtrackCuts::kAny);
284   trkQualityCut->SetITSmaxChi2(36);
285   trkQualityCut->AddStatusFlag(AliESDtrack::kTPCin   , kTRUE);//already in defaults 2011
286   trkQualityCut->AddStatusFlag(AliESDtrack::kTPCrefit, kTRUE);//already in defaults 2011
287   trkQualityCut->AddStatusFlag(AliESDtrack::kITSrefit, kTRUE);//already in defaults 2011
288
289   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kFilterBitCustom) {
290     trkQualityCut->SetCheckOnlyFilterBit(kTRUE);
291   } 
292   
293   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kStdLooserDCAXY){
294     trkQualityCut->SetDCARmax(2.4);
295   } else {
296     trkQualityCut->SetDCARPtFormula("0.0105+0.0350/pt^1.1");
297   }
298   
299   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kStdLooserDCAZ){
300     trkQualityCut->SetDCAZmax(3.2);
301   } else {
302     trkQualityCut->SetDCAZmax(2.0); 
303   }
304   
305   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kStdCrossedRows60){
306     trkQualityCut->SetMinNCrossedRowsTPC(60, kTRUE);
307   }
308   
309   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kStdCrossedRows80){
310     trkQualityCut->SetMinNCrossedRowsTPC(80, kTRUE);
311   }
312   
313   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kStdRowsToCls075){
314     trkQualityCut->SetMinNCrossedRowsOverFindableClsTPC(0.75, kTRUE);
315   }
316   
317   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kStdRowsToCls085){
318     trkQualityCut->SetMinNCrossedRowsOverFindableClsTPC(0.85, kTRUE);
319   }
320   
321   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kStdCls70){
322     trkQualityCut->SetAODTestFilterBit(10);
323     trkQualityCut->SetTPCminNClusters(70);
324   }
325   
326   if (customQualityCutsID==AliRsnCutSetDaughterParticle::kStdChi2TPCCls35){
327     trkQualityCut->SetTPCmaxChi2(3.5);
328   }
329   
330   trkQualityCut->SetPtRange(0.15, 20.0);
331   trkQualityCut->SetEtaRange(-0.8, 0.8);
332   
333   Printf(Form("::::: SetCustomQualityCut:: using custom track quality cuts #%i",customQualityCutsID));
334   trkQualityCut->Print();
335   return kTRUE;
336 }