1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
20 #include <TInterpreter.h>
25 #include "AliAnalysisTaskSE.h"
26 #include "AliAnalysisManager.h"
27 #include "AliAnalysisCuts.h"
28 #include "AliAnalysisDataSlot.h"
29 #include "AliAnalysisDataContainer.h"
31 #include "AliESDEvent.h"
32 #include "AliESDfriend.h"
34 #include "AliAODEvent.h"
35 #include "AliAODHeader.h"
36 #include "AliAODVZERO.h"
37 #include "AliTOFHeader.h"
38 #include "AliAODTracklets.h"
39 #include "AliAODCaloCells.h"
40 #include "AliAODCaloTrigger.h"
41 #include "AliAODMCParticle.h"
42 #include "AliVEvent.h"
43 #include "AliAODHandler.h"
44 #include "AliAODInputHandler.h"
45 #include "AliMCEventHandler.h"
46 #include "AliInputEventHandler.h"
47 #include "AliMultiInputEventHandler.h"
48 #include "AliESDInputHandler.h"
49 #include "AliMCEvent.h"
52 #include "AliAODDimuon.h"
55 ClassImp(AliAnalysisTaskSE)
57 ////////////////////////////////////////////////////////////////////////
58 AliAODHeader* AliAnalysisTaskSE::fgAODHeader = NULL;
59 AliTOFHeader* AliAnalysisTaskSE::fgTOFHeader = NULL;
60 AliAODVZERO* AliAnalysisTaskSE::fgAODVZERO = NULL;
61 TClonesArray* AliAnalysisTaskSE::fgAODTracks = NULL;
62 TClonesArray* AliAnalysisTaskSE::fgAODVertices = NULL;
63 TClonesArray* AliAnalysisTaskSE::fgAODV0s = NULL;
64 TClonesArray* AliAnalysisTaskSE::fgAODPMDClusters = NULL;
65 TClonesArray* AliAnalysisTaskSE::fgAODJets = NULL;
66 TClonesArray* AliAnalysisTaskSE::fgAODFMDClusters = NULL;
67 TClonesArray* AliAnalysisTaskSE::fgAODCaloClusters = NULL;
68 TClonesArray* AliAnalysisTaskSE::fgAODMCParticles = NULL;
69 AliAODTracklets* AliAnalysisTaskSE::fgAODTracklets = NULL;
70 AliAODCaloCells* AliAnalysisTaskSE::fgAODEmcalCells = NULL;
71 AliAODCaloCells* AliAnalysisTaskSE::fgAODPhosCells = NULL;
72 AliAODCaloTrigger* AliAnalysisTaskSE::fgAODEMCALTrigger = NULL;
73 AliAODCaloTrigger* AliAnalysisTaskSE::fgAODPHOSTrigger = NULL;
74 TClonesArray* AliAnalysisTaskSE::fgAODDimuons = NULL;
75 TClonesArray* AliAnalysisTaskSE::fgAODHmpidRings = NULL;
77 AliAnalysisTaskSE::AliAnalysisTaskSE():
87 fCurrentRunNumber(-1),
89 fOfflineTriggerMask(0),
90 fMultiInputHandler(0),
93 // Default constructor
96 AliAnalysisTaskSE::AliAnalysisTaskSE(const char* name):
97 AliAnalysisTask(name, "AnalysisTaskSE"),
106 fCurrentRunNumber(-1),
108 fOfflineTriggerMask(0),
109 fMultiInputHandler(0),
112 // Default constructor
113 DefineInput (0, TChain::Class());
114 DefineOutput(0, TTree::Class());
117 AliAnalysisTaskSE::AliAnalysisTaskSE(const AliAnalysisTaskSE& obj):
118 AliAnalysisTask(obj),
127 fCurrentRunNumber(-1),
129 fOfflineTriggerMask(0),
130 fMultiInputHandler(obj.fMultiInputHandler),
131 fMCEventHandler(obj.fMCEventHandler)
136 fInputEvent = obj.fInputEvent;
137 fESDfriend = obj.fESDfriend;
138 fInputHandler = obj.fInputHandler;
139 fOutputAOD = obj.fOutputAOD;
140 fMCEvent = obj.fMCEvent;
142 fCurrentRunNumber = obj.fCurrentRunNumber;
143 fHistosQA = obj.fHistosQA;
148 AliAnalysisTaskSE& AliAnalysisTaskSE::operator=(const AliAnalysisTaskSE& other)
151 if(&other == this) return *this;
152 AliAnalysisTask::operator=(other);
154 AliAnalysisTask::operator=(other);
155 fDebug = other.fDebug;
156 fEntry = other.fEntry;
157 fInputEvent = other.fInputEvent;
158 fESDfriend = other.fESDfriend;
159 fInputHandler = other.fInputHandler;
160 fOutputAOD = other.fOutputAOD;
161 fMCEvent = other.fMCEvent;
162 fTreeA = other.fTreeA;
163 fCurrentRunNumber = other.fCurrentRunNumber;
164 fHistosQA = other.fHistosQA;
165 fOfflineTriggerMask = other.fOfflineTriggerMask;
166 fMultiInputHandler = other.fMultiInputHandler;
167 fMCEventHandler = other.fMCEventHandler;
171 //______________________________________________________________________________
172 void AliAnalysisTaskSE::ConnectInputData(Option_t* /*option*/)
174 // Connect the input data
175 if (fDebug > 1) printf("AnalysisTaskSE::ConnectInputData() \n");
177 // Connect input handlers (multi input handler is handled)
178 ConnectMultiHandler();
180 if (fInputHandler && fInputHandler->GetTree()) {
181 if (fInputHandler->GetTree()->GetBranch("ESDfriend."))
182 fESDfriend = ((AliESDInputHandler*)fInputHandler)->GetESDfriend();
184 fInputEvent = fInputHandler->GetEvent();
185 } else if( fMCEvent ) {
186 AliWarning("No Input Event Handler connected, only MC Truth Event Handler") ;
188 AliError("No Input Event Handler connected") ;
191 // Disconnect multi handler
192 DisconnectMultiHandler();
195 void AliAnalysisTaskSE::CreateOutputObjects()
197 // Create the output container
200 if (fDebug > 1) printf("AnalysisTaskSE::CreateOutPutData() \n");
202 AliAODHandler* handler = dynamic_cast<AliAODHandler*>
203 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
205 Bool_t merging = kFALSE;
206 AliAODInputHandler* aodIH = static_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
208 if (aodIH->GetMergeEvents()) merging = kTRUE;
211 // Check if AOD replication has been required
213 fOutputAOD = handler->GetAOD();
214 fTreeA = handler->GetTree();
215 if (fOutputAOD && !(handler->IsStandard())) {
216 if ((handler->NeedsHeaderReplication() || merging) && !(fgAODHeader))
218 if (fDebug > 1) AliInfo("Replicating header");
219 fgAODHeader = new AliAODHeader;
220 handler->AddBranch("AliAODHeader", &fgAODHeader);
222 if ((handler->NeedsTOFHeaderReplication() || merging) && !(fgTOFHeader))
224 if (fDebug > 1) AliInfo("Replicating TOFheader");
225 fgTOFHeader = new AliTOFHeader;
226 handler->AddBranch("AliTOFHeader", &fgTOFHeader);
228 if ((handler->NeedsVZEROReplication() || merging) && !(fgAODVZERO))
230 if (fDebug > 1) AliInfo("Replicating VZERO");
231 fgAODVZERO = new AliAODVZERO;
232 handler->AddBranch("AliAODVZERO", &fgAODVZERO);
235 if ((handler->NeedsTracksBranchReplication() || merging) && !(fgAODTracks))
237 if (fDebug > 1) AliInfo("Replicating track branch\n");
238 fgAODTracks = new TClonesArray("AliAODTrack",500);
239 fgAODTracks->SetName("tracks");
240 handler->AddBranch("TClonesArray", &fgAODTracks);
242 if ((handler->NeedsVerticesBranchReplication() || merging) && !(fgAODVertices))
244 if (fDebug > 1) AliInfo("Replicating vertices branch\n");
245 fgAODVertices = new TClonesArray("AliAODVertex",500);
246 fgAODVertices->SetName("vertices");
247 handler->AddBranch("TClonesArray", &fgAODVertices);
249 if ((handler->NeedsV0sBranchReplication()) && !(fgAODV0s))
251 if (fDebug > 1) AliInfo("Replicating V0s branch\n");
252 fgAODV0s = new TClonesArray("AliAODv0",500);
253 fgAODV0s->SetName("v0s");
254 handler->AddBranch("TClonesArray", &fgAODV0s);
256 if ((handler->NeedsTrackletsBranchReplication()) && !(fgAODTracklets))
258 if (fDebug > 1) AliInfo("Replicating Tracklets branch\n");
259 fgAODTracklets = new AliAODTracklets("tracklets","tracklets");
260 handler->AddBranch("AliAODTracklets", &fgAODTracklets);
262 if ((handler->NeedsPMDClustersBranchReplication()) && !(fgAODPMDClusters))
264 if (fDebug > 1) AliInfo("Replicating PMDClusters branch\n");
265 fgAODPMDClusters = new TClonesArray("AliAODPmdCluster",500);
266 fgAODPMDClusters->SetName("pmdClusters");
267 handler->AddBranch("TClonesArray", &fgAODPMDClusters);
269 if ((handler->NeedsJetsBranchReplication() || merging) && !(fgAODJets))
271 if (fDebug > 1) AliInfo("Replicating Jets branch\n");
272 fgAODJets = new TClonesArray("AliAODJet",500);
273 fgAODJets->SetName("jets");
274 handler->AddBranch("TClonesArray", &fgAODJets);
276 if ((handler->NeedsFMDClustersBranchReplication()) && !(fgAODFMDClusters))
278 AliInfo("Replicating FMDClusters branch\n");
279 fgAODFMDClusters = new TClonesArray("AliAODFmdCluster",500);
280 fgAODFMDClusters->SetName("fmdClusters");
281 handler->AddBranch("TClonesArray", &fgAODFMDClusters);
283 if ((handler->NeedsCaloClustersBranchReplication() || merging) && !(fgAODCaloClusters))
285 if (fDebug > 1) AliInfo("Replicating CaloClusters branch\n");
286 fgAODCaloClusters = new TClonesArray("AliAODCaloCluster",500);
287 fgAODCaloClusters->SetName("caloClusters");
288 handler->AddBranch("TClonesArray", &fgAODCaloClusters);
290 fgAODEmcalCells = new AliAODCaloCells("emcalCells","emcalCells",AliVCaloCells::kEMCALCell);
291 handler->AddBranch("AliAODCaloCells", &fgAODEmcalCells);
293 fgAODPhosCells = new AliAODCaloCells("phosCells","phosCells",AliVCaloCells::kPHOSCell);
294 handler->AddBranch("AliAODCaloCells", &fgAODPhosCells);
296 if ((handler->NeedsCaloTriggerBranchReplication() || merging) && !(fgAODEMCALTrigger))
298 if (fDebug > 1) AliInfo("Replicating EMCAL Calo Trigger branches\n");
299 fgAODEMCALTrigger = new AliAODCaloTrigger("emcalTrigger","emcalTrigger");
300 handler->AddBranch("AliAODCaloTrigger", &fgAODEMCALTrigger);
302 if ((handler->NeedsCaloTriggerBranchReplication() || merging) && !(fgAODPHOSTrigger))
304 if (fDebug > 1) AliInfo("Replicating PHOS Calo Trigger branches\n");
305 fgAODPHOSTrigger = new AliAODCaloTrigger("phosTrigger","phosTrigger");
306 handler->AddBranch("AliAODCaloTrigger", &fgAODPHOSTrigger);
308 if ((handler->NeedsMCParticlesBranchReplication() || merging) && !(fgAODMCParticles))
310 if (fDebug > 1) AliInfo("Replicating MCParticles branch\n");
311 fgAODMCParticles = new TClonesArray("AliAODMCParticle",500);
312 fgAODMCParticles->SetName("mcparticles");
313 handler->AddBranch("TClonesArray", &fgAODMCParticles);
315 if ((handler->NeedsDimuonsBranchReplication() || merging) && !(fgAODDimuons))
317 if (fDebug > 1) AliInfo("Replicating dimuon branch\n");
318 fgAODDimuons = new TClonesArray("AliAODDimuon",0);
319 fgAODDimuons->SetName("dimuons");
320 handler->AddBranch("TClonesArray", &fgAODDimuons);
322 if ((handler->NeedsHMPIDBranchReplication() || merging) && !(fgAODHmpidRings))
324 if (fDebug > 1) AliInfo("Replicating HMPID branch\n");
325 fgAODHmpidRings = new TClonesArray("AliAODHMPIDrings",0);
326 fgAODHmpidRings->SetName("hmpidRings");
327 handler->AddBranch("TClonesArray", &fgAODHmpidRings);
331 // cache the pointerd in the AODEvent
332 fOutputAOD->GetStdContent();
335 ConnectMultiHandler();
336 UserCreateOutputObjects();
337 DisconnectMultiHandler();
340 void AliAnalysisTaskSE::Exec(Option_t* option)
343 // Exec analysis of one event
345 ConnectMultiHandler();
346 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
347 if (mgr->GetDebugLevel() > 1) {
348 if (!mgr->GetTopTasks()->FindObject(this))
349 printf(" -> Executing sub-task %s\n", GetName());
352 printf("Task is active %5d\n", IsActive());
354 if (fDebug > 1) AliInfo("AliAnalysisTaskSE::Exec() \n");
356 AliAODHandler* handler = dynamic_cast<AliAODHandler*>(mgr->GetOutputEventHandler());
358 AliAODInputHandler* aodH = dynamic_cast<AliAODInputHandler*>(fInputHandler);
360 // Was event selected ? If no event selection mechanism, the event SHOULD be selected (AG)
361 UInt_t isSelected = AliVEvent::kAny;
362 if( fInputHandler && (fInputHandler->GetEventSelection() || aodH)) {
363 // Get the actual offline trigger mask for the event and AND it with the
364 // requested mask. If no mask requested select by default the event.
365 if (fOfflineTriggerMask)
366 isSelected = fOfflineTriggerMask & fInputHandler->IsEventSelected();
368 // Functionality below moved in the filter tasks (AG)
369 // if (handler) handler->SetFillAOD(isSelected);
371 if( fInputHandler ) {
372 fEntry = fInputHandler->GetReadEntry();
373 fESDfriend = ((AliESDInputHandler*)fInputHandler)->GetESDfriend();
377 // Notify the change of run number
378 if (InputEvent() && (InputEvent()->GetRunNumber() != fCurrentRunNumber)) {
379 fCurrentRunNumber = InputEvent()->GetRunNumber();
381 } else if( fMCEvent )
382 fEntry = fMCEvent->Header()->GetEvent();
384 if ( !((Entry()-1)%100) && fDebug > 0)
385 AliInfo(Form("%s ----> Processing event # %lld", CurrentFileName(), Entry()));
387 if (aodH) fMCEvent = aodH->MCEvent();
389 if (handler && aodH) {
390 Bool_t merging = aodH->GetMergeEvents();
391 // Do not analyze merged events if last embedded file has less events than normal event,
392 // skip analysis after last embeded event
394 if(aodH->GetReadEntry() + aodH->GetMergeOffset() >= aodH->GetTreeToMerge()->GetEntriesFast()){
395 // Do I need to add the lines before the return?
396 // Added protection in case the derived task is not an AOD producer.
397 AliAnalysisDataSlot *out0 = GetOutputSlot(0);
398 if (out0 && out0->IsConnected()) PostData(0, fTreeA);
399 DisconnectMultiHandler();
404 AliAODEvent* aod = dynamic_cast<AliAODEvent*>(InputEvent());
406 if (aod && !(handler->IsStandard()) && !(handler->AODIsReplicated())) {
407 if ((handler->NeedsHeaderReplication() || merging) && (fgAODHeader))
409 // copy the contents by assigment
410 *fgAODHeader = *(aod->GetHeader());
412 if ((handler->NeedsTOFHeaderReplication() || merging) && (fgTOFHeader))
414 if (aod->GetTOFHeader()) *fgTOFHeader = *(aod->GetTOFHeader());
416 if ((handler->NeedsVZEROReplication() || merging) && (fgAODVZERO) && aod->GetVZEROData())
418 *fgAODVZERO = *(aod->GetVZEROData());
421 if ((handler->NeedsTracksBranchReplication() || (merging && aodH->GetMergeTracks())) && (fgAODTracks))
423 TClonesArray* tracks = aod->GetTracks();
424 new (fgAODTracks) TClonesArray(*tracks);
426 if ((handler->NeedsVerticesBranchReplication() || merging) && (fgAODVertices))
428 TClonesArray* vertices = aod->GetVertices();
429 new (fgAODVertices) TClonesArray(*vertices);
431 if ((handler->NeedsV0sBranchReplication()) && (fgAODV0s))
433 TClonesArray* v0s = aod->GetV0s();
434 new (fgAODV0s) TClonesArray(*v0s);
436 if ((handler->NeedsTrackletsBranchReplication()) && (fgAODTracklets))
438 *fgAODTracklets = *aod->GetTracklets();
440 if ((handler->NeedsPMDClustersBranchReplication()) && (fgAODPMDClusters))
442 TClonesArray* pmdClusters = aod->GetPmdClusters();
443 new (fgAODPMDClusters) TClonesArray(*pmdClusters);
445 if ((handler->NeedsJetsBranchReplication() || (merging &&aodH->GetMergeTracks())) && (fgAODJets))
447 TClonesArray* jets = aod->GetJets();
448 new (fgAODJets) TClonesArray(*jets);
450 if ((handler->NeedsFMDClustersBranchReplication()) && (fgAODFMDClusters))
452 TClonesArray* fmdClusters = aod->GetFmdClusters();
453 new (fgAODFMDClusters) TClonesArray(*fmdClusters);
455 if ((handler->NeedsCaloClustersBranchReplication() ||
456 (merging && (aodH->GetMergeEMCALClusters() || aodH->GetMergePHOSClusters())))
457 && (fgAODCaloClusters))
459 TClonesArray* caloClusters = aod->GetCaloClusters();
460 new (fgAODCaloClusters) TClonesArray(*caloClusters);
463 if ((handler->NeedsMCParticlesBranchReplication() || merging) && (fgAODMCParticles))
465 TClonesArray* mcParticles = (TClonesArray*) (aod->FindListObject("mcparticles"));
467 new (fgAODMCParticles) TClonesArray(*mcParticles);
470 if ((handler->NeedsDimuonsBranchReplication() || (merging && aodH->GetMergeTracks())) && (fgAODDimuons))
472 fgAODDimuons->Clear();
473 TClonesArray& dimuons = *fgAODDimuons;
474 TClonesArray& tracksnew = *fgAODTracks;
476 Int_t nMuonTrack[100];
477 for(Int_t imuon = 0; imuon < 100; imuon++) nMuonTrack[imuon] = 0;
479 for(Int_t ii=0; ii < fgAODTracks->GetEntries(); ii++){
480 AliAODTrack *track = (AliAODTrack*) fgAODTracks->At(ii);
481 if(track->IsMuonTrack()) {
482 nMuonTrack[nMuons]= ii;
488 for(Int_t i = 0; i < nMuons; i++){
489 Int_t index0 = nMuonTrack[i];
490 for(Int_t j = i+1; j < nMuons; j++){
491 Int_t index1 = nMuonTrack[j];
492 tracksnew.At(index0)->ResetBit(kIsReferenced);
493 tracksnew.At(index0)->SetUniqueID(0);
494 tracksnew.At(index1)->ResetBit(kIsReferenced);
495 tracksnew.At(index1)->SetUniqueID(0);
496 new(dimuons[jDimuons++]) AliAODDimuon(tracksnew.At(index0),tracksnew.At(index1));
501 if ((handler->NeedsHMPIDBranchReplication()) && (fgAODHmpidRings))
503 TClonesArray* hmpidRings = aod->GetHMPIDrings();
504 new (fgAODHmpidRings) TClonesArray(*hmpidRings);
509 // Additional merging if needed
514 TClonesArray* mcparticles = (TClonesArray*) ((aodH->GetEventToMerge())->FindListObject("mcparticles"));
516 Int_t npart = mcparticles->GetEntries();
517 nc = fgAODMCParticles->GetEntries();
519 for (Int_t i = 0; i < npart; i++) {
520 AliAODMCParticle* particle = (AliAODMCParticle*) mcparticles->At(i);
521 new((*fgAODMCParticles)[nc++]) AliAODMCParticle(*particle);
526 TClonesArray* tracks = aodH->GetEventToMerge()->GetTracks();
527 if(tracks && aodH->GetMergeTracks()){
528 Int_t ntr = tracks->GetEntries();
529 nc = fgAODTracks->GetEntries();
530 for (Int_t i = 0; i < ntr; i++) {
531 AliAODTrack* track = (AliAODTrack*) tracks->At(i);
532 AliAODTrack* newtrack = new((*fgAODTracks)[nc++]) AliAODTrack(*track);
533 newtrack->SetLabel(newtrack->GetLabel() + fgAODMCParticles->GetEntries());
536 for (Int_t i = 0; i < nc; i++)
538 AliAODTrack* track = (AliAODTrack*) fgAODTracks->At(i);
539 track->ResetBit(kIsReferenced);
540 track->SetUniqueID(0);
545 TClonesArray* clusters = aodH->GetEventToMerge()->GetCaloClusters();
546 if( clusters && (aodH->GetMergeEMCALClusters() || aodH->GetMergePHOSClusters())) {
547 Int_t ncl = clusters->GetEntries();
548 nc = fgAODCaloClusters->GetEntries();
549 for (Int_t i = 0; i < ncl; i++) {
550 AliAODCaloCluster* cluster = (AliAODCaloCluster*) clusters->At(i);
551 if(cluster->IsEMCAL() && !aodH->GetMergeEMCALClusters() ) continue;
552 if(cluster->IsPHOS() && !aodH->GetMergePHOSClusters() ) continue;
553 new((*fgAODCaloClusters)[nc++]) AliAODCaloCluster(*cluster);
558 //*fgAODEmcalCells = *(aod->GetEMCALCells()); // This will be valid after 10.Mar.2011.
559 if(aodH->GetMergeEMCALCells())
561 AliAODCaloCells* copycells = aod->GetEMCALCells();
562 fgAODEmcalCells->CreateContainer(copycells->GetNumberOfCells());
563 nc = copycells->GetNumberOfCells();
565 while( nc-- ){ fgAODEmcalCells->SetCell(nc,copycells->GetCellNumber(nc),copycells->GetAmplitude(nc),
566 copycells->GetTime(nc),copycells->GetMCLabel(nc),copycells->GetEFraction(nc)); }
568 AliAODCaloCells* cellsA = aodH->GetEventToMerge()->GetEMCALCells();
571 Int_t ncells = cellsA->GetNumberOfCells();
572 nc = fgAODEmcalCells->GetNumberOfCells();
574 for (Int_t i = 0; i < ncells; i++)
576 Int_t cn = cellsA->GetCellNumber(i);
577 Int_t pos = fgAODEmcalCells->GetCellPosition(cn);
581 Double_t amp = cellsA->GetAmplitude(i) + fgAODEmcalCells->GetAmplitude(pos);
583 //Check if it is MC, depending on that assing the mc lable, time and e fraction
584 // Double_t time = 0;
587 if(cellsA->GetMCLabel(i) >= 0 && fgAODEmcalCells->GetMCLabel(pos) < 0)
589 mclabel = cellsA->GetMCLabel(i) ;
590 // time = fgAODEmcalCells->GetTime(pos) ; // Time from data
591 if(amp > 0) efrac = cellsA->GetAmplitude(i) / amp;
593 else if(fgAODEmcalCells->GetMCLabel(pos) >= 0 && cellsA->GetMCLabel(i) < 0)
595 mclabel = fgAODEmcalCells->GetMCLabel(pos) ;
596 // time = cellsA->GetTime(i) ; // Time from data
597 if(amp > 0) efrac = fgAODEmcalCells->GetAmplitude(pos) / amp;
600 { // take all from input
601 mclabel = cellsA->GetMCLabel(i) ;
602 // time = cellsA->GetTime(i) ;
603 if(amp > 0) efrac = cellsA->GetAmplitude(i) / amp;
606 fgAODEmcalCells->SetCell(pos, cn, amp,cellsA->GetTime(i),mclabel,efrac);
610 AliAODCaloCells* copycells1 = new AliAODCaloCells(*fgAODEmcalCells);
611 fgAODEmcalCells->CreateContainer(nc+1);
612 Int_t nn = copycells1->GetNumberOfCells();
614 while( nn-- ){ fgAODEmcalCells->SetCell(nn,copycells1->GetCellNumber(nn),copycells1->GetAmplitude(nn),
615 copycells1->GetTime(nn),copycells1->GetMCLabel(nn),copycells1->GetEFraction(nn)); }
617 fgAODEmcalCells->SetCell(nc++,cn,cellsA->GetAmplitude(i),cellsA->GetTime(i), cellsA->GetMCLabel(i),1.);
622 fgAODEmcalCells->Sort();
624 } // merge emcal cells
628 //*fgAODPhosCells = *(aod->GetPHOSCells()); // This will be valid after 10.Mar.2011.
629 if(aodH->GetMergePHOSCells())
631 AliAODCaloCells* copycells = aod->GetPHOSCells();
632 fgAODPhosCells->CreateContainer(copycells->GetNumberOfCells());
633 nc = copycells->GetNumberOfCells();
635 while( nc-- ){ fgAODPhosCells->SetCell(nc,copycells->GetCellNumber(nc),copycells->GetAmplitude(nc),
636 copycells->GetTime(nc),copycells->GetMCLabel(nc),copycells->GetEFraction(nc)); }
638 AliAODCaloCells* cellsP = aodH->GetEventToMerge()->GetPHOSCells();
641 Int_t ncellsP = cellsP->GetNumberOfCells();
642 nc = fgAODPhosCells->GetNumberOfCells();
644 for (Int_t i = 0; i < ncellsP; i++)
646 Int_t cn = cellsP->GetCellNumber(i);
647 Int_t pos = fgAODPhosCells->GetCellPosition(cn);
651 Double_t amp = cellsP->GetAmplitude(i) + fgAODPhosCells->GetAmplitude(pos);
653 //Check if it is MC, depending on that assing the mc lable, time and e fraction
654 // Double_t time = 0;
657 if(cellsP->GetMCLabel(i) >= 0 && fgAODPhosCells->GetMCLabel(pos) < 0)
659 mclabel = cellsP->GetMCLabel(i) ;
660 // time = fgAODPhosCells->GetTime(pos) ; // Time from data
661 if(amp > 0) efrac = cellsP->GetAmplitude(i) / amp;
663 else if(fgAODPhosCells->GetMCLabel(pos) >= 0 && cellsP->GetMCLabel(i) < 0)
665 mclabel = fgAODPhosCells->GetMCLabel(pos) ;
666 // time = cellsP->GetTime(i) ; // Time from data
667 if(amp > 0) efrac = fgAODPhosCells->GetAmplitude(pos) / amp;
670 { // take all from input
671 mclabel = cellsP->GetMCLabel(i) ;
672 // time = cellsP->GetTime(i) ;
673 if(amp > 0) efrac = cellsP->GetAmplitude(i) / amp;
676 fgAODPhosCells->SetCell(pos, cn, amp,cellsP->GetTime(i),mclabel,efrac);
680 AliAODCaloCells* copycells1 = new AliAODCaloCells(*fgAODPhosCells);
681 fgAODPhosCells->CreateContainer(nc+1);
682 Int_t nn = copycells1->GetNumberOfCells();
684 while( nn-- ){ fgAODPhosCells->SetCell(nn,copycells1->GetCellNumber(nn),copycells1->GetAmplitude(nn),
685 copycells1->GetTime(nn),copycells1->GetMCLabel(nn),copycells1->GetEFraction(nn)); }
687 fgAODPhosCells->SetCell(nc++,cn,cellsP->GetAmplitude(i),cellsP->GetTime(i), cellsP->GetMCLabel(i),1.);
692 fgAODPhosCells->Sort();
694 } // Merge PHOS Cells
696 if (aodH->GetMergeEMCALTrigger() && aod->GetCaloTrigger("EMCAL"))
698 Int_t tsEMCAL[48][64], px, py, ts;
699 Float_t foEMCAL[48][64], am;
701 for (Int_t i = 0; i < 48; i++) for (Int_t j = 0; j < 64; j++)
707 AliAODCaloTrigger& trg0 = *(aod->GetCaloTrigger("EMCAL"));
711 trg0.GetPosition(px, py);
713 if (px > -1 && py > -1)
715 trg0.GetL1TimeSum(ts);
716 if (ts > -1) tsEMCAL[px][py] += ts;
718 trg0.GetAmplitude(am);
719 if (am > -1) foEMCAL[px][py] += am;
723 AliAODCaloTrigger& trg1 = *((aodH->GetEventToMerge())->GetCaloTrigger("EMCAL"));
728 trg1.GetPosition(px, py);
730 if (px > -1 && py > -1)
732 trg1.GetL1TimeSum(ts);
733 if (ts > -1) tsEMCAL[px][py] += ts;
735 trg1.GetAmplitude(am);
736 if (am > -1) foEMCAL[px][py] += am;
741 for (Int_t i = 0; i < 48; i++)
742 for (Int_t j = 0; j < 64; j++)
743 if (tsEMCAL[i][j] || foEMCAL[i][j]) nEntries++;
745 fgAODEMCALTrigger->Allocate(nEntries);
746 Int_t timesL0[10]; for (int i = 0; i < 10; i++) timesL0[i] = -1;
748 for (Int_t i = 0; i < 48; i++)
749 for (Int_t j = 0; j < 64; j++)
750 if (tsEMCAL[i][j] || foEMCAL[i][j])
751 fgAODEMCALTrigger->Add(i, j, foEMCAL[i][j], -1., timesL0, 0, tsEMCAL[i][j], 0);
754 if (aodH->GetMergePHOSTrigger())
756 // To be implemented by PHOS
760 handler->SetAODIsReplicated();
765 // Call the user analysis
766 if (isSelected) UserExec(option);
768 // Added protection in case the derived task is not an AOD producer.
769 AliAnalysisDataSlot *out0 = GetOutputSlot(0);
770 if (out0 && out0->IsConnected()) PostData(0, fTreeA);
772 DisconnectMultiHandler();
775 const char* AliAnalysisTaskSE::CurrentFileName()
777 // Returns the current file name
779 return fInputHandler->GetTree()->GetCurrentFile()->GetName();
781 return ((AliMCEventHandler*) ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler()))->TreeK()->GetCurrentFile()->GetName();
785 void AliAnalysisTaskSE::AddAODBranch(const char* cname, void* addobj, const char *fname)
787 // Add a new branch to the aod tree
788 AliAODHandler* handler = dynamic_cast<AliAODHandler*>
789 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
791 handler->AddBranch(cname, addobj, fname);
795 Bool_t AliAnalysisTaskSE::IsStandardAOD() const
797 // Check if the output AOD handler is configured for standard or delta AOD.
798 // Users should first check that AODEvent() returns non-null.
799 AliAODHandler* handler = dynamic_cast<AliAODHandler*>
800 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
802 Error("IsStandardAOD", "No AOD handler. Please use AODEvent() to check this first");
805 return handler->IsStandard();
808 Bool_t AliAnalysisTaskSE::Notify()
810 return (UserNotify());
813 const AliEventTag *AliAnalysisTaskSE::EventTag() const
815 // Returns tag for the current event, if any. The return value should always be checked by the user.
816 if (!fInputHandler) {
817 Error("EventTag", "Input handler not yet available. Call this in UserExec");
820 return fInputHandler->GetEventTag();
823 void AliAnalysisTaskSE::LoadBranches() const
825 // Load all branches declared in fBranchNames data member of the parent class.
826 // Should be called in UserExec.
827 if (!fInputHandler) {
828 Error("LoadBranches", "Input handler not available yet. Call this in UserExec");
831 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
832 if (mgr->GetAutoBranchLoading()) return;
833 TString taskbranches;
834 GetBranches(fInputHandler->GetDataType(), taskbranches);
835 if (taskbranches.IsNull()) return;
836 TObjArray *arr = taskbranches.Tokenize(",");
839 while ((obj=next())) mgr->LoadBranch(obj->GetName());
844 //_________________________________________________________________________________________________
845 void AliAnalysisTaskSE::ConnectMultiHandler()
848 // Connect MultiHandler
850 fInputHandler = (AliInputEventHandler *)((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
851 fMultiInputHandler = dynamic_cast<AliMultiInputEventHandler *>(fInputHandler);
852 if (fMultiInputHandler) {
853 fInputHandler = dynamic_cast<AliInputEventHandler *>(fMultiInputHandler->GetFirstInputEventHandler());
854 fMCEventHandler = dynamic_cast<AliInputEventHandler *>(fMultiInputHandler->GetFirstMCEventHandler());
856 fMCEventHandler = dynamic_cast<AliInputEventHandler *>((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
858 if (fMCEventHandler) fMCEvent = fMCEventHandler->MCEvent();
861 //_________________________________________________________________________________________________
862 void AliAnalysisTaskSE::DisconnectMultiHandler()
865 // Disconnect MultiHandler
867 if (fMultiInputHandler) fInputHandler = fMultiInputHandler;