]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliAnalysisTaskDxHFECorrelation.cxx
Add HFE in the include path for par files (Matthias)
[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 "AliAnalysisManager.h"
33 #include "AliLog.h"
34 #include "AliESDInputHandler.h"
35 #include "AliAnalysisDataSlot.h"
36 #include "AliAnalysisDataContainer.h"
37 #include "AliAnalysisManager.h"
38 #include "AliVertexerTracks.h"
39 #include "AliAODHandler.h"
40 #include "AliInputEventHandler.h"
41 #include "AliAODEvent.h"
42 #include "AliAODVertex.h"
43 #include "AliAODTrack.h"
44 #include "AliAODMCHeader.h"
45 #include "AliAODMCParticle.h"
46 #include "AliAODRecoDecayHF2Prong.h"
47 #include "AliAODRecoCascadeHF.h"
48 #include "AliRDHFCutsD0toKpi.h"
49 #include "AliPID.h"
50 #include "AliPIDResponse.h"
51 #include "AliHFEcontainer.h"
52 #include "AliHFEpid.h"
53 #include "AliHFEpidBase.h"
54 #include "AliHFEcuts.h"
55 #include "AliHFEtools.h"
56 #include "TObject.h"
57 #include "TChain.h"
58 #include "TSystem.h"
59 #include "TFile.h"
60 #include <memory>
61
62 using namespace std;
63
64 /// ROOT macro for the implementation of ROOT specific class methods
65 ClassImp(AliAnalysisTaskDxHFECorrelation)
66
67 AliAnalysisTaskDxHFECorrelation::AliAnalysisTaskDxHFECorrelation(const char* opt)
68   : AliAnalysisTaskSE("AliAnalysisTaskDxHFECorrelation")
69   , fOutput(0)
70   , fOption(opt)
71   , fCorrelation(NULL)
72   , fD0s(NULL)
73   , fElectrons(NULL)
74   , fCutsD0(NULL)
75   , fCutsHFE(NULL)
76   , fPID(NULL)
77   , fFillOnlyD0D0bar(0)
78   , fUseMC(kFALSE)
79 {
80   // constructor
81   //
82   //
83   DefineSlots();
84   fPID = new AliHFEpid("hfePid");
85
86 }
87
88 int AliAnalysisTaskDxHFECorrelation::DefineSlots()
89 {
90   // define the data slots
91   DefineInput(0, TChain::Class());
92   DefineOutput(1, TList::Class());
93   DefineOutput(2,AliRDHFCutsD0toKpi::Class());
94   return 0;
95 }
96
97 AliAnalysisTaskDxHFECorrelation::~AliAnalysisTaskDxHFECorrelation()
98 {
99   // destructor
100   //
101   //
102
103   // histograms are in the output list and deleted when the output
104   // list is deleted by the TSelector dtor
105
106   if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
107     delete fOutput;
108     fOutput = 0;
109   }
110   if (fD0s) delete fD0s;
111   fD0s=NULL;
112   if (fElectrons) delete fElectrons;
113   fElectrons=NULL;
114   if (fCorrelation) delete fCorrelation;
115   fCorrelation=NULL;
116   if (fCutsD0) delete fCutsD0;
117   fCutsD0=NULL;
118   if (fCutsHFE) delete fCutsHFE;
119   fCutsHFE=NULL;
120   if(fPID) delete fPID;
121   fPID=NULL;
122
123
124 }
125
126 void AliAnalysisTaskDxHFECorrelation::UserCreateOutputObjects()
127 {
128   // create result objects and add to output list
129
130   //Initialize PID for electron selection
131   if(!fPID->GetNumberOfPIDdetectors()) { 
132     fPID->AddDetector("TOF",0);
133     fPID->AddDetector("TPC",1);
134   }
135   fPID->InitializePID();
136
137   fOutput = new TList;
138   fOutput->SetOwner();
139
140   // setting up for D0s
141   TString selectionD0Options;
142   switch (fFillOnlyD0D0bar) {
143   case 1: selectionD0Options+="FillOnlyD0 "; break;
144   case 2: selectionD0Options+="FillOnlyD0bar "; break;
145   default: selectionD0Options+="FillD0D0bar ";
146   }
147
148   if(fUseMC) fD0s=new AliDxHFEParticleSelectionMCD0(selectionD0Options);
149   else fD0s=new AliDxHFEParticleSelectionD0(selectionD0Options);
150   fD0s->SetCuts(fCutsD0);
151   fD0s->Init();
152
153   //Electrons
154   if(fUseMC) fElectrons=new AliDxHFEParticleSelectionMCEl;
155   else fElectrons=new AliDxHFEParticleSelectionEl;
156   fElectrons->SetCuts(fPID, AliDxHFEParticleSelectionEl::kCutPID);
157   fElectrons->SetCuts(fCutsHFE, AliDxHFEParticleSelectionEl::kCutHFE);
158   fElectrons->Init();
159
160   //Correlation
161   fCorrelation=new AliDxHFECorrelation;
162   fCorrelation->SetCuts(dynamic_cast<AliRDHFCutsD0toKpi*>(fCutsD0));
163   fCorrelation->Init();
164
165   // Fix for merging:
166   // Retrieving the individual objects created
167   // and storing them instead of fD0s, fElectrons etc.. 
168   TList *list =(TList*)fD0s->GetControlObjects();
169   TObject *obj=NULL;
170
171   TIter next(list);
172   while((obj = next())){
173     fOutput->Add(obj);
174   }
175
176   list=(TList*)fCorrelation->GetControlObjects();
177   next=TIter(list);
178   while((obj= next())){
179     fOutput->Add(obj);
180   }
181
182   list=(TList*)fElectrons->GetControlObjects();
183   next=TIter(list);
184   while((obj = next()))
185     fOutput->Add(obj);
186
187   if (fCutsD0) {
188     AliRDHFCutsD0toKpi* cuts=dynamic_cast<AliRDHFCutsD0toKpi*>(fCutsD0);
189     if (!cuts) {
190       AliFatal(Form("cut object %s is of incorrect type %s, expecting AliRDHFCutsD0toKpi", fCutsD0->GetName(), fCutsD0->ClassName()));
191       return;
192     }
193   }
194
195   // TODO: why copy? cleanup?
196   AliRDHFCutsD0toKpi* copyfCuts=new AliRDHFCutsD0toKpi(dynamic_cast<AliRDHFCutsD0toKpi&>(*fCutsD0));
197   const char* nameoutput=GetOutputSlot(2)->GetContainer()->GetName();
198   copyfCuts->SetName(nameoutput);
199
200   // all tasks must post data once for all outputs
201   PostData(1, fOutput);
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 }