]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliAnalysisTaskReducedTree.cxx
o add Reduced Event
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliAnalysisTaskReducedTree.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 ///////////////////////////////////////////////////////////////////////////
17 //                                                                       //
18 //    Analysis task for creating a reduced data tree                     //
19 //                                                                       //
20 ///////////////////////////////////////////////////////////////////////////
21
22
23 #include <TChain.h>
24 #include <TH1D.h>
25 #include <TFile.h>
26
27 #include <AliCFContainer.h>
28 #include <AliInputEventHandler.h>
29 #include <AliESDInputHandler.h>
30 #include <AliAODInputHandler.h>
31 #include <AliAnalysisManager.h>
32 #include <AliVEvent.h>
33 #include <AliESDEvent.h>
34 #include <AliAODEvent.h>
35 #include <AliAODTrack.h>
36 #include <AliTriggerAnalysis.h>
37 #include <AliESDtrackCuts.h>
38 #include <AliVZDC.h>
39 #include <AliESDv0.h>
40 #include <AliESDv0Cuts.h>
41 #include <AliVCluster.h>
42 #include "AliDielectron.h"
43 #include "AliDielectronHistos.h"
44 #include "AliDielectronMC.h"
45 #include "AliDielectronVarManager.h"
46 #include "AliFlowTrackCuts.h"
47 #include "AliFlowBayesianPID.h"
48
49 #include "AliReducedEvent.h"
50 #include "AliAnalysisTaskReducedTree.h"
51
52 ClassImp(AliAnalysisTaskReducedTree)
53
54
55 //_________________________________________________________________________________
56 AliAnalysisTaskReducedTree::AliAnalysisTaskReducedTree() :
57   AliAnalysisTaskSE(),
58   fListDielectron(),
59   fListHistos(),
60   fSelectPhysics(kFALSE),
61   fTriggerMask(AliVEvent::kAny),
62   fRejectPileup(kFALSE),
63   fFillTrackInfo(kTRUE),
64   fFillDielectronInfo(kTRUE),
65   fFillV0Info(kTRUE),
66   fFillCaloClusterInfo(kTRUE),
67   fFillFriendInfo(kTRUE),
68   fEventFilter(0x0),
69   fTrackFilter(0x0),
70   fFlowTrackFilter(0x0),
71   fK0sCuts(0x0),
72   fLambdaCuts(0x0),
73   fK0sPionCuts(0x0),
74   fLambdaProtonCuts(0x0),
75   fLambdaPionCuts(0x0),
76   fK0sMassRange(),
77   fLambdaMassRange(),
78   fV0Histos(0x0),
79   fTreeFile(0x0),
80   fTree(0x0),
81   fFriendTreeFile(0x0),
82   fFriendTree(0x0),
83   fReducedEvent(0x0),
84   fReducedEventFriend(0x0)
85 {
86   //
87   // Constructor
88   //
89 }
90
91 //_________________________________________________________________________________
92 AliAnalysisTaskReducedTree::AliAnalysisTaskReducedTree(const char *name) :
93   AliAnalysisTaskSE(name),
94   fListDielectron(),
95   fListHistos(),
96   fSelectPhysics(kFALSE),
97   fTriggerMask(AliVEvent::kAny),
98   fRejectPileup(kFALSE),
99   fFillTrackInfo(kTRUE),
100   fFillDielectronInfo(kTRUE),
101   fFillV0Info(kTRUE),
102   fFillCaloClusterInfo(kTRUE),
103   fFillFriendInfo(kTRUE),
104   fEventFilter(0x0),
105   fTrackFilter(0x0),
106   fFlowTrackFilter(0x0),
107   fK0sCuts(0x0),
108   fLambdaCuts(0x0),
109   fK0sPionCuts(0x0),
110   fLambdaProtonCuts(0x0),
111   fLambdaPionCuts(0x0),
112   fK0sMassRange(),
113   fLambdaMassRange(),
114   fV0Histos(0x0),
115   fTreeFile(0x0),
116   fTree(0x0),
117   fFriendTreeFile(0x0),
118   fFriendTree(0x0),
119   fReducedEvent(0x0),
120   fReducedEventFriend(0x0)
121 {
122   //
123   // Constructor
124   //
125   fK0sMassRange[0] = 0.4; fK0sMassRange[1] = 0.6;
126   fLambdaMassRange[0] = 1.08; fLambdaMassRange[1] = 1.15;
127   
128   DefineInput(0,TChain::Class());
129   DefineOutput(1, TList::Class());   // QA histograms
130   DefineOutput(2, TTree::Class());   // reduced information tree
131   if(fFillFriendInfo) DefineOutput(3, TTree::Class());   // reduced information tree with friends
132   //DefineOutput(2, TTree::Class());   // reduced information tree with friends
133   
134   fListHistos.SetName("QAhistograms");
135   fListDielectron.SetOwner();
136   fListHistos.SetOwner(kFALSE);
137 }
138
139
140 //_________________________________________________________________________________
141 void AliAnalysisTaskReducedTree::UserCreateOutputObjects()
142 {
143   //
144   // Add all histogram manager histogram lists to the output TList
145   //
146
147   if (!fListHistos.IsEmpty() || fTree || fFriendTree) return; //already initialised
148
149   TIter nextDie(&fListDielectron);
150   AliDielectron *die=0;
151   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
152     die->Init();
153     if (die->GetHistogramList()) fListHistos.Add(const_cast<THashList*>(die->GetHistogramList()));
154   }  
155   if(fV0Histos) fListHistos.Add(const_cast<THashList*>(fV0Histos->GetHistogramList()));
156
157   if(fFillFriendInfo) {
158     fFriendTree = new TTree("DstFriendTree","Reduced ESD information");
159     fReducedEventFriend = new AliReducedEventFriend();
160     fFriendTree->Branch("Event",&fReducedEventFriend,16000,99);
161   }
162   
163   //fTreeFile = new TFile("dstTree.root", "RECREATE");
164   fTree = new TTree("DstTree","Reduced ESD information");
165   fReducedEvent = new AliReducedEvent("DstEvent");
166   fTree->Branch("Event",&fReducedEvent,16000,99);
167     
168   PostData(1, &fListHistos);
169   PostData(2, fTree);
170   if(fFillFriendInfo) PostData(3, fFriendTree);
171   //PostData(2, fFriendTree);
172 }
173
174 //_________________________________________________________________________________
175 void AliAnalysisTaskReducedTree::UserExec(Option_t *option)
176 {
177   //
178   // Main loop. Called for every event
179   //  
180   option = option;
181   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
182   Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
183   Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODInputHandler::Class();
184   
185   AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
186   if (!inputHandler) return;
187   
188   if ( inputHandler->GetPIDResponse() ){
189     AliDielectronVarManager::SetPIDResponse( inputHandler->GetPIDResponse() );
190   } else {
191     AliFatal("This task needs the PID response attached to the input event handler!");
192   }
193   
194   // Was event selected ?
195   UInt_t isSelected = AliVEvent::kAny;
196   if(fSelectPhysics && inputHandler){
197     if((isESD && inputHandler->GetEventSelection()) || isAOD){
198       isSelected = inputHandler->IsEventSelected();
199       isSelected&=fTriggerMask;
200     }
201   }
202
203   if(isSelected==0) {
204     return;
205   }
206
207   // fill event histograms before event filter
208   Double_t values[AliDielectronVarManager::kNMaxValues]={0};
209   AliDielectronVarManager::Fill(InputEvent(),values);
210   
211   TIter nextDie(&fListDielectron);
212   AliDielectron *die=0;
213   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
214     AliDielectronHistos *h=die->GetHistoManager();
215     if (h){
216       if (h->GetHistogramList()->FindObject("Event_noCuts"))
217         h->FillClass("Event_noCuts",AliDielectronVarManager::kNMaxValues,values);
218     }
219   }
220   nextDie.Reset();
221
222   //event filter
223   if (fEventFilter) {
224     if (!fEventFilter->IsSelected(InputEvent())) return;
225   }
226   
227   //pileup
228   if (fRejectPileup){
229     if (InputEvent()->IsPileupFromSPD(3,0.8,3.,2.,5.)) return;
230   }
231
232   //bz for AliKF
233   Double_t bz = InputEvent()->GetMagneticField();
234   AliKFParticle::SetField( bz );
235
236   //Process event in all AliDielectron instances
237   fReducedEvent->ClearEvent();
238   if(fFillFriendInfo) fReducedEventFriend->ClearEvent();
239   FillEventInfo();
240   if(fFillV0Info) FillV0PairInfo();
241   
242   Short_t idie=0;
243   if(fFillDielectronInfo) {
244     while((die=static_cast<AliDielectron*>(nextDie()))){
245       die->Process(InputEvent());
246       FillDielectronPairInfo(die, idie);
247       ++idie;
248     }
249   }
250   nextDie.Reset();
251   
252   if(fFillTrackInfo) FillTrackInfo();
253   if(fFillFriendInfo) FillFriendEventInfo();              // Q-vector calculation
254   
255   fTree->Fill();
256   if(fFillFriendInfo) fFriendTree->Fill();
257       
258   // if there are candidate pairs, add the information to the reduced tree
259   PostData(1, &fListHistos);
260   PostData(2, fTree);
261   if(fFillFriendInfo) PostData(3, fFriendTree);
262   //PostData(2, fFriendTree);
263 }
264
265
266 //_________________________________________________________________________________
267 void AliAnalysisTaskReducedTree::FillEventInfo() 
268 {
269   //
270   // fill reduced event information
271   //
272   AliVEvent* event = InputEvent();
273   // Was event selected ?
274   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
275   Bool_t isESD = (event->IsA()==AliESDEvent::Class());
276   Bool_t isAOD = (event->IsA()==AliAODEvent::Class());
277   
278   AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
279   UInt_t isSelected = AliVEvent::kAny;
280   if(inputHandler){
281     if((isESD && inputHandler->GetEventSelection()) || isAOD){
282       isSelected = inputHandler->IsEventSelected();
283       isSelected&=fTriggerMask;
284     }
285   }
286
287   Double_t values[AliDielectronVarManager::kNMaxValues];
288   AliDielectronVarManager::Fill(event, values);
289   
290   fReducedEvent->fRunNo       = event->GetRunNumber();
291   fReducedEvent->fBC          = event->GetBunchCrossNumber();
292   fReducedEvent->fTriggerMask = event->GetTriggerMask();
293   fReducedEvent->fIsPhysicsSelection = (isSelected!=0 ? kTRUE : kFALSE);
294   AliVVertex* eventVtx = 0x0;
295   if(isESD) eventVtx = const_cast<AliESDVertex*>((static_cast<AliESDEvent*>(event))->GetPrimaryVertexTracks());
296   if(isAOD) eventVtx = const_cast<AliAODVertex*>((static_cast<AliAODEvent*>(event))->GetPrimaryVertex());
297   if(eventVtx) {
298     fReducedEvent->fVtx[0] = (isESD ? ((AliESDVertex*)eventVtx)->GetXv() : ((AliAODVertex*)eventVtx)->GetX());
299     fReducedEvent->fVtx[1] = (isESD ? ((AliESDVertex*)eventVtx)->GetYv() : ((AliAODVertex*)eventVtx)->GetY());
300     fReducedEvent->fVtx[2] = (isESD ? ((AliESDVertex*)eventVtx)->GetZv() : ((AliAODVertex*)eventVtx)->GetZ());
301     fReducedEvent->fNVtxContributors = eventVtx->GetNContributors();
302   }
303   if(isESD) {
304     eventVtx = const_cast<AliESDVertex*>((static_cast<AliESDEvent*>(event))->GetPrimaryVertexTPC());
305     if(eventVtx) {
306       fReducedEvent->fVtxTPC[0] = ((AliESDVertex*)eventVtx)->GetXv();
307       fReducedEvent->fVtxTPC[1] = ((AliESDVertex*)eventVtx)->GetYv();
308       fReducedEvent->fVtxTPC[2] = ((AliESDVertex*)eventVtx)->GetZv();
309       fReducedEvent->fNVtxTPCContributors = eventVtx->GetNContributors();
310     }
311   }
312   
313   AliCentrality *centrality = event->GetCentrality();
314   if(centrality) {
315     fReducedEvent->fCentrality[0] = centrality->GetCentralityPercentile("V0M");
316     fReducedEvent->fCentrality[1] = centrality->GetCentralityPercentile("CL1");
317     fReducedEvent->fCentrality[2] = centrality->GetCentralityPercentile("TRK");
318     fReducedEvent->fCentrality[3] = centrality->GetCentralityPercentile("ZEMvsZDC");    
319     fReducedEvent->fCentQuality   = centrality->GetQuality();
320   }
321   
322   fReducedEvent->fNtracks[0] = event->GetNumberOfTracks();
323   fReducedEvent->fSPDntracklets = values[AliDielectronVarManager::kNaccTrckltsEsd10Corr];
324
325   AliVVZERO* vzero = event->GetVZEROData();
326   for(Int_t i=0;i<64;++i) 
327     fReducedEvent->fVZEROMult[i] = vzero->GetMultiplicity(i);  
328
329   AliESDZDC* zdc = (isESD ? (static_cast<AliESDEvent*>(event))->GetESDZDC() : 0x0);
330   if(zdc) {
331     for(Int_t i=0; i<4; ++i)  fReducedEvent->fZDCnEnergy[i]   = zdc->GetZN1TowerEnergy()[i];
332     for(Int_t i=4; i<8; ++i)  fReducedEvent->fZDCnEnergy[i]   = zdc->GetZN2TowerEnergy()[i-5];
333   }
334   
335   // EMCAL/PHOS clusters
336   if(fFillCaloClusterInfo) FillCaloClusters();
337   
338   // TODO FMD multiplicities
339   
340 }
341
342
343 //_________________________________________________________________________________
344 void AliAnalysisTaskReducedTree::FillCaloClusters() {
345   //
346   // Fill info about the calorimeter clusters
347   //
348   AliVEvent* event = InputEvent();
349   Int_t nclusters = event->GetNumberOfCaloClusters();
350
351   fReducedEvent->fNCaloClusters = 0;
352   for(Int_t iclus=0; iclus<nclusters; ++iclus) {
353     AliVCluster* cluster = event->GetCaloCluster(iclus);
354     
355     TClonesArray& clusters = *(fReducedEvent->fCaloClusters);
356     AliReducedCaloCluster *reducedCluster=new(clusters[fReducedEvent->fNCaloClusters]) AliReducedCaloCluster();
357     
358     reducedCluster->fType    = (cluster->IsEMCAL() ? AliReducedCaloCluster::kEMCAL : AliReducedCaloCluster::kPHOS);
359     reducedCluster->fEnergy  = cluster->E();
360     reducedCluster->fTrackDx = cluster->GetTrackDx();
361     reducedCluster->fTrackDz = cluster->GetTrackDz();
362     fReducedEvent->fNCaloClusters += 1;
363   }  // end loop over clusters
364 }
365
366
367 //_________________________________________________________________________________
368 void AliAnalysisTaskReducedTree::FillFriendEventInfo() {
369   //
370   // Fill event info into the friend tree
371   //
372   // Add here calculated Q-vector components from all detectors
373   for(Int_t idet=0; idet<AliReducedEventFriend::kNdetectors; ++idet) {
374     fReducedEvent->GetQvector(fReducedEventFriend->fQvector[idet], idet);
375     for(Int_t ih=0; ih<fgkNMaxHarmonics; ++ih)
376       fReducedEventFriend->fEventPlaneStatus[idet][ih] = AliReducedEventFriend::kRaw;
377   }
378 }
379
380
381 //_________________________________________________________________________________
382 void AliAnalysisTaskReducedTree::FillTrackInfo() 
383 {
384   //
385   // fill reduced track information
386   //
387   AliVEvent* event = InputEvent();
388   Bool_t isESD = (event->IsA()==AliESDEvent::Class());
389   Bool_t isAOD = (event->IsA()==AliAODEvent::Class());
390
391   Int_t ntracks=event->GetNumberOfTracks();
392   for(Int_t itrack=0; itrack<ntracks; ++itrack){
393     AliVParticle *particle=event->GetTrack(itrack);
394     //apply track cuts
395     if(fTrackFilter && !fTrackFilter->IsSelected(particle)) continue;
396     
397     TClonesArray& tracks = *(fReducedEvent->fTracks);
398     AliReducedTrack *reducedParticle=new(tracks[fReducedEvent->fNtracks[1]]) AliReducedTrack();
399         
400     Double_t values[AliDielectronVarManager::kNMaxValues];
401     AliDielectronVarManager::Fill(particle, values);
402     reducedParticle->fStatus        = (ULong_t)values[AliDielectronVarManager::kTrackStatus];
403     reducedParticle->fGlobalPhi     = values[AliDielectronVarManager::kPhi];
404     reducedParticle->fGlobalPt      = values[AliDielectronVarManager::kPt]*values[AliDielectronVarManager::kCharge];
405     reducedParticle->fGlobalEta     = values[AliDielectronVarManager::kEta];    
406     reducedParticle->fMomentumInner = values[AliDielectronVarManager::kPIn];
407     reducedParticle->fDCA[0]        = values[AliDielectronVarManager::kImpactParXY];
408     reducedParticle->fDCA[1]        = values[AliDielectronVarManager::kImpactParZ];
409     
410     reducedParticle->fITSclusterMap = values[AliDielectronVarManager::kITSclusterMap];
411     reducedParticle->fITSsignal     = values[AliDielectronVarManager::kITSsignal];
412     
413     reducedParticle->fTPCNcls      = (UChar_t)values[AliDielectronVarManager::kNclsTPC];
414     reducedParticle->fTPCNclsF     = (UChar_t)values[AliDielectronVarManager::kNFclsTPC];
415     reducedParticle->fTPCNclsIter1 = (UChar_t)values[AliDielectronVarManager::kNclsTPCiter1];
416     reducedParticle->fTPCsignal    = values[AliDielectronVarManager::kTPCsignal];
417     reducedParticle->fTPCnSig[0]   = values[AliDielectronVarManager::kTPCnSigmaEle];
418     reducedParticle->fTPCnSig[1]   = values[AliDielectronVarManager::kTPCnSigmaPio];
419     reducedParticle->fTPCnSig[2]   = values[AliDielectronVarManager::kTPCnSigmaKao];
420     reducedParticle->fTPCnSig[3]   = values[AliDielectronVarManager::kTPCnSigmaPro];
421     
422     reducedParticle->fTOFbeta      = values[AliDielectronVarManager::kTOFbeta];
423     reducedParticle->fTOFnSig[0]   = values[AliDielectronVarManager::kTOFnSigmaEle];
424     reducedParticle->fTOFnSig[1]   = values[AliDielectronVarManager::kTOFnSigmaPio];
425     reducedParticle->fTOFnSig[2]   = values[AliDielectronVarManager::kTOFnSigmaKao];
426     reducedParticle->fTOFnSig[3]   = values[AliDielectronVarManager::kTOFnSigmaPro];
427
428     reducedParticle->fTRDpid[0]    = values[AliDielectronVarManager::kTRDprobEle];
429     reducedParticle->fTRDpid[1]    = values[AliDielectronVarManager::kTRDprobPio];
430     
431     if(fFlowTrackFilter) {
432       // switch on the first bit if this particle should be used for the event plane
433       if(fFlowTrackFilter->IsSelected(particle)) reducedParticle->fFlags |= (1<<0);
434     }
435     
436     if(isESD){
437       AliESDtrack *track=static_cast<AliESDtrack*>(particle);
438       reducedParticle->fTrackId          = (UShort_t)track->GetID();
439       reducedParticle->fTPCCrossedRows   = (UChar_t)track->GetTPCCrossedRows();
440       reducedParticle->fTPCClusterMap    = EncodeTPCClusterMap(track);
441       const AliExternalTrackParam* tpcInner = track->GetTPCInnerParam();
442       reducedParticle->fTPCPhi        = (tpcInner ? tpcInner->Phi() : 0.0);
443       reducedParticle->fTPCPt         = (tpcInner ? tpcInner->Pt() : 0.0);
444       reducedParticle->fTPCEta        = (tpcInner ? tpcInner->Eta() : 0.0);
445       reducedParticle->fTRDntracklets[0] = track->GetTRDntracklets();
446       reducedParticle->fTRDntracklets[1] = track->GetTRDntrackletsPID();
447       if(track->IsEMCAL()) reducedParticle->fCaloClusterId = track->GetEMCALcluster();
448       if(track->IsPHOS()) reducedParticle->fCaloClusterId = track->GetPHOScluster();
449     }
450     if(isAOD) {
451       AliAODTrack *track=static_cast<AliAODTrack*>(particle);
452       reducedParticle->fTrackId = track->GetID(); 
453       if(track->IsEMCAL()) reducedParticle->fCaloClusterId = track->GetEMCALcluster();
454       if(track->IsPHOS()) reducedParticle->fCaloClusterId = track->GetPHOScluster();
455     }
456
457     fReducedEvent->fNtracks[1] += 1;
458   }
459 }
460
461
462 //_________________________________________________________________________________
463 void AliAnalysisTaskReducedTree::FillDielectronPairInfo(AliDielectron* die, Short_t iDie) 
464 {
465   //
466   // fill reduced pair information
467   //
468   Bool_t hasMC=AliDielectronMC::Instance()->HasMC();
469
470   for(Int_t iType=0; iType<3; ++iType) {
471     
472     const TObjArray* array = die->GetPairArray(iType);
473     if(!array || array->GetEntriesFast()==0) continue;
474     
475     for(Int_t iCandidate=0; iCandidate<array->GetEntriesFast(); ++iCandidate) {
476       AliDielectronPair* pair = (AliDielectronPair*)array->At(iCandidate);
477       Double_t values[AliDielectronVarManager::kNMaxValues];
478       AliDielectronVarManager::Fill(pair, values);
479       
480       TClonesArray& tracks = *(fReducedEvent->fCandidates);
481       AliReducedPair *reducedParticle= 
482          new (tracks[fReducedEvent->fNV0candidates[1]+fReducedEvent->fNDielectronCandidates]) AliReducedPair();
483       // !!! hardcoded flag for dielectron id 
484       reducedParticle->fCandidateId  = (iDie==0 ? AliReducedPair::kJpsiToEE : AliReducedPair::kPhiToKK);
485       reducedParticle->fPairType     = values[AliDielectronVarManager::kPairType];
486       reducedParticle->fLegIds[0]    = (UShort_t)(static_cast<AliVTrack*>(pair->GetFirstDaughter()))->GetID();
487       reducedParticle->fLegIds[1]    = (UShort_t)(static_cast<AliVTrack*>(pair->GetSecondDaughter()))->GetID();
488       reducedParticle->fMass[0]      = values[AliDielectronVarManager::kM];
489       reducedParticle->fMass[1]      = -999.;
490       reducedParticle->fMass[2]      = -999.;
491       reducedParticle->fPhi          = values[AliDielectronVarManager::kPhi];  // in the [-pi,pi] interval
492       if(reducedParticle->fPhi<0.0) reducedParticle->fPhi = 2.0*TMath::Pi() + reducedParticle->fPhi;  // converted to [0,2pi]
493       reducedParticle->fPt           = values[AliDielectronVarManager::kPt];
494       reducedParticle->fEta          = values[AliDielectronVarManager::kEta];
495       reducedParticle->fLxy          = values[AliDielectronVarManager::kPseudoProperTime];
496       reducedParticle->fLxyErr       = values[AliDielectronVarManager::kPseudoProperTimeErr];
497       reducedParticle->fOpeningAngle = values[AliDielectronVarManager::kOpeningAngle];     
498      
499       reducedParticle->fMCid         = 0;
500       if(hasMC) {
501         AliDielectronMC::Instance()->ConnectMCEvent();
502         const TObjArray* mcSignals = die->GetMCSignals();
503         for(Int_t iSig=0; iSig<mcSignals->GetEntries(); ++iSig) {
504           if(iSig>31) break;
505           AliDielectronMC *mc=AliDielectronMC::Instance();
506           if(mc->IsMCTruth(pair, (AliDielectronSignalMC*)mcSignals->At(iSig))) {
507             reducedParticle->fMCid = reducedParticle->fMCid | (1<<iSig);
508           }
509         }
510       }   // end if has MC
511       fReducedEvent->fNDielectronCandidates += 1;
512     }    // end loop over candidates
513   }    // end loop over pair type
514 }
515
516
517 //_________________________________________________________________________________
518 void AliAnalysisTaskReducedTree::FillV0PairInfo() 
519 {
520   //
521   // fill reduced pair information
522   //
523   AliESDEvent* esd = (AliESDEvent*)InputEvent();
524   const AliESDVertex *primaryVertex = esd->GetPrimaryVertex();
525   AliKFVertex primaryVertexKF(*primaryVertex);
526   
527   fReducedEvent->fNV0candidates[0] = InputEvent()->GetNumberOfV0s();
528   
529   Double_t valuesPos[AliDielectronVarManager::kNMaxValues];
530   Double_t valuesNeg[AliDielectronVarManager::kNMaxValues];
531    
532   for(Int_t iV0=0; iV0<InputEvent()->GetNumberOfV0s(); ++iV0) {   // loop over V0s
533     AliESDv0 *v0 = esd->GetV0(iV0);
534        
535     AliESDtrack* legPos = esd->GetTrack(v0->GetPindex());
536     AliESDtrack* legNeg = esd->GetTrack(v0->GetNindex());
537  
538     if(legPos->GetSign() == legNeg->GetSign()) {
539       continue;
540     }
541
542     Bool_t v0ChargesAreCorrect = (legPos->GetSign()==+1 ? kTRUE : kFALSE);
543     legPos = (!v0ChargesAreCorrect ? esd->GetTrack(v0->GetNindex()) : legPos);
544     legNeg = (!v0ChargesAreCorrect ? esd->GetTrack(v0->GetPindex()) : legNeg); 
545     
546     // apply the K0s filter
547     Bool_t goodK0s = kTRUE;
548     if(fK0sPionCuts && (!fK0sPionCuts->IsSelected(legPos) || !fK0sPionCuts->IsSelected(legNeg))) goodK0s = kFALSE;
549     if(fK0sCuts) {
550       TList k0sList;
551       k0sList.Add(v0);
552       k0sList.Add(legPos); k0sList.Add(legNeg);
553       k0sList.Add(const_cast<AliESDVertex*>(primaryVertex));
554       if(!fK0sCuts->IsSelected(&k0sList)) goodK0s = kFALSE;
555     }
556     
557     // apply the lambda filter
558     Bool_t goodLambda = kTRUE;
559     if(fLambdaProtonCuts && !fLambdaProtonCuts->IsSelected(legPos)) goodLambda = kFALSE;
560     if(fLambdaPionCuts && !fLambdaPionCuts->IsSelected(legNeg)) goodLambda = kFALSE;
561     Bool_t goodALambda = kTRUE;
562     if(fLambdaProtonCuts && !fLambdaProtonCuts->IsSelected(legNeg)) goodALambda = kFALSE;
563     if(fLambdaPionCuts && !fLambdaPionCuts->IsSelected(legPos)) goodALambda = kFALSE;
564     if(fLambdaCuts) {
565       TList lambdaList;
566       lambdaList.Add(v0);
567       lambdaList.Add(legPos); lambdaList.Add(legNeg);
568       lambdaList.Add(const_cast<AliESDVertex*>(primaryVertex));
569       if(!fLambdaCuts->IsSelected(&lambdaList)) {goodLambda = kFALSE; goodALambda = kFALSE;}
570     }
571     
572     if(!(goodK0s || goodLambda || goodALambda)) continue;
573     
574     // Fill the V0 information into the tree for 3 hypothesis: K0s, Lambda, Anti-Lambda
575     AliReducedPair* k0sReducedPair     = FillV0PairInfo(v0, AliReducedPair::kK0sToPiPi,     legPos, legNeg, &primaryVertexKF, v0ChargesAreCorrect);
576     AliReducedPair* lambdaReducedPair  = FillV0PairInfo(v0, AliReducedPair::kLambda0ToPPi,  legPos, legNeg, &primaryVertexKF, v0ChargesAreCorrect);
577     AliReducedPair* alambdaReducedPair = FillV0PairInfo(v0, AliReducedPair::kALambda0ToPPi, legPos, legNeg, &primaryVertexKF, v0ChargesAreCorrect);
578
579     if(goodK0s && k0sReducedPair->fMass[0]>fK0sMassRange[0] && k0sReducedPair->fMass[0]<fK0sMassRange[1]) {
580       TClonesArray& tracks = *(fReducedEvent->fCandidates);
581       AliReducedPair *goodK0sPair = new (tracks[fReducedEvent->fNV0candidates[1]]) AliReducedPair(*k0sReducedPair);
582       goodK0sPair->fMass[0] = k0sReducedPair->fMass[0];
583       goodK0sPair->fMass[1] = lambdaReducedPair->fMass[0];
584       goodK0sPair->fMass[2] = alambdaReducedPair->fMass[0];
585       fReducedEvent->fNV0candidates[1] += 1;
586     } else {goodK0s=kFALSE;}
587     if(goodLambda && lambdaReducedPair->fMass[0]>fLambdaMassRange[0] && lambdaReducedPair->fMass[0]<fLambdaMassRange[1]) {
588       TClonesArray& tracks = *(fReducedEvent->fCandidates);
589       AliReducedPair *goodLambdaPair = new (tracks[fReducedEvent->fNV0candidates[1]]) AliReducedPair(*lambdaReducedPair);
590       fReducedEvent->fNV0candidates[1] += 1;
591       goodLambdaPair->fMass[0] = k0sReducedPair->fMass[0];
592       goodLambdaPair->fMass[1] = lambdaReducedPair->fMass[0];
593       goodLambdaPair->fMass[2] = alambdaReducedPair->fMass[0];
594     } else {goodLambda=kFALSE;}
595     if(goodALambda && alambdaReducedPair->fMass[0]>fLambdaMassRange[0] && alambdaReducedPair->fMass[0]<fLambdaMassRange[1]) {
596       TClonesArray& tracks = *(fReducedEvent->fCandidates);
597       AliReducedPair *goodALambdaPair = new (tracks[fReducedEvent->fNV0candidates[1]]) AliReducedPair(*alambdaReducedPair);
598       fReducedEvent->fNV0candidates[1] += 1;  
599       goodALambdaPair->fMass[0] = k0sReducedPair->fMass[0];
600       goodALambdaPair->fMass[1] = lambdaReducedPair->fMass[0];
601       goodALambdaPair->fMass[2] = alambdaReducedPair->fMass[0];
602     } else {goodALambda = kFALSE;}
603     delete k0sReducedPair;
604     delete lambdaReducedPair;
605     delete alambdaReducedPair;
606
607     if(!(goodK0s || goodLambda || goodALambda)) continue;
608     
609     //  Fill histograms and the CF container
610     AliDielectronVarManager::Fill(legPos, valuesPos);
611     AliDielectronVarManager::Fill(legNeg, valuesNeg);
612     
613     if(fV0Histos && fV0Histos->GetHistogramList()->FindObject("V0Track_Pos"))
614       fV0Histos->FillClass("V0Track_Pos", AliDielectronVarManager::kNMaxValues, valuesPos);
615     if(fV0Histos && fV0Histos->GetHistogramList()->FindObject("V0Track_Neg"))
616       fV0Histos->FillClass("V0Track_Neg", AliDielectronVarManager::kNMaxValues, valuesNeg);    
617   }   // end loop over V0s
618 }
619
620
621 //_________________________________________________________________________________
622 AliReducedPair* AliAnalysisTaskReducedTree::FillV0PairInfo(AliESDv0* v0, Int_t id, 
623                                                     AliESDtrack* legPos, AliESDtrack* legNeg,
624                                                     AliKFVertex* vtxKF, Bool_t chargesAreCorrect) {
625   //
626   // Create a reduced V0 object and fill it
627   //
628   AliReducedPair* reducedPair=new AliReducedPair();  
629   reducedPair->fCandidateId = id;
630   reducedPair->fPairType    = v0->GetOnFlyStatus();    // on the fly status
631   //reducedPair->fOnTheFly    = v0->GetOnFlyStatus();
632   reducedPair->fLegIds[0]   = legPos->GetID();
633   reducedPair->fLegIds[1]   = legNeg->GetID();
634   if(!reducedPair->fPairType) {    // offline
635     UInt_t pidPos = AliPID::kPion;
636     if(id==AliReducedPair::kLambda0ToPPi) pidPos = AliPID::kProton;
637     UInt_t pidNeg = AliPID::kPion;
638     if(id==AliReducedPair::kALambda0ToPPi) pidNeg = AliPID::kProton;
639     reducedPair->fMass[0]      = v0->GetEffMass(pidPos, pidNeg);
640     reducedPair->fPhi          = v0->Phi();
641     if(reducedPair->fPhi<0.0) reducedPair->fPhi = 2.0*TMath::Pi() + reducedPair->fPhi;  // converted to [0,2pi]
642     reducedPair->fPt           = v0->Pt();
643     reducedPair->fEta          = v0->Eta();
644     reducedPair->fLxy          = 0.0;           //TODO
645     reducedPair->fOpeningAngle = 0.0;
646   }
647   else {
648     const AliExternalTrackParam *negHelix=v0->GetParamN();
649     const AliExternalTrackParam *posHelix=v0->GetParamP();
650     if(!chargesAreCorrect) {
651       negHelix = v0->GetParamP();
652       posHelix = v0->GetParamN();
653     }
654     AliKFParticle negKF(*(negHelix),(id==AliReducedPair::kALambda0ToPPi ? -2212 : -211));
655     AliKFParticle posKF(*(posHelix),(id==AliReducedPair::kLambda0ToPPi ? +2212 : +211));
656     AliKFParticle v0Refit;
657     v0Refit += negKF;
658     v0Refit += posKF;
659     Double_t massFit=0.0, massErrFit=0.0;
660     v0Refit.GetMass(massFit,massErrFit);
661     reducedPair->fMass[0] = massFit;
662     reducedPair->fPhi          = v0Refit.GetPhi();
663     if(reducedPair->fPhi<0.0) reducedPair->fPhi = 2.0*TMath::Pi() + reducedPair->fPhi;  // converted to [0,2pi]
664     reducedPair->fPt           = v0Refit.GetPt();
665     reducedPair->fEta          = v0Refit.GetEta();
666     reducedPair->fLxy          = v0Refit.GetPseudoProperDecayTime(*vtxKF, massFit);
667     reducedPair->fOpeningAngle = negKF.GetAngle(posKF);
668   }
669   return reducedPair;
670 }
671
672
673 //_________________________________________________________________________________
674 UChar_t AliAnalysisTaskReducedTree::EncodeTPCClusterMap(AliESDtrack* track) {
675   //
676   // Encode the TPC cluster map into an UChar_t
677   // Divide the 159 bits from the bit map into 8 groups of adiacent clusters
678   // For each group enable its corresponding bit if in that group there are more clusters compared to
679   // a threshold.
680   //
681   const UChar_t threshold=5;
682   TBits tpcClusterMap = track->GetTPCClusterMap();
683   UChar_t map=0;
684   UChar_t n=0;
685   UChar_t j=0;
686   for(UChar_t i=0; i<8; ++i) {
687     n=0;
688     for(j=i*20; j<(i+1)*20 && j<159; ++j) n+=tpcClusterMap.TestBitNumber(j);
689     if(n>=threshold) map |= (1<<i);
690   }
691   return map;
692 }
693
694
695 //_________________________________________________________________________________
696 void AliAnalysisTaskReducedTree::FinishTaskOutput()
697 {
698   //
699   // Finish Task 
700   //
701   //fTreeFile->Write();
702   //fTreeFile->Close();
703 }