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