]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/macros/train/LHC2010-7TeV-phi/fixed_rapidity/AddRsnEfficiency.C
2ca86d38cc927ceabc3a0800dd7a0f2c1a84b684
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / macros / train / LHC2010-7TeV-phi / fixed_rapidity / 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, const char *path)
12 {
13   // load useful macros
14   gROOT->LoadMacro(Form("%s/QualityCutsITS.C", path));
15   gROOT->LoadMacro(Form("%s/QualityCutsTPC.C", path));
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   task[0]->SelectCollisionCandidates();
25   task[1]->SelectCollisionCandidates();
26
27   // pair definition: 
28   // phi --> K+ K-
29   AliRsnPairDef *pairPhi = new AliRsnPairDef(AliPID::kKaon, '+', AliPID::kKaon, '-', 333, 1.019455);
30   
31   // axis definition
32   // 0) invariant mass
33   // 1) transverse momentum
34   // 2) rapidity
35   // 3) multiplicity
36   Double_t     mult[]   = {0., 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12., 13., 14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25., 30., 35., 40., 50., 60., 70., 80., 90., 100., 120., 140., 160., 180., 200., 1E+8};
37   Int_t        nmult    = sizeof(mult) / sizeof(mult[0]);
38   AliRsnValue *axisIM   = new AliRsnValue("IM"  , AliRsnValue::kPairInvMass     , 0.9, 1.4, 0.001);
39   AliRsnValue *axisPt   = new AliRsnValue("PT"  , AliRsnValue::kPairPt          , 0.0, 5.0, 0.100);
40   AliRsnValue *axisY    = new AliRsnValue("Y"   , AliRsnValue::kPairY           ,-0.5, 0.5, 0.500);
41   AliRsnValue *axisMult = new AliRsnValue("Mult", AliRsnValue::kEventMultESDCuts, nmult, mult);
42   
43   // initialize the support object: AliESDtrackCuts
44   // configured using the standard values
45   AliESDtrackCuts *cuts = new AliESDtrackCuts(QualityCutsTPC());
46   axisMult->SetSupportObject(cuts);
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   cutVertex->SetCheckPileUp(kTRUE);
53   
54   // define standard 2010 track quality/PID cuts:
55   // - first  index: [0] = no PID, [1] = PID
56   // - second index: [0] = no ITS, [1] = ITS
57   AliRsnCutESD2010 *cuts2010[2][2];
58   cuts2010[0][0] = new AliRsnCutESD2010("cutESD2010nopidNoSA");
59   cuts2010[0][1] = new AliRsnCutESD2010("cutESD2010nopidSA");
60   cuts2010[1][0] = new AliRsnCutESD2010("cutESD2010pidNoSA");
61   cuts2010[1][1] = new AliRsnCutESD2010("cutESD2010pidSA");
62   // define Bethe-Bloch parameters (only for MC, since this computes efficiency)
63   Double_t bbPar[5] = {2.15898 / 50.0, 1.75295E1, 3.40030E-9, 1.96178, 3.91720};
64   // since both indexes are 0/1, the boolean settings are done according to them, for clarity
65   for (Int_t ipid = 0; ipid < 2; ipid++)
66   {
67     for (Int_t iits = 0; iits < 2; iits++)
68     {
69       // all work with MC here
70       cuts2010[ipid][iits]->SetMC(kTRUE);
71       
72       // PID reference is kaons
73       cuts2010[ipid][iits]->SetPID(AliPID::kKaon);
74       
75       // all use global tracks
76       cuts2010[ipid][iits]->SetUseITSTPC(kTRUE);
77       
78       // other flags, depend on indexes
79       cuts2010[ipid][iits]->SetUseITSSA((Bool_t)iits);
80       cuts2010[ipid][iits]->SetCheckITS((Bool_t)ipid);
81       cuts2010[ipid][iits]->SetCheckTPC((Bool_t)ipid);
82       cuts2010[ipid][iits]->SetCheckTOF((Bool_t)ipid);
83       
84       // basic quality settings
85       cuts2010[ipid][iits]->CopyCutsTPC(QualityCutsTPC());
86       cuts2010[ipid][iits]->CopyCutsITS(QualityCutsITS());
87       
88       // (unused for No PID) set the ITS PID-related variables
89       cuts2010[ipid][iits]->SetITSband(3.0);
90       
91       // (unused for No PID) set the TPC PID-related variables
92       cuts2010[ipid][iits]->SetTPCrange(3.0, 5.0);
93       cuts2010[ipid][iits]->SetTPCpLimit(0.35);
94       cuts2010[ipid][iits]->GetESDpid()->GetTPCResponse().SetBetheBlochParameters(bbPar[0], bbPar[1], bbPar[2], bbPar[3], bbPar[4]);
95       
96       // (unused for No PID) set the TOF PID-related variables
97       cuts2010[ipid][iits]->SetTOFrange(-3.0, 3.0);
98     }
99   }
100
101   // define cut on dip angle:
102   AliRsnCutValue *cutDip = new AliRsnCutValue("cutDip", AliRsnValue::kPairDipAngle, 0.02, 1.01);
103   
104   // define a common path for the output file
105   Char_t commonPath[500];
106   sprintf(commonPath, "%s", AliAnalysisManager::GetCommonFileName());
107   
108   // add all the steps
109   // two-folded loop on the two tasks, where one contains the ITS-SA and the other doesn't
110   for (Int_t itask = 0; itask < 2; itask++)
111   {
112     // add pair definition, to choose the checked resonance
113     task[itask]->AddPairDef(pairPhi);
114
115     // add the output histogram axis
116     task[itask]->AddAxis(axisIM);
117     task[itask]->AddAxis(axisPt);
118     task[itask]->AddAxis(axisY);
119     task[itask]->AddAxis(axisMult);
120     
121     // add the cut on primary vertex
122     task[itask]->GetEventCuts()->AddCut(cutVertex);
123     task[itask]->GetEventCuts()->SetCutScheme(cutVertex->GetName());
124
125     //
126     // *** STEP 0 - All resonances which decay in the specified pair
127     //
128     // This step does not need any kind of definition, since
129     // its requirement is automatically checked during execution,
130     // but to avoid segfaults, it is better to initialize a cut manager.
131     //
132     AliRsnCutManager *mgr_step0 = new AliRsnCutManager("mc_step0", "");
133     
134     //
135     // *** STEP 1 - All resonances which decay into tracked particles
136     //
137     // This step does not need any kind of definition, since
138     // its requirement is automatically checked during execution,
139     // but to avoid segfaults, it is better to initialize a cut manager.
140     //
141     AliRsnCutManager *mgr_step1 = new AliRsnCutManager("esd_step0", "");
142
143     //
144     // *** STEP 2 - Reconstruction & track quality
145     //
146     // Define a cut on track quality, disabling the PID cuts (first index = [0])
147     //
148     AliRsnCutManager *mgr_step2 = new AliRsnCutManager("esd_step2", "");
149     AliRsnCutSet     *set_step2 = mgr_step2->GetCommonDaughterCuts();
150     
151     set_step2->AddCut(cuts2010[0][itask]);
152     set_step2->SetCutScheme(cuts2010[0][itask]->GetName());
153     
154     //
155     // *** STEP 3 - PID
156     //
157     // Define a cut on track quality, enabling the PID cuts (first index = [1])
158     //
159     AliRsnCutManager *mgr_step3 = new AliRsnCutManager("esd_step3", "");
160     AliRsnCutSet     *set_step3 = mgr_step3->GetCommonDaughterCuts();
161     
162     set_step3->AddCut(cuts2010[1][itask]);
163     set_step3->SetCutScheme(cuts2010[1][itask]->GetName());
164
165     //
166     // *** STEP 4 - Dip angle
167     //
168     // Add a cut on the pair dip angle
169     //
170     AliRsnCutManager *mgr_step4 = new AliRsnCutManager("esd_step4", "");
171     AliRsnCutSet     *set_step4 = mgr_step4->GetMotherCuts();
172     
173     set_step4->AddCut(cutDip);
174     set_step4->SetCutScheme(Form("%s", cutDip->GetName()));
175     
176     // add all steps to the task:
177     // - first step computed on MC
178     // - all other steps computed on reconstruction
179     task[itask]->AddStepMC (mgr_step0);
180     task[itask]->AddStepESD(mgr_step1);
181     task[itask]->AddStepESD(mgr_step2);
182     task[itask]->AddStepESD(mgr_step3);
183     task[itask]->AddStepESD(mgr_step4);
184     
185     // add the task to the manager and connect to input
186     mgr->AddTask(task[itask]);
187     mgr->ConnectInput(task[itask], 0, mgr->GetCommonInputContainer());
188     
189     // create paths for the output in the common file
190     TString infoname(task[itask]->GetName());
191     TString histname(task[itask]->GetName());
192     infoname.ReplaceAll("TaskEff", "Info");
193     histname.ReplaceAll("TaskEff", "Hist");
194     AliAnalysisDataContainer *outputInfo = mgr->CreateContainer(infoname.Data(), TList::Class(), AliAnalysisManager::kOutputContainer, commonPath);
195     AliAnalysisDataContainer *outputHist = mgr->CreateContainer(histname.Data(), TList::Class(), AliAnalysisManager::kOutputContainer, commonPath);
196     mgr->ConnectOutput(task[itask], 1, outputInfo);
197     mgr->ConnectOutput(task[itask], 2, outputHist);
198   }
199
200   return kTRUE;
201 }