]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliAnalysisTaskDxHFECorrelation.cxx
7cedf9b99e9cc9c00844ca4c5fa157fe14a21c3a
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliAnalysisTaskDxHFECorrelation.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   AliAnalysisTaskDxHFECorrelation.cxx
21 /// @author Sedat Altinpinar, Hege Erdal, Matthias Richter
22 /// @date   2012-03-19
23 /// @brief  AnalysisTask D0 - HFE correlation
24 ///
25
26 #include "AliAnalysisTaskDxHFECorrelation.h"
27 #include "AliDxHFECorrelation.h"
28 #include "AliDxHFEParticleSelectionD0.h"
29 #include "AliDxHFEParticleSelectionMCD0.h"
30 #include "AliDxHFEParticleSelectionEl.h"
31 #include "AliDxHFEParticleSelectionMCEl.h"
32 #include "AliDxHFEParticleSelection.h"
33 #include "AliAnalysisManager.h"
34 #include "AliLog.h"
35 #include "AliESDInputHandler.h"
36 #include "AliAnalysisDataSlot.h"
37 #include "AliAnalysisDataContainer.h"
38 #include "AliAnalysisManager.h"
39 #include "AliVertexerTracks.h"
40 #include "AliAODHandler.h"
41 #include "AliInputEventHandler.h"
42 #include "AliAODEvent.h"
43 #include "AliAODVertex.h"
44 #include "AliAODTrack.h"
45 #include "AliAODMCHeader.h"
46 #include "AliAODMCParticle.h"
47 #include "AliAODRecoDecayHF2Prong.h"
48 #include "AliAODRecoCascadeHF.h"
49 #include "AliRDHFCutsD0toKpi.h"
50 #include "AliPID.h"
51 #include "AliPIDResponse.h"
52 #include "AliHFEcontainer.h"
53 #include "AliHFEpid.h"
54 #include "AliHFEpidBase.h"
55 #include "AliHFEcuts.h"
56 #include "AliHFEtools.h"
57 #include "TObject.h"
58 #include "TChain.h"
59 #include "TSystem.h"
60 #include "TFile.h"
61 #include <memory>
62
63 using namespace std;
64
65 /// ROOT macro for the implementation of ROOT specific class methods
66 ClassImp(AliAnalysisTaskDxHFECorrelation)
67
68 AliAnalysisTaskDxHFECorrelation::AliAnalysisTaskDxHFECorrelation(const char* opt)
69   : AliAnalysisTaskSE("AliAnalysisTaskDxHFECorrelation")
70   , fOutput(0)
71   , fOption(opt)
72   , fCorrelation(NULL)
73   , fD0s(NULL)
74   , fElectrons(NULL)
75   , fCutsD0(NULL)
76   , fCutsHFE(NULL)
77   , fPID(NULL)
78   , fFillOnlyD0D0bar(0)
79   , fUseMC(kFALSE)
80 {
81   // constructor
82   //
83   //
84   DefineSlots();
85   fPID = new AliHFEpid("hfePid");
86
87 }
88
89 int AliAnalysisTaskDxHFECorrelation::DefineSlots()
90 {
91   // define the data slots
92   DefineInput(0, TChain::Class());
93   DefineOutput(1, TList::Class());
94   DefineOutput(2,AliRDHFCutsD0toKpi::Class());
95   return 0;
96 }
97
98 AliAnalysisTaskDxHFECorrelation::~AliAnalysisTaskDxHFECorrelation()
99 {
100   // destructor
101   //
102   //
103
104   // histograms are in the output list and deleted when the output
105   // list is deleted by the TSelector dtor
106
107   if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
108     delete fOutput;
109     fOutput = 0;
110   }
111   if (fD0s) delete fD0s;
112   fD0s=NULL;
113   if (fElectrons) delete fElectrons;
114   fElectrons=NULL;
115   if (fCorrelation) delete fCorrelation;
116   fCorrelation=NULL;
117   // external object, do not delete
118   fCutsD0=NULL;
119   // external object, do not delete
120   fCutsHFE=NULL;
121   if(fPID) delete fPID;
122   fPID=NULL;
123
124
125 }
126
127 void AliAnalysisTaskDxHFECorrelation::UserCreateOutputObjects()
128 {
129   // create result objects and add to output list
130
131   //Initialize PID for electron selection
132   if(!fPID->GetNumberOfPIDdetectors()) { 
133     fPID->AddDetector("TOF",0);
134     fPID->AddDetector("TPC",1);
135   }
136   fPID->InitializePID();
137
138   fOutput = new TList;
139   fOutput->SetOwner();
140
141   // setting up for D0s
142   TString selectionD0Options;
143   switch (fFillOnlyD0D0bar) {
144   case 1: selectionD0Options+="FillOnlyD0 "; break;
145   case 2: selectionD0Options+="FillOnlyD0bar "; break;
146   default: selectionD0Options+="FillD0D0bar ";
147   }
148
149   if(fUseMC) fD0s=new AliDxHFEParticleSelectionMCD0(selectionD0Options);
150   else fD0s=new AliDxHFEParticleSelectionD0(selectionD0Options);
151   fD0s->SetCuts(fCutsD0);
152   fD0s->Init();
153
154   //Electrons
155   if(fUseMC) fElectrons=new AliDxHFEParticleSelectionMCEl;
156   else fElectrons=new AliDxHFEParticleSelectionEl;
157   fElectrons->SetCuts(fPID, AliDxHFEParticleSelectionEl::kCutPID);
158   fElectrons->SetCuts(fCutsHFE, AliDxHFEParticleSelectionEl::kCutHFE);
159   fElectrons->Init();
160
161   //Correlation
162   fCorrelation=new AliDxHFECorrelation;
163   fCorrelation->SetCuts(dynamic_cast<AliRDHFCutsD0toKpi*>(fCutsD0));
164   // TODO: check if we can get rid of the mc flag in the correlation analysis class
165   fCorrelation->SetUseMC(fUseMC);
166   fCorrelation->Init();
167
168   // Fix for merging:
169   // Retrieving the individual objects created
170   // and storing them instead of fD0s, fElectrons etc.. 
171   TList *list =(TList*)fD0s->GetControlObjects();
172   TObject *obj=NULL;
173
174   TIter next(list);
175   while((obj = next())){
176     fOutput->Add(obj);
177   }
178
179   list=(TList*)fCorrelation->GetControlObjects();
180   next=TIter(list);
181   while((obj= next())){
182     fOutput->Add(obj);
183   }
184
185   list=(TList*)fElectrons->GetControlObjects();
186   next=TIter(list);
187   while((obj = next()))
188     fOutput->Add(obj);
189
190   PostData(1, fOutput);
191   if (fCutsD0) {
192     AliRDHFCutsD0toKpi* cuts=dynamic_cast<AliRDHFCutsD0toKpi*>(fCutsD0);
193     if (!cuts) {
194       AliFatal(Form("cut object %s is of incorrect type %s, expecting AliRDHFCutsD0toKpi", fCutsD0->GetName(), fCutsD0->ClassName()));
195       return;
196     }
197     // TODO: why copy? cleanup?
198     AliRDHFCutsD0toKpi* copyfCuts=new AliRDHFCutsD0toKpi(*cuts);
199     const char* nameoutput=GetOutputSlot(2)->GetContainer()->GetName();
200     copyfCuts->SetName(nameoutput);
201
202     PostData(2,copyfCuts);
203   }
204 }
205
206 void AliAnalysisTaskDxHFECorrelation::UserExec(Option_t* /*option*/)
207 {
208   // process the event
209   TObject* pInput=InputEvent();
210   if (!pInput) {
211     AliError("failed to get input");
212     return;
213   }
214   AliVEvent *pEvent = dynamic_cast<AliVEvent*>(pInput);
215   TClonesArray *inputArray=0;
216
217   if(!pEvent && AODEvent() && IsStandardAOD()) { //Not sure if this is needed.. Keep it for now. 
218     // In case there is an AOD handler writing a standard AOD, use the AOD 
219     // event in memory rather than the input (ESD) event.    
220     pEvent = dynamic_cast<AliAODEvent*> (AODEvent());
221     // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
222     // have to taken from the AOD event hold by the AliAODExtension
223     AliAODHandler* aodHandler = (AliAODHandler*) 
224       ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
225     
226     if(aodHandler->GetExtensions()) {
227       AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
228       AliAODEvent* aodFromExt = ext->GetAOD();
229       inputArray=(TClonesArray*)aodFromExt->GetList()->FindObject("D0toKpi");
230     }
231   } else if(pEvent) {
232     inputArray=(TClonesArray*)pEvent->GetList()->FindObject("D0toKpi");
233   }
234   if(!inputArray || !pEvent) {
235     AliError("Input branch not found!\n");
236     return;
237   }
238   // fix for temporary bug in ESDfilter
239   // the AODs with null vertex pointer didn't pass the PhysSel
240   if(!pEvent->GetPrimaryVertex() || TMath::Abs(pEvent->GetMagneticField())<0.001){
241     AliDebug(2,"Rejected at GetPrimaryvertex");
242     return;
243   }
244
245   fCorrelation->HistogramEventProperties(AliDxHFECorrelation::kEventsAll);
246   AliRDHFCuts* cutsd0=dynamic_cast<AliRDHFCuts*>(fCutsD0);
247   if (!cutsd0) return; // Fatal thrown already in initialization
248
249   if(!cutsd0->IsEventSelected(pEvent)) {
250     AliDebug(2,"rejected at IsEventSelected");
251     return;
252   }
253
254   if(!fPID->IsInitialized()){ 
255     // Initialize PID with the given run number
256     AliWarning("PID not initialised, get from Run no");
257     fPID->InitializePID(pEvent->GetRunNumber());
258   }
259
260   AliPIDResponse *pidResponse = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler())->GetPIDResponse();
261   if(!pidResponse){
262     AliDebug(1, "Using default PID Response");
263     pidResponse = AliHFEtools::GetDefaultPID(kFALSE, fInputEvent->IsA() == AliAODEvent::Class()); 
264   }
265  
266   fPID->SetPIDResponse(pidResponse);
267
268
269   Int_t nInD0toKpi = inputArray->GetEntriesFast();
270
271   fCorrelation->HistogramEventProperties(AliDxHFECorrelation::kEventsSel);
272
273   //TODO: pSelectedD0s will now contain both D0 and D0bar. Probably add argument to Select
274   // which determines which ones to return. 
275   std::auto_ptr<TObjArray> pSelectedD0s(fD0s->Select(inputArray,pEvent));
276   if(! pSelectedD0s.get()) {
277     return;
278   }
279   Int_t nD0Selected = pSelectedD0s->GetEntriesFast();
280
281   // No need to go further if no D0s are found. Not sure if this is the best way though..
282   // At the moment this means no control histos can be implemented in AliDxHFECorrelation
283   if(nD0Selected==0){
284     //AliInfo("No D0s found in this event");
285     return;
286   }
287   fCorrelation->HistogramEventProperties(AliDxHFECorrelation::kEventsD0);
288
289
290   std::auto_ptr<TObjArray> pSelectedElectrons(fElectrons->Select(pEvent));
291
292   // TODO: use the array of selected track for something, right now
293   // only the control histograms of the selection class are filled
294   // note: the pointer is deleted automatically once the scope is left
295   // if the array should be published, the auto pointer must be released
296   // first, however some other cleanup will be necessary in that case
297   // probably a clone with a reduced AliVParticle implementation is
298   // appropriate.
299
300
301   if(! pSelectedElectrons.get()) {
302     return;
303   }
304
305   Int_t nElSelected =  pSelectedElectrons->GetEntriesFast();
306
307   // No need to go further if no electrons are found. Not sure if this is the best way though..
308   // At the moment this means no control histos can be implemented in AliDxHFECorrelation
309   if(nElSelected==0){
310     //AliInfo("No electrons found in this event");
311     return;
312   }
313
314   AliDebug(4,Form("Number of D0->Kpi Start: %d , End: %d    Electrons Selected: %d\n", nInD0toKpi, nD0Selected, nElSelected));
315
316   fCorrelation->HistogramEventProperties(AliDxHFECorrelation::kEventsD0e);
317
318   //This is only called if there are electrons and D0s present.
319   int iResult=fCorrelation->Fill(pSelectedD0s.get(), pSelectedElectrons.get());
320
321   if (iResult<0) {
322     AliError(Form("%s processing failed with error %d", fCorrelation->GetName(), iResult));
323   }
324
325   PostData(1, fOutput);
326 }
327
328 void AliAnalysisTaskDxHFECorrelation::FinishTaskOutput()
329 {
330   // end of the processing
331 }
332
333 void AliAnalysisTaskDxHFECorrelation::Terminate(Option_t *)
334 {
335   // last action on the client
336   fOutput = dynamic_cast<TList*> (GetOutputData(1));
337   if (!fOutput) {
338     AliFatal("failed to get output container");
339     return;
340   }
341
342 }