]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/hfe/macros/configs/pPb/ConfigHFEnpepPb.C
updtaes for pPb analysis
[u/mrichter/AliRoot.git] / PWGHF / hfe / macros / configs / pPb / ConfigHFEnpepPb.C
1 AliAnalysisTaskHFE* ConfigHFEnpepPb(Bool_t useMC, TString appendix,\r
2                 UChar_t TPCcl=70, UChar_t TPCclPID = 80, \r
3                 UChar_t ITScl=3, Double_t DCAxy=1000., Double_t DCAz=1000., \r
4                 Double_t* tpcdEdxcutlow=NULL, Double_t* tpcdEdxcuthigh=NULL, \r
5                 Double_t TOFs=3., Int_t TOFmis=0, \r
6                 Int_t itshitpixel = 0, Int_t icent, \r
7                 Double_t assETA=0.8, Int_t assITS=2, \r
8                 Int_t assTPCcl=100, Int_t assTPCPIDcl=80, \r
9                 Double_t assDCAr=1.0, Double_t assDCAz=2.0, \r
10                 Double_t *assTPCSminus=NULL, Double_t *assTPCSplus=NULL)\r
11 {\r
12   Bool_t kAnalyseTaggedTracks = kTRUE;\r
13  \r
14   //***************************************//\r
15   //        Setting up the HFE cuts        //\r
16   //***************************************//\r
17   \r
18   AliHFEcuts *hfecuts = new AliHFEcuts(appendix,"HFE cuts pPb");\r
19   //hfecuts->SetQAOn();\r
20   hfecuts->CreateStandardCuts();\r
21   hfecuts->SetMinNClustersTPC(TPCcl);\r
22   hfecuts->SetMinNClustersTPCPID(TPCclPID);\r
23   hfecuts->SetMinNClustersITS(ITScl);\r
24   hfecuts->SetMinRatioTPCclusters(0.6);\r
25   hfecuts->SetTPCmodes(AliHFEextraCuts::kFound, AliHFEextraCuts::kFoundOverFindable);\r
26   hfecuts->SetCutITSpixel(itshitpixel);\r
27   hfecuts->SetCheckITSLayerStatus(kFALSE);\r
28   hfecuts->SetEtaRange(0.8);\r
29   hfecuts->SetMaxImpactParam(DCAxy,DCAz);\r
30   hfecuts->SetUseMixedVertex(kTRUE);\r
31   hfecuts->SetVertexRange(10.);\r
32 \r
33   // TOF settings:\r
34   Int_t usetof=0;\r
35   Bool_t kTOFmis=kFALSE;\r
36   if (TOFs>0.){\r
37     usetof = 1;\r
38     printf("CONFIGURATION FILE: TOF is used \n");\r
39     hfecuts->SetTOFPIDStep(kTRUE);\r
40     printf("CONFIGURATION FILE: TOF PID step is requested !!!! \n");\r
41     if (TOFmis>0){\r
42       kTOFmis = kTRUE;\r
43       printf("CONFIGURATION FILE: TOF mismatch rejection is set ON \n");\r
44     }\r
45   }\r
46 \r
47   //***************************************//\r
48   //        Setting up the task            //\r
49   //***************************************//\r
50 \r
51   AliAnalysisTaskHFE *task = new AliAnalysisTaskHFE(Form("HFEtask%s",appendix.Data()));\r
52   printf("task %p\n", task);\r
53   task->SetpPbAnalysis();\r
54   task->SetHFECuts(hfecuts);\r
55   task->SetRemovePileUp(kFALSE);\r
56   task->GetPIDQAManager()->SetHighResolutionHistos();\r
57 \r
58   // Setttings for pPb\r
59   task    -> SetRemoveFirstEventInChunk();\r
60   hfecuts -> SetUseCorrelationVertex();\r
61   hfecuts -> SetSPDVtxResolutionCut();\r
62 \r
63   //***************************************//\r
64   //          Variable manager             //\r
65   //***************************************//\r
66 \r
67   // Define Variables\r
68   Double_t ptbinning[36] = {0., 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1., 1.1, 1.2, 1.3, 1.4, 1.5, 1.75, 2., 2.25, 2.5, 2.75, 3., 3.5, 4., 4.5, 5., 5.5, 6., 7., 8., 10., 12., 14., 16., 18., 20.};\r
69   Double_t etabinning[17] = {-0.8, -0.7, -0.6, -0.5, -0.4, -0.3, -0.2, -0.1, 0., 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8};\r
70 \r
71   Int_t sizept=(sizeof(ptbinning)/sizeof(double))-1;\r
72   Int_t sizeeta=(sizeof(etabinning)/sizeof(double))-1;\r
73 \r
74   AliHFEvarManager *vm = task->GetVarManager();\r
75   vm->AddVariable("pt", sizept, ptbinning);\r
76   vm->AddVariable("eta", sizeeta, -0.8,0.8);\r
77   vm->AddVariable("phi",21, -0, 2*TMath::Pi());\r
78   vm->AddVariable("charge");\r
79   vm->AddVariable("source");\r
80   vm->AddVariable("centrality");\r
81 \r
82   // Determine the centrality estimator\r
83   task->SetCentralityEstimator("V0A");\r
84   if (icent == 2) task->SetCentralityEstimator("V0M");\r
85   else if (icent == 3) task->SetCentralityEstimator("CL1");\r
86   else if (icent == 4) task->SetCentralityEstimator("ZNA");\r
87 \r
88   // For the moment, remove the part dedicated to the background subtraction.\r
89   // It will be implemented in a different way, reading it from a root file.\r
90  \r
91 \r
92   //***************************************//\r
93   //          Configure the PID            //\r
94   //***************************************//\r
95 \r
96   // Define PID\r
97   AliHFEpid *pid = task->GetPID();\r
98   if(useMC) pid->SetHasMCData(kTRUE);\r
99 \r
100   if (usetof){\r
101     pid->AddDetector("TOF", 0);\r
102     pid->AddDetector("TPC", 1);\r
103   } else {\r
104     pid->AddDetector("TPC", 0);\r
105   }\r
106   \r
107   // Configure TPC PID\r
108   if(!useMC){\r
109     Double_t paramsTPCdEdxcutlow[12] ={0.0, 0.0, 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0};\r
110     if(tpcdEdxcutlow) memcpy(paramsTPCdEdxcutlow,tpcdEdxcutlow,sizeof(paramsTPCdEdxcutlow));\r
111 \r
112     Double_t paramsTPCdEdxcuthigh[12] ={3.0, 3.0, 3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0};\r
113     if(tpcdEdxcuthigh) memcpy(paramsTPCdEdxcuthigh,tpcdEdxcuthigh,sizeof(paramsTPCdEdxcuthigh));\r
114     \r
115     char *cutmodel;\r
116     cutmodel="pol0";\r
117     \r
118     for(Int_t a=0;a<11;a++)\r
119     {\r
120       //   cout << a << " " << paramsTPCdEdxcut[a] << endl;\r
121       Double_t tpcparamlow[1]={paramsTPCdEdxcutlow[a]};\r
122       Float_t tpcparamhigh=paramsTPCdEdxcuthigh[a];\r
123       pid->ConfigureTPCcentralityCut(a,cutmodel,tpcparamlow,tpcparamhigh);\r
124     }\r
125   }\r
126 \r
127   // Configure TOF PID\r
128   if (usetof){\r
129     pid->ConfigureTOF(TOFs);\r
130     AliHFEpidTOF *tofpid = pid->GetDetPID(AliHFEpid::kTOFpid);\r
131     if (kTOFmis){\r
132       tofpid->SetRejectTOFmismatch();\r
133     }\r
134   }\r
135 \r
136   // To make different upper TOF cut to see contamination effect\r
137   // The below two lines should be removed after this check\r
138   //AliHFEpidTOF *tofpid = pid->GetDetPID(AliHFEpid::kTOFpid);\r
139   //if(TOFs<3.) tofpid->SetTOFnSigmaBand(-3,TOFs); //only to check the assymmetric tof cut\r
140   \r
141   //***************************************//\r
142   //       Configure NPE plugin            //\r
143   //***************************************//\r
144 \r
145   AliHFENonPhotonicElectron *backe = new AliHFENonPhotonicElectron(Form("HFEBackGroundSubtractionPID2%s",appendix.Data()),"Background subtraction");  //appendix\r
146     //Setting the Cuts for the Associated electron-pool\r
147   AliHFEcuts *hfeBackgroundCuts = new AliHFEcuts(Form("HFEBackSub%s",appendix.Data()),"Background sub Cuts");\r
148   hfeBackgroundCuts->SetEtaRange(assETA);\r
149   hfeBackgroundCuts->SetPtRange(0.1,1e10);\r
150 \r
151   hfeBackgroundCuts->SetMaxChi2perClusterTPC(4);\r
152   hfeBackgroundCuts->SetMinNClustersITS(assITS);\r
153   hfeBackgroundCuts->SetMinNClustersTPC(assTPCcl);\r
154   hfeBackgroundCuts->SetMinNClustersTPCPID(assTPCPIDcl);\r
155   hfeBackgroundCuts->SetQAOn();                         // QA\r
156 \r
157   AliHFEpid *pidbackground = backe->GetPIDBackground();\r
158   if(useMC) pidbackground->SetHasMCData(kTRUE);\r
159   pidbackground->AddDetector("TPC", 0);\r
160   Double_t paramsTPCdEdxcutlowAssoc[12] ={-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0};\r
161   if(assTPCSminus) memcpy(paramsTPCdEdxcutlowAssoc,assTPCSminus,sizeof(paramsTPCdEdxcutlowAssoc));\r
162 \r
163   Double_t paramsTPCdEdxcuthighAssoc[12] ={3.0, 3.0, 3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0};\r
164   if(assTPCSplus) memcpy(paramsTPCdEdxcuthighAssoc,assTPCSplus,sizeof(paramsTPCdEdxcuthighAssoc));\r
165     \r
166   char *cutmodelAssoc;\r
167   cutmodelAssoc="pol0";\r
168   for(Int_t a=0;a<11;a++)\r
169   {\r
170     //   cout << a << " " << paramsTPCdEdxcut[a] << endl;\r
171     Double_t tpcparamlow[1]={paramsTPCdEdxcutlowAssoc[a]};\r
172     Float_t tpcparamhigh=paramsTPCdEdxcuthighAssoc[a];\r
173     pidbackground->ConfigureTPCcentralityCut(a,cutmodelAssoc,tpcparamlow,tpcparamhigh);\r
174   }\r
175   backe->GetPIDBackgroundQAManager()->SetHighResolutionHistos();\r
176   backe->SetHFEBackgroundCuts(hfeBackgroundCuts);\r
177 \r
178   task->SetHFEBackgroundSubtraction(backe);\r
179 \r
180   //***************************************//\r
181   //          V0 tagged tracks             //\r
182   //***************************************//\r
183 \r
184   if(kAnalyseTaggedTracks){\r
185     AliHFEcuts *v0trackCuts = new AliHFEcuts("V0trackCuts", "Track Cuts for tagged track Analysis");\r
186     v0trackCuts->CreateStandardCuts();\r
187     v0trackCuts->SetMinNClustersTPC(TPCcl);\r
188     v0trackCuts->SetMinRatioTPCclusters(0.6);\r
189     v0trackCuts->SetTPCmodes(AliHFEextraCuts::kFound, AliHFEextraCuts::kFoundOverFindable);\r
190     v0trackCuts->SetMinNClustersITS(1);\r
191     v0trackCuts->SetCutITSpixel(AliHFEextraCuts::kFirst);\r
192     v0trackCuts->SetCheckITSLayerStatus(kFALSE);\r
193     v0trackCuts->UnsetVertexRequirement();\r
194     //hfecuts->SetSigmaToVertex(10);\r
195     v0trackCuts->SetTOFPIDStep(kTRUE);\r
196     v0trackCuts->SetQAOn();\r
197 \r
198     task->SwitchOnPlugin(AliAnalysisTaskHFE::kTaggedTrackAnalysis);\r
199     task->SetTaggedTrackCuts(v0trackCuts);\r
200     task->SetCleanTaggedTrack(kTRUE);\r
201   }\r
202 \r
203   // QA\r
204   printf("task %p\n", task);\r
205   task->SetQAOn(AliAnalysisTaskHFE::kPIDqa);\r
206   task->SetQAOn(AliAnalysisTaskHFE::kMCqa);\r
207   task->SwitchOnPlugin(AliAnalysisTaskHFE::kNonPhotonicElectron);\r
208   task->SwitchOnPlugin(AliAnalysisTaskHFE::kDEstep);\r
209 \r
210   printf("*************************************\n");\r
211   printf("Configuring standard Task:\n");\r
212   task->PrintStatus();\r
213   pid->PrintStatus();\r
214   printf("*************************************\n");\r
215   return task;\r
216 }\r