]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTaskSE.cxx
add histograms to check the mass of selected mesons
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisTaskSE.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 /* $Id$ */
17  
18 #include <TROOT.h>
19 #include <TSystem.h>
20 #include <TInterpreter.h>
21 #include <TChain.h>
22 #include <TFile.h>
23 #include <TList.h>
24
25 #include "AliAnalysisTaskSE.h"
26 #include "AliAnalysisManager.h"
27 #include "AliAnalysisCuts.h"
28 #include "AliAnalysisDataSlot.h"
29 #include "AliAnalysisDataContainer.h"
30
31 #include "AliESDEvent.h"
32 #include "AliESDfriend.h"
33 #include "AliESD.h"
34 #include "AliAODEvent.h"
35 #include "AliAODHeader.h"
36 #include "AliAODTracklets.h"
37 #include "AliAODCaloCells.h"
38 #include "AliAODCaloTrigger.h"
39 #include "AliAODMCParticle.h"
40 #include "AliVEvent.h"
41 #include "AliAODHandler.h"
42 #include "AliAODInputHandler.h"
43 #include "AliMCEventHandler.h"
44 #include "AliInputEventHandler.h"
45 #include "AliMultiInputEventHandler.h"
46 #include "AliESDInputHandler.h"
47 #include "AliMCEvent.h"
48 #include "AliStack.h"
49 #include "AliLog.h"
50 #include "AliAODDimuon.h"
51
52
53 ClassImp(AliAnalysisTaskSE)
54
55 ////////////////////////////////////////////////////////////////////////
56 AliAODHeader*    AliAnalysisTaskSE::fgAODHeader         = NULL;
57 TClonesArray*    AliAnalysisTaskSE::fgAODTracks         = NULL;
58 TClonesArray*    AliAnalysisTaskSE::fgAODVertices       = NULL;
59 TClonesArray*    AliAnalysisTaskSE::fgAODV0s            = NULL;
60 TClonesArray*    AliAnalysisTaskSE::fgAODPMDClusters    = NULL;
61 TClonesArray*    AliAnalysisTaskSE::fgAODJets           = NULL;
62 TClonesArray*    AliAnalysisTaskSE::fgAODFMDClusters    = NULL;
63 TClonesArray*    AliAnalysisTaskSE::fgAODCaloClusters   = NULL;
64 TClonesArray*    AliAnalysisTaskSE::fgAODMCParticles    = NULL;
65 AliAODTracklets* AliAnalysisTaskSE::fgAODTracklets      = NULL;
66 AliAODCaloCells* AliAnalysisTaskSE::fgAODEmcalCells     = NULL;
67 AliAODCaloCells* AliAnalysisTaskSE::fgAODPhosCells      = NULL;
68 AliAODCaloTrigger* AliAnalysisTaskSE::fgAODEMCALTrigger = NULL;
69 AliAODCaloTrigger* AliAnalysisTaskSE::fgAODPHOSTrigger  = NULL;
70 TClonesArray*    AliAnalysisTaskSE::fgAODDimuons        = NULL;
71 TClonesArray*    AliAnalysisTaskSE::fgAODHmpidRings     = NULL;
72
73 AliAnalysisTaskSE::AliAnalysisTaskSE():
74     AliAnalysisTask(),
75     fDebug(0),
76     fEntry(0),
77     fInputEvent(0x0),
78     fESDfriend(0x0),
79     fInputHandler(0x0),
80     fOutputAOD(0x0),
81     fMCEvent(0x0),
82     fTreeA(0x0),
83     fCurrentRunNumber(-1),
84     fHistosQA(0x0),
85     fOfflineTriggerMask(0),
86     fMultiInputHandler(0),
87     fMCEventHandler(0)
88 {
89   // Default constructor
90 }
91
92 AliAnalysisTaskSE::AliAnalysisTaskSE(const char* name):
93     AliAnalysisTask(name, "AnalysisTaskSE"),
94     fDebug(0),
95     fEntry(0),
96     fInputEvent(0x0),
97     fESDfriend(0x0),
98     fInputHandler(0x0),
99     fOutputAOD(0x0),
100     fMCEvent(0x0),
101     fTreeA(0x0),
102     fCurrentRunNumber(-1),
103     fHistosQA(0x0),
104     fOfflineTriggerMask(0),
105     fMultiInputHandler(0),
106     fMCEventHandler(0)
107 {
108   // Default constructor
109     DefineInput (0, TChain::Class());
110     DefineOutput(0,  TTree::Class());
111 }
112
113 AliAnalysisTaskSE::AliAnalysisTaskSE(const AliAnalysisTaskSE& obj):
114     AliAnalysisTask(obj),
115     fDebug(0),
116     fEntry(0),
117     fInputEvent(0x0),
118     fESDfriend(0x0),
119     fInputHandler(0x0),
120     fOutputAOD(0x0),
121     fMCEvent(0x0),
122     fTreeA(0x0),
123     fCurrentRunNumber(-1),
124     fHistosQA(0x0),
125     fOfflineTriggerMask(0),
126     fMultiInputHandler(obj.fMultiInputHandler),
127     fMCEventHandler(obj.fMCEventHandler)
128 {
129 // Copy constructor
130     fDebug            = obj.fDebug;
131     fEntry            = obj.fEntry;
132     fInputEvent       = obj.fInputEvent;
133     fESDfriend        = obj.fESDfriend;
134     fInputHandler     = obj.fInputHandler;
135     fOutputAOD        = obj.fOutputAOD;
136     fMCEvent          = obj.fMCEvent;
137     fTreeA            = obj.fTreeA;    
138     fCurrentRunNumber = obj.fCurrentRunNumber;
139     fHistosQA         = obj.fHistosQA;
140
141 }
142
143
144 AliAnalysisTaskSE& AliAnalysisTaskSE::operator=(const AliAnalysisTaskSE& other)
145 {
146 // Assignment
147   if(&other == this) return *this;
148   AliAnalysisTask::operator=(other);
149
150     AliAnalysisTask::operator=(other);
151     fDebug            = other.fDebug;
152     fEntry            = other.fEntry;
153     fInputEvent       = other.fInputEvent;
154     fESDfriend        = other.fESDfriend;
155     fInputHandler     = other.fInputHandler;
156     fOutputAOD        = other.fOutputAOD;
157     fMCEvent          = other.fMCEvent;
158     fTreeA            = other.fTreeA;    
159     fCurrentRunNumber = other.fCurrentRunNumber;
160     fHistosQA         = other.fHistosQA;
161     fOfflineTriggerMask = other.fOfflineTriggerMask;
162     fMultiInputHandler  = other.fMultiInputHandler;
163     fMCEventHandler     = other.fMCEventHandler;
164     return *this;
165 }
166
167 //______________________________________________________________________________
168 void AliAnalysisTaskSE::ConnectInputData(Option_t* /*option*/)
169 {
170 // Connect the input data
171     if (fDebug > 1) printf("AnalysisTaskSE::ConnectInputData() \n");
172
173    // Connect input handlers (multi input handler is handled)
174     ConnectMultiHandler();
175     
176     if (fInputHandler) {
177         if ((fInputHandler->GetTree())->GetBranch("ESDfriend."))
178             fESDfriend = ((AliESDInputHandler*)fInputHandler)->GetESDfriend();
179
180         fInputEvent = fInputHandler->GetEvent();
181     } else if( fMCEvent ) {
182          AliWarning("No Input Event Handler connected, only MC Truth Event Handler") ; 
183     } else {
184          AliError("No Input Event Handler connected") ; 
185          return ; 
186     }
187     // Disconnect multi handler
188     DisconnectMultiHandler();
189 }
190
191 void AliAnalysisTaskSE::CreateOutputObjects()
192 {
193 // Create the output container
194 //
195 //  Default AOD
196     if (fDebug > 1) printf("AnalysisTaskSE::CreateOutPutData() \n");
197
198     AliAODHandler* handler = dynamic_cast<AliAODHandler*> 
199          ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
200     
201     Bool_t merging = kFALSE;
202     AliAODInputHandler* aodIH = static_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
203     if (aodIH) {
204         if (aodIH->GetMergeEvents()) merging = kTRUE;
205     }
206
207     // Check if AOD replication has been required
208     if (handler) {
209         fOutputAOD   = handler->GetAOD();
210         fTreeA = handler->GetTree();
211         if (fOutputAOD && !(handler->IsStandard())) {
212             if ((handler->NeedsHeaderReplication() || merging) && !(fgAODHeader)) 
213                 {
214                  if (fDebug > 1) AliInfo("Replicating header");
215                  fgAODHeader = new AliAODHeader;
216                  handler->AddBranch("AliAODHeader", &fgAODHeader);
217                 }
218             if ((handler->NeedsTracksBranchReplication() || merging) && !(fgAODTracks))      
219             {   
220                 if (fDebug > 1) AliInfo("Replicating track branch\n");
221                 fgAODTracks = new TClonesArray("AliAODTrack",500);
222                 fgAODTracks->SetName("tracks");
223                 handler->AddBranch("TClonesArray", &fgAODTracks);
224             }    
225             if ((handler->NeedsVerticesBranchReplication() || merging) && !(fgAODVertices))
226             {
227                 if (fDebug > 1) AliInfo("Replicating vertices branch\n");
228                 fgAODVertices = new TClonesArray("AliAODVertex",500);
229                 fgAODVertices->SetName("vertices");
230                 handler->AddBranch("TClonesArray", &fgAODVertices);
231             }   
232             if ((handler->NeedsV0sBranchReplication()) && !(fgAODV0s))    
233             {   
234                 if (fDebug > 1) AliInfo("Replicating V0s branch\n");
235                 fgAODV0s = new TClonesArray("AliAODv0",500);
236                 fgAODV0s->SetName("v0s");
237                 handler->AddBranch("TClonesArray", &fgAODV0s);
238             }
239             if ((handler->NeedsTrackletsBranchReplication()) && !(fgAODTracklets))        
240             {   
241                 if (fDebug > 1) AliInfo("Replicating Tracklets branch\n");
242                 fgAODTracklets = new AliAODTracklets("tracklets","tracklets");
243                 handler->AddBranch("AliAODTracklets", &fgAODTracklets);
244             }
245             if ((handler->NeedsPMDClustersBranchReplication()) && !(fgAODPMDClusters))    
246             {   
247                 if (fDebug > 1) AliInfo("Replicating PMDClusters branch\n");
248                 fgAODPMDClusters = new TClonesArray("AliAODPmdCluster",500);
249                 fgAODPMDClusters->SetName("pmdClusters");
250                 handler->AddBranch("TClonesArray", &fgAODPMDClusters);
251             }
252             if ((handler->NeedsJetsBranchReplication() || merging) && !(fgAODJets))       
253             {   
254                 if (fDebug > 1) AliInfo("Replicating Jets branch\n");
255                 fgAODJets = new TClonesArray("AliAODJet",500);
256                 fgAODJets->SetName("jets");
257                 handler->AddBranch("TClonesArray", &fgAODJets);
258             }
259             if ((handler->NeedsFMDClustersBranchReplication()) && !(fgAODFMDClusters))    
260             {   
261                 AliInfo("Replicating FMDClusters branch\n");
262                 fgAODFMDClusters = new TClonesArray("AliAODFmdCluster",500);
263                 fgAODFMDClusters->SetName("fmdClusters");
264                 handler->AddBranch("TClonesArray", &fgAODFMDClusters);
265             }
266             if ((handler->NeedsCaloClustersBranchReplication() || merging) && !(fgAODCaloClusters))       
267             {   
268                 if (fDebug > 1) AliInfo("Replicating CaloClusters branch\n");
269                 fgAODCaloClusters = new TClonesArray("AliAODCaloCluster",500);
270                 fgAODCaloClusters->SetName("caloClusters");
271                 handler->AddBranch("TClonesArray", &fgAODCaloClusters);
272
273                 fgAODEmcalCells = new AliAODCaloCells("emcalCells","emcalCells",AliVCaloCells::kEMCALCell);
274                 handler->AddBranch("AliAODCaloCells", &fgAODEmcalCells);
275
276                 fgAODPhosCells = new AliAODCaloCells("phosCells","phosCells",AliVCaloCells::kPHOSCell);
277                 handler->AddBranch("AliAODCaloCells", &fgAODPhosCells);
278                 }
279             if ((handler->NeedsCaloTriggerBranchReplication() || merging) && !(fgAODEMCALTrigger))        
280             {   
281                 if (fDebug > 1) AliInfo("Replicating EMCAL Calo Trigger branches\n");
282                 fgAODEMCALTrigger = new AliAODCaloTrigger("emcalTrigger","emcalTrigger");
283                 handler->AddBranch("AliAODCaloTrigger", &fgAODEMCALTrigger);
284                 }
285                 if ((handler->NeedsCaloTriggerBranchReplication() || merging) && !(fgAODPHOSTrigger))     
286                 {   
287                 if (fDebug > 1) AliInfo("Replicating PHOS Calo Trigger branches\n");
288                 fgAODPHOSTrigger = new AliAODCaloTrigger("phosTrigger","phosTrigger");
289                 handler->AddBranch("AliAODCaloTrigger", &fgAODPHOSTrigger);
290             }
291             if ((handler->NeedsMCParticlesBranchReplication() || merging) && !(fgAODMCParticles))         
292             {   
293                 if (fDebug > 1) AliInfo("Replicating MCParticles branch\n");
294                 fgAODMCParticles = new TClonesArray("AliAODMCParticle",500);
295                 fgAODMCParticles->SetName("mcparticles");
296                 handler->AddBranch("TClonesArray", &fgAODMCParticles);
297             }
298             if ((handler->NeedsDimuonsBranchReplication() || merging) && !(fgAODDimuons))      
299             {   
300                 if (fDebug > 1) AliInfo("Replicating dimuon branch\n");
301                 fgAODDimuons = new TClonesArray("AliAODDimuon",0);
302                 fgAODDimuons->SetName("dimuons");
303                 handler->AddBranch("TClonesArray", &fgAODDimuons);
304             }    
305             if ((handler->NeedsHMPIDBranchReplication() || merging) && !(fgAODHmpidRings))      
306             {   
307                 if (fDebug > 1) AliInfo("Replicating HMPID branch\n");
308                 fgAODHmpidRings = new TClonesArray("AliAODHMPIDrings",0);
309                 fgAODHmpidRings->SetName("hmpidRings");
310                 handler->AddBranch("TClonesArray", &fgAODHmpidRings);
311             }
312                 
313
314             // cache the pointerd in the AODEvent
315             fOutputAOD->GetStdContent();
316         }
317     }
318     ConnectMultiHandler();
319     UserCreateOutputObjects();
320     DisconnectMultiHandler();
321 }
322
323 void AliAnalysisTaskSE::Exec(Option_t* option)
324 {
325 //
326 // Exec analysis of one event
327
328     ConnectMultiHandler();
329
330     if ( fDebug >= 10)
331       printf("Task is active %5d\n", IsActive());
332     
333     if (fDebug > 1) AliInfo("AliAnalysisTaskSE::Exec() \n");
334 //
335     AliAODHandler* handler = dynamic_cast<AliAODHandler*> 
336         ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
337
338     AliAODInputHandler* aodH = dynamic_cast<AliAODInputHandler*>(fInputHandler);
339 //
340 // Was event selected ? If no event selection mechanism, the event SHOULD be selected (AG)
341     UInt_t isSelected = AliVEvent::kAny;
342     if( fInputHandler && (fInputHandler->GetEventSelection() || aodH)) {
343       // Get the actual offline trigger mask for the event and AND it with the
344       // requested mask. If no mask requested select by default the event.
345       if (fOfflineTriggerMask)
346          isSelected = fOfflineTriggerMask & fInputHandler->IsEventSelected();
347     }
348 //  Functionality below moved in the filter tasks (AG)
349 //    if (handler) handler->SetFillAOD(isSelected);
350
351     if( fInputHandler ) {
352         fEntry = fInputHandler->GetReadEntry();
353         fESDfriend = ((AliESDInputHandler*)fInputHandler)->GetESDfriend();
354     }
355     
356
357 // Notify the change of run number
358     if (InputEvent() && (InputEvent()->GetRunNumber() != fCurrentRunNumber)) {
359         fCurrentRunNumber = InputEvent()->GetRunNumber();
360         NotifyRun();
361     }    
362            
363     else if( fMCEvent )
364        fEntry = fMCEvent->Header()->GetEvent(); 
365     if ( !((Entry()-1)%100) && fDebug > 0) 
366          AliInfo(Form("%s ----> Processing event # %lld", CurrentFileName(), Entry()));
367
368     
369     
370
371     if (handler && aodH) {
372         fMCEvent = aodH->MCEvent();
373         Bool_t merging = aodH->GetMergeEvents();
374       
375   // Do not analyze merged events if last embedded file has less events than normal event, 
376   // skip analysis after last embeded event 
377   if(merging){
378     if(aodH->GetReadEntry() + aodH->GetMergeOffset() >= aodH->GetTreeToMerge()->GetEntriesFast()){
379       //printf("Skip Entry %lld, Offset %d, Tree Entries %d\n",aodH->GetReadEntry(),aodH->GetMergeOffset(), aodH->GetTreeToMerge()->GetEntries());
380           
381       // Do I need to add the lines before the return?
382       // Added protection in case the derived task is not an AOD producer.
383       AliAnalysisDataSlot *out0 = GetOutputSlot(0);
384       if (out0 && out0->IsConnected()) PostData(0, fTreeA);    
385           
386       DisconnectMultiHandler();
387           
388       return;
389     }
390     //else   printf("MERGE Entry %lld, Offset %d, Tree Entries %d\n",aodH->GetReadEntry(),aodH->GetMergeOffset(), aodH->GetTreeToMerge()->GetEntries());
391   }
392       
393         AliAODEvent* aod = dynamic_cast<AliAODEvent*>(InputEvent());
394
395         if (aod && !(handler->IsStandard()) && !(handler->AODIsReplicated())) {
396             if ((handler->NeedsHeaderReplication() || merging) && (fgAODHeader))
397             {
398               // copy the contents by assigment
399               *fgAODHeader =  *(aod->GetHeader());
400             }
401             if ((handler->NeedsTracksBranchReplication() || (merging &&  aodH->GetMergeTracks())) && (fgAODTracks))
402             {
403                 TClonesArray* tracks = aod->GetTracks();
404                 new (fgAODTracks) TClonesArray(*tracks);
405             }
406             if ((handler->NeedsVerticesBranchReplication() || merging) && (fgAODVertices))
407             {
408                 TClonesArray* vertices = aod->GetVertices();
409                 new (fgAODVertices) TClonesArray(*vertices);
410             }
411             if ((handler->NeedsV0sBranchReplication()) && (fgAODV0s))
412             {
413                 TClonesArray* v0s = aod->GetV0s();
414                 new (fgAODV0s) TClonesArray(*v0s);
415             }
416             if ((handler->NeedsTrackletsBranchReplication()) && (fgAODTracklets))
417             {
418               *fgAODTracklets = *aod->GetTracklets();
419             }
420             if ((handler->NeedsPMDClustersBranchReplication()) && (fgAODPMDClusters))
421             {
422                 TClonesArray* pmdClusters = aod->GetPmdClusters();
423                 new (fgAODPMDClusters) TClonesArray(*pmdClusters);
424             }
425             if ((handler->NeedsJetsBranchReplication() || (merging &&aodH->GetMergeTracks())) && (fgAODJets))
426             {
427                 TClonesArray* jets = aod->GetJets();
428                 new (fgAODJets) TClonesArray(*jets);
429             }
430             if ((handler->NeedsFMDClustersBranchReplication()) && (fgAODFMDClusters))
431             {
432                 TClonesArray* fmdClusters = aod->GetFmdClusters();
433                 new (fgAODFMDClusters) TClonesArray(*fmdClusters);
434             }
435             if ((handler->NeedsCaloClustersBranchReplication() || 
436                  (merging && (aodH->GetMergeEMCALClusters() || aodH->GetMergePHOSClusters()))) 
437                 && (fgAODCaloClusters))
438             {
439                 TClonesArray* caloClusters = aod->GetCaloClusters();
440                 new (fgAODCaloClusters) TClonesArray(*caloClusters);
441             }
442
443             if ((handler->NeedsMCParticlesBranchReplication() || merging) && (fgAODMCParticles))
444             {
445                 TClonesArray* mcParticles = (TClonesArray*) (aod->FindListObject("mcparticles"));
446                 if( mcParticles )
447                   new (fgAODMCParticles) TClonesArray(*mcParticles);
448             }
449             
450             if ((handler->NeedsDimuonsBranchReplication() || (merging && aodH->GetMergeTracks())) && (fgAODDimuons))
451             {
452                 fgAODDimuons->Clear();
453                 TClonesArray& dimuons = *fgAODDimuons;
454                 TClonesArray& tracksnew = *fgAODTracks;
455                 
456                 Int_t nMuonTrack[100]; 
457                 for(Int_t imuon = 0; imuon < 100; imuon++) nMuonTrack[imuon] = 0;
458                 Int_t nMuons=0;
459                 for(Int_t ii=0; ii < fgAODTracks->GetEntries(); ii++){
460                     AliAODTrack *track = (AliAODTrack*) fgAODTracks->At(ii);
461                     if(track->IsMuonTrack()) {
462                         nMuonTrack[nMuons]= ii;
463                         nMuons++;
464                     }  
465                 }
466                 Int_t jDimuons=0;
467                 if(nMuons >= 2){
468                     for(Int_t i = 0; i < nMuons; i++){
469                         Int_t index0 = nMuonTrack[i];
470                         for(Int_t j = i+1; j < nMuons; j++){
471                             Int_t index1 = nMuonTrack[j];
472                             tracksnew.At(index0)->ResetBit(kIsReferenced);
473                             tracksnew.At(index0)->SetUniqueID(0); 
474                             tracksnew.At(index1)->ResetBit(kIsReferenced);
475                             tracksnew.At(index1)->SetUniqueID(0);
476                             new(dimuons[jDimuons++]) AliAODDimuon(tracksnew.At(index0),tracksnew.At(index1));
477                         }
478                     }    
479                 }
480             }
481             if ((handler->NeedsHMPIDBranchReplication()) && (fgAODHmpidRings))
482             {
483                 TClonesArray* hmpidRings = aod->GetHMPIDrings();
484                 new (fgAODHmpidRings) TClonesArray(*hmpidRings);
485             }
486             
487             
488             
489             // Additional merging if needed
490             if (merging) {
491               Int_t nc;
492
493                 // mcParticles
494                 TClonesArray* mcparticles = (TClonesArray*) ((aodH->GetEventToMerge())->FindListObject("mcparticles"));
495                 if( mcparticles ){
496                   Int_t npart = mcparticles->GetEntries();
497                   nc = fgAODMCParticles->GetEntries();
498                   
499                   for (Int_t i = 0; i < npart; i++) {
500                     AliAODMCParticle* particle = (AliAODMCParticle*) mcparticles->At(i);
501                     new((*fgAODMCParticles)[nc++]) AliAODMCParticle(*particle);
502                   }
503                 }
504
505                 // tracks
506                 TClonesArray* tracks = aodH->GetEventToMerge()->GetTracks();
507                 if(tracks && aodH->GetMergeTracks()){
508                   Int_t ntr = tracks->GetEntries();
509                   nc  = fgAODTracks->GetEntries();      
510                   for (Int_t i = 0; i < ntr; i++) {
511                     AliAODTrack*    track = (AliAODTrack*) tracks->At(i);
512                     AliAODTrack* newtrack = new((*fgAODTracks)[nc++]) AliAODTrack(*track);
513                     newtrack->SetLabel(newtrack->GetLabel() + fgAODMCParticles->GetEntries());
514                   }
515                   
516                   for (Int_t i = 0; i < nc; i++) 
517                     {
518                       AliAODTrack* track = (AliAODTrack*) fgAODTracks->At(i);
519                       track->ResetBit(kIsReferenced);
520                       track->SetUniqueID(0);
521                     }
522                 }
523                 
524                 // clusters
525                 TClonesArray* clusters = aodH->GetEventToMerge()->GetCaloClusters();
526                 if( clusters  && (aodH->GetMergeEMCALClusters() || aodH->GetMergePHOSClusters())) {
527                   Int_t ncl  = clusters->GetEntries();
528                   nc         =  fgAODCaloClusters->GetEntries();
529                   for (Int_t i = 0; i < ncl; i++) {
530                     AliAODCaloCluster*    cluster = (AliAODCaloCluster*) clusters->At(i);
531                     if(cluster->IsEMCAL() && !aodH->GetMergeEMCALClusters() ) continue;
532                     if(cluster->IsPHOS()  && !aodH->GetMergePHOSClusters()  ) continue;   
533                     new((*fgAODCaloClusters)[nc++]) AliAODCaloCluster(*cluster);
534                   }
535                 }
536
537                 // EMCAL cells
538                 //*fgAODEmcalCells =  *(aod->GetEMCALCells()); // This will be valid after 10.Mar.2011.
539                 if(aodH->GetMergeEMCALCells()) 
540                   {
541                     AliAODCaloCells* copycells = aod->GetEMCALCells();
542                     fgAODEmcalCells->CreateContainer(copycells->GetNumberOfCells());
543                     nc  = copycells->GetNumberOfCells();
544                     
545                     while( nc-- ){ fgAODEmcalCells->SetCell(nc,copycells->GetCellNumber(nc),copycells->GetAmplitude(nc),
546                                                             copycells->GetTime(nc),copycells->GetMCLabel(nc),copycells->GetEFraction(nc)); }
547                     
548                     AliAODCaloCells* cellsA = aodH->GetEventToMerge()->GetEMCALCells();
549                     if( cellsA )
550                       {
551                         Int_t ncells  = cellsA->GetNumberOfCells();
552                         nc = fgAODEmcalCells->GetNumberOfCells();
553                         
554                         for (Int_t i  = 0; i < ncells; i++) 
555                           {
556                             Int_t cn  = cellsA->GetCellNumber(i);
557                             Int_t pos = fgAODEmcalCells->GetCellPosition(cn);
558                             
559                             if (pos >= 0) 
560                               {
561                                 Double_t amp = cellsA->GetAmplitude(i) + fgAODEmcalCells->GetAmplitude(pos);
562                                 
563                                 //Check if it is MC, depending on that assing the mc lable, time and e fraction
564                                 Double_t time    = 0;
565                                 Int_t    mclabel =-1;
566                                 Double_t efrac   = 0;
567                                 if(cellsA->GetMCLabel(i) >= 0 && fgAODEmcalCells->GetMCLabel(i) < 0)
568                                   {
569                                     mclabel = cellsA->GetMCLabel(i) ;
570                                     time    = fgAODEmcalCells->GetTime(i) ; // Time from data
571                                     if(amp > 0) efrac = cellsA->GetAmplitude(i) / amp;
572                                   }
573                                 else if(fgAODEmcalCells->GetMCLabel(i) >= 0 &&  cellsA->GetMCLabel(i) < 0)
574                                   {
575                                     mclabel = fgAODEmcalCells->GetMCLabel(i) ;
576                                     time    = cellsA->GetTime(i) ; // Time from data
577                                     if(amp > 0) efrac = fgAODEmcalCells->GetAmplitude(i) / amp;              
578                                   }
579                                 else 
580                                   { // take all from input
581                                     mclabel = cellsA->GetMCLabel(i) ;
582                                     time    = cellsA->GetTime(i) ; 
583                                     if(amp > 0) efrac = cellsA->GetAmplitude(i) / amp;  
584                                   }
585                                 
586                                 fgAODEmcalCells->SetCell(pos, cn, amp,cellsA->GetTime(i),mclabel,efrac);
587                                 
588                               } else 
589                               {
590                                 AliAODCaloCells* copycells1 = new AliAODCaloCells(*fgAODEmcalCells);
591                                 fgAODEmcalCells->CreateContainer(nc+1);
592                                 Int_t nn = copycells1->GetNumberOfCells();
593                                 
594                                 while( nn-- ){ fgAODEmcalCells->SetCell(nn,copycells1->GetCellNumber(nn),copycells1->GetAmplitude(nn),
595                                                                         copycells1->GetTime(nn),copycells1->GetMCLabel(nn),0.); }
596                                 
597                                 fgAODEmcalCells->SetCell(nc++,cn,cellsA->GetAmplitude(i),cellsA->GetTime(i), cellsA->GetMCLabel(i),0.);
598                                 
599                                 delete copycells1;
600                               }
601                           }
602                         fgAODEmcalCells->Sort();
603                       }
604                   } // merge emcal cells
605                 
606                 
607                 // PHOS cells
608                 //*fgAODPhosCells =  *(aod->GetPHOSCells()); // This will be valid after 10.Mar.2011.
609                 if(aodH->GetMergePHOSCells()) 
610                   {
611                     AliAODCaloCells* copycells = aod->GetPHOSCells();
612                     fgAODPhosCells->CreateContainer(copycells->GetNumberOfCells());
613                     nc  = copycells->GetNumberOfCells();
614                     
615                     while( nc-- ){ fgAODPhosCells->SetCell(nc,copycells->GetCellNumber(nc),copycells->GetAmplitude(nc),
616                                                            copycells->GetTime(nc),copycells->GetMCLabel(nc),copycells->GetEFraction(nc)); }
617                     
618                     AliAODCaloCells* cellsP = aodH->GetEventToMerge()->GetPHOSCells();
619                     if( cellsP )
620                       {
621                         Int_t ncellsP  = cellsP->GetNumberOfCells();
622                         nc = fgAODPhosCells->GetNumberOfCells();
623                         
624                         for (Int_t i  = 0; i < ncellsP; i++) 
625                           {
626                             Int_t cn  = cellsP->GetCellNumber(i);
627                             Int_t pos = fgAODPhosCells->GetCellPosition(cn);
628                             
629                             if (pos >= 0) 
630                               {
631                                 Double_t amp = cellsP->GetAmplitude(i) + fgAODPhosCells->GetAmplitude(pos);
632                                 
633                                 //Check if it is MC, depending on that assing the mc lable, time and e fraction
634                                 Double_t time    = 0;
635                                 Int_t    mclabel =-1;
636                                 Double_t    efrac   = 0;
637                                 if(cellsP->GetMCLabel(i) >= 0 && fgAODPhosCells->GetMCLabel(i) < 0)
638                                   {
639                                     mclabel = cellsP->GetMCLabel(i) ;
640                                     time    = fgAODPhosCells->GetTime(i) ; // Time from data
641                                     if(amp > 0) efrac = cellsP->GetAmplitude(i) / amp;
642                                   }
643                                 else if(fgAODPhosCells->GetMCLabel(i) >= 0 &&  cellsP->GetMCLabel(i) < 0)
644                                   {
645                                     mclabel = fgAODPhosCells->GetMCLabel(i) ;
646                                     time    = cellsP->GetTime(i) ; // Time from data
647                                     if(amp > 0) efrac = fgAODPhosCells->GetAmplitude(i) / amp;              
648                                   }
649                                 else 
650                                   { // take all from input
651                                     mclabel = cellsP->GetMCLabel(i) ;
652                                     time    = cellsP->GetTime(i) ; 
653                                     if(amp > 0) efrac = cellsP->GetAmplitude(i) / amp;  
654                                   }
655                                 
656                                 fgAODPhosCells->SetCell(pos, cn, amp,cellsP->GetTime(i),mclabel,efrac);                
657                                 
658                               } else 
659                               {
660                                 AliAODCaloCells* copycells1 = new AliAODCaloCells(*fgAODPhosCells);
661                                 fgAODPhosCells->CreateContainer(nc+1);
662                                 Int_t nn = copycells1->GetNumberOfCells();
663                                 
664                                 while( nn-- ){ fgAODPhosCells->SetCell(nn,copycells1->GetCellNumber(nn),copycells1->GetAmplitude(nn), 
665                                                                        copycells1->GetTime(nn),copycells1->GetMCLabel(nn),0.); }
666                                 
667                                 fgAODPhosCells->SetCell(nc++,cn,cellsP->GetAmplitude(i),cellsP->GetTime(i), cellsP->GetMCLabel(i),0.);
668                                 
669                                 delete copycells1;
670                               }
671                           }
672                         fgAODPhosCells->Sort();
673                       }
674                   } // Merge PHOS Cells
675                 
676                 if (aodH->GetMergeEMCALTrigger() && aod->GetCaloTrigger("EMCAL")) 
677                 {
678                         Int_t   tsEMCAL[48][64], px, py, ts;
679                         Float_t foEMCAL[48][64], am;
680                         
681                         for (Int_t i = 0; i < 48; i++) for (Int_t j = 0; j < 64; j++) 
682                         {
683                                 tsEMCAL[i][j] = 0;
684                                 foEMCAL[i][j] = 0.;
685                         }
686       
687       AliAODCaloTrigger& trg0 = *(aod->GetCaloTrigger("EMCAL"));
688       trg0.Reset();
689       while (trg0.Next())
690       {
691         trg0.GetPosition(px, py);
692         
693         if (px > -1 && py > -1) 
694         {
695           trg0.GetL1TimeSum(ts);
696           if (ts > -1) tsEMCAL[px][py] += ts;
697           
698           trg0.GetAmplitude(am);
699           if (am > -1) foEMCAL[px][py] += am;
700         }
701       }
702       
703       AliAODCaloTrigger& trg1 = *((aodH->GetEventToMerge())->GetCaloTrigger("EMCAL"));
704       
705       trg1.Reset();
706       while (trg1.Next())
707       {
708         trg1.GetPosition(px, py);
709         
710         if (px > -1 && py > -1) 
711         {
712           trg1.GetL1TimeSum(ts);
713           if (ts > -1) tsEMCAL[px][py] += ts;
714           
715           trg1.GetAmplitude(am);
716           if (am > -1) foEMCAL[px][py] += am;
717         }
718       }
719       
720       int nEntries = 0;
721       for (Int_t i = 0; i < 48; i++) 
722         for (Int_t j = 0; j < 64; j++) 
723           if (tsEMCAL[i][j] || foEMCAL[i][j]) nEntries++;
724       
725       fgAODEMCALTrigger->Allocate(nEntries);
726       Int_t timesL0[10]; for (int i = 0; i < 10; i++) timesL0[i] = -1;
727       
728       for (Int_t i = 0; i < 48; i++) 
729         for (Int_t j = 0; j < 64; j++) 
730           if (tsEMCAL[i][j] || foEMCAL[i][j]) 
731             fgAODEMCALTrigger->Add(i, j, foEMCAL[i][j], -1., timesL0, 0, tsEMCAL[i][j], 0);
732     }
733                 
734                 if (aodH->GetMergePHOSTrigger()) 
735                 {
736                         // To be implemented by PHOS
737                 }
738             } // merging
739             
740             handler->SetAODIsReplicated();
741         }
742     }
743
744
745 // Call the user analysis    
746     if (!fMCEventHandler) {
747         if (isSelected) 
748             UserExec(option);
749     } else {
750         if (isSelected && (fMCEventHandler->InitOk())) 
751             UserExec(option);
752     }
753     
754 // Added protection in case the derived task is not an AOD producer.
755     AliAnalysisDataSlot *out0 = GetOutputSlot(0);
756     if (out0 && out0->IsConnected()) PostData(0, fTreeA);    
757
758     DisconnectMultiHandler();
759 }
760
761 const char* AliAnalysisTaskSE::CurrentFileName()
762 {
763 // Returns the current file name    
764     if( fInputHandler )
765       return fInputHandler->GetTree()->GetCurrentFile()->GetName();
766     else if( fMCEvent )
767       return ((AliMCEventHandler*) ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler()))->TreeK()->GetCurrentFile()->GetName();
768     else return "";
769 }
770
771 void AliAnalysisTaskSE::AddAODBranch(const char* cname, void* addobj, const char *fname)
772 {
773     // Add a new branch to the aod tree
774     AliAODHandler* handler = dynamic_cast<AliAODHandler*> 
775         ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
776     if (handler) {
777         handler->AddBranch(cname, addobj, fname);
778     }
779 }
780
781 Bool_t AliAnalysisTaskSE::IsStandardAOD() const
782 {
783 // Check if the output AOD handler is configured for standard or delta AOD.
784 // Users should first check that AODEvent() returns non-null.
785     AliAODHandler* handler = dynamic_cast<AliAODHandler*> 
786          ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
787     if (!handler) {
788        Error("IsStandardAOD", "No AOD handler. Please use AODEvent() to check this first");
789        return kTRUE;
790     }
791     return handler->IsStandard();   
792 }
793
794 Bool_t AliAnalysisTaskSE::Notify()
795 {
796     return (UserNotify());
797 }
798
799 const AliEventTag *AliAnalysisTaskSE::EventTag() const
800 {
801 // Returns tag for the current event, if any. The return value should always be checked by the user.
802    if (!fInputHandler) {
803       Error("EventTag", "Input handler not yet available. Call this in UserExec");
804       return NULL;
805    }
806    return fInputHandler->GetEventTag();
807 }
808
809 void AliAnalysisTaskSE::LoadBranches() const
810 {
811 // Load all branches declared in fBranchNames data member of the parent class.
812 // Should be called in UserExec.
813   if (!fInputHandler) {
814      Error("LoadBranches", "Input handler not available yet. Call this in UserExec");
815      return;
816   }
817   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
818   if (mgr->GetAutoBranchLoading()) return;
819   TString taskbranches;
820   GetBranches(fInputHandler->GetDataType(), taskbranches);
821   if (taskbranches.IsNull()) return;
822   TObjArray *arr = taskbranches.Tokenize(",");
823   TIter next(arr);
824   TObject *obj;
825   while ((obj=next())) mgr->LoadBranch(obj->GetName());
826 }
827
828
829 //_________________________________________________________________________________________________
830 void AliAnalysisTaskSE::ConnectMultiHandler()
831 {
832    //
833    // Connect MultiHandler
834    //
835    fInputHandler = (AliInputEventHandler *)((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
836    fMultiInputHandler = dynamic_cast<AliMultiInputEventHandler *>(fInputHandler);
837    if (fMultiInputHandler) {
838       fInputHandler = dynamic_cast<AliInputEventHandler *>(fMultiInputHandler->GetFirstInputEventHandler());
839       fMCEventHandler = dynamic_cast<AliMCEventHandler *>(fMultiInputHandler->GetFirstMCEventHandler());
840    } else { 
841       fMCEventHandler = dynamic_cast<AliMCEventHandler *>((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
842    }
843    if (fMCEventHandler) fMCEvent = fMCEventHandler->MCEvent();
844 }
845
846 //_________________________________________________________________________________________________
847 void AliAnalysisTaskSE::DisconnectMultiHandler()
848 {
849    //
850    // Disconnect MultiHandler
851    //
852    if (fMultiInputHandler) fInputHandler = fMultiInputHandler;
853 }