]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliAnalysisTaskDxHFEParticleSelection.cxx
#97492 Request to: patch AliSimulation; port to Release; make tag on release; for...
[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 "AliAnalysisManager.h"
31 #include "AliAnalysisCuts.h"
32 #include "AliLog.h"
33 #include "TH1F.h"
34 #include "AliESDInputHandler.h"
35 #include "AliAODHandler.h"
36 #include "AliAODRecoDecayHF2Prong.h"
37 #include "AliRDHFCutsD0toKpi.h"
38 #include "TChain.h"
39 #include "TSystem.h"
40 #include "TFile.h"
41 #include "TObjArray.h"
42 #include <memory>
43
44 using namespace std;
45
46 /// ROOT macro for the implementation of ROOT specific class methods
47 ClassImp(AliAnalysisTaskDxHFEParticleSelection)
48
49 AliAnalysisTaskDxHFEParticleSelection::AliAnalysisTaskDxHFEParticleSelection(const char* opt)
50   : AliAnalysisTaskSE("AliAnalysisTaskDxHFEParticleSelection")
51   , fOutput(0)
52   , fOption(opt)
53   , fCuts(NULL)
54   , fSelector(NULL)
55   , fUseMC(kFALSE)
56   , fFillOnlyD0D0bar(0)
57 {
58   // constructor
59   //
60   //
61
62   DefineSlots();
63 }
64
65 int AliAnalysisTaskDxHFEParticleSelection::DefineSlots()
66 {
67   // define the data slots
68   DefineInput(0, TChain::Class());
69   DefineOutput(1, TList::Class());
70   return 0;
71 }
72
73 AliAnalysisTaskDxHFEParticleSelection::~AliAnalysisTaskDxHFEParticleSelection()
74 {
75   // destructor
76   //
77   //
78
79   // histograms are in the output list and deleted when the output
80   // list is deleted by the TSelector dtor
81
82   if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
83     delete fOutput;
84     fOutput = 0;
85   }
86
87   if (fSelector) {
88     fSelector->Clear();
89     delete fSelector;
90     fSelector=NULL;
91   }
92
93   if (fCuts) {
94     delete fCuts;
95     fCuts=NULL;
96   }
97 }
98
99 void AliAnalysisTaskDxHFEParticleSelection::UserCreateOutputObjects()
100 {
101   // create result objects and add to output list
102
103   fOutput = new TList;
104   fOutput->SetOwner();
105
106   // setting up for D0s
107   TString selectionD0Options;
108   switch (fFillOnlyD0D0bar) {
109   case 1: selectionD0Options+="FillOnlyD0 "; break;
110   case 2: selectionD0Options+="FillOnlyD0bar "; break;
111   default: selectionD0Options+="FillD0D0bar ";
112   }
113
114   if(fUseMC) fSelector=new AliDxHFEParticleSelectionMCD0(selectionD0Options);
115   else fSelector=new AliDxHFEParticleSelectionD0(selectionD0Options);
116   fSelector->SetCuts(fCuts);
117   fSelector->Init();
118
119   // Retrieving the list containing histos and THnSparse
120   // and storing them instead of fSelector
121   // Fix to be able to merge
122   TList *list =(TList*)fSelector->GetControlObjects();
123   TObject *obj=NULL;
124
125   TIter next(list);
126   while((obj = next())){
127     fOutput->Add(obj);
128   }
129
130   // all tasks must post data once for all outputs
131   PostData(1, fOutput);
132 }
133
134 void AliAnalysisTaskDxHFEParticleSelection::UserExec(Option_t* /*option*/)
135 {
136   // process the event
137
138   // TODO: implement correct input, this is likely not to be the
139   // ESD
140   TObject* pInput=InputEvent();
141   if (!pInput) {
142     AliError("failed to get input");
143     return;
144   }
145
146   // check if input is an ESD
147   AliVEvent *pEvent = dynamic_cast<AliVEvent*>(pInput);
148   TClonesArray *inputArray=0;
149
150   if(!pEvent && AODEvent() && IsStandardAOD()) { //Not sure if this is needed.. Keep it for now. 
151     // In case there is an AOD handler writing a standard AOD, use the AOD 
152     // event in memory rather than the input (ESD) event.    
153     pEvent = AODEvent();
154     // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
155     // have to taken from the AOD event hold by the AliAODExtension
156     AliAODHandler* aodHandler = (AliAODHandler*) 
157       ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
158     
159     if(aodHandler->GetExtensions()) {
160       AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
161       AliAODEvent* aodFromExt = ext->GetAOD();
162       inputArray=(TClonesArray*)aodFromExt->GetList()->FindObject("D0toKpi");
163     }
164   } else if(pEvent) {
165     inputArray=(TClonesArray*)pEvent->GetList()->FindObject("D0toKpi");
166   }
167   if(!inputArray || !pEvent) {
168     AliError("Input branch not found!\n");
169     return;
170   }
171   // fix for temporary bug in ESDfilter
172   // the AODs with null vertex pointer didn't pass the PhysSel
173   if(!pEvent->GetPrimaryVertex() || TMath::Abs(pEvent->GetMagneticField())<0.001){
174     AliDebug(2,"Rejected at GetPrimaryvertex");
175     return;
176   }
177
178   fSelector->HistogramEventProperties(AliDxHFEParticleSelection::kEventsAll);
179   AliRDHFCuts* cuts=dynamic_cast<AliRDHFCuts*>(fCuts);
180   if (!cuts) return; // Fatal thrown already in initialization
181
182   if(!cuts->IsEventSelected(pEvent)) {
183     AliDebug(2,"rejected at IsEventSelected");
184     return;
185   }
186
187   //  Int_t nInD0toKpi = inputArray->GetEntriesFast();
188
189   fSelector->HistogramEventProperties(AliDxHFEParticleSelection::kEventsSel);
190              
191   std::auto_ptr<TObjArray> selectedTracks(fSelector->Select(inputArray,pEvent));
192   // TODO: use the array of selected track for something, right now
193   // only the control histograms of the selection class are filled
194   // note: the pointer is deleted automatically once the scope is left
195   // if the array should be published, the auto pointer must be released
196   // first, however some other cleanup will be necessary in that case
197   // probably a clone with a reduced AliVParticle implementation is
198   // appropriate.
199
200   if(! selectedTracks.get()) {
201     cout << "No selected D0s in this event" << endl; 
202     return;
203   }
204
205   //Test to see if I have read in D0s and retrieved them after selection
206   Int_t nD0Selected = selectedTracks->GetEntriesFast();
207
208   fSelector->HistogramEventProperties(AliDxHFEParticleSelection::kEventsWithParticle);
209
210   for(Int_t iD0toKpi = 0; iD0toKpi < nD0Selected; iD0toKpi++) {
211     AliAODRecoDecayHF2Prong *particle = (AliAODRecoDecayHF2Prong*)selectedTracks->UncheckedAt(iD0toKpi);
212     if (!particle) continue;
213   }
214
215   PostData(1, fOutput);
216 }
217
218 void AliAnalysisTaskDxHFEParticleSelection::FinishTaskOutput()
219 {
220   // end of the processing
221 }
222
223 void AliAnalysisTaskDxHFEParticleSelection::Terminate(Option_t *)
224 {
225   // last action on the client
226   fOutput = dynamic_cast<TList*> (GetOutputData(1));
227   if (!fOutput) {
228     AliFatal("failed to get output container");
229     return;
230   }
231
232 }