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 **************************************************************************/
16 //_________________________________________________________________________
17 // This analysis provides a new list of clusters to be used in other analysis
19 // Author: Gustavo Conesa Balbastre,
20 // Adapted from analysis class from Deepa Thomas
23 //_________________________________________________________________________
27 #include <TRefArray.h>
28 #include <TClonesArray.h>
30 #include <TGeoManager.h>
32 #include <TInterpreter.h>
35 // --- AliRoot Analysis Steering
36 #include "AliAnalysisTask.h"
37 #include "AliAnalysisManager.h"
38 #include "AliESDEvent.h"
39 #include "AliGeomManager.h"
40 #include "AliVCaloCells.h"
41 #include "AliAODCaloCluster.h"
42 #include "AliCDBManager.h"
43 #include "AliCDBStorage.h"
44 #include "AliCDBEntry.h"
46 #include "AliVEventHandler.h"
47 #include "AliAODInputHandler.h"
48 #include "AliOADBContainer.h"
51 #include "AliEMCALAfterBurnerUF.h"
52 #include "AliEMCALGeometry.h"
53 #include "AliEMCALClusterizerNxN.h"
54 #include "AliEMCALClusterizerv1.h"
55 #include "AliEMCALClusterizerv2.h"
56 #include "AliEMCALRecPoint.h"
57 #include "AliEMCALDigit.h"
58 #include "AliCaloCalibPedestal.h"
59 #include "AliEMCALCalibData.h"
61 #include "AliAnalysisTaskEMCALClusterize.h"
63 ClassImp(AliAnalysisTaskEMCALClusterize)
65 //______________________________________________________________________________
66 AliAnalysisTaskEMCALClusterize::AliAnalysisTaskEMCALClusterize(const char *name)
67 : AliAnalysisTaskSE(name)
69 , fGeom(0), fGeomName("")
70 , fGeomMatrixSet(kFALSE), fLoadGeomMatrices(kFALSE)
71 , fCalibData(0), fPedestalData(0)
72 , fOCDBpath(""), fAccessOCDB(kFALSE)
73 , fDigitsArr(0), fClusterArr(0), fCaloClusterArr(0)
74 , fRecParam(0), fClusterizer(0)
75 , fUnfolder(0), fJustUnfold(kFALSE)
76 , fOutputAODBranch(0), fOutputAODBranchName(""), fOutputAODBranchSet(0)
77 , fFillAODFile(kFALSE), fFillAODHeader(0)
78 , fFillAODCaloCells(0), fRun(-1)
79 , fRecoUtils(0), fConfigName("")
80 , fCellLabels(), fCellSecondLabels(), fCellTime()
81 , fCellMatchdEta(), fCellMatchdPhi()
82 , fMaxEvent(0), fDoTrackMatching(kFALSE)
83 , fSelectCell(kFALSE), fSelectCellMinE(0), fSelectCellMinFrac(0)
84 , fRemoveLEDEvents(kTRUE),fRemoveExoticEvents(kFALSE)
85 , fImportGeometryFromFile(kFALSE), fImportGeometryFilePath("")
86 , fOADBSet(kFALSE), fAccessOADB(kTRUE), fOADBFilePath("")
87 , fCentralityClass("")
91 for(Int_t i = 0; i < 12; i++) fGeomMatrix[i] = 0;
92 for(Int_t j = 0; j < 24*48*11; j++)
95 fCellSecondLabels[j] = -1;
97 fCellMatchdEta[j] = -999;
98 fCellMatchdPhi[j] = -999;
101 fCentralityBin[0] = fCentralityBin[1]=-1;
105 //______________________________________________________________
106 AliAnalysisTaskEMCALClusterize::AliAnalysisTaskEMCALClusterize()
107 : AliAnalysisTaskSE("DefaultAnalysis_AliAnalysisTaskEMCALClusterize")
109 , fGeom(0), fGeomName("")
110 , fGeomMatrixSet(kFALSE), fLoadGeomMatrices(kFALSE)
111 , fCalibData(0), fPedestalData(0)
112 , fOCDBpath(""), fAccessOCDB(kFALSE)
113 , fDigitsArr(0), fClusterArr(0), fCaloClusterArr(0)
114 , fRecParam(0), fClusterizer(0)
115 , fUnfolder(0), fJustUnfold(kFALSE)
116 , fOutputAODBranch(0), fOutputAODBranchName(""), fOutputAODBranchSet(0)
117 , fFillAODFile(kFALSE), fFillAODHeader(0)
118 , fFillAODCaloCells(0), fRun(-1)
119 , fRecoUtils(0), fConfigName("")
120 , fCellLabels(), fCellSecondLabels(), fCellTime()
121 , fCellMatchdEta(), fCellMatchdPhi()
122 , fMaxEvent(0), fDoTrackMatching(kFALSE)
123 , fSelectCell(kFALSE), fSelectCellMinE(0), fSelectCellMinFrac(0)
124 , fRemoveLEDEvents(kTRUE), fRemoveExoticEvents(kFALSE)
125 , fImportGeometryFromFile(kFALSE), fImportGeometryFilePath("")
126 , fOADBSet(kFALSE), fAccessOADB(kTRUE), fOADBFilePath("")
127 , fCentralityClass("")
131 for(Int_t i = 0; i < 12; i++) fGeomMatrix[i] = 0;
132 for(Int_t j = 0; j < 24*48*11; j++)
135 fCellSecondLabels[j] = -1;
137 fCellMatchdEta[j] = -999;
138 fCellMatchdPhi[j] = -999;
141 fCentralityBin[0] = fCentralityBin[1]=-1;
146 //_______________________________________________________________
147 AliAnalysisTaskEMCALClusterize::~AliAnalysisTaskEMCALClusterize()
153 fDigitsArr->Clear("C");
159 fClusterArr->Delete();
165 fCaloClusterArr->Delete();
166 delete fCaloClusterArr;
169 if(fClusterizer) delete fClusterizer;
170 if(fUnfolder) delete fUnfolder;
171 if(fRecoUtils) delete fRecoUtils;
175 //_______________________________________________
176 void AliAnalysisTaskEMCALClusterize::AccessOADB()
178 // Set the AODB calibration, bad channels etc. parameters at least once
179 // alignment matrices from OADB done in SetGeometryMatrices
182 if(fOADBSet) return ;
184 Int_t runnumber = InputEvent()->GetRunNumber() ;
185 TString pass = GetPass();
187 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Get AODB parameters from EMCAL in %s for run %d, and <%s> \n",fOADBFilePath.Data(),runnumber,pass.Data());
189 Int_t nSM = fGeom->GetNumberOfSuperModules();
192 if(fRecoUtils->IsBadChannelsRemovalSwitchedOn())
194 AliOADBContainer *contBC=new AliOADBContainer("");
195 contBC->InitFromFile(Form("%s/EMCALBadChannels.root",fOADBFilePath.Data()),"AliEMCALBadChannels");
197 TObjArray *arrayBC=(TObjArray*)contBC->GetObject(runnumber);
201 fRecoUtils->SwitchOnDistToBadChannelRecalculation();
202 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Remove EMCAL bad cells \n");
204 for (Int_t i=0; i<nSM; ++i)
206 TH2I *hbm = fRecoUtils->GetEMCALChannelStatusMap(i);
211 hbm=(TH2I*)arrayBC->FindObject(Form("EMCALBadChannelMap_Mod%d",i));
215 AliError(Form("Can not get EMCALBadChannelMap_Mod%d",i));
219 hbm->SetDirectory(0);
220 fRecoUtils->SetEMCALChannelStatusMap(i,hbm);
223 } else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT remove EMCAL bad channels\n"); // run array
226 // Energy Recalibration
227 if(fRecoUtils->IsRecalibrationOn())
229 AliOADBContainer *contRF=new AliOADBContainer("");
231 contRF->InitFromFile(Form("%s/EMCALRecalib.root",fOADBFilePath.Data()),"AliEMCALRecalib");
233 TObjArray *recal=(TObjArray*)contRF->GetObject(runnumber);
237 TObjArray *recalpass=(TObjArray*)recal->FindObject(pass);
241 TObjArray *recalib=(TObjArray*)recalpass->FindObject("Recalib");
245 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Recalibrate EMCAL \n");
246 for (Int_t i=0; i<nSM; ++i)
248 TH2F *h = fRecoUtils->GetEMCALChannelRecalibrationFactors(i);
253 h = (TH2F*)recalib->FindObject(Form("EMCALRecalFactors_SM%d",i));
257 AliError(Form("Could not load EMCALRecalFactors_SM%d",i));
263 fRecoUtils->SetEMCALChannelRecalibrationFactors(i,h);
265 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params object array \n"); // array ok
266 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params for pass\n"); // array pass ok
267 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params for run\n"); // run number array ok
269 } // Recalibration on
271 // Energy Recalibration, apply on top of previous calibration factors
272 if(fRecoUtils->IsRunDepRecalibrationOn())
274 AliOADBContainer *contRFTD=new AliOADBContainer("");
276 contRFTD->InitFromFile(Form("%s/EMCALTemperatureCorrCalib.root",fOADBFilePath.Data()),"AliEMCALRunDepTempCalibCorrections");
278 TH1S *htd=(TH1S*)contRFTD->GetObject(runnumber);
282 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Recalibrate (Temperature) EMCAL \n");
284 for (Int_t ism=0; ism<nSM; ++ism)
286 for (Int_t icol=0; icol<48; ++icol)
288 for (Int_t irow=0; irow<24; ++irow)
290 Float_t factor = fRecoUtils->GetEMCALChannelRecalibrationFactor(ism,icol,irow);
292 Int_t absID = fGeom->GetAbsCellIdFromCellIndexes(ism, irow, icol); // original calibration factor
293 factor *= htd->GetBinContent(absID) / 10000. ; // correction dependent on T
294 //printf("\t ism %d, icol %d, irow %d,absID %d, corrA %2.3f, corrB %2.3f, corrAB %2.3f\n",ism, icol, irow, absID,
295 // GetEMCALChannelRecalibrationFactor(ism,icol,irow) , htd->GetBinContent(absID) / 10000., factor);
296 fRecoUtils->SetEMCALChannelRecalibrationFactor(ism,icol,irow,factor);
300 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL with T variations, no params TH1 \n");
301 } // Run by Run T calibration
303 // Time Recalibration
304 if(fRecoUtils->IsTimeRecalibrationOn())
306 AliOADBContainer *contTRF=new AliOADBContainer("");
308 contTRF->InitFromFile(Form("%s/EMCALTimeCalib.root",fOADBFilePath.Data()),"AliEMCALTimeCalib");
310 TObjArray *trecal=(TObjArray*)contTRF->GetObject(runnumber);
314 TString passTmp = pass;
315 if(pass!="pass1" && pass!="pass2") passTmp = "pass2"; // TEMPORARY FIX FOR LHC11a analysis
316 TObjArray *trecalpass=(TObjArray*)trecal->FindObject(passTmp);
320 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Time Recalibrate EMCAL \n");
321 for (Int_t ibc = 0; ibc < 4; ++ibc)
323 TH1F *h = fRecoUtils->GetEMCALChannelTimeRecalibrationFactors(ibc);
328 h = (TH1F*)trecalpass->FindObject(Form("hAllTimeAvBC%d",ibc));
332 AliError(Form("Could not load hAllTimeAvBC%d",ibc));
338 fRecoUtils->SetEMCALChannelTimeRecalibrationFactors(ibc,h);
339 } // bunch crossing loop
340 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate time EMCAL, no params for pass\n"); // array pass ok
341 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate time EMCAL, no params for run\n"); // run number array ok
343 } // Time recalibration on
345 // Parameters already set once, so do not it again
350 //_________________________________________________
351 Bool_t AliAnalysisTaskEMCALClusterize::AccessOCDB()
353 //Access to OCDB stuff
355 fEvent = InputEvent();
358 Warning("AccessOCDB","Event not available!!!");
362 if (fEvent->GetRunNumber()==fRun)
364 fRun = fEvent->GetRunNumber();
366 if(DebugLevel() > 1 )
367 printf("AliAnalysisTaksEMCALClusterize::AccessOCDB() - Begin");
369 AliCDBManager *cdb = AliCDBManager::Instance();
372 if (fOCDBpath.Length())
374 cdb->SetDefaultStorage(fOCDBpath.Data());
375 printf("AliAnalysisTaksEMCALClusterize::AccessOCDB() - Default storage %s",fOCDBpath.Data());
378 cdb->SetRun(fEvent->GetRunNumber());
381 // EMCAL from RAW OCDB
382 if (fOCDBpath.Contains("alien:"))
384 cdb->SetSpecificStorage("EMCAL/Calib/Data","alien://Folder=/alice/data/2010/OCDB");
385 cdb->SetSpecificStorage("EMCAL/Calib/Pedestals","alien://Folder=/alice/data/2010/OCDB");
388 TString path = cdb->GetDefaultStorage()->GetBaseFolder();
392 //Get calibration parameters
395 AliCDBEntry *entry = (AliCDBEntry*)
396 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
398 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
402 AliFatal("Calibration parameters not found in CDB!");
404 //Get calibration parameters
407 AliCDBEntry *entry = (AliCDBEntry*)
408 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
410 if (entry) fPedestalData = (AliCaloCalibPedestal*) entry->GetObject();
414 AliFatal("Dead map not found in CDB!");
419 //_____________________________________________________
420 void AliAnalysisTaskEMCALClusterize::CheckAndGetEvent()
422 // Get the input event, it can depend in embedded events what you want to get
423 // Also check if the quality of the event is good if not reject it
427 AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
428 Int_t eventN = Entry();
429 if(aodIH) eventN = aodIH->GetReadEntry();
431 if (eventN > fMaxEvent)
434 //printf("Clusterizer --- Event %d-- \n",eventN);
436 //Check if input event are embedded events
437 //If so, take output event
438 if (aodIH && aodIH->GetMergeEvents())
442 if(!aodIH->GetMergeEMCALCells())
443 AliFatal("Events merged but not EMCAL cells, check analysis settings!");
447 printf("AliAnalysisTaksEMCALClusterize::UserExec() - Use embedded events\n");
449 printf("\t InputEvent N Clusters %d, N Cells %d\n",InputEvent()->GetNumberOfCaloClusters(),
450 InputEvent()->GetEMCALCells()->GetNumberOfCells());
452 printf("\t MergedEvent N Clusters %d, N Cells %d\n",aodIH->GetEventToMerge()->GetNumberOfCaloClusters(),
453 aodIH->GetEventToMerge()->GetEMCALCells()->GetNumberOfCells());
455 for (Int_t icl=0; icl < aodIH->GetEventToMerge()->GetNumberOfCaloClusters(); icl++)
457 AliAODCaloCluster *sigCluster = aodIH->GetEventToMerge()->GetCaloCluster(icl);
458 if(sigCluster->IsEMCAL()) printf("\t \t Signal cluster: i %d, E %f\n",icl,sigCluster->E());
461 printf("\t OutputEvent N Clusters %d, N Cells %d\n", AODEvent()->GetNumberOfCaloClusters(),
462 AODEvent()->GetEMCALCells()->GetNumberOfCells());
467 fEvent = InputEvent();
468 if(fFillAODCaloCells) FillAODCaloCells();
469 if(fFillAODHeader) FillAODHeader();
474 Error("UserExec","Event not available");
478 //-------------------------------------------------------------------------------------
479 // Reject events if LED was firing, use only for LHC11a data
480 // Reject event if triggered by exotic cell and remove exotic cells if not triggered
481 //-------------------------------------------------------------------------------------
483 if( IsLEDEvent( InputEvent()->GetRunNumber() ) ) { fEvent = 0x0 ; return ; }
485 if( IsExoticEvent() ) { fEvent = 0x0 ; return ; }
487 //-------------------------------------------------------------------------------------
488 // Set the cluster array in the event (output or input)
489 //-------------------------------------------------------------------------------------
493 //Magic line to write events to AOD filem put after event rejection
494 AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
496 else if( !fOutputAODBranchSet )
498 // Create array and put it in the input event, if output AOD not selected, only once
499 InputEvent()->AddObject(fOutputAODBranch);
500 fOutputAODBranchSet = kTRUE;
501 printf("AliAnalysisTaskEMCALClusterize::UserExec() - Add AOD branch <%s> to input event\n",fOutputAODBranchName.Data());
506 //____________________________________________________
507 void AliAnalysisTaskEMCALClusterize::ClusterizeCells()
509 // Recluster calocells, transform them into digits,
510 // feed the clusterizer with them and get new list of clusters
512 //In case of MC, first loop on the clusters and fill MC label to array
513 Int_t nClusters = fEvent->GetNumberOfCaloClusters();
514 Int_t nClustersOrg = 0;
516 AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
517 if(aodIH && aodIH->GetEventToMerge()) //Embedding
518 nClusters = aodIH->GetEventToMerge()->GetNumberOfCaloClusters(); //Get clusters directly from embedded signal
520 for (Int_t i = 0; i < nClusters; i++)
522 AliVCluster *clus = 0;
523 if(aodIH && aodIH->GetEventToMerge()) //Embedding
524 clus = aodIH->GetEventToMerge()->GetCaloCluster(i); //Get clusters directly from embedded signal
526 clus = fEvent->GetCaloCluster(i);
532 Int_t label = clus->GetLabel();
534 //printf("Org cluster E %f, Time %e, Id = ", clus->E(), clus->GetTOF() );
535 if (clus->GetNLabels()>=2) label2 = clus->GetLabelAt(1) ;
536 UShort_t * index = clus->GetCellsAbsId() ;
537 for(Int_t icell=0; icell < clus->GetNCells(); icell++ )
539 fCellLabels[index[icell]] = label;
540 fCellSecondLabels[index[icell]] = label2;
541 fCellTime[index[icell]] = clus->GetTOF();
542 fCellMatchdEta[index[icell]] = clus->GetTrackDz();
543 fCellMatchdPhi[index[icell]] = clus->GetTrackDx();
544 //printf(" %d,", index[icell] );
551 // Transform CaloCells into Digits
558 AliVCaloCells *cells = fEvent->GetEMCALCells();
560 Int_t bc = InputEvent()->GetBunchCrossNumber();
561 for (Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
563 // Get cell values, recalibrate and not include bad channels found in analysis, nor cells with too low energy, nor exotic cell
564 id = cells->GetCellNumber(icell);
565 Bool_t accept = fRecoUtils->AcceptCalibrateCell(id,bc,amp,time,cells);
567 // Do not include cells with too low energy, nor exotic cell
568 if(amp < fRecParam->GetMinECut() ) accept = kFALSE;
570 // In case of AOD analysis cell time is 0, approximate replacing by time of the cluster the digit belongs.
573 time = fCellTime[id];
574 //printf("cell %d time org %f - ",id, time*1.e9);
575 fRecoUtils->RecalibrateCellTime(id,bc,time);
576 //printf("recal %f\n",time*1.e9);
579 if( accept && fRecoUtils->IsExoticCell(id,cells,bc))
586 fCellLabels[id] =-1; //reset the entry in the array for next event
587 fCellSecondLabels[id]=-1; //reset the entry in the array for next event
589 fCellMatchdEta[id] =-999;
590 fCellMatchdPhi[id] =-999;
591 if( DebugLevel() > 2 )
592 printf("AliAnalysisTaksEMCALClusterize::ClusterizeCells() - Remove channel absId %d, index %d of %d, amp %f, time %f\n",
593 id,icell, cells->GetNumberOfCells(), amp, time*1.e9);
597 //Create the digit, put a fake primary deposited energy to trick the clusterizer when checking the most likely primary
598 new((*fDigitsArr)[idigit]) AliEMCALDigit( fCellLabels[id], fCellLabels[id],id, amp, time,AliEMCALDigit::kHG,idigit, 0, 0, 1);
600 fCellLabels[id] =-1; //reset the entry in the array for next event
607 //-------------------------------------------------------------------------------------
608 //Do the clusterization
609 //-------------------------------------------------------------------------------------
611 fClusterizer->Digits2Clusters("");
613 //-------------------------------------------------------------------------------------
614 //Transform the recpoints into AliVClusters
615 //-------------------------------------------------------------------------------------
617 RecPoints2Clusters();
621 printf("AliAnalysisTaksEMCALClusterize::UserExec() - No array with CaloClusters, input RecPoints entries %d\n",fClusterArr->GetEntriesFast());
625 if( DebugLevel() > 0 )
627 printf("AliAnalysisTaksEMCALClusterize::ClusterizeCells() - N clusters: before recluster %d, after recluster %d\n",nClustersOrg, fCaloClusterArr->GetEntriesFast());
629 if(fCaloClusterArr->GetEntriesFast() != fClusterArr->GetEntriesFast())
631 printf("\t Some RecRoints not transformed into CaloClusters (clusterizer %d, unfold %d): Input entries %d - Output entries %d - %d (not fast)\n",
632 fRecParam->GetClusterizerFlag(),fRecParam->GetUnfold(),
633 fClusterArr->GetEntriesFast(), fCaloClusterArr->GetEntriesFast(), fCaloClusterArr->GetEntries());
637 //Reset the array with second labels for this event
638 memset(fCellSecondLabels, -1, sizeof(fCellSecondLabels));
642 //_____________________________________________________
643 void AliAnalysisTaskEMCALClusterize::ClusterUnfolding()
645 // Take the event clusters and unfold them
647 AliVCaloCells *cells = fEvent->GetEMCALCells();
648 Double_t cellAmplitude = 0;
649 Double_t cellTime = 0;
650 Short_t cellNumber = 0;
651 Short_t cellMCLabel = 0;
652 Double_t cellEFrac = 0;
653 Int_t nClustersOrg = 0;
655 // Fill the array with the EMCAL clusters, copy them
656 for (Int_t i = 0; i < fEvent->GetNumberOfCaloClusters(); i++)
658 AliVCluster *clus = fEvent->GetCaloCluster(i);
661 //recalibrate/remove bad channels/etc if requested
662 if(fRecoUtils->ClusterContainsBadChannel(fGeom,clus->GetCellsAbsId(), clus->GetNCells()))
667 if(fRecoUtils->IsRecalibrationOn())
670 fRecoUtils->RecalibrateClusterEnergy(fGeom, clus, cells);
673 for (Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
675 if (cells->GetCell(icell, cellNumber, cellAmplitude, cellTime, cellMCLabel, cellEFrac) != kTRUE)
678 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
679 fGeom->GetCellIndex(cellNumber,imod,iTower,iIphi,iIeta);
680 fGeom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
682 //Do not include bad channels found in analysis?
683 if( fRecoUtils->IsBadChannelsRemovalSwitchedOn() &&
684 fRecoUtils->GetEMCALChannelStatus(imod, ieta, iphi))
686 fCellLabels[cellNumber] =-1; //reset the entry in the array for next event
687 fCellSecondLabels[cellNumber]=-1; //reset the entry in the array for next event
688 fCellTime[cellNumber] = 0.;
689 fCellMatchdEta[cellNumber] =-999;
690 fCellMatchdPhi[cellNumber] =-999;
694 cells->SetCell(icell, cellNumber, cellAmplitude*fRecoUtils->GetEMCALChannelRecalibrationFactor(imod,ieta,iphi),cellTime);
699 //Cast to ESD or AOD, needed to create the cluster array
700 AliESDCaloCluster * esdCluster = dynamic_cast<AliESDCaloCluster*> (clus);
701 AliAODCaloCluster * aodCluster = dynamic_cast<AliAODCaloCluster*> (clus);
705 fCaloClusterArr->Add( new AliESDCaloCluster(*esdCluster) );
709 fCaloClusterArr->Add( new AliAODCaloCluster(*aodCluster) );
712 Warning("UserExec()"," - Wrong CaloCluster type?");
719 fUnfolder->UnfoldClusters(fCaloClusterArr, cells);
726 //_____________________________________________________
727 void AliAnalysisTaskEMCALClusterize::FillAODCaloCells()
729 // Put calo cells in standard branch
730 AliVCaloCells &eventEMcells = *(fEvent->GetEMCALCells());
731 Int_t nEMcell = eventEMcells.GetNumberOfCells() ;
733 AliAODCaloCells &aodEMcells = *(AODEvent()->GetEMCALCells());
734 aodEMcells.CreateContainer(nEMcell);
735 aodEMcells.SetType(AliVCaloCells::kEMCALCell);
736 Double_t calibFactor = 1.;
737 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
739 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
740 fGeom->GetCellIndex(eventEMcells.GetCellNumber(iCell),imod,iTower,iIphi,iIeta);
741 fGeom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
743 if(fRecoUtils->IsRecalibrationOn())
745 calibFactor = fRecoUtils->GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
748 if(!fRecoUtils->GetEMCALChannelStatus(imod, ieta, iphi))
749 { //Channel is not declared as bad
750 aodEMcells.SetCell(iCell,eventEMcells.GetCellNumber(iCell),eventEMcells.GetAmplitude(iCell)*calibFactor,
751 eventEMcells.GetTime(iCell),eventEMcells.GetMCLabel(iCell),eventEMcells.GetEFraction(iCell));
755 aodEMcells.SetCell(iCell,eventEMcells.GetCellNumber(iCell),0,-1,-1,0);
762 //__________________________________________________
763 void AliAnalysisTaskEMCALClusterize::FillAODHeader()
765 //Put event header information in standard AOD branch
767 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (fEvent);
768 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (fEvent);
772 for (Int_t i = 0; i < 6; i++) covVtx[i] = 0.;
774 AliAODHeader* header = AODEvent()->GetHeader();
775 header->SetRunNumber(fEvent->GetRunNumber());
779 TTree* tree = fInputHandler->GetTree();
782 TFile* file = tree->GetCurrentFile();
783 if (file) header->SetESDFileName(file->GetName());
786 else if (aodevent) header->SetESDFileName(aodevent->GetHeader()->GetESDFileName());
788 header->SetBunchCrossNumber(fEvent->GetBunchCrossNumber());
789 header->SetOrbitNumber(fEvent->GetOrbitNumber());
790 header->SetPeriodNumber(fEvent->GetPeriodNumber());
791 header->SetEventType(fEvent->GetEventType());
794 if(fEvent->GetCentrality())
795 header->SetCentrality(new AliCentrality(*(fEvent->GetCentrality())));
797 header->SetCentrality(0);
800 header->SetOfflineTrigger(fInputHandler->IsEventSelected()); // propagate the decision of the physics selection
801 if (esdevent) header->SetFiredTriggerClasses(esdevent->GetFiredTriggerClasses());
802 else if (aodevent) header->SetFiredTriggerClasses(aodevent->GetFiredTriggerClasses());
804 header->SetTriggerMask(fEvent->GetTriggerMask());
805 header->SetTriggerCluster(fEvent->GetTriggerCluster());
809 header->SetL0TriggerInputs(esdevent->GetHeader()->GetL0TriggerInputs());
810 header->SetL1TriggerInputs(esdevent->GetHeader()->GetL1TriggerInputs());
811 header->SetL2TriggerInputs(esdevent->GetHeader()->GetL2TriggerInputs());
815 header->SetL0TriggerInputs(aodevent->GetHeader()->GetL0TriggerInputs());
816 header->SetL1TriggerInputs(aodevent->GetHeader()->GetL1TriggerInputs());
817 header->SetL2TriggerInputs(aodevent->GetHeader()->GetL2TriggerInputs());
820 header->SetMagneticField(fEvent->GetMagneticField());
821 //header->SetMuonMagFieldScale(esdevent->GetCurrentDip()/6000.);
823 header->SetZDCN1Energy(fEvent->GetZDCN1Energy());
824 header->SetZDCP1Energy(fEvent->GetZDCP1Energy());
825 header->SetZDCN2Energy(fEvent->GetZDCN2Energy());
826 header->SetZDCP2Energy(fEvent->GetZDCP2Energy());
827 header->SetZDCEMEnergy(fEvent->GetZDCEMEnergy(0),fEvent->GetZDCEMEnergy(1));
829 Float_t diamxy[2]={fEvent->GetDiamondX(),fEvent->GetDiamondY()};
831 fEvent->GetDiamondCovXY(diamcov);
832 header->SetDiamond(diamxy,diamcov);
833 if (esdevent) header->SetDiamondZ(esdevent->GetDiamondZ(),esdevent->GetSigma2DiamondZ());
834 else if (aodevent) header->SetDiamondZ(aodevent->GetDiamondZ(),aodevent->GetSigma2DiamondZ());
837 Int_t nVertices = 1 ;/* = prim. vtx*/;
838 Int_t nCaloClus = fEvent->GetNumberOfCaloClusters();
840 AODEvent()->ResetStd(0, nVertices, 0, 0, 0, nCaloClus, 0, 0);
842 // Access to the AOD container of vertices
843 TClonesArray &vertices = *(AODEvent()->GetVertices());
846 // Add primary vertex. The primary tracks will be defined
847 // after the loops on the composite objects (V0, cascades, kinks)
848 fEvent->GetPrimaryVertex()->GetXYZ(pos);
852 esdevent->GetPrimaryVertex()->GetCovMatrix(covVtx);
853 chi = esdevent->GetPrimaryVertex()->GetChi2toNDF();
857 aodevent->GetPrimaryVertex()->GetCovMatrix(covVtx);
858 chi = aodevent->GetPrimaryVertex()->GetChi2perNDF();//Different from ESD?
861 AliAODVertex * primary = new(vertices[jVertices++])
862 AliAODVertex(pos, covVtx, chi, NULL, -1, AliAODVertex::kPrimary);
863 primary->SetName(fEvent->GetPrimaryVertex()->GetName());
864 primary->SetTitle(fEvent->GetPrimaryVertex()->GetTitle());
868 //___________________________________________________________
869 void AliAnalysisTaskEMCALClusterize::FillCaloClusterInEvent()
871 // Get the CaloClusters array, do some final calculations
872 // and put the clusters in the output or input event
873 // as a separate branch
875 //First recalculate track-matching for the new clusters
878 fRecoUtils->FindMatches(fEvent,fCaloClusterArr,fGeom);
880 //Put the new clusters in the AOD list
882 Int_t kNumberOfCaloClusters = fCaloClusterArr->GetEntriesFast();
884 for(Int_t i = 0; i < kNumberOfCaloClusters; i++)
886 AliAODCaloCluster *newCluster = (AliAODCaloCluster *) fCaloClusterArr->At(i);
888 newCluster->SetID(i);
890 // Correct cluster energy non linearity
891 newCluster->SetE(fRecoUtils->CorrectClusterEnergyLinearity(newCluster));
896 Int_t trackIndex = fRecoUtils->GetMatchedTrackIndex(i);
899 newCluster->AddTrackMatched(fEvent->GetTrack(trackIndex));
901 printf("AliAnalysisTaksEMCALClusterize::UserExec() - Matched Track index %d to new cluster %d \n",trackIndex,i);
904 Float_t dR = 999., dZ = 999.;
905 fRecoUtils->GetMatchedResiduals(newCluster->GetID(),dR,dZ);
906 newCluster->SetTrackDistance(dR,dZ);
910 {// Assign previously assigned matched track in reco, very very rough
911 Int_t absId0 = newCluster->GetCellsAbsId()[0]; // Assign match of first cell in cluster
912 newCluster->SetTrackDistance(fCellMatchdPhi[absId0],fCellMatchdEta[absId0]);
915 //printf("New cluster E %f, Time %e, Id = ", newCluster->E(), newCluster->GetTOF() );
916 //for(Int_t icell=0; icell < newCluster->GetNCells(); icell++ ) printf(" %d,", newCluster->GetCellsAbsId() [icell] );
919 // Calculate distance to bad channel for new cluster. Make sure you give the list of bad channels.
920 fRecoUtils->RecalculateClusterDistanceToBadChannel(fGeom, fEvent->GetEMCALCells(), newCluster);
922 new((*fOutputAODBranch)[i]) AliAODCaloCluster(*newCluster);
924 if(DebugLevel() > 1 )
925 printf("AliAnalysisTaksEMCALClusterize::UserExec() - New cluster %d of %d, energy %f\n",newCluster->GetID(), kNumberOfCaloClusters, newCluster->E());
929 fOutputAODBranch->Expand(kNumberOfCaloClusters); // resize TObjArray to 'remove' slots
934 //_______________________________________________
935 TString AliAnalysisTaskEMCALClusterize::GetPass()
937 // Get passx from filename.
939 if (!AliAnalysisManager::GetAnalysisManager()->GetTree())
941 AliError("AliAnalysisTaskEMCALClusterize::GetPass() - Pointer to tree = 0, returning null\n");
945 if (!AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile())
947 AliError("AliAnalysisTaskEMCALClusterize::GetPass() - Null pointer input file, returning null\n");
951 TString pass(AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile()->GetName());
952 if (pass.Contains("ass1")) return TString("pass1");
953 else if (pass.Contains("ass2")) return TString("pass2");
954 else if (pass.Contains("ass3")) return TString("pass3");
955 else if (pass.Contains("ass4")) return TString("pass4");
956 else if (pass.Contains("ass5")) return TString("pass5");
958 // No condition fullfilled, give a default value
959 printf("AliAnalysisTaskEMCALClusterize::GetPass() - Pass number string not found \n");
964 //_________________________________________
965 void AliAnalysisTaskEMCALClusterize::Init()
967 //Init analysis with configuration macro if available
970 if(fOADBFilePath == "") fOADBFilePath = "$ALICE_ROOT/OADB/EMCAL" ;
972 fBranchNames = "ESD:AliESDHeader.,EMCALCells.";
974 if(!fRecParam) fRecParam = new AliEMCALRecParam;
975 if(!fRecoUtils) fRecoUtils = new AliEMCALRecoUtils();
977 if(fMaxEvent <= 0) fMaxEvent = 1000000000;
978 if(fSelectCellMinE <= 0) fSelectCellMinE = 0.005;
979 if(fSelectCellMinFrac <= 0) fSelectCellMinFrac = 0.001;
982 if(fCentralityClass == "") fCentralityClass = "V0M";
984 if (fOCDBpath == "") fOCDBpath = "raw://" ;
985 if (fOutputAODBranchName == "") fOutputAODBranchName = "newEMCALClusters" ;
987 if(gROOT->LoadMacro(fConfigName) >=0)
989 printf("AliAnalysisTaksEMCALClusterize::Init() - Configure analysis with %s\n",fConfigName.Data());
990 AliAnalysisTaskEMCALClusterize *clus = (AliAnalysisTaskEMCALClusterize*)gInterpreter->ProcessLine("ConfigEMCALClusterize()");
991 fGeomName = clus->fGeomName;
992 fLoadGeomMatrices = clus->fLoadGeomMatrices;
993 fOCDBpath = clus->fOCDBpath;
994 fAccessOCDB = clus->fAccessOCDB;
995 fRecParam = clus->fRecParam;
996 fJustUnfold = clus->fJustUnfold;
997 fFillAODFile = clus->fFillAODFile;
998 fRecoUtils = clus->fRecoUtils;
999 fConfigName = clus->fConfigName;
1000 fMaxEvent = clus->fMaxEvent;
1001 fDoTrackMatching = clus->fDoTrackMatching;
1002 fOutputAODBranchName = clus->fOutputAODBranchName;
1003 for(Int_t i = 0; i < 12; i++) fGeomMatrix[i] = clus->fGeomMatrix[i] ;
1004 fCentralityClass = clus->fCentralityClass;
1005 fCentralityBin[0] = clus->fCentralityBin[0];
1006 fCentralityBin[1] = clus->fCentralityBin[1];
1009 // Init geometry, I do not like much to do it like this ...
1010 if(fImportGeometryFromFile && !gGeoManager)
1012 if (fImportGeometryFilePath == "") fImportGeometryFilePath = "$ALICE_ROOT/OADB/EMCAL/geometry_2011.root" ; // "$ALICE_ROOT/EVE/alice-data/default_geo.root"
1013 printf("AliAnalysisTaskEMCALClusterize::Init() - Import %s\n",fImportGeometryFilePath.Data());
1014 TGeoManager::Import(fImportGeometryFilePath) ;
1019 //_______________________________________________________
1020 void AliAnalysisTaskEMCALClusterize::InitClusterization()
1022 //Select clusterization/unfolding algorithm and set all the needed parameters
1026 // init the unfolding afterburner
1028 fUnfolder = new AliEMCALAfterBurnerUF(fRecParam->GetW0(),fRecParam->GetLocMaxCut(),fRecParam->GetMinECut());
1032 //First init the clusterizer
1033 delete fClusterizer;
1034 if (fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv1)
1035 fClusterizer = new AliEMCALClusterizerv1 (fGeom, fCalibData, fPedestalData);
1036 else if(fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv2)
1037 fClusterizer = new AliEMCALClusterizerv2(fGeom, fCalibData, fPedestalData);
1038 else if(fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerNxN)
1040 fClusterizer = new AliEMCALClusterizerNxN(fGeom, fCalibData, fPedestalData);
1041 fClusterizer->SetNRowDiff(fRecParam->GetNRowDiff());
1042 fClusterizer->SetNColDiff(fRecParam->GetNColDiff());
1046 AliFatal(Form("Clusterizer < %d > not available", fRecParam->GetClusterizerFlag()));
1049 //Now set the parameters
1050 fClusterizer->SetECAClusteringThreshold( fRecParam->GetClusteringThreshold() );
1051 fClusterizer->SetECALogWeight ( fRecParam->GetW0() );
1052 fClusterizer->SetMinECut ( fRecParam->GetMinECut() );
1053 fClusterizer->SetUnfolding ( fRecParam->GetUnfold() );
1054 fClusterizer->SetECALocalMaxCut ( fRecParam->GetLocMaxCut() );
1055 fClusterizer->SetTimeCut ( fRecParam->GetTimeCut() );
1056 fClusterizer->SetTimeMin ( fRecParam->GetTimeMin() );
1057 fClusterizer->SetTimeMax ( fRecParam->GetTimeMax() );
1058 fClusterizer->SetInputCalibrated ( kTRUE );
1059 fClusterizer->SetJustClusters ( kTRUE );
1061 // Initialize the cluster rec points and digits arrays and get them.
1062 fClusterizer->SetOutput(0);
1063 fClusterArr = const_cast<TObjArray *>(fClusterizer->GetRecPoints());
1064 fDigitsArr = fClusterizer->GetDigits();
1066 //In case of unfolding after clusterization is requested, set the corresponding parameters
1067 if(fRecParam->GetUnfold())
1070 for (i = 0; i < 8; i++)
1072 fClusterizer->SetSSPars(i, fRecParam->GetSSPars(i));
1073 }//end of loop over parameters
1075 for (i = 0; i < 3; i++)
1077 fClusterizer->SetPar5 (i, fRecParam->GetPar5(i));
1078 fClusterizer->SetPar6 (i, fRecParam->GetPar6(i));
1079 }//end of loop over parameters
1081 fClusterizer->InitClusterUnfolding();
1086 //_________________________________________________
1087 void AliAnalysisTaskEMCALClusterize::InitGeometry()
1089 // Init geometry and set the geometry matrix, for the first event, skip the rest
1090 // Also set once the run dependent calibrations
1092 if(fGeomMatrixSet) return;
1094 Int_t runnumber = InputEvent()->GetRunNumber() ;
1099 if (runnumber < 140000) fGeomName = "EMCAL_FIRSTYEARV1";
1100 else if(runnumber < 171000) fGeomName = "EMCAL_COMPLETEV1";
1101 else fGeomName = "EMCAL_COMPLETE12SMV1";
1102 printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - Set EMCAL geometry name to <%s> for run %d\n",fGeomName.Data(),runnumber);
1105 fGeom = AliEMCALGeometry::GetInstance(fGeomName);
1109 printf("AliAnalysisTaskEMCALClusterize::InitGeometry(run=%d)",runnumber);
1110 if (!gGeoManager) printf(" - Careful!, gGeoManager not loaded, load misalign matrices");
1113 } // geometry pointer did not exist before
1115 if(fLoadGeomMatrices)
1117 printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - Load user defined EMCAL geometry matrices\n");
1119 // OADB if available
1120 AliOADBContainer emcGeoMat("AliEMCALgeo");
1121 emcGeoMat.InitFromFile(Form("%s/EMCALlocal2master.root",fOADBFilePath.Data()),"AliEMCALgeo");
1122 TObjArray *matEMCAL=(TObjArray*)emcGeoMat.GetObject(runnumber,"EmcalMatrices");
1125 for(Int_t mod=0; mod < (fGeom->GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
1128 if (!fGeomMatrix[mod]) // Get it from OADB
1131 printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - EMCAL matrices SM %d, %p\n",
1132 mod,((TGeoHMatrix*) matEMCAL->At(mod)));
1133 //((TGeoHMatrix*) matEMCAL->At(mod))->Print();
1135 fGeomMatrix[mod] = (TGeoHMatrix*) matEMCAL->At(mod) ;
1138 if(fGeomMatrix[mod])
1140 if(DebugLevel() > 1)
1141 fGeomMatrix[mod]->Print();
1143 fGeom->SetMisalMatrix(fGeomMatrix[mod],mod) ;
1146 fGeomMatrixSet=kTRUE;
1150 else if(!gGeoManager)
1152 printf("AliAnalysisTaksEMCALClusterize::InitGeometry() - Get geo matrices from data");
1153 //Still not implemented in AOD, just a workaround to be able to work at least with ESDs
1154 if(!strcmp(fEvent->GetName(),"AliAODEvent"))
1156 if(DebugLevel() > 1)
1157 Warning("UserExec","Use ideal geometry, values geometry matrix not kept in AODs.");
1161 if(DebugLevel() > 1)
1162 printf("AliAnalysisTaksEMCALClusterize::InitGeometry() - AliAnalysisTaskEMCALClusterize Load Misaligned matrices.");
1164 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(fEvent) ;
1168 Error("InitGeometry"," - This event does not contain ESDs?");
1169 if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);
1173 for(Int_t mod=0; mod < (fGeom->GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
1175 if(DebugLevel() > 1)
1176 esd->GetEMCALMatrix(mod)->Print();
1178 if(esd->GetEMCALMatrix(mod)) fGeom->SetMisalMatrix(esd->GetEMCALMatrix(mod),mod) ;
1182 fGeomMatrixSet=kTRUE;
1185 }//Load matrices from Data
1189 //____________________________________________________
1190 Bool_t AliAnalysisTaskEMCALClusterize::IsExoticEvent()
1193 // Check if event is exotic, get an exotic cell and compare with triggered patch
1194 // If there is a match remove event ... to be completed, filled with something provisional
1196 if(!fRemoveExoticEvents) return kFALSE;
1199 AliVCaloCells * cells = fEvent->GetEMCALCells();
1200 Float_t totCellE = 0;
1201 Int_t bc = InputEvent()->GetBunchCrossNumber();
1202 for(Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
1206 Double_t tcell = 0 ;
1208 Int_t absID = cells->GetCellNumber(icell);
1209 Bool_t accept = fRecoUtils->AcceptCalibrateCell(absID,bc,ecell,tcell,cells);
1210 if(accept && !fRecoUtils->IsExoticCell(absID,cells,bc)) totCellE += ecell;
1213 // TString triggerclasses = "";
1214 // if(esdevent) triggerclasses = esdevent ->GetFiredTriggerClasses();
1215 // else triggerclasses = ((AliAODEvent*)event)->GetFiredTriggerClasses();
1217 // printf("AliAnalysisTaskEMCALClusterize - reject event %d with cluster - reject event with ncells in SM3 %d and SM4 %d\n",(Int_t)Entry(),ncellsSM3, ncellsSM4);
1218 // if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);;
1222 //printf("TotE cell %f\n",totCellE);
1223 if(totCellE < 1) return kTRUE;
1229 //________________________________________________________________
1230 Bool_t AliAnalysisTaskEMCALClusterize::IsLEDEvent(const Int_t run)
1232 //Check if event is LED
1234 if(!fRemoveLEDEvents) return kFALSE;
1236 //check events of LHC11a period
1237 if(run < 146858 || run > 146860) return kFALSE ;
1239 // Count number of cells with energy larger than 0.1 in SM3, cut on this number
1240 Int_t ncellsSM3 = 0;
1241 AliVCaloCells * cells = fEvent->GetEMCALCells();
1242 for(Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
1244 if(cells->GetAmplitude(icell) > 0.1 && cells->GetCellNumber(icell)/(24*48)==3) ncellsSM3++;
1247 TString triggerclasses = "";
1249 AliESDEvent *esdevent = dynamic_cast<AliESDEvent*>(fEvent);
1250 if(esdevent) triggerclasses = esdevent ->GetFiredTriggerClasses();
1251 else triggerclasses = ((AliAODEvent*)fEvent)->GetFiredTriggerClasses();
1253 Int_t ncellcut = 21;
1254 if(triggerclasses.Contains("EMC")) ncellcut = 35;
1256 if( ncellsSM3 >= ncellcut)
1258 printf("AliAnalysisTaksEMCALClusterize::IsLEDEvent() - reject event %d with ncells in SM3 %d\n",(Int_t)Entry(),ncellsSM3);
1259 if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);;
1267 //_______________________________________________________
1268 void AliAnalysisTaskEMCALClusterize::RecPoints2Clusters()
1270 // Restore clusters from recPoints
1271 // Cluster energy, global position, cells and their amplitude
1272 // fractions are restored
1275 for(Int_t i = 0; i < fClusterArr->GetEntriesFast(); i++)
1277 AliEMCALRecPoint *recPoint = (AliEMCALRecPoint*) fClusterArr->At(i);
1279 const Int_t ncells = recPoint->GetMultiplicity();
1280 Int_t ncellsTrue = 0;
1282 if(recPoint->GetEnergy() < fRecParam->GetClusteringThreshold()) continue;
1284 // cells and their amplitude fractions
1285 UShort_t absIds[ncells];
1286 Double32_t ratios[ncells];
1288 //For later check embedding
1289 AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
1291 Float_t clusterE = 0;
1292 for (Int_t c = 0; c < ncells; c++)
1294 AliEMCALDigit *digit = (AliEMCALDigit*) fDigitsArr->At(recPoint->GetDigitsList()[c]);
1296 absIds[ncellsTrue] = digit->GetId();
1297 ratios[ncellsTrue] = recPoint->GetEnergiesList()[c]/digit->GetAmplitude();
1299 // In case of unfolding, remove digits with unfolded energy too low
1302 if (recPoint->GetEnergiesList()[c] < fSelectCellMinE || ratios[ncellsTrue] < fSelectCellMinFrac)
1304 if(DebugLevel() > 1)
1306 printf("AliAnalysisTaksEMCALClusterize::RecPoints2Clusters() - Too small energy in cell of cluster: cluster cell %f, digit %f\n",
1307 recPoint->GetEnergiesList()[c],digit->GetAmplitude());
1315 //Recalculate cluster energy and number of cluster cells in case any of the cells was rejected
1317 clusterE +=recPoint->GetEnergiesList()[c];
1319 // In case of embedding, fill ratio with amount of signal,
1320 if (aodIH && aodIH->GetMergeEvents())
1322 //AliVCaloCells* inEMCALCells = InputEvent()->GetEMCALCells();
1323 AliVCaloCells* meEMCALCells = aodIH->GetEventToMerge()->GetEMCALCells();
1324 AliVCaloCells* ouEMCALCells = AODEvent()->GetEMCALCells();
1326 Float_t sigAmplitude = meEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1327 //Float_t bkgAmplitude = inEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1328 Float_t sumAmplitude = ouEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1329 //printf("\t AbsID %d, amplitude : bkg %f, sigAmplitude %f, summed %f - %f\n",absIds[ncellsTrue], bkgAmplitude, sigAmplitude, sumAmplitude, digit->GetAmplitude());
1331 if(sumAmplitude > 0) ratios[ncellsTrue] = sigAmplitude/sumAmplitude;
1332 //printf("\t \t ratio %f\n",ratios[ncellsTrue]);
1336 }// cluster cell loop
1340 if (DebugLevel() > 1)
1341 printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Skipping cluster with no cells avobe threshold E = %f, ncells %d\n",
1342 recPoint->GetEnergy(), ncells);
1346 //if(ncellsTrue != ncells) printf("Old E %f, ncells %d; New E %f, ncells %d\n",recPoint->GetEnergy(),ncells,clusterE,ncellsTrue);
1348 if(clusterE < fRecParam->GetClusteringThreshold())
1351 printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Remove cluster with energy below seed threshold %f\n",clusterE);
1358 // calculate new cluster position
1360 recPoint->EvalGlobalPosition(fRecParam->GetW0(), fDigitsArr);
1361 recPoint->GetGlobalPosition(gpos);
1364 // create a new cluster
1366 (*fCaloClusterArr)[j] = new AliAODCaloCluster() ;
1367 AliAODCaloCluster *clus = dynamic_cast<AliAODCaloCluster *>( fCaloClusterArr->At(j) ) ;
1369 clus->SetType(AliVCluster::kEMCALClusterv1);
1370 clus->SetE(clusterE);
1371 clus->SetPosition(g);
1372 clus->SetNCells(ncellsTrue);
1373 clus->SetCellsAbsId(absIds);
1374 clus->SetCellsAmplitudeFraction(ratios);
1375 clus->SetChi2(-1); //not yet implemented
1376 clus->SetTOF(recPoint->GetTime()) ; //time-of-flight
1377 clus->SetNExMax(recPoint->GetNExMax()); //number of local maxima
1378 clus->SetDistanceToBadChannel(recPoint->GetDistanceToBadTower());
1380 if(ncells == ncellsTrue)
1382 Float_t elipAxis[2];
1383 recPoint->GetElipsAxis(elipAxis);
1384 clus->SetM02(elipAxis[0]*elipAxis[0]) ;
1385 clus->SetM20(elipAxis[1]*elipAxis[1]) ;
1386 clus->SetDispersion(recPoint->GetDispersion());
1388 else if(fSelectCell)
1390 // In case some cells rejected, in unfolding case, recalculate
1391 // shower shape parameters and position
1392 if(DebugLevel() > 1)
1393 printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Cells removed from cluster (ncells %d, ncellsTrue %d), recalculate Shower Shape\n",ncells,ncellsTrue);
1395 AliVCaloCells* cells = 0x0;
1396 if (aodIH && aodIH->GetMergeEvents()) cells = AODEvent() ->GetEMCALCells();
1397 else cells = InputEvent()->GetEMCALCells();
1399 fRecoUtils->RecalculateClusterShowerShapeParameters(fGeom,cells,clus);
1400 fRecoUtils->RecalculateClusterPID(clus);
1401 fRecoUtils->RecalculateClusterPosition(fGeom,cells,clus);
1406 Int_t parentMult = 0;
1407 Int_t *parentList = recPoint->GetParents(parentMult);
1408 clus->SetLabel(parentList, parentMult);
1410 //Write the second major contributor to each MC cluster.
1412 for ( Int_t iLoopCell = 0 ; iLoopCell < clus->GetNCells() ; iLoopCell++ )
1415 Int_t idCell = clus->GetCellAbsId(iLoopCell) ;
1418 iNewLabel = 1 ; //iNewLabel makes sure we don't write twice the same label.
1419 for ( UInt_t iLoopLabels = 0 ; iLoopLabels < clus->GetNLabels() ; iLoopLabels++ )
1421 if ( fCellSecondLabels[idCell] == -1 ) iNewLabel = 0; // -1 is never a good second label.
1422 if ( fCellSecondLabels[idCell] == clus->GetLabelAt(iLoopLabels) ) iNewLabel = 0;
1426 Int_t * newLabelArray = new Int_t[clus->GetNLabels()+1] ;
1427 for ( UInt_t iLoopNewLabels = 0 ; iLoopNewLabels < clus->GetNLabels() ; iLoopNewLabels++ )
1429 newLabelArray[iLoopNewLabels] = clus->GetLabelAt(iLoopNewLabels) ;
1432 newLabelArray[clus->GetNLabels()] = fCellSecondLabels[idCell] ;
1433 clus->SetLabel(newLabelArray,clus->GetNLabels()+1) ;
1434 delete [] newLabelArray;
1436 }//positive cell number
1444 //____________________________________________________________
1445 void AliAnalysisTaskEMCALClusterize::UserCreateOutputObjects()
1447 // Init geometry, create list of output clusters
1450 fOutputAODBranch = new TClonesArray("AliAODCaloCluster", 0);
1452 if(fOutputAODBranchName.Length()==0)
1454 fOutputAODBranchName = "newEMCALClustersArray";
1455 printf("Cluster branch name not set, set it to newEMCALClustersArray \n");
1458 fOutputAODBranch->SetName(fOutputAODBranchName);
1462 //fOutputAODBranch = new TClonesArray("AliAODCaloCluster", 0);
1465 //fOutputAODBranch->SetOwner(kFALSE);
1467 AddAODBranch("TClonesArray", &fOutputAODBranch);
1472 //_______________________________________________________
1473 void AliAnalysisTaskEMCALClusterize::UserExec(Option_t *)
1475 // Do clusterization event by event, execute different steps
1476 // 1) Do some checks on the kind of events (ESD, AOD) or if some filtering is needed, initializations
1477 // load things and clear arrays
1478 // 2) Clusterize a) just unfolding existing clusters (fJustUnfold)
1479 // b) recluster cells
1480 // + convert cells into digits (calibrating them)
1481 // + recluster digits into recPoints with the chosen clusterizer (V1, V1+Unfold,V2, NxN)
1482 // with methods in AliEMCALClusterizer
1483 // + transform recPoints into CaloClusters
1484 // 3) Do some final calculations in the found clusters (track-matching) and put them in an AOD branch
1490 //Remove the contents of AOD branch output list set in the previous event
1491 fOutputAODBranch->Clear("C");
1495 //Check if there is a centrality value, PbPb analysis, and if a centrality bin selection is requested
1496 //If we need a centrality bin, we select only those events in the corresponding bin.
1497 if( GetCentrality() && fCentralityBin[0] >= 0 && fCentralityBin[1] >= 0 )
1499 Float_t cen = GetEventCentrality();
1500 if(cen > fCentralityBin[1] || cen < fCentralityBin[0]) return ; //reject events out of bin.
1503 // intermediate array with new clusters : init the array only once or clear from previous event
1504 if(!fCaloClusterArr) fCaloClusterArr = new TObjArray(10000);
1505 else fCaloClusterArr->Delete();//Clear("C"); it leaks?
1507 //Get the event, do some checks and settings
1508 CheckAndGetEvent() ;
1512 if(DebugLevel() > 0 ) printf("AliAnalysisTaksEMCALClusterize::UserExec() - Skip Event %d", (Int_t) Entry());
1516 //Init pointers, geometry, clusterizer, ocdb, aodb
1518 InitGeometry(); // only once, must be done before OADB, geo OADB accessed here
1520 if(fAccessOCDB) AccessOCDB();
1521 if(fAccessOADB) AccessOADB(); // only once
1523 InitClusterization();
1529 if (fJustUnfold) ClusterUnfolding();
1530 else ClusterizeCells() ;
1535 FillCaloClusterInEvent();