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