]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliAnalysisTaskDxHFEParticleSelection.cxx
Add sigma2 jet shape and fill constituent info. for subtracted jets
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliAnalysisTaskDxHFEParticleSelection.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE Project            * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  Sedat Altinpinar <Sedat.Altinpinar@cern.ch>           *
9 //*                  Hege Erdal       <hege.erdal@gmail.com>               *
10 //*                                                                        *
11 //* Permission to use, copy, modify and distribute this software and its   *
12 //* documentation strictly for non-commercial purposes is hereby granted   *
13 //* without fee, provided that the above copyright notice appears in all   *
14 //* copies and that both the copyright notice and this permission notice   *
15 //* appear in the supporting documentation. The authors make no claims     *
16 //* about the suitability of this software for any purpose. It is          *
17 //* provided "as is" without express or implied warranty.                  *
18 //**************************************************************************
19
20 /// @file   AliAnalysisTaskDxHFEParticleSelection.cxx
21 /// @author Sedat Altinpinar, Hege Erdal, Matthias Richter
22 /// @date   2012-03-19
23 /// @brief  AnalysisTask electron selection for D0 - HFE correlation
24 ///
25
26 #include "AliAnalysisTaskDxHFEParticleSelection.h"
27 #include "AliDxHFEParticleSelection.h"
28 #include "AliDxHFEParticleSelectionD0.h"
29 #include "AliDxHFEParticleSelectionMCD0.h"
30 #include "AliDxHFEParticleSelectionEl.h"
31 #include "AliDxHFEParticleSelectionMCEl.h"
32 #include "AliAnalysisManager.h"
33 #include "AliAnalysisCuts.h"
34 #include "AliLog.h"
35 #include "TH1F.h"
36 #include "AliESDInputHandler.h"
37 #include "AliAODHandler.h"
38 #include "AliAODRecoDecayHF2Prong.h"
39 #include "AliRDHFCutsD0toKpi.h"
40 #include "AliHFEtools.h"
41 #include "TChain.h"
42 #include "TSystem.h"
43 #include "TObjArray.h"
44 #include "AliAODMCParticle.h"
45 #include "TFile.h"
46 #include "TList.h"
47 #include "TObjArray.h"
48 #include <memory>
49
50 using namespace std;
51
52 /// ROOT macro for the implementation of ROOT specific class methods
53 ClassImp(AliAnalysisTaskDxHFEParticleSelection)
54
55 AliAnalysisTaskDxHFEParticleSelection::AliAnalysisTaskDxHFEParticleSelection(const char* opt)
56   : AliAnalysisTaskSE("AliAnalysisTaskDxHFEParticleSelection")
57   , fOutput(0)
58   , fOption(opt)
59   , fCutList(NULL)
60   , fCutsD0(NULL)
61   , fSelector(NULL)
62   , fUseMC(kFALSE)
63   , fFillOnlyD0D0bar(0)
64   , fSelectedTracks(NULL)
65   , fMCArray(NULL)
66   , fParticleType(kD0)
67   , fSystem(0)
68   , fUseKine(kFALSE)
69 {
70   // constructor
71   //
72   //
73
74   DefineSlots();
75 }
76
77 int AliAnalysisTaskDxHFEParticleSelection::DefineSlots()
78 {
79   // define the data slots
80   DefineInput(0, TChain::Class());
81   DefineOutput(1, TList::Class());
82   DefineOutput(2, TList::Class());
83   return 0;
84 }
85
86 AliAnalysisTaskDxHFEParticleSelection::~AliAnalysisTaskDxHFEParticleSelection()
87 {
88   // destructor
89   //
90   //
91
92   // histograms are in the output list and deleted when the output
93   // list is deleted by the TSelector dtor
94
95   if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
96     delete fOutput;
97     fOutput = 0;
98   }
99
100   if (fSelector) {
101     fSelector->Clear();
102     delete fSelector;
103     fSelector=NULL;
104   }
105   if(fSelectedTracks){
106     delete fSelectedTracks;
107     fSelectedTracks=NULL;
108   }
109   if(fMCArray) delete fMCArray;
110   fMCArray=NULL;
111   // external object, do not delete
112   fCutsD0=NULL;
113
114 }
115
116 void AliAnalysisTaskDxHFEParticleSelection::UserCreateOutputObjects()
117 {
118   // create result objects and add to output list
119
120   Int_t iResult=0;
121
122   fOutput = new TList;
123   fOutput->SetOwner();
124
125   ParseArguments(fOption.Data());
126
127   // REVIEW: this has only effect if the list is deleted, it should
128   // be done where the list is created, in any case, all objects of the
129   // list are presumably external objects, they should not be deleted in
130   // this class
131   fCutList->SetOwner(); // NOt sure if needed
132
133   // setting up for D0s
134   if(fParticleType==kD0){
135
136     if(fUseMC) fSelector=new AliDxHFEParticleSelectionMCD0(fOption);
137     else fSelector=new AliDxHFEParticleSelectionD0(fOption);
138     fSelector->SetCuts(fCutList,AliDxHFEParticleSelectionD0::kCutList);
139     iResult=fSelector->Init();
140
141     TObject *obj=NULL;
142     TIter iter(fCutList);
143     while((obj = iter())){
144       fCutsD0=dynamic_cast<AliRDHFCuts*>(obj);
145       if (!fCutsD0) {
146         AliError(Form("Cut object is not of required type AliRDHFCuts but %s", obj->ClassName()));
147       }
148     }
149   }
150   
151
152   // Setting up for electrons
153   if(fParticleType==kElectron){
154     if(fUseMC) fSelector=new AliDxHFEParticleSelectionMCEl(fOption);
155     else fSelector=new AliDxHFEParticleSelectionEl(fOption);
156     fSelector->SetCuts(fCutList, AliDxHFEParticleSelectionEl::kCutList);
157     iResult=fSelector->Init();
158
159     // If running on electrons, for now use RDHFCuts for event selection. 
160     // Set up default cut object:
161     AliRDHFCutsD0toKpi* cuts=new AliRDHFCutsD0toKpi();
162     cuts->SetStandardCutsPP2010();
163     fCutsD0=cuts;
164   }
165
166   if (iResult<0) {
167     AliFatal(Form("initialization of worker class instance fElectrons failed with error %d", iResult));
168   }
169
170
171   // Retrieving the list containing histos and THnSparse
172   // and storing them instead of fSelector
173   // Fix to be able to merge
174   TList *list =(TList*)fSelector->GetControlObjects();
175   TObject *obj=NULL;
176
177   TIter next(list);
178   while((obj = next())){
179     fOutput->Add(obj);
180   }
181
182   // all tasks must post data once for all outputs
183   PostData(1, fOutput);
184   PostData(2, fCutList);
185 }
186
187 void AliAnalysisTaskDxHFEParticleSelection::UserExec(Option_t* /*option*/)
188 {
189   // process the event
190
191   // TODO: implement correct input, this is likely not to be the
192   // ESD
193   TObject* pInput=InputEvent();
194   if (!pInput) {
195     AliError("failed to get input");
196     return;
197   }
198
199   // check if input is an ESD
200   AliVEvent *pEvent = dynamic_cast<AliVEvent*>(pInput);
201   TClonesArray *inputArray=0;
202
203   if(!pEvent && AODEvent() && IsStandardAOD()) { //Not sure if this is needed.. Keep it for now. 
204     // In case there is an AOD handler writing a standard AOD, use the AOD 
205     // event in memory rather than the input (ESD) event.    
206     pEvent = AODEvent();
207     // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
208     // have to taken from the AOD event hold by the AliAODExtension
209     AliAODHandler* aodHandler = (AliAODHandler*) 
210       ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
211     
212     if(aodHandler->GetExtensions()) {
213       AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
214       AliAODEvent* aodFromExt = ext->GetAOD();
215       inputArray=(TClonesArray*)aodFromExt->GetList()->FindObject("D0toKpi");
216     }
217   } else if(pEvent) {
218     inputArray=(TClonesArray*)pEvent->GetList()->FindObject("D0toKpi");
219   }
220   if(!inputArray || !pEvent) {
221     AliError("Input branch not found!\n");
222     return;
223   }
224   // fix for temporary bug in ESDfilter
225   // the AODs with null vertex pointer didn't pass the PhysSel
226   if(!pEvent->GetPrimaryVertex() || TMath::Abs(pEvent->GetMagneticField())<0.001){
227     AliDebug(2,"Rejected at GetPrimaryvertex");
228     return;
229   }
230
231   fSelector->HistogramEventProperties(AliDxHFEParticleSelection::kEventsAll);
232
233   //At the moment: Use AliRDHFCuts for event selection
234   AliRDHFCuts* cutsd0=dynamic_cast<AliRDHFCuts*>(fCutsD0);
235   if (!cutsd0) return; // Fatal thrown already in initialization
236
237   if(!cutsd0->IsEventSelected(pEvent)) {
238     AliDebug(2,"rejected at IsEventSelected");
239     return;
240   }
241
242   fSelector->HistogramEventProperties(AliDxHFEParticleSelection::kEventsSel);
243
244   if(fUseMC && fUseKine){
245     fMCArray = dynamic_cast<TObjArray*>(pEvent->FindListObject(AliAODMCParticle::StdBranchName()));
246     if(!fMCArray){
247       AliError("Array of MC particles not found");
248       return;
249     }
250   }
251
252   if (fSelectedTracks) delete fSelectedTracks;
253
254   if(fParticleType==kElectron){
255     // Gets the PID response from the analysis manager
256     AliPIDResponse *pidResponse = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler())->GetPIDResponse();
257     if(!pidResponse){
258       // TODO: consider issuing fatal instead of debug in case pidresponse not available
259       AliDebug(1, "Using default PID Response");
260       pidResponse = AliHFEtools::GetDefaultPID(kFALSE, fInputEvent->IsA() == AliAODEvent::Class()); 
261     }
262   
263     // Fetching the PID objects from the list, checks if the objects are AliHFEpids
264     // If so, checks if they are initialized and also sets the pidresponse
265     TObject *obj=NULL;
266     TIter next(fCutList);
267     while((obj = next())){
268       if(strcmp(obj->ClassName(),"AliHFEpid")==0){
269         AliHFEpid* pidObj=(AliHFEpid*)obj;
270         if(!pidObj->IsInitialized()){
271           AliWarning("PID not initialised, get from Run no");
272           pidObj->InitializePID(pEvent->GetRunNumber());
273         }
274         pidObj->SetPIDResponse(pidResponse);
275       }
276     }
277
278     // Also sends the pidresponse to the particle selection class for electron
279     fSelector->SetPIDResponse(pidResponse); 
280
281     if(fUseKine) fSelectedTracks=(fSelector->Select(fMCArray, pEvent));
282     else fSelectedTracks=(fSelector->Select(pEvent));
283   }
284   else{
285     if(fUseKine)  fSelectedTracks=(fSelector->Select(fMCArray,pEvent));
286     else fSelectedTracks=(fSelector->Select(inputArray,pEvent));
287   }
288
289   Int_t nSelected = fSelectedTracks->GetEntriesFast();
290   if(nSelected>0)
291     fSelector->HistogramEventProperties(AliDxHFEParticleSelection::kEventsWithParticle);
292
293   PostData(1, fOutput);
294 }
295
296 int AliAnalysisTaskDxHFEParticleSelection::ParseArguments(const char* arguments)
297 {
298   // parse arguments and set internal flags
299   TString strArguments(arguments);
300   auto_ptr<TObjArray> tokens(strArguments.Tokenize(" "));
301   if (!tokens.get()) return 0;
302
303   TIter next(tokens.get());
304   TObject* token;
305   while ((token=next())) {
306     TString argument=token->GetName();
307    
308     if (argument.BeginsWith("mc")) {
309       fUseMC=true;
310       continue;
311     }
312     if (argument.BeginsWith("system=")) {
313       argument.ReplaceAll("system=", "");
314       if (argument.CompareTo("pp")==0) {fSystem=0; }
315       else if (argument.CompareTo("Pb-Pb")==0){ fSystem=1;}
316       else if (argument.CompareTo("p-Pb")==0){ fSystem=2;}
317       else {
318         AliWarning(Form("can not set collision system, unknown parameter '%s'", argument.Data()));
319         // TODO: check what makes sense
320         fSystem=0;
321       }
322       continue;
323     }
324     if (argument.BeginsWith("usekine") || argument.BeginsWith("kine")) {
325       fUseKine=true;
326       AliInfo("Running on MC stack");
327       continue;
328     }
329     if (argument.BeginsWith("fillD0scheme=")){
330       argument.ReplaceAll("fillD0scheme=", "");
331       if (argument.CompareTo("both")==0){ fFillOnlyD0D0bar=0;}
332       else if (argument.CompareTo("D0")==0){ fFillOnlyD0D0bar=1;}
333       else if (argument.CompareTo("D0bar")==0){ fFillOnlyD0D0bar=2;}
334       else {
335         AliWarning(Form("can not set D0 filling scheme, unknown parameter '%s'", argument.Data()));
336         fFillOnlyD0D0bar=0;
337       }
338       continue;
339     }
340     if (argument.BeginsWith("particle=")){
341       argument.ReplaceAll("particle=","");     
342       if (argument.CompareTo("D0")==0){ fParticleType=kD0; AliInfo("Running over D0s ");}
343       else if (argument.CompareTo("electron")==0){ fParticleType=kElectron; AliInfo("Running over Electrons ");}
344       else{
345         AliWarning(Form("can not set particle, unknown parameter '%s'. Particle set to D0", argument.Data()));
346         fParticleType=kD0;
347       }
348       continue;
349     }
350     AliWarning(Form("unknown argument '%s'", argument.Data()));
351       
352   }
353
354   return 0;
355 }
356
357 void AliAnalysisTaskDxHFEParticleSelection::FinishTaskOutput()
358 {
359   // end of the processing
360 }
361
362 void AliAnalysisTaskDxHFEParticleSelection::Terminate(Option_t *)
363 {
364   // last action on the client
365   fOutput = dynamic_cast<TList*> (GetOutputData(1));
366   if (!fOutput) {
367     AliFatal("failed to get output container");
368     return;
369   }
370
371 }