]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/macros/train/LHC2010-7TeV-phi/RsnConfig.C
Update of analysis macros
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / macros / train / LHC2010-7TeV-phi / RsnConfig.C
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,
36   const char *path,
37   Int_t       multMin = 0,
38   Int_t       multMax = 0
39 )
40 {
41   // load useful macros
42   gROOT->LoadMacro(Form("%s/QualityCutsITS.C", path));
43   gROOT->LoadMacro(Form("%s/QualityCutsTPC.C", path));
44   
45   // interpret the useful information from second argument
46   TString opt(options);
47   Bool_t isSim  = opt.Contains("sim");
48   Bool_t isData = opt.Contains("data");
49   if (!isSim && !isData)
50   {
51     Error("RsnConfig", "Required to know if working on data or MC");
52     return kFALSE;
53   }
54   
55   // interpret the specific info from third argument
56   // which should be fixed in the various calls to this function
57   TString conf(config);
58   Bool_t addPID    = conf.Contains("pid");
59   Bool_t addITSSA  = conf.Contains("its");
60   Bool_t addDipCut = conf.Contains("dip");
61       
62   // generate a common suffix depending on chosen options
63   TString suffix;
64   if (addPID)    suffix += "_pid";
65   if (addITSSA)  suffix += "_its";
66   if (addDipCut) suffix += "_dip";
67   Info("RsnConfig", "=== Specific configuration: %s ====================================================", config);
68   Info("RsnConfig", "=== suffix used           : %s ====================================================", suffix.Data());
69
70   // retrieve analysis manager & task
71   AliAnalysisManager *mgr  = AliAnalysisManager::GetAnalysisManager();
72   AliRsnAnalysisSE   *task = (AliRsnAnalysisSE*)mgr->GetTask(taskName);
73
74   // for safety, return if no task is passed
75   if (!task)
76   {
77     Error("RsnConfig2010PhiFcn", "Task not found");
78     return kFALSE;
79   }
80   
81   //
82   // -- Setup event cuts (added directly to task) ---------------------------------------------------
83   //
84   
85   // define a common cut on primary vertex, which also checks pile-up
86   AliRsnCutPrimaryVertex *cutVertex  = new AliRsnCutPrimaryVertex("cutVertex", 10.0, 0, kFALSE);
87   cutVertex->SetCheckPileUp(kTRUE);
88   task->GetEventCuts()->AddCut(cutVertex);
89   
90   // if at least one of last two arguments is not zero, 
91   // add a multiplicity cut using those arguments as range limits
92   if (multMin > 0 || multMax > 0)
93   {
94     ::Info("RsnConfig.C", "Adding multiplicity cut: %d --> %d", multMin, multMax);
95     AliRsnCutValue *cutMult = new AliRsnCutValue(Form("cutMult_%d-%d", multMin, multMax), AliRsnValue::kEventMultESDCuts, (Double_t)multMin, (Double_t)multMax);
96     
97     // initialize the support object: AliESDtrackCuts
98     // configured using the standard values
99     AliESDtrackCuts *cuts = new AliESDtrackCuts(QualityCutsTPC());
100     cutMult->GetValueObj()->SetSupportObject(cuts);
101     
102     // add the cut and set the cut string
103     task->GetEventCuts()->AddCut(cutMult);
104     task->GetEventCuts()->SetCutScheme(Form("cutVertex&%s", cutMult->GetName()));
105   }
106   else
107   {
108     // if no mult cut is added, only primary vertex cut is used
109     task->GetEventCuts()->SetCutScheme("cutVertex");
110   }
111
112   //
113   // -- Setup pairs ---------------------------------------------------------------------------------
114   //
115
116   // decay channels
117   AliRsnPairDef *pairDefPM = new AliRsnPairDef(AliPID::kKaon, '+', AliPID::kKaon, '-', 333, 1.019455);
118   AliRsnPairDef *pairDefPP = new AliRsnPairDef(AliPID::kKaon, '+', AliPID::kKaon, '+', 333, 1.019455);
119   AliRsnPairDef *pairDefMM = new AliRsnPairDef(AliPID::kKaon, '-', AliPID::kKaon, '-', 333, 1.019455);
120
121   // computation objects
122   AliRsnPairFunctions *pairPM = new AliRsnPairFunctions(Form("PairPM%s", suffix.Data()), pairDefPM);
123   AliRsnPairFunctions *truePM = new AliRsnPairFunctions(Form("TruePM%s", suffix.Data()), pairDefPM);
124   AliRsnPairFunctions *pairPP = new AliRsnPairFunctions(Form("PairPP%s", suffix.Data()), pairDefPP);
125   AliRsnPairFunctions *pairMM = new AliRsnPairFunctions(Form("PairMM%s", suffix.Data()), pairDefMM);
126
127   //
128   // -- Setup cuts ----------------------------------------------------------------------------------
129   //
130
131   // track cut -----------------------
132   // --> global cuts for 2010 analysis
133   // --> most options are set to right values by default
134   // --> second argument in constructor tells if we are working in simulation or not
135   AliRsnCutESD2010 *cuts2010 = new AliRsnCutESD2010(Form("cuts2010%s", suffix.Data()), isSim);
136   // --> set the reference particle for PID
137   cuts2010->SetPID(AliPID::kKaon);
138   // --> include or not the ITS standalone (TPC is always in)
139   cuts2010->SetUseITSTPC(kTRUE);
140   cuts2010->SetUseITSSA (addITSSA);
141   // --> set the quality cuts using the general macro and using the 'Copy()' method in AliESDtrackCuts
142   cuts2010->CopyCutsTPC(QualityCutsTPC());
143   cuts2010->CopyCutsITS(QualityCutsITS());
144   // --> set values for PID flags, depending on the choice expressed in the options
145   cuts2010->SetCheckITS (addPID);
146   cuts2010->SetCheckTPC (addPID);
147   cuts2010->SetCheckTOF (addPID);
148   // --> set the ITS PID-related variables
149   cuts2010->SetITSband(3.0);
150   // --> set the TPC PID-related variables
151   Double_t bbPar[5];
152   if (isSim)
153   {
154     bbPar[0] = 2.15898 / 50.0;
155     bbPar[1] = 1.75295E1;
156     bbPar[2] = 3.40030E-9;
157     bbPar[3] = 1.96178;
158     bbPar[4] = 3.91720;
159   }
160   else
161   {
162     bbPar[0] = 1.41543 / 50.0;
163     bbPar[1] = 2.63394E1;
164     bbPar[2] = 5.0411E-11;
165     bbPar[3] = 2.12543;
166     bbPar[4] = 4.88663;
167   }
168   cuts2010->SetTPCrange(3.0, 5.0);
169   cuts2010->SetTPCpLimit(0.35);
170   cuts2010->GetESDpid()->GetTPCResponse().SetBetheBlochParameters(bbPar[0], bbPar[1], bbPar[2], bbPar[3], bbPar[4]);
171   // --> set the TOF PID-related variables
172   cuts2010->SetTOFrange(-3.0, 3.0);
173   
174   // pair cut ----------------------------------------
175   // --> dip angle between daughters: (it is a cosine)
176   AliRsnCutValue *cutDip = new AliRsnCutValue("cutDip", AliRsnValue::kPairDipAngle, 0.02, 1.01);
177
178   // setup cut set for tracks------------------------------------------------------------
179   // --> in this case, only common cuts are applied, depending if working with ESD or AOD
180   // --> these cuts are added always
181   pairPM->GetCutManager()->GetCommonDaughterCuts()->AddCut(cuts2010);
182   truePM->GetCutManager()->GetCommonDaughterCuts()->AddCut(cuts2010);
183   pairPP->GetCutManager()->GetCommonDaughterCuts()->AddCut(cuts2010);
184   pairMM->GetCutManager()->GetCommonDaughterCuts()->AddCut(cuts2010);
185   
186   pairPM->GetCutManager()->GetCommonDaughterCuts()->SetCutScheme(cuts2010->GetName());
187   truePM->GetCutManager()->GetCommonDaughterCuts()->SetCutScheme(cuts2010->GetName());
188   pairPP->GetCutManager()->GetCommonDaughterCuts()->SetCutScheme(cuts2010->GetName());
189   pairMM->GetCutManager()->GetCommonDaughterCuts()->SetCutScheme(cuts2010->GetName());
190   
191   // cut set for pairs---------------------
192   // --> add dip angle cut only if required
193   if (addDipCut)
194   {
195     pairPM->GetCutManager()->GetMotherCuts()->AddCut(cutDip);
196     truePM->GetCutManager()->GetMotherCuts()->AddCut(cutDip);
197     pairPP->GetCutManager()->GetMotherCuts()->AddCut(cutDip);
198     pairMM->GetCutManager()->GetMotherCuts()->AddCut(cutDip);
199     
200     pairPM->GetCutManager()->GetMotherCuts()->SetCutScheme(cutDip->GetName());
201     truePM->GetCutManager()->GetMotherCuts()->SetCutScheme(cutDip->GetName());
202     pairPP->GetCutManager()->GetMotherCuts()->SetCutScheme(cutDip->GetName());
203     pairMM->GetCutManager()->GetMotherCuts()->SetCutScheme(cutDip->GetName());
204   }
205   
206   // set additional option for true pairs
207   truePM->SetOnlyTrue  (kTRUE);
208   truePM->SetCheckDecay(kTRUE);
209
210   //
211   // -- Setup functions -----------------------------------------------------------------------------
212   //
213
214   // axis definition
215   // 0) invariant mass
216   // 1) transverse momentum
217   // 2) rapidity
218   AliRsnValue *axisIM = new AliRsnValue("IM", AliRsnValue::kPairInvMass,  0.9,  1.3, 0.001);
219   AliRsnValue *axisPt = new AliRsnValue("PT", AliRsnValue::kPairPt     ,  0.0, 10.0, 0.100);
220   AliRsnValue *axisY  = new AliRsnValue("Y" , AliRsnValue::kPairY      , -1.1,  1.1, 0.100);
221
222   // create function and add axes
223   AliRsnFunction *fcn = new AliRsnFunction;
224   if ( !fcn->AddAxis(axisIM  ) ) return kFALSE;
225   if ( !fcn->AddAxis(axisPt  ) ) return kFALSE;
226   if ( !fcn->AddAxis(axisY   ) ) return kFALSE;
227
228   // add functions to pairs
229   pairPM->AddFunction(fcn);
230   truePM->AddFunction(fcn);
231   pairPP->AddFunction(fcn);
232   pairMM->AddFunction(fcn);
233
234   //
235   // -- Conclusion ----------------------------------------------------------------------------------
236   //
237
238   // add all created AliRsnPair objects to the AliRsnAnalysisManager in the task
239   task->GetAnalysisManager()->Add(pairPM);
240   task->GetAnalysisManager()->Add(pairPP);
241   task->GetAnalysisManager()->Add(pairMM);
242   if (isSim) task->GetAnalysisManager()->Add(truePM);
243
244   return kTRUE;
245 }