]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/macros/train/LHC2010-7TeV-phi/RsnConfig.C
Added cut for checking multiplicity and done some adaptments for AOD analysis with...
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / macros / train / LHC2010-7TeV-phi / RsnConfig.C
1 //
2 // This function configures the entire task for all resonances the user is interested in.
3 // This is done by creating all configuration objects which are defined in the package.
4 //
5 // Generally speaking, one has to define the following objects for each resonance:
6 //
7 //  1 - an AliRsnPairDef to define the resonance decay channel to be studied
8 //  2 - an AliRsnPair{Ntuple|Functions} where the output is stored
9 //  3 - one or more AliRsnCut objects to define track selections
10 //      which will have then to be organized into AliRsnCutSet objects
11 //  4 - an AliRsnCutManager to include all cuts to be applied (see point 3)
12 //  5 - definitions to build the TNtuple or histograms which are returned
13 //
14 // The return value is used to know if the configuration was successful
15 //
16 Bool_t RsnConfig
17 (
18   const char *taskName, 
19   const char *options,
20   const char *config
21 )
22 {
23   // load useful macros
24   gROOT->LoadMacro("$(ALICE_ROOT)/PWG2/RESONANCES/macros/train/LHC2010-7TeV-phi/ConfigESDCutsITS.C");
25   gROOT->LoadMacro("$(ALICE_ROOT)/PWG2/RESONANCES/macros/train/LHC2010-7TeV-phi/ConfigESDCutsTPC.C");
26   
27   // interpret the useful information from second argument
28   TString opt(options);
29   Bool_t isSim  = opt.Contains("sim");
30   Bool_t isData = opt.Contains("data");
31   Bool_t isESD  = opt.Contains("ESD");
32   Bool_t isAOD  = opt.Contains("AOD");
33   
34   // interpret the specific info from third argument
35   // which should be fixed in the various calls to this function
36   TString opt(config);
37   Bool_t realPID     = opt.Contains("realistic");
38   Bool_t perfPID     = opt.Contains("perfect");
39   Bool_t addITSSA    = opt.Contains("its");
40   Bool_t dipAngleCut = opt.Contains("dip");
41   Int_t  typePID     = 0;
42   if (realPID) typePID = 1;
43   else if (perfPID) typePID = 2;
44       
45   // info
46   const Char_t *pid[3] = {"nopid", "realistic", "perfect"};
47   Info("RsnConfig2010PhiFcn", "=== Specific configuration: %s ===", config);
48   Info("RsnConfig2010PhiFcn", "--> PID           : %s", pid[typePID]);
49   Info("RsnConfig2010PhiFcn", "--> ITS standalone: %s", (addITSSA ? "INCLUDED" : "EXCLUDED"));
50   Info("RsnConfig2010PhiFcn", "--> dip-angle cut : %s", (dipAngleCut ? "INCLUDED" : "EXCLUDED"));
51   
52   // generate a common suffix depending on chosen options
53   TString suffix(pid[typePID]);
54   if (addITSSA) suffix += "_sa"; else suffix += "_nosa";
55   if (dipAngleCut) suffix += "_dip";
56   Info("RsnConfig2010PhiFcn", "--> suffix used   : %s", suffix.Data());
57
58   // retrieve analysis manager & task
59   AliAnalysisManager *mgr  = AliAnalysisManager::GetAnalysisManager();
60   AliRsnAnalysisSE   *task = (AliRsnAnalysisSE*)mgr->GetTask(taskName);
61
62   // for safety, return if no task is passed
63   if (!task)
64   {
65     Error("RsnConfig2010PhiFcn", "Task not found");
66     return kFALSE;
67   }
68
69   //
70   // -- Setup pairs ---------------------------------------------------------------------------------
71   //
72
73   // decay channels
74   AliRsnPairDef *pairDefPM = new AliRsnPairDef(AliPID::kKaon, '+', AliPID::kKaon, '-', 333, 1.019455);
75   AliRsnPairDef *pairDefPP = new AliRsnPairDef(AliPID::kKaon, '+', AliPID::kKaon, '+', 333, 1.019455);
76   AliRsnPairDef *pairDefMM = new AliRsnPairDef(AliPID::kKaon, '-', AliPID::kKaon, '-', 333, 1.019455);
77
78   // computation objects
79   AliRsnPairFunctions *pairPM = new AliRsnPairFunctions(Form("PairPM_%s", suffix.Data()), pairDefPM);
80   AliRsnPairFunctions *truePM = new AliRsnPairFunctions(Form("TruePM_%s", suffix.Data()), pairDefPM);
81   AliRsnPairFunctions *pairPP = new AliRsnPairFunctions(Form("PairPP_%s", suffix.Data()), pairDefPP);
82   AliRsnPairFunctions *pairMM = new AliRsnPairFunctions(Form("PairMM_%s", suffix.Data()), pairDefMM);
83
84   //
85   // -- Setup cuts ----------------------------------------------------------------------------------
86   //
87
88   // track cut -----------------------
89   // --> global cuts for 2010 analysis
90   // --> most options are set to right values by default
91   AliRsnCutESD2010 *cuts2010_esd = new AliRsnCutESD2010(Form("cuts2010_esd_%s", suffix.Data()));
92   AliRsnCutAOD2010 *cuts2010_aod = new AliRsnCutAOD2010(Form("cuts2010_aod_%s", suffix.Data()));
93   // ----> set the flag for sim/data management (which sets some other options)
94   cuts2010_esd->SetMC(isSim);
95   cuts2010_aod->SetMC(isSim);
96   // ----> include or not the ITS standalone (TPC is always in)
97   cuts2010_esd->SetUseGlobal(kTRUE);
98   cuts2010_esd->SetUseITSSA (addITSSA);
99   cuts2010_aod->SetUseGlobal(kTRUE);
100   cuts2010_aod->SetUseITSSA (addITSSA);
101   // ----> require to check PID or not, depending on the label
102   if (realPID)
103   {
104     // if doing realistic PID, it must be activated
105     cuts2010_esd->SetCheckITS (kTRUE);
106     cuts2010_esd->SetCheckTPC (kTRUE);
107     cuts2010_esd->SetCheckTOF (kTRUE);
108     cuts2010_aod->SetCheckITS (kTRUE);
109     cuts2010_aod->SetCheckTPC (kTRUE);
110     cuts2010_aod->SetCheckTOF (kTRUE);
111   }
112   else
113   {
114     // otherwise (both for no pid and perfect PID)
115     // the PID cuts are deactivated
116     cuts2010_esd->SetCheckITS (kFALSE);
117     cuts2010_esd->SetCheckTPC (kFALSE);
118     cuts2010_esd->SetCheckTOF (kFALSE);
119     cuts2010_aod->SetCheckITS (kFALSE);
120     cuts2010_aod->SetCheckTPC (kFALSE);
121     cuts2010_aod->SetCheckTOF (kFALSE);
122   }
123   // ----> set all other defaults
124   ConfigESDCutsTPC(cuts2010_esd->GetCutsTPC());
125   ConfigESDCutsITS(cuts2010_esd->GetCutsITS());
126   
127   // track cut -----------------------------
128   // --> perfect PID for check of PID issues
129   AliRsnCutPID *cutPID = new AliRsnCutPID("cutPID", AliPID::kKaon, 0.0, kTRUE);
130   
131   // pair cut ----------------------
132   // --> dip angle between daughters
133   AliRsnCutStd *cutDip = new AliRsnCutStd("cutDip", AliRsnCut::kMother, AliRsnCutStd::kDipAngle, 0.02, 1.1);
134
135   // cut set for tracks------------------------
136   // --> only common cuts for tracks are needed
137   // --> standard 2010 cuts are applied always
138   TString       cutSchemeTrack;
139   AliRsnCutSet *cutSetDaughterCommon = new AliRsnCutSet("commonDaughterCuts", AliRsnCut::kDaughter);
140   if (isESD)
141   {
142     cutSetDaughterCommon->AddCut(cuts2010_esd);
143     cutSchemeTrack += cuts2010_esd->GetName();
144   }
145   else if (isAOD)
146   {
147     cutSetDaughterCommon->AddCut(cuts2010_aod);
148     cutSchemeTrack += cuts2010_aod->GetName();
149   }
150   else
151   {
152     Error("Required ESD or AOD");
153     return kFALSE;
154   }
155   if (perfPID)
156   {
157     cutSetDaughterCommon->AddCut(cutPID);
158     cutSchemeTrack += "&cutPID";
159   }
160   cutSetDaughterCommon->SetCutScheme(cutSchemeTrack.Data());
161   
162   // cut set for pairs---------------------------------------
163   // --> add dip angle cut (but then include only if required)
164   AliRsnCutSet *cutSetPair = new AliRsnCutSet("cutsPair", AliRsnCut::kMother);
165   cutSetPair->AddCut(cutDip);
166   cutSetPair->SetCutScheme(cutDip->GetName());
167   
168   // configure cut managers -------------------
169   // --> track cuts are always defined, so add them by default
170   pairPM->GetCutManager()->SetCommonDaughterCuts(cutSetDaughterCommon);
171   truePM->GetCutManager()->SetCommonDaughterCuts(cutSetDaughterCommon);
172   pairPP->GetCutManager()->SetCommonDaughterCuts(cutSetDaughterCommon);
173   pairMM->GetCutManager()->SetCommonDaughterCuts(cutSetDaughterCommon);
174   // --> pair cut is added only when dip angle cut is required
175   if (dipAngleCut)
176   {
177     pairPM->GetCutManager()->SetMotherCuts(cutSetPair);
178     truePM->GetCutManager()->SetMotherCuts(cutSetPair);
179     pairPP->GetCutManager()->SetMotherCuts(cutSetPair);
180     pairMM->GetCutManager()->SetMotherCuts(cutSetPair);
181   }
182   
183   // set additional option for true pairs when needed
184   truePM->SetOnlyTrue  (kTRUE);
185   truePM->SetCheckDecay(kTRUE);
186
187   //
188   // -- Setup functions -----------------------------------------------------------------------------
189   //
190
191   // function axes
192   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};
193   Double_t y   [] = {-0.8, -0.7, -0.6, 0.6, 0.7, 0.8};
194   Double_t mult[] = {0.0, 5.0, 9.0, 14.0, 22.0, 1000.0};
195   Int_t    npt    = sizeof(pt  ) / sizeof(pt  [0]);
196   Int_t    ny     = sizeof(y   ) / sizeof(y   [0]);
197   Int_t    nmult  = sizeof(mult) / sizeof(mult[0]);
198   AliRsnValue *axisIM   = new AliRsnValue("IM"  , AliRsnValue::kPairInvMass     , 1000 , 0.9,  1.9);
199   AliRsnValue *axisPt   = new AliRsnValue("PT"  , AliRsnValue::kPairPt          , npt  , pt);
200   AliRsnValue *axisY    = new AliRsnValue("Y"   , AliRsnValue::kPairY           , ny   , y);
201   //AliRsnValue *axisMult = new AliRsnValue("Mult", AliRsnValue::kEventMultESDcuts, nmult, mult);
202   AliRsnValue *axisMult = new AliRsnValue("Mult", AliRsnValue::kEventMultESDcuts, 100, 0, 100);
203   ConfigESDCutsTPC(axisMult->GetCuts());
204
205   // create function and add axes
206   AliRsnFunction *fcnImPtY = new AliRsnFunction;
207   //fcnImPtY->AddAxis(axisIM);
208   //fcnImPtY->AddAxis(axisPt);
209   //fcnImPtY->AddAxis(axisY);
210   fcnImPtY->AddAxis(axisMult);
211
212   // add functions to pairs
213   pairPM->AddFunction(fcnImPtY);
214   truePM->AddFunction(fcnImPtY);
215   pairPP->AddFunction(fcnImPtY);
216   pairMM->AddFunction(fcnImPtY);
217
218   //
219   // -- Conclusion ----------------------------------------------------------------------------------
220   //
221
222   // add all created AliRsnPair objects to the AliRsnAnalysisManager in the task
223   task->GetAnalysisManager()->Add(pairPM);
224   task->GetAnalysisManager()->Add(pairPP);
225   task->GetAnalysisManager()->Add(pairMM);
226   if (isSim) task->GetAnalysisManager()->Add(truePM);
227
228   return kTRUE;
229 }