]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/macros/train/LHC2010-7TeV-phi/AddRsnEfficiency.C
Block of updates on RSN package:
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / macros / train / LHC2010-7TeV-phi / AddRsnEfficiency.C
1 //
2 // This macro add an analysis task for computing efficiency.
3 // It will have as output an AliCFContainer with several steps:
4 //
5 //  0) all resonances in MC which decay in the pair specified
6 //  1) subset of (0) whose daughters are in acceptance
7 //  2) subset of (1) whose daughters satisfy quality track cuts (covariance, chi square && nTPCclusters)
8 //  3) subset of (2) whose daughters satisfy primary track cuts (nsigma to vertex, no kink daughters)
9 //  4) subset of (3) whose daughters satisty the BB TPC compatibility cut at 3 sigma
10 //
11 Bool_t AddRsnEfficiency(const char *dataLabel)
12 {
13   // load useful macros
14   gROOT->LoadMacro("$(ALICE_ROOT)/PWG2/RESONANCES/macros/train/LHC2010-7TeV-phi/ConfigESDCutsITS.C");
15   gROOT->LoadMacro("$(ALICE_ROOT)/PWG2/RESONANCES/macros/train/LHC2010-7TeV-phi/ConfigESDCutsTPC.C");
16   
17   // retrieve analysis manager
18   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
19   
20   // create task
21   AliRsnAnalysisEffSE *task[2];
22   task[0] = new AliRsnAnalysisEffSE("RsnTaskEffNoSA");
23   task[1] = new AliRsnAnalysisEffSE("RsnTaskEffSA");
24
25   // pair definition: 
26   // phi --> K+ K-
27   AliRsnPairDef *pairPhi = new AliRsnPairDef(AliPID::kKaon, '+', AliPID::kKaon, '-', 333, 1.019455);
28   
29   // axis definition
30   // 0) invariant mass
31   // 1) transverse momentum
32   // 2) rapidity
33   // 3) multiplicity
34   Double_t pt  [] = {0.0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6, 1.7, 1.8, 1.9, 2.0, 2.2, 2.4, 2.6, 2.8, 3.0, 3.2, 3.4, 3.6, 3.8, 4.0, 4.5, 5.0, 5.5, 6.0, 7.0, 8.0, 9.0, 10.0};
35   Double_t y   [] = {-1.0, -0.9, -0.8, -0.7, -0.6, -0.5, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0};
36   Double_t mult[] = {0.0, 5.0, 9.0, 14.0, 22.0, 100000.0};
37   Int_t    npt    = sizeof(pt  ) / sizeof(pt  [0]);
38   Int_t    ny     = sizeof(y   ) / sizeof(y   [0]);
39   Int_t    nmult  = sizeof(mult) / sizeof(mult[0]);
40   AliRsnValue *axisIM   = new AliRsnValue("IM"  , AliRsnValue::kPairInvMass     , 500  , 0.9,  1.4);
41   //AliRsnValue *axisPt   = new AliRsnValue("PT"  , AliRsnValue::kPairPt          , npt  , pt);
42   //AliRsnValue *axisY    = new AliRsnValue("Y"   , AliRsnValue::kPairY           , ny   , y);
43   //AliRsnValue *axisMult = new AliRsnValue("Mult", AliRsnValue::kEventMultESDcuts, nmult, mult);
44   AliRsnValue *axisPt   = new AliRsnValue("PT"  , AliRsnValue::kPairPt          , 100,  0.0, 10.0);
45   AliRsnValue *axisY    = new AliRsnValue("Y"   , AliRsnValue::kPairY           ,  20, -1.0,  1.0);
46   ConfigESDCutsTPC(axisMult->GetCuts());
47   
48   // define cuts for event selection:
49   // this will determine the filling of bins in the "info" histograms
50   // and should be computed as additional correction factor in efficiency
51   AliRsnCutPrimaryVertex *cutVertex = new AliRsnCutPrimaryVertex("cutVertex", 10.0, 0, kFALSE);
52   
53   // define standard 2010 track quality/PID cuts:
54   // - first  index: [0] = no PID, [1] = PID
55   // - second index: [0] = no ITS, [1] = ITS
56   AliRsnCutESD2010 *cuts2010[2][2];
57   cuts2010[0][0] = new AliRsnCutESD2010("cutESD2010nopidNoSA");
58   cuts2010[0][1] = new AliRsnCutESD2010("cutESD2010nopidSA");
59   cuts2010[1][0] = new AliRsnCutESD2010("cutESD2010pidNoSA");
60   cuts2010[1][1] = new AliRsnCutESD2010("cutESD2010pidSA");
61   // since both indexes are 0/1, the boolean settings
62   // are done according to them, for clarity
63   for (Int_t ipid = 0; ipid < 2; ipid++)
64   {
65     for (Int_t iits = 0; iits < 2; iits++)
66     {
67       // all work with MC here
68       cuts2010[ipid][iits]->SetMC(kTRUE);
69       
70       // all use global tracks
71       cuts2010[ipid][iits]->SetUseGlobal(kTRUE);
72       
73       // other flags, depend on indexes
74       cuts2010[ipid][iits]->SetUseITSSA((Bool_t)iits);
75       cuts2010[ipid][iits]->SetCheckITS((Bool_t)ipid);
76       cuts2010[ipid][iits]->SetCheckTPC((Bool_t)ipid);
77       cuts2010[ipid][iits]->SetCheckTOF((Bool_t)ipid);
78       
79       // basic quality settings
80       ConfigESDCutsITS(cuts2010[ipid][iits]->GetCutsITS());
81       ConfigESDCutsTPC(cuts2010[ipid][iits]->GetCutsTPC());
82     }
83   }
84
85   // define cut on dip angle:
86   AliRsnCutStd *cutDip = new AliRsnCutStd("cutDip", AliRsnCut::kMother, AliRsnCutStd::kDipAngle, 0.0, 0.04);
87   
88   // define a common path for the output file
89   Char_t commonPath[500];
90   sprintf(commonPath, "%s", AliAnalysisManager::GetCommonFileName());
91   
92   // add all the steps
93   // two-folded loop on the two tasks, where one contains the ITS-SA and the other doesn't
94   for (Int_t itask = 0; itask < 2; itask++)
95   {
96     // add pair definition, to choose the checked resonance
97     task[itask]->AddPairDef(pairPhi);
98
99     // add the output histogram axis
100     task[itask]->AddAxis(axisIM);
101     task[itask]->AddAxis(axisPt);
102     task[itask]->AddAxis(axisY);
103     task[itask]->AddAxis(axisMult);
104     
105     // add the cut only when working on ESD, not on MC only
106     task[itask]->GetEventCuts()->AddCut(cutVertex);
107     task[itask]->GetEventCuts()->SetCutScheme("cutVertex");
108
109     //
110     // *** STEP 0 - All resonances which decay in the specified pair
111     //
112     // This step does not need any kind of definition, since
113     // its requirement is automatically checked during execution,
114     // but to avoid segfaults, it is better to initialize a cut manager.
115     //
116     AliRsnCutManager *mgr_step0 = new AliRsnCutManager("mc_step0", "");
117     
118     //
119     // *** STEP 1 - All resonances which decay into tracked particles
120     //
121     // This step does not need any kind of definition, since
122     // its requirement is automatically checked during execution,
123     // but to avoid segfaults, it is better to initialize a cut manager.
124     //
125     AliRsnCutManager *mgr_step1 = new AliRsnCutManager("reco_step0", "");
126
127     //
128     // *** STEP 2 - Reconstruction & track quality
129     //
130       
131     AliRsnCutSet     *set_step2 = new AliRsnCutSet("cuts_step2", AliRsnCut::kDaughter);
132     AliRsnCutManager *mgr_step2 = new AliRsnCutManager("esd_step2", "");
133     
134     set_step2->AddCut(cuts2010[0][itask]);
135     set_step2->SetCutScheme(cuts2010[0][itask]->GetName());
136     mgr_step2->SetCommonDaughterCuts(set_step2);
137     
138     //
139     // *** STEP 3 - PID
140     //
141     
142     AliRsnCutSet     *set_step3 = new AliRsnCutSet("cuts_step3", AliRsnCut::kDaughter);
143     AliRsnCutManager *mgr_step3 = new AliRsnCutManager("esd_step3", "");
144     
145     set_step3->AddCut(cuts2010[1][itask]);
146     set_step3->SetCutScheme(cuts2010[1][itask]->GetName());
147     mgr_step3->SetCommonDaughterCuts(set_step3);
148
149     //
150     // *** STEP 4 - Dip angle
151     //
152     
153     AliRsnCutSet     *set_step4 = new AliRsnCutSet("cuts_step4", AliRsnCut::kMother);
154     AliRsnCutManager *mgr_step4 = new AliRsnCutManager("esd_step4", "");
155     
156     set_step4->AddCut(cutDip);
157     set_step4->SetCutScheme(Form("!%s", cutDip->GetName()));
158     mgr_step4->SetMotherCuts(set_step4);
159     
160     // add all steps to the task
161     task[itask]->AddStepMC (mgr_step0);
162     task[itask]->AddStepESD(mgr_step1);
163     task[itask]->AddStepESD(mgr_step2);
164     task[itask]->AddStepESD(mgr_step3);
165     task[itask]->AddStepESD(mgr_step4);
166     
167     // add the task to the manager and connect to input
168     mgr->AddTask(task[itask]);
169     mgr->ConnectInput(task[itask], 0, mgr->GetCommonInputContainer());
170     
171     // create paths for the output in the common file
172     TString infoname(task[itask]->GetName());
173     TString histname(task[itask]->GetName());
174     infoname.ReplaceAll("TaskEff", "Info");
175     histname.ReplaceAll("TaskEff", "Hist");
176     AliAnalysisDataContainer *outputInfo = mgr->CreateContainer(infoname.Data(), TList::Class(), AliAnalysisManager::kOutputContainer, commonPath);
177     AliAnalysisDataContainer *outputHist = mgr->CreateContainer(histname.Data(), TList::Class(), AliAnalysisManager::kOutputContainer, commonPath);
178     mgr->ConnectOutput(task[itask], 1, outputInfo);
179     mgr->ConnectOutput(task[itask], 2, outputHist);
180   }
181
182   return kTRUE;
183 }