]>
Commit | Line | Data |
---|---|---|
b0d3d6b9 | 1 | /* |
2 | #include <TROOT.h> | |
3 | #include <TString.h> | |
4 | #include <AliAnalysisManager.h> | |
5 | #include <AliRsnAnalysisSE.h> | |
6 | #include <AliRsnCutESD2010.h> | |
7 | #include <AliRsnCutValue.h> | |
8 | #include <AliRsnPairFunctions.h> | |
9 | #include <AliRsnFunction.h> | |
10 | #include <AliRsnCutPrimaryVertex.h> | |
11 | ||
12 | #include "config/QualityCutsITS.C" | |
13 | #include "config/QualityCutsTPC.C" | |
14 | */ | |
15 | ||
16 | // | |
17 | // This function configures the entire task for all resonances the user is interested in. | |
18 | // This is done by creating all configuration objects which are defined in the package. | |
19 | // | |
20 | // Generally speaking, one has to define the following objects for each resonance: | |
21 | // | |
22 | // 1 - an AliRsnPairDef to define the resonance decay channel to be studied | |
23 | // 2 - an AliRsnPair{Ntuple|Functions} where the output is stored | |
24 | // 3 - one or more AliRsnCut objects to define track selections | |
25 | // which will have then to be organized into AliRsnCutSet objects | |
26 | // 4 - an AliRsnCutManager to include all cuts to be applied (see point 3) | |
27 | // 5 - definitions to build the TNtuple or histograms which are returned | |
28 | // | |
29 | // The return value is used to know if the configuration was successful | |
30 | // | |
31 | Bool_t RsnConfig | |
32 | ( | |
33 | const char *taskName, | |
34 | const char *options, | |
35 | const char *config, | |
f073c807 | 36 | const char *path |
b0d3d6b9 | 37 | ) |
38 | { | |
39 | // load useful macros | |
40 | gROOT->LoadMacro(Form("%s/QualityCutsITS.C", path)); | |
41 | gROOT->LoadMacro(Form("%s/QualityCutsTPC.C", path)); | |
42 | ||
43 | // interpret the useful information from second argument | |
44 | TString opt(options); | |
45 | Bool_t isSim = opt.Contains("sim"); | |
46 | Bool_t isData = opt.Contains("data"); | |
47 | if (!isSim && !isData) | |
48 | { | |
49 | Error("RsnConfig", "Required to know if working on data or MC"); | |
50 | return kFALSE; | |
51 | } | |
52 | ||
53 | // interpret the specific info from third argument | |
54 | // which should be fixed in the various calls to this function | |
55 | TString conf(config); | |
56 | Bool_t addPID = conf.Contains("pid"); | |
57 | Bool_t addITSSA = conf.Contains("its"); | |
58 | Bool_t addDipCut = conf.Contains("dip"); | |
59 | ||
60 | // generate a common suffix depending on chosen options | |
61 | TString suffix; | |
62 | if (addPID) suffix += "_pid"; | |
63 | if (addITSSA) suffix += "_its"; | |
64 | if (addDipCut) suffix += "_dip"; | |
65 | Info("RsnConfig", "=== Specific configuration: %s ====================================================", config); | |
66 | Info("RsnConfig", "=== suffix used : %s ====================================================", suffix.Data()); | |
67 | ||
68 | // retrieve analysis manager & task | |
69 | AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager(); | |
70 | AliRsnAnalysisSE *task = (AliRsnAnalysisSE*)mgr->GetTask(taskName); | |
71 | ||
72 | // for safety, return if no task is passed | |
73 | if (!task) | |
74 | { | |
75 | Error("RsnConfig2010PhiFcn", "Task not found"); | |
76 | return kFALSE; | |
77 | } | |
78 | ||
79 | // | |
80 | // -- Setup event cuts (added directly to task) --------------------------------------------------- | |
81 | // | |
82 | ||
83 | // define a common cut on primary vertex, which also checks pile-up | |
84 | AliRsnCutPrimaryVertex *cutVertex = new AliRsnCutPrimaryVertex("cutVertex", 10.0, 0, kFALSE); | |
85 | cutVertex->SetCheckPileUp(kTRUE); | |
86 | task->GetEventCuts()->AddCut(cutVertex); | |
87 | task->GetEventCuts()->SetCutScheme("cutVertex"); | |
88 | ||
89 | // | |
90 | // -- Setup pairs --------------------------------------------------------------------------------- | |
91 | // | |
92 | ||
93 | // decay channels | |
94 | AliRsnPairDef *pairDefPM = new AliRsnPairDef(AliPID::kKaon, '+', AliPID::kKaon, '-', 333, 1.019455); | |
95 | AliRsnPairDef *pairDefPP = new AliRsnPairDef(AliPID::kKaon, '+', AliPID::kKaon, '+', 333, 1.019455); | |
96 | AliRsnPairDef *pairDefMM = new AliRsnPairDef(AliPID::kKaon, '-', AliPID::kKaon, '-', 333, 1.019455); | |
97 | ||
98 | // computation objects | |
99 | AliRsnPairFunctions *pairPM = new AliRsnPairFunctions(Form("PairPM%s", suffix.Data()), pairDefPM); | |
100 | AliRsnPairFunctions *truePM = new AliRsnPairFunctions(Form("TruePM%s", suffix.Data()), pairDefPM); | |
101 | AliRsnPairFunctions *pairPP = new AliRsnPairFunctions(Form("PairPP%s", suffix.Data()), pairDefPP); | |
102 | AliRsnPairFunctions *pairMM = new AliRsnPairFunctions(Form("PairMM%s", suffix.Data()), pairDefMM); | |
103 | ||
104 | // | |
105 | // -- Setup cuts ---------------------------------------------------------------------------------- | |
106 | // | |
107 | ||
108 | // track cut ----------------------- | |
109 | // --> global cuts for 2010 analysis | |
110 | // --> most options are set to right values by default | |
111 | // --> second argument in constructor tells if we are working in simulation or not | |
112 | AliRsnCutESD2010 *cuts2010 = new AliRsnCutESD2010(Form("cuts2010%s", suffix.Data()), isSim); | |
113 | // --> set the reference particle for PID | |
114 | cuts2010->SetPID(AliPID::kKaon); | |
115 | // --> include or not the ITS standalone (TPC is always in) | |
116 | cuts2010->SetUseITSTPC(kTRUE); | |
117 | cuts2010->SetUseITSSA (addITSSA); | |
118 | // --> set the quality cuts using the general macro and using the 'Copy()' method in AliESDtrackCuts | |
119 | cuts2010->CopyCutsTPC(QualityCutsTPC()); | |
120 | cuts2010->CopyCutsITS(QualityCutsITS()); | |
121 | // --> set values for PID flags, depending on the choice expressed in the options | |
122 | cuts2010->SetCheckITS (addPID); | |
123 | cuts2010->SetCheckTPC (addPID); | |
124 | cuts2010->SetCheckTOF (addPID); | |
125 | // --> set the ITS PID-related variables | |
126 | cuts2010->SetITSband(3.0); | |
127 | // --> set the TPC PID-related variables | |
128 | Double_t bbPar[5]; | |
129 | if (isSim) | |
130 | { | |
131 | bbPar[0] = 2.15898 / 50.0; | |
132 | bbPar[1] = 1.75295E1; | |
133 | bbPar[2] = 3.40030E-9; | |
134 | bbPar[3] = 1.96178; | |
135 | bbPar[4] = 3.91720; | |
136 | } | |
137 | else | |
138 | { | |
139 | bbPar[0] = 1.41543 / 50.0; | |
140 | bbPar[1] = 2.63394E1; | |
141 | bbPar[2] = 5.0411E-11; | |
142 | bbPar[3] = 2.12543; | |
143 | bbPar[4] = 4.88663; | |
144 | } | |
145 | cuts2010->SetTPCrange(3.0, 5.0); | |
146 | cuts2010->SetTPCpLimit(0.35); | |
147 | cuts2010->GetESDpid()->GetTPCResponse().SetBetheBlochParameters(bbPar[0], bbPar[1], bbPar[2], bbPar[3], bbPar[4]); | |
148 | // --> set the TOF PID-related variables | |
149 | cuts2010->SetTOFrange(-3.0, 3.0); | |
150 | ||
151 | // pair cut ---------------------------------------- | |
152 | // --> dip angle between daughters: (it is a cosine) | |
b0d3d6b9 | 153 | AliRsnCutValue *cutY = new AliRsnCutValue("cutY" , AliRsnValue::kPairY , -0.5 , 0.5 ); |
154 | cutY->GetValueObj()->SetSupportObject(pairDefPM); | |
155 | ||
156 | // setup cut set for tracks------------------------------------------------------------ | |
157 | // --> in this case, only common cuts are applied, depending if working with ESD or AOD | |
158 | // --> these cuts are added always | |
159 | pairPM->GetCutManager()->GetCommonDaughterCuts()->AddCut(cuts2010); | |
160 | truePM->GetCutManager()->GetCommonDaughterCuts()->AddCut(cuts2010); | |
161 | pairPP->GetCutManager()->GetCommonDaughterCuts()->AddCut(cuts2010); | |
162 | pairMM->GetCutManager()->GetCommonDaughterCuts()->AddCut(cuts2010); | |
163 | ||
164 | pairPM->GetCutManager()->GetCommonDaughterCuts()->SetCutScheme(cuts2010->GetName()); | |
165 | truePM->GetCutManager()->GetCommonDaughterCuts()->SetCutScheme(cuts2010->GetName()); | |
166 | pairPP->GetCutManager()->GetCommonDaughterCuts()->SetCutScheme(cuts2010->GetName()); | |
167 | pairMM->GetCutManager()->GetCommonDaughterCuts()->SetCutScheme(cuts2010->GetName()); | |
168 | ||
169 | // setup cut set for pairs--------------- | |
170 | // --> add rapidity range cut | |
171 | TString scheme(cutY->GetName()); | |
172 | pairPM->GetCutManager()->GetMotherCuts()->AddCut(cutY); | |
173 | truePM->GetCutManager()->GetMotherCuts()->AddCut(cutY); | |
174 | pairPP->GetCutManager()->GetMotherCuts()->AddCut(cutY); | |
175 | pairMM->GetCutManager()->GetMotherCuts()->AddCut(cutY); | |
f073c807 | 176 | |
177 | pairPM->GetCutManager()->GetMotherCuts()->SetCutScheme(cutY->GetName()); | |
178 | truePM->GetCutManager()->GetMotherCuts()->SetCutScheme(cutY->GetName()); | |
179 | pairPP->GetCutManager()->GetMotherCuts()->SetCutScheme(cutY->GetName()); | |
180 | pairMM->GetCutManager()->GetMotherCuts()->SetCutScheme(cutY->GetName()); | |
b0d3d6b9 | 181 | |
182 | // set additional option for true pairs | |
183 | truePM->SetOnlyTrue (kTRUE); | |
184 | truePM->SetCheckDecay(kTRUE); | |
185 | ||
186 | // | |
187 | // -- Setup functions ----------------------------------------------------------------------------- | |
188 | // | |
189 | ||
190 | // axis definition | |
191 | // 0) invariant mass | |
192 | // 1) transverse momentum | |
193 | // 2) multiplicity | |
194 | 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}; | |
195 | Int_t nmult = sizeof(mult) / sizeof(mult[0]); | |
196 | AliRsnValue *axisIM = new AliRsnValue("IM", AliRsnValue::kPairInvMass , 0.9, 1.4, 0.001); | |
197 | AliRsnValue *axisPt = new AliRsnValue("PT", AliRsnValue::kPairPt , 0.0, 5.0, 0.100); | |
198 | AliRsnValue *axisMult = new AliRsnValue("M" , AliRsnValue::kEventMultESDCuts, nmult, mult); | |
199 | ||
200 | // initialize the support object: AliESDtrackCuts | |
201 | // configured using the standard values | |
202 | AliESDtrackCuts *cuts = new AliESDtrackCuts(QualityCutsTPC()); | |
203 | axisMult->SetSupportObject(cuts); | |
204 | ||
205 | // create function and add axes | |
206 | AliRsnFunction *fcn = new AliRsnFunction; | |
207 | if ( !fcn->AddAxis(axisIM ) ) return kFALSE; | |
208 | if ( !fcn->AddAxis(axisPt ) ) return kFALSE; | |
209 | if ( !fcn->AddAxis(axisMult) ) return kFALSE; | |
210 | ||
211 | // add functions to pairs | |
212 | pairPM->AddFunction(fcn); | |
213 | truePM->AddFunction(fcn); | |
214 | pairPP->AddFunction(fcn); | |
215 | pairMM->AddFunction(fcn); | |
216 | ||
217 | // | |
218 | // -- Conclusion ---------------------------------------------------------------------------------- | |
219 | // | |
220 | ||
221 | // add all created AliRsnPair objects to the AliRsnAnalysisManager in the task | |
222 | task->GetAnalysisManager()->Add(pairPM); | |
223 | task->GetAnalysisManager()->Add(pairPP); | |
224 | task->GetAnalysisManager()->Add(pairMM); | |
225 | if (isSim) task->GetAnalysisManager()->Add(truePM); | |
226 | ||
227 | return kTRUE; | |
228 | } |