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"
49 #include "AliAODMCParticle.h"
52 #include "AliEMCALAfterBurnerUF.h"
53 #include "AliEMCALGeometry.h"
54 #include "AliEMCALClusterizerNxN.h"
55 #include "AliEMCALClusterizerv1.h"
56 #include "AliEMCALClusterizerv2.h"
57 #include "AliEMCALRecPoint.h"
58 #include "AliEMCALDigit.h"
59 #include "AliCaloCalibPedestal.h"
60 #include "AliEMCALCalibData.h"
62 #include "AliAnalysisTaskEMCALClusterize.h"
64 ClassImp(AliAnalysisTaskEMCALClusterize)
66 //______________________________________________________________________________
67 AliAnalysisTaskEMCALClusterize::AliAnalysisTaskEMCALClusterize(const char *name)
68 : AliAnalysisTaskSE(name)
70 , fGeom(0), fGeomName("")
71 , fGeomMatrixSet(kFALSE), fLoadGeomMatrices(kFALSE)
72 , fCalibData(0), fPedestalData(0)
73 , fOCDBpath(""), fAccessOCDB(kFALSE)
74 , fDigitsArr(0), fClusterArr(0), fCaloClusterArr(0)
75 , fRecParam(0), fClusterizer(0)
76 , fUnfolder(0), fJustUnfold(kFALSE)
77 , fOutputAODBranch(0), fOutputAODBranchName(""), fOutputAODBranchSet(0)
78 , fFillAODFile(kFALSE), fFillAODHeader(0)
79 , fFillAODCaloCells(0), fRun(-1)
80 , fRecoUtils(0), fConfigName("")
82 , fCellLabels(), fCellSecondLabels(), fCellTime()
83 , fCellMatchdEta(), fCellMatchdPhi()
84 , fRecalibrateWithClusterTime(0)
85 , fMaxEvent(0), fDoTrackMatching(kFALSE)
86 , fSelectCell(kFALSE), fSelectCellMinE(0), fSelectCellMinFrac(0)
87 , fRemoveLEDEvents(kTRUE),fRemoveExoticEvents(kFALSE)
88 , fImportGeometryFromFile(kFALSE), fImportGeometryFilePath("")
89 , fOADBSet(kFALSE), fAccessOADB(kTRUE), fOADBFilePath("")
90 , fCentralityClass(""), fSelectEMCALEvent(0)
91 , fEMCALEnergyCut(0.), fEMCALNcellsCut (0)
92 , fSetCellMCLabelFromCluster(0)
93 , fRemapMCLabelForAODs(0)
97 for(Int_t i = 0; i < 12; i++) fGeomMatrix[i] = 0;
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 , fOrgClusterCellId()
121 , fCellLabels(), fCellSecondLabels(), fCellTime()
122 , fCellMatchdEta(), fCellMatchdPhi()
123 , fRecalibrateWithClusterTime(0)
124 , fMaxEvent(0), fDoTrackMatching(kFALSE)
125 , fSelectCell(kFALSE), fSelectCellMinE(0), fSelectCellMinFrac(0)
126 , fRemoveLEDEvents(kTRUE), fRemoveExoticEvents(kFALSE)
127 , fImportGeometryFromFile(kFALSE), fImportGeometryFilePath("")
128 , fOADBSet(kFALSE), fAccessOADB(kTRUE), fOADBFilePath("")
129 , fCentralityClass(""), fSelectEMCALEvent(0)
130 , fEMCALEnergyCut(0.), fEMCALNcellsCut (0)
131 , fSetCellMCLabelFromCluster(0)
132 , fRemapMCLabelForAODs(0)
136 for(Int_t i = 0; i < 12; i++) fGeomMatrix[i] = 0;
140 fCentralityBin[0] = fCentralityBin[1]=-1;
145 //_______________________________________________________________
146 AliAnalysisTaskEMCALClusterize::~AliAnalysisTaskEMCALClusterize()
152 fDigitsArr->Clear("C");
158 fClusterArr->Delete();
164 fCaloClusterArr->Delete();
165 delete fCaloClusterArr;
168 if(fClusterizer) delete fClusterizer;
169 if(fUnfolder) delete fUnfolder;
170 if(fRecoUtils) delete fRecoUtils;
174 //_______________________________________________________
175 Bool_t AliAnalysisTaskEMCALClusterize::AcceptEventEMCAL()
177 // Accept event given there is a EMCAL cluster with enough energy, and not noisy, exotic
179 if(!fSelectEMCALEvent) return kTRUE; // accept
181 if(fEMCALEnergyCut <= 0) return kTRUE; // accept
183 Int_t nCluster = InputEvent() -> GetNumberOfCaloClusters();
184 AliVCaloCells * caloCell = InputEvent() -> GetEMCALCells();
185 Int_t bc = InputEvent() -> GetBunchCrossNumber();
187 for(Int_t icalo = 0; icalo < nCluster; icalo++)
189 AliVCluster *clus = (AliVCluster*) (InputEvent()->GetCaloCluster(icalo));
191 if( ( clus->IsEMCAL() ) && ( clus->GetNCells() > fEMCALNcellsCut ) && ( clus->E() > fEMCALEnergyCut ) &&
192 fRecoUtils->IsGoodCluster(clus,fGeom,caloCell,bc))
196 printf("AliAnalysisTaskEMCALClusterize::AcceptEventEMCAL() - Accept : E %2.2f > %2.2f, nCells %d > %d \n",
197 clus->E(), fEMCALEnergyCut, clus->GetNCells(), fEMCALNcellsCut);
205 printf("AliAnalysisTaskEMCALClusterize::AcceptEventEMCAL() - Reject \n");
211 //_______________________________________________
212 void AliAnalysisTaskEMCALClusterize::AccessOADB()
214 // Set the AODB calibration, bad channels etc. parameters at least once
215 // alignment matrices from OADB done in SetGeometryMatrices
218 if(fOADBSet) return ;
220 Int_t runnumber = InputEvent()->GetRunNumber() ;
221 TString pass = GetPass();
223 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Get AODB parameters from EMCAL in %s for run %d, and <%s> \n",fOADBFilePath.Data(),runnumber,pass.Data());
225 Int_t nSM = fGeom->GetNumberOfSuperModules();
228 if(fRecoUtils->IsBadChannelsRemovalSwitchedOn())
230 AliOADBContainer *contBC=new AliOADBContainer("");
231 contBC->InitFromFile(Form("%s/EMCALBadChannels.root",fOADBFilePath.Data()),"AliEMCALBadChannels");
233 TObjArray *arrayBC=(TObjArray*)contBC->GetObject(runnumber);
237 fRecoUtils->SwitchOnDistToBadChannelRecalculation();
238 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Remove EMCAL bad cells \n");
240 for (Int_t i=0; i<nSM; ++i)
242 TH2I *hbm = fRecoUtils->GetEMCALChannelStatusMap(i);
247 hbm=(TH2I*)arrayBC->FindObject(Form("EMCALBadChannelMap_Mod%d",i));
251 AliError(Form("Can not get EMCALBadChannelMap_Mod%d",i));
255 hbm->SetDirectory(0);
256 fRecoUtils->SetEMCALChannelStatusMap(i,hbm);
259 } else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT remove EMCAL bad channels\n"); // run array
262 // Energy Recalibration
263 if(fRecoUtils->IsRecalibrationOn())
265 AliOADBContainer *contRF=new AliOADBContainer("");
267 contRF->InitFromFile(Form("%s/EMCALRecalib.root",fOADBFilePath.Data()),"AliEMCALRecalib");
269 TObjArray *recal=(TObjArray*)contRF->GetObject(runnumber);
273 TObjArray *recalpass=(TObjArray*)recal->FindObject(pass);
277 TObjArray *recalib=(TObjArray*)recalpass->FindObject("Recalib");
281 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Recalibrate EMCAL \n");
282 for (Int_t i=0; i<nSM; ++i)
284 TH2F *h = fRecoUtils->GetEMCALChannelRecalibrationFactors(i);
289 h = (TH2F*)recalib->FindObject(Form("EMCALRecalFactors_SM%d",i));
293 AliError(Form("Could not load EMCALRecalFactors_SM%d",i));
299 fRecoUtils->SetEMCALChannelRecalibrationFactors(i,h);
301 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params object array \n"); // array ok
302 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params for pass\n"); // array pass ok
303 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params for run\n"); // run number array ok
305 } // Recalibration on
307 // Energy Recalibration, apply on top of previous calibration factors
308 if(fRecoUtils->IsRunDepRecalibrationOn())
310 AliOADBContainer *contRFTD=new AliOADBContainer("");
312 contRFTD->InitFromFile(Form("%s/EMCALTemperatureCorrCalib.root",fOADBFilePath.Data()),"AliEMCALRunDepTempCalibCorrections");
314 TH1S *htd=(TH1S*)contRFTD->GetObject(runnumber);
316 //If it did not exist for this run, get closes one
319 AliWarning(Form("No TemperatureCorrCalib Objects for run: %d",runnumber));
320 // let's get the closest runnumber instead then..
323 Int_t maxEntry = contRFTD->GetNumberOfEntries();
325 while ( (ic < maxEntry) && (contRFTD->UpperLimit(ic) < runnumber) ) {
330 Int_t closest = lower;
331 if ( (ic<maxEntry) &&
332 (contRFTD->LowerLimit(ic)-runnumber) < (runnumber - contRFTD->UpperLimit(lower)) ) {
336 AliWarning(Form("TemperatureCorrCalib Objects found closest id %d from run: %d", closest, contRFTD->LowerLimit(closest)));
337 htd = (TH1S*) contRFTD->GetObjectByIndex(closest);
342 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Recalibrate (Temperature) EMCAL \n");
344 for (Int_t ism=0; ism<nSM; ++ism)
346 for (Int_t icol=0; icol<48; ++icol)
348 for (Int_t irow=0; irow<24; ++irow)
350 Float_t factor = fRecoUtils->GetEMCALChannelRecalibrationFactor(ism,icol,irow);
352 Int_t absID = fGeom->GetAbsCellIdFromCellIndexes(ism, irow, icol); // original calibration factor
353 factor *= htd->GetBinContent(absID) / 10000. ; // correction dependent on T
354 //printf("\t ism %d, icol %d, irow %d,absID %d, corrA %2.3f, corrB %2.3f, corrAB %2.3f\n",ism, icol, irow, absID,
355 // GetEMCALChannelRecalibrationFactor(ism,icol,irow) , htd->GetBinContent(absID) / 10000., factor);
356 fRecoUtils->SetEMCALChannelRecalibrationFactor(ism,icol,irow,factor);
360 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL with T variations, no params TH1 \n");
361 } // Run by Run T calibration
363 // Time Recalibration
364 if(fRecoUtils->IsTimeRecalibrationOn())
366 AliOADBContainer *contTRF=new AliOADBContainer("");
368 contTRF->InitFromFile(Form("%s/EMCALTimeCalib.root",fOADBFilePath.Data()),"AliEMCALTimeCalib");
370 TObjArray *trecal=(TObjArray*)contTRF->GetObject(runnumber);
374 TString passTmp = pass;
375 if(pass!="pass1" && pass!="pass2") passTmp = "pass2"; // TEMPORARY FIX FOR LHC11a analysis
376 TObjArray *trecalpass=(TObjArray*)trecal->FindObject(passTmp);
380 printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Time Recalibrate EMCAL \n");
381 for (Int_t ibc = 0; ibc < 4; ++ibc)
383 TH1F *h = fRecoUtils->GetEMCALChannelTimeRecalibrationFactors(ibc);
388 h = (TH1F*)trecalpass->FindObject(Form("hAllTimeAvBC%d",ibc));
392 AliError(Form("Could not load hAllTimeAvBC%d",ibc));
398 fRecoUtils->SetEMCALChannelTimeRecalibrationFactors(ibc,h);
399 } // bunch crossing loop
400 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate time EMCAL, no params for pass\n"); // array pass ok
401 }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate time EMCAL, no params for run\n"); // run number array ok
403 } // Time recalibration on
405 // Parameters already set once, so do not it again
410 //_________________________________________________
411 Bool_t AliAnalysisTaskEMCALClusterize::AccessOCDB()
413 //Access to OCDB stuff
415 fEvent = InputEvent();
418 Warning("AccessOCDB","Event not available!!!");
422 if (fEvent->GetRunNumber()==fRun)
424 fRun = fEvent->GetRunNumber();
426 if(DebugLevel() > 1 )
427 printf("AliAnalysisTaksEMCALClusterize::AccessOCDB() - Begin");
429 AliCDBManager *cdb = AliCDBManager::Instance();
432 if (fOCDBpath.Length())
434 cdb->SetDefaultStorage(fOCDBpath.Data());
435 printf("AliAnalysisTaksEMCALClusterize::AccessOCDB() - Default storage %s",fOCDBpath.Data());
438 cdb->SetRun(fEvent->GetRunNumber());
441 // EMCAL from RAW OCDB
442 if (fOCDBpath.Contains("alien:"))
444 cdb->SetSpecificStorage("EMCAL/Calib/Data","alien://Folder=/alice/data/2010/OCDB");
445 cdb->SetSpecificStorage("EMCAL/Calib/Pedestals","alien://Folder=/alice/data/2010/OCDB");
448 TString path = cdb->GetDefaultStorage()->GetBaseFolder();
452 //Get calibration parameters
455 AliCDBEntry *entry = (AliCDBEntry*)
456 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
458 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
462 AliFatal("Calibration parameters not found in CDB!");
464 //Get calibration parameters
467 AliCDBEntry *entry = (AliCDBEntry*)
468 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
470 if (entry) fPedestalData = (AliCaloCalibPedestal*) entry->GetObject();
474 AliFatal("Dead map not found in CDB!");
479 //_____________________________________________________
480 void AliAnalysisTaskEMCALClusterize::CheckAndGetEvent()
482 // Get the input event, it can depend in embedded events what you want to get
483 // Also check if the quality of the event is good if not reject it
487 //Process events if there is a high energy cluster
488 if(!AcceptEventEMCAL()) return ;
490 AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
491 Int_t eventN = Entry();
492 if(aodIH) eventN = aodIH->GetReadEntry();
494 if (eventN > fMaxEvent)
497 //printf("Clusterizer --- Event %d-- \n",eventN);
499 //Check if input event are embedded events
500 //If so, take output event
501 if (aodIH && aodIH->GetMergeEvents())
505 if(!aodIH->GetMergeEMCALCells())
506 AliFatal("Events merged but not EMCAL cells, check analysis settings!");
510 printf("AliAnalysisTaksEMCALClusterize::UserExec() - Use embedded events\n");
512 printf("\t InputEvent N Clusters %d, N Cells %d\n",InputEvent()->GetNumberOfCaloClusters(),
513 InputEvent()->GetEMCALCells()->GetNumberOfCells());
515 printf("\t MergedEvent N Clusters %d, N Cells %d\n",aodIH->GetEventToMerge()->GetNumberOfCaloClusters(),
516 aodIH->GetEventToMerge()->GetEMCALCells()->GetNumberOfCells());
518 for (Int_t icl=0; icl < aodIH->GetEventToMerge()->GetNumberOfCaloClusters(); icl++)
520 AliAODCaloCluster *sigCluster = aodIH->GetEventToMerge()->GetCaloCluster(icl);
521 if(sigCluster->IsEMCAL()) printf("\t \t Signal cluster: i %d, E %f\n",icl,sigCluster->E());
524 printf("\t OutputEvent N Clusters %d, N Cells %d\n", AODEvent()->GetNumberOfCaloClusters(),
525 AODEvent()->GetEMCALCells()->GetNumberOfCells());
530 fEvent = InputEvent();
531 if(fFillAODCaloCells) FillAODCaloCells();
532 if(fFillAODHeader) FillAODHeader();
537 Error("UserExec","Event not available");
541 //-------------------------------------------------------------------------------------
542 // Reject events if LED was firing, use only for LHC11a data
543 // Reject event if triggered by exotic cell and remove exotic cells if not triggered
544 //-------------------------------------------------------------------------------------
546 if( IsLEDEvent( InputEvent()->GetRunNumber() ) ) { fEvent = 0x0 ; return ; }
548 if( IsExoticEvent() ) { fEvent = 0x0 ; return ; }
550 //-------------------------------------------------------------------------------------
551 // Set the cluster array in the event (output or input)
552 //-------------------------------------------------------------------------------------
556 //Magic line to write events to AOD filem put after event rejection
557 AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
559 else if( !fOutputAODBranchSet )
561 // Create array and put it in the input event, if output AOD not selected, only once
562 InputEvent()->AddObject(fOutputAODBranch);
563 fOutputAODBranchSet = kTRUE;
564 printf("AliAnalysisTaskEMCALClusterize::UserExec() - Add AOD branch <%s> to input event\n",fOutputAODBranchName.Data());
569 //____________________________________________________
570 void AliAnalysisTaskEMCALClusterize::ClusterizeCells()
572 // Recluster calocells, transform them into digits,
573 // feed the clusterizer with them and get new list of clusters
575 //In case of MC, first loop on the clusters and fill MC label to array
576 Int_t nClusters = fEvent->GetNumberOfCaloClusters();
577 Int_t nClustersOrg = 0;
579 AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
580 if(aodIH && aodIH->GetEventToMerge()) //Embedding
581 nClusters = aodIH->GetEventToMerge()->GetNumberOfCaloClusters(); //Get clusters directly from embedded signal
585 for (Int_t i = 0; i < nClusters; i++)
587 AliVCluster *clus = 0;
588 if(aodIH && aodIH->GetEventToMerge()) //Embedding
589 clus = aodIH->GetEventToMerge()->GetCaloCluster(i); //Get clusters directly from embedded signal
591 clus = fEvent->GetCaloCluster(i);
597 Int_t label = clus->GetLabel();
599 //printf("Org cluster E %f, Time %e, Index = %d, ID %d, MC label %d\n", clus->E(), clus->GetTOF(),i, clus->GetID(),label );
600 //printf("Original list of labels from old cluster : \n");
601 //for(Int_t imc = 0; imc < clus->GetNLabels(); imc++) printf("\t Label %d\n",clus->GetLabelAt(imc));
603 if (clus->GetNLabels()>=2) label2 = clus->GetLabelAt(1) ;
604 UShort_t * index = clus->GetCellsAbsId() ;
605 for(Int_t icell=0; icell < clus->GetNCells(); icell++ )
607 //printf("\t cell %d, MC label %d\n",index[icell],fEvent->GetEMCALCells()->GetCellMCLabel(index[icell]));
608 fOrgClusterCellId[index[icell]] = i;
609 fCellLabels[index[icell]] = label;
610 fCellSecondLabels[index[icell]] = label2;
611 fCellTime[index[icell]] = clus->GetTOF();
612 fCellMatchdEta[index[icell]] = clus->GetTrackDz();
613 fCellMatchdPhi[index[icell]] = clus->GetTrackDx();
614 //printf(" %d,", index[icell] );
621 // Transform CaloCells into Digits
628 AliVCaloCells *cells = fEvent->GetEMCALCells();
630 Int_t bc = InputEvent()->GetBunchCrossNumber();
631 for (Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
633 // Get cell values, recalibrate and not include bad channels found in analysis, nor cells with too low energy, nor exotic cell
634 id = cells->GetCellNumber(icell);
635 Bool_t accept = fRecoUtils->AcceptCalibrateCell(id,bc,amp,time,cells);
637 // Do not include cells with too low energy, nor exotic cell
638 if( amp < fRecParam->GetMinECut() ||
639 time > fRecParam->GetTimeMax() ||
640 time < fRecParam->GetTimeMin() ) accept = kFALSE;
642 // In case of old AOD analysis cell time is -1 s, approximate replacing by time of the cluster the digit belongs.
643 if (fRecalibrateWithClusterTime)
645 time = fCellTime[id];
646 //printf("cell %d time org %f - ",id, time*1.e9);
647 fRecoUtils->RecalibrateCellTime(id,bc,time);
648 //printf("recal %f\n",time*1.e9);
652 if (accept && fRecoUtils->IsExoticCell(id,cells,bc))
657 if( DebugLevel() > 2 )
658 printf("AliAnalysisTaksEMCALClusterize::ClusterizeCells() - Remove channel absId %d, index %d of %d, amp %f, time %f\n",
659 id,icell, cells->GetNumberOfCells(), amp, time*1.e9);
663 Int_t mcLabel = cells->GetMCLabel(icell);
664 //printf("AliAnalysisTaksEMCALClusterize::ClusterizeCells() - cell %d, mc label %d\n",id,mcLabel);
666 //if(fCellLabels[id]!=mcLabel)printf("mcLabel %d - %d\n",mcLabel,fCellLabels[id]);
667 if ( fSetCellMCLabelFromCluster == 1 ) mcLabel = fCellLabels[id]; // Older aliroot MC productions
668 else if( fSetCellMCLabelFromCluster == 0 && fRemapMCLabelForAODs) RemapMCLabelForAODs(mcLabel);
669 else mcLabel = -1; // found later
671 //printf("\t new label %d\n",mcLabel);
673 // Create the digit, put a fake primary deposited energy to trick the clusterizer
674 // when checking the most likely primary
676 new((*fDigitsArr)[idigit]) AliEMCALDigit( mcLabel, mcLabel, id, amp, time,AliEMCALDigit::kHG,idigit, 0, 0, amp);
677 // Last parameter should be MC deposited energy, since it is not available, add just the cell amplitude so that
678 // we give more weight to the MC label of the cell with highest energy in the cluster
685 //-------------------------------------------------------------------------------------
686 //Do the clusterization
687 //-------------------------------------------------------------------------------------
689 fClusterizer->Digits2Clusters("");
691 //-------------------------------------------------------------------------------------
692 //Transform the recpoints into AliVClusters
693 //-------------------------------------------------------------------------------------
695 RecPoints2Clusters();
699 printf("AliAnalysisTaksEMCALClusterize::UserExec() - No array with CaloClusters, input RecPoints entries %d\n",fClusterArr->GetEntriesFast());
703 if( DebugLevel() > 0 )
705 printf("AliAnalysisTaksEMCALClusterize::ClusterizeCells() - N clusters: before recluster %d, after recluster %d\n",nClustersOrg, fCaloClusterArr->GetEntriesFast());
707 if(fCaloClusterArr->GetEntriesFast() != fClusterArr->GetEntriesFast())
709 printf("\t Some RecRoints not transformed into CaloClusters (clusterizer %d, unfold %d): Input entries %d - Output entries %d - %d (not fast)\n",
710 fRecParam->GetClusterizerFlag(),fRecParam->GetUnfold(),
711 fClusterArr->GetEntriesFast(), fCaloClusterArr->GetEntriesFast(), fCaloClusterArr->GetEntries());
716 //_____________________________________________________
717 void AliAnalysisTaskEMCALClusterize::ClusterUnfolding()
719 // Take the event clusters and unfold them
721 AliVCaloCells *cells = fEvent->GetEMCALCells();
722 Double_t cellAmplitude = 0;
723 Double_t cellTime = 0;
724 Short_t cellNumber = 0;
725 Int_t cellMCLabel = 0;
726 Double_t cellEFrac = 0;
727 Int_t nClustersOrg = 0;
729 // Fill the array with the EMCAL clusters, copy them
730 for (Int_t i = 0; i < fEvent->GetNumberOfCaloClusters(); i++)
732 AliVCluster *clus = fEvent->GetCaloCluster(i);
735 //recalibrate/remove bad channels/etc if requested
736 if(fRecoUtils->ClusterContainsBadChannel(fGeom,clus->GetCellsAbsId(), clus->GetNCells()))
741 if(fRecoUtils->IsRecalibrationOn())
744 fRecoUtils->RecalibrateClusterEnergy(fGeom, clus, cells);
747 for (Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
749 if (cells->GetCell(icell, cellNumber, cellAmplitude, cellTime, cellMCLabel, cellEFrac) != kTRUE)
752 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
753 fGeom->GetCellIndex(cellNumber,imod,iTower,iIphi,iIeta);
754 fGeom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
756 //Do not include bad channels found in analysis?
757 if( fRecoUtils->IsBadChannelsRemovalSwitchedOn() &&
758 fRecoUtils->GetEMCALChannelStatus(imod, ieta, iphi))
761 cells->SetCell(icell, cellNumber, cellAmplitude*fRecoUtils->GetEMCALChannelRecalibrationFactor(imod,ieta,iphi),cellTime);
766 //Cast to ESD or AOD, needed to create the cluster array
767 AliESDCaloCluster * esdCluster = dynamic_cast<AliESDCaloCluster*> (clus);
768 AliAODCaloCluster * aodCluster = dynamic_cast<AliAODCaloCluster*> (clus);
772 fCaloClusterArr->Add( new AliESDCaloCluster(*esdCluster) );
776 fCaloClusterArr->Add( new AliAODCaloCluster(*aodCluster) );
779 Warning("UserExec()"," - Wrong CaloCluster type?");
786 fUnfolder->UnfoldClusters(fCaloClusterArr, cells);
793 //_____________________________________________________
794 void AliAnalysisTaskEMCALClusterize::FillAODCaloCells()
796 // Put calo cells in standard branch
797 AliVCaloCells &eventEMcells = *(fEvent->GetEMCALCells());
798 Int_t nEMcell = eventEMcells.GetNumberOfCells() ;
800 AliAODCaloCells &aodEMcells = *(AODEvent()->GetEMCALCells());
801 aodEMcells.CreateContainer(nEMcell);
802 aodEMcells.SetType(AliVCaloCells::kEMCALCell);
803 Double_t calibFactor = 1.;
804 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
806 Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1;
807 fGeom->GetCellIndex(eventEMcells.GetCellNumber(iCell),imod,iTower,iIphi,iIeta);
808 fGeom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);
810 if(fRecoUtils->IsRecalibrationOn())
812 calibFactor = fRecoUtils->GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
815 if(!fRecoUtils->GetEMCALChannelStatus(imod, ieta, iphi))
816 { //Channel is not declared as bad
817 aodEMcells.SetCell(iCell,eventEMcells.GetCellNumber(iCell),eventEMcells.GetAmplitude(iCell)*calibFactor,
818 eventEMcells.GetTime(iCell),eventEMcells.GetMCLabel(iCell),eventEMcells.GetEFraction(iCell));
822 aodEMcells.SetCell(iCell,eventEMcells.GetCellNumber(iCell),0,-1,-1,0);
829 //__________________________________________________
830 void AliAnalysisTaskEMCALClusterize::FillAODHeader()
832 //Put event header information in standard AOD branch
834 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (fEvent);
835 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (fEvent);
839 for (Int_t i = 0; i < 6; i++) covVtx[i] = 0.;
841 AliAODHeader* header = AODEvent()->GetHeader();
842 header->SetRunNumber(fEvent->GetRunNumber());
846 TTree* tree = fInputHandler->GetTree();
849 TFile* file = tree->GetCurrentFile();
850 if (file) header->SetESDFileName(file->GetName());
853 else if (aodevent) header->SetESDFileName(aodevent->GetHeader()->GetESDFileName());
855 header->SetBunchCrossNumber(fEvent->GetBunchCrossNumber());
856 header->SetOrbitNumber(fEvent->GetOrbitNumber());
857 header->SetPeriodNumber(fEvent->GetPeriodNumber());
858 header->SetEventType(fEvent->GetEventType());
861 if(fEvent->GetCentrality())
862 header->SetCentrality(new AliCentrality(*(fEvent->GetCentrality())));
864 header->SetCentrality(0);
867 header->SetOfflineTrigger(fInputHandler->IsEventSelected()); // propagate the decision of the physics selection
868 if (esdevent) header->SetFiredTriggerClasses(esdevent->GetFiredTriggerClasses());
869 else if (aodevent) header->SetFiredTriggerClasses(aodevent->GetFiredTriggerClasses());
871 header->SetTriggerMask(fEvent->GetTriggerMask());
872 header->SetTriggerCluster(fEvent->GetTriggerCluster());
876 header->SetL0TriggerInputs(esdevent->GetHeader()->GetL0TriggerInputs());
877 header->SetL1TriggerInputs(esdevent->GetHeader()->GetL1TriggerInputs());
878 header->SetL2TriggerInputs(esdevent->GetHeader()->GetL2TriggerInputs());
882 header->SetL0TriggerInputs(aodevent->GetHeader()->GetL0TriggerInputs());
883 header->SetL1TriggerInputs(aodevent->GetHeader()->GetL1TriggerInputs());
884 header->SetL2TriggerInputs(aodevent->GetHeader()->GetL2TriggerInputs());
887 header->SetMagneticField(fEvent->GetMagneticField());
888 //header->SetMuonMagFieldScale(esdevent->GetCurrentDip()/6000.);
890 header->SetZDCN1Energy(fEvent->GetZDCN1Energy());
891 header->SetZDCP1Energy(fEvent->GetZDCP1Energy());
892 header->SetZDCN2Energy(fEvent->GetZDCN2Energy());
893 header->SetZDCP2Energy(fEvent->GetZDCP2Energy());
894 header->SetZDCEMEnergy(fEvent->GetZDCEMEnergy(0),fEvent->GetZDCEMEnergy(1));
896 Float_t diamxy[2]={fEvent->GetDiamondX(),fEvent->GetDiamondY()};
898 fEvent->GetDiamondCovXY(diamcov);
899 header->SetDiamond(diamxy,diamcov);
900 if (esdevent) header->SetDiamondZ(esdevent->GetDiamondZ(),esdevent->GetSigma2DiamondZ());
901 else if (aodevent) header->SetDiamondZ(aodevent->GetDiamondZ(),aodevent->GetSigma2DiamondZ());
904 Int_t nVertices = 1 ;/* = prim. vtx*/;
905 Int_t nCaloClus = fEvent->GetNumberOfCaloClusters();
907 AODEvent()->ResetStd(0, nVertices, 0, 0, 0, nCaloClus, 0, 0);
909 // Access to the AOD container of vertices
910 TClonesArray &vertices = *(AODEvent()->GetVertices());
913 // Add primary vertex. The primary tracks will be defined
914 // after the loops on the composite objects (V0, cascades, kinks)
915 fEvent->GetPrimaryVertex()->GetXYZ(pos);
919 esdevent->GetPrimaryVertex()->GetCovMatrix(covVtx);
920 chi = esdevent->GetPrimaryVertex()->GetChi2toNDF();
924 aodevent->GetPrimaryVertex()->GetCovMatrix(covVtx);
925 chi = aodevent->GetPrimaryVertex()->GetChi2perNDF();//Different from ESD?
928 AliAODVertex * primary = new(vertices[jVertices++])
929 AliAODVertex(pos, covVtx, chi, NULL, -1, AliAODVertex::kPrimary);
930 primary->SetName(fEvent->GetPrimaryVertex()->GetName());
931 primary->SetTitle(fEvent->GetPrimaryVertex()->GetTitle());
935 //___________________________________________________________
936 void AliAnalysisTaskEMCALClusterize::FillCaloClusterInEvent()
938 // Get the CaloClusters array, do some final calculations
939 // and put the clusters in the output or input event
940 // as a separate branch
942 //First recalculate track-matching for the new clusters
945 fRecoUtils->FindMatches(fEvent,fCaloClusterArr,fGeom);
947 //Put the new clusters in the AOD list
949 Int_t kNumberOfCaloClusters = fCaloClusterArr->GetEntriesFast();
951 for(Int_t i = 0; i < kNumberOfCaloClusters; i++)
953 AliAODCaloCluster *newCluster = (AliAODCaloCluster *) fCaloClusterArr->At(i);
955 newCluster->SetID(i);
957 // Correct cluster energy non linearity
958 newCluster->SetE(fRecoUtils->CorrectClusterEnergyLinearity(newCluster));
963 Int_t trackIndex = fRecoUtils->GetMatchedTrackIndex(i);
966 newCluster->AddTrackMatched(fEvent->GetTrack(trackIndex));
968 printf("AliAnalysisTaksEMCALClusterize::UserExec() - Matched Track index %d to new cluster %d \n",trackIndex,i);
971 Float_t dR = 999., dZ = 999.;
972 fRecoUtils->GetMatchedResiduals(newCluster->GetID(),dR,dZ);
973 newCluster->SetTrackDistance(dR,dZ);
977 {// Assign previously assigned matched track in reco, very very rough
978 Int_t absId0 = newCluster->GetCellsAbsId()[0]; // Assign match of first cell in cluster
979 newCluster->SetTrackDistance(fCellMatchdPhi[absId0],fCellMatchdEta[absId0]);
982 //printf("New cluster E %f, Time %e, Id = ", newCluster->E(), newCluster->GetTOF() );
983 //for(Int_t icell=0; icell < newCluster->GetNCells(); icell++ ) printf(" %d,", newCluster->GetCellsAbsId() [icell] );
986 // Calculate distance to bad channel for new cluster. Make sure you give the list of bad channels.
987 fRecoUtils->RecalculateClusterDistanceToBadChannel(fGeom, fEvent->GetEMCALCells(), newCluster);
989 new((*fOutputAODBranch)[i]) AliAODCaloCluster(*newCluster);
991 if(DebugLevel() > 1 )
992 printf("AliAnalysisTaksEMCALClusterize::UserExec() - New cluster %d of %d, energy %f\n",newCluster->GetID(), kNumberOfCaloClusters, newCluster->E());
996 fOutputAODBranch->Expand(kNumberOfCaloClusters); // resize TObjArray to 'remove' slots
1001 //_______________________________________________
1002 TString AliAnalysisTaskEMCALClusterize::GetPass()
1004 // Get passx from filename.
1006 if (!AliAnalysisManager::GetAnalysisManager()->GetTree())
1008 AliError("AliAnalysisTaskEMCALClusterize::GetPass() - Pointer to tree = 0, returning null\n");
1012 if (!AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile())
1014 AliError("AliAnalysisTaskEMCALClusterize::GetPass() - Null pointer input file, returning null\n");
1018 TString pass(AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile()->GetName());
1019 if (pass.Contains("ass1")) return TString("pass1");
1020 else if (pass.Contains("ass2")) return TString("pass2");
1021 else if (pass.Contains("ass3")) return TString("pass3");
1022 else if (pass.Contains("ass4")) return TString("pass4");
1023 else if (pass.Contains("ass5")) return TString("pass5");
1025 // No condition fullfilled, give a default value
1026 printf("AliAnalysisTaskEMCALClusterize::GetPass() - Pass number string not found \n");
1031 //_________________________________________
1032 void AliAnalysisTaskEMCALClusterize::Init()
1034 //Init analysis with configuration macro if available
1037 if(fOADBFilePath == "") fOADBFilePath = "$ALICE_ROOT/OADB/EMCAL" ;
1039 fBranchNames = "ESD:AliESDHeader.,EMCALCells.";
1041 if(!fRecParam) fRecParam = new AliEMCALRecParam;
1042 if(!fRecoUtils) fRecoUtils = new AliEMCALRecoUtils();
1044 if(fMaxEvent <= 0) fMaxEvent = 1000000000;
1045 if(fSelectCellMinE <= 0) fSelectCellMinE = 0.005;
1046 if(fSelectCellMinFrac <= 0) fSelectCellMinFrac = 0.001;
1049 if(fCentralityClass == "") fCentralityClass = "V0M";
1051 if (fOCDBpath == "") fOCDBpath = "raw://" ;
1052 if (fOutputAODBranchName == "") fOutputAODBranchName = "newEMCALClusters" ;
1054 if(gROOT->LoadMacro(fConfigName) >=0)
1056 printf("AliAnalysisTaksEMCALClusterize::Init() - Configure analysis with %s\n",fConfigName.Data());
1057 AliAnalysisTaskEMCALClusterize *clus = (AliAnalysisTaskEMCALClusterize*)gInterpreter->ProcessLine("ConfigEMCALClusterize()");
1058 fGeomName = clus->fGeomName;
1059 fLoadGeomMatrices = clus->fLoadGeomMatrices;
1060 fOCDBpath = clus->fOCDBpath;
1061 fAccessOCDB = clus->fAccessOCDB;
1062 fRecParam = clus->fRecParam;
1063 fJustUnfold = clus->fJustUnfold;
1064 fFillAODFile = clus->fFillAODFile;
1065 fRecoUtils = clus->fRecoUtils;
1066 fConfigName = clus->fConfigName;
1067 fMaxEvent = clus->fMaxEvent;
1068 fDoTrackMatching = clus->fDoTrackMatching;
1069 fOutputAODBranchName = clus->fOutputAODBranchName;
1070 for(Int_t i = 0; i < 12; i++) fGeomMatrix[i] = clus->fGeomMatrix[i] ;
1071 fCentralityClass = clus->fCentralityClass;
1072 fCentralityBin[0] = clus->fCentralityBin[0];
1073 fCentralityBin[1] = clus->fCentralityBin[1];
1076 // Init geometry, I do not like much to do it like this ...
1077 if(fImportGeometryFromFile && !gGeoManager)
1079 if (fImportGeometryFilePath == "") fImportGeometryFilePath = "$ALICE_ROOT/OADB/EMCAL/geometry_2011.root" ; // "$ALICE_ROOT/EVE/alice-data/default_geo.root"
1080 printf("AliAnalysisTaskEMCALClusterize::Init() - Import %s\n",fImportGeometryFilePath.Data());
1081 TGeoManager::Import(fImportGeometryFilePath) ;
1086 //_______________________________________________________
1087 void AliAnalysisTaskEMCALClusterize::InitClusterization()
1089 //Select clusterization/unfolding algorithm and set all the needed parameters
1093 // init the unfolding afterburner
1095 fUnfolder = new AliEMCALAfterBurnerUF(fRecParam->GetW0(),fRecParam->GetLocMaxCut(),fRecParam->GetMinECut());
1099 //First init the clusterizer
1100 delete fClusterizer;
1101 if (fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv1)
1102 fClusterizer = new AliEMCALClusterizerv1 (fGeom, fCalibData, fPedestalData);
1103 else if(fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv2)
1104 fClusterizer = new AliEMCALClusterizerv2(fGeom, fCalibData, fPedestalData);
1105 else if(fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerNxN)
1107 fClusterizer = new AliEMCALClusterizerNxN(fGeom, fCalibData, fPedestalData);
1108 fClusterizer->SetNRowDiff(fRecParam->GetNRowDiff());
1109 fClusterizer->SetNColDiff(fRecParam->GetNColDiff());
1113 AliFatal(Form("Clusterizer < %d > not available", fRecParam->GetClusterizerFlag()));
1116 //Now set the parameters
1117 fClusterizer->SetECAClusteringThreshold( fRecParam->GetClusteringThreshold() );
1118 fClusterizer->SetECALogWeight ( fRecParam->GetW0() );
1119 fClusterizer->SetMinECut ( fRecParam->GetMinECut() );
1120 fClusterizer->SetUnfolding ( fRecParam->GetUnfold() );
1121 fClusterizer->SetECALocalMaxCut ( fRecParam->GetLocMaxCut() );
1122 fClusterizer->SetTimeCut ( fRecParam->GetTimeCut() );
1123 fClusterizer->SetTimeMin ( fRecParam->GetTimeMin() );
1124 fClusterizer->SetTimeMax ( fRecParam->GetTimeMax() );
1125 fClusterizer->SetInputCalibrated ( kTRUE );
1126 fClusterizer->SetJustClusters ( kTRUE );
1128 // Initialize the cluster rec points and digits arrays and get them.
1129 fClusterizer->SetOutput(0);
1130 fClusterArr = const_cast<TObjArray *>(fClusterizer->GetRecPoints());
1131 fDigitsArr = fClusterizer->GetDigits();
1133 //In case of unfolding after clusterization is requested, set the corresponding parameters
1134 if(fRecParam->GetUnfold())
1137 for (i = 0; i < 8; i++)
1139 fClusterizer->SetSSPars(i, fRecParam->GetSSPars(i));
1140 }//end of loop over parameters
1142 for (i = 0; i < 3; i++)
1144 fClusterizer->SetPar5 (i, fRecParam->GetPar5(i));
1145 fClusterizer->SetPar6 (i, fRecParam->GetPar6(i));
1146 }//end of loop over parameters
1148 fClusterizer->InitClusterUnfolding();
1153 //_________________________________________________
1154 void AliAnalysisTaskEMCALClusterize::InitGeometry()
1156 // Init geometry and set the geometry matrix, for the first event, skip the rest
1157 // Also set once the run dependent calibrations
1159 if(fGeomMatrixSet) return;
1161 Int_t runnumber = InputEvent()->GetRunNumber() ;
1166 if (runnumber < 140000) fGeomName = "EMCAL_FIRSTYEARV1";
1167 else if(runnumber < 171000) fGeomName = "EMCAL_COMPLETEV1";
1168 else fGeomName = "EMCAL_COMPLETE12SMV1";
1169 printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - Set EMCAL geometry name to <%s> for run %d\n",fGeomName.Data(),runnumber);
1172 fGeom = AliEMCALGeometry::GetInstance(fGeomName);
1176 printf("AliAnalysisTaskEMCALClusterize::InitGeometry(run=%d)",runnumber);
1177 if (!gGeoManager) printf(" - Careful!, gGeoManager not loaded, load misalign matrices");
1180 } // geometry pointer did not exist before
1182 if(fLoadGeomMatrices)
1184 printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - Load user defined EMCAL geometry matrices\n");
1186 // OADB if available
1187 AliOADBContainer emcGeoMat("AliEMCALgeo");
1188 emcGeoMat.InitFromFile(Form("%s/EMCALlocal2master.root",fOADBFilePath.Data()),"AliEMCALgeo");
1189 TObjArray *matEMCAL=(TObjArray*)emcGeoMat.GetObject(runnumber,"EmcalMatrices");
1192 for(Int_t mod=0; mod < (fGeom->GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
1195 if (!fGeomMatrix[mod]) // Get it from OADB
1198 printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - EMCAL matrices SM %d, %p\n",
1199 mod,((TGeoHMatrix*) matEMCAL->At(mod)));
1200 //((TGeoHMatrix*) matEMCAL->At(mod))->Print();
1202 fGeomMatrix[mod] = (TGeoHMatrix*) matEMCAL->At(mod) ;
1205 if(fGeomMatrix[mod])
1207 if(DebugLevel() > 1)
1208 fGeomMatrix[mod]->Print();
1210 fGeom->SetMisalMatrix(fGeomMatrix[mod],mod) ;
1213 fGeomMatrixSet=kTRUE;
1217 else if(!gGeoManager)
1219 printf("AliAnalysisTaksEMCALClusterize::InitGeometry() - Get geo matrices from data");
1220 //Still not implemented in AOD, just a workaround to be able to work at least with ESDs
1221 if(!strcmp(fEvent->GetName(),"AliAODEvent"))
1223 if(DebugLevel() > 1)
1224 Warning("UserExec","Use ideal geometry, values geometry matrix not kept in AODs.");
1228 if(DebugLevel() > 1)
1229 printf("AliAnalysisTaksEMCALClusterize::InitGeometry() - AliAnalysisTaskEMCALClusterize Load Misaligned matrices.");
1231 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(fEvent) ;
1235 Error("InitGeometry"," - This event does not contain ESDs?");
1236 if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);
1240 for(Int_t mod=0; mod < (fGeom->GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
1242 if(DebugLevel() > 1)
1243 esd->GetEMCALMatrix(mod)->Print();
1245 if(esd->GetEMCALMatrix(mod)) fGeom->SetMisalMatrix(esd->GetEMCALMatrix(mod),mod) ;
1249 fGeomMatrixSet=kTRUE;
1252 }//Load matrices from Data
1256 //____________________________________________________
1257 Bool_t AliAnalysisTaskEMCALClusterize::IsExoticEvent()
1260 // Check if event is exotic, get an exotic cell and compare with triggered patch
1261 // If there is a match remove event ... to be completed, filled with something provisional
1263 if(!fRemoveExoticEvents) return kFALSE;
1266 AliVCaloCells * cells = fEvent->GetEMCALCells();
1267 Float_t totCellE = 0;
1268 Int_t bc = InputEvent()->GetBunchCrossNumber();
1269 for(Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
1273 Double_t tcell = 0 ;
1275 Int_t absID = cells->GetCellNumber(icell);
1276 Bool_t accept = fRecoUtils->AcceptCalibrateCell(absID,bc,ecell,tcell,cells);
1277 if(accept && !fRecoUtils->IsExoticCell(absID,cells,bc)) totCellE += ecell;
1280 // TString triggerclasses = "";
1281 // if(esdevent) triggerclasses = esdevent ->GetFiredTriggerClasses();
1282 // else triggerclasses = ((AliAODEvent*)event)->GetFiredTriggerClasses();
1284 // printf("AliAnalysisTaskEMCALClusterize - reject event %d with cluster - reject event with ncells in SM3 %d and SM4 %d\n",(Int_t)Entry(),ncellsSM3, ncellsSM4);
1285 // if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);;
1289 //printf("TotE cell %f\n",totCellE);
1290 if(totCellE < 1) return kTRUE;
1296 //________________________________________________________________
1297 Bool_t AliAnalysisTaskEMCALClusterize::IsLEDEvent(const Int_t run)
1299 //Check if event is LED
1301 if(!fRemoveLEDEvents) return kFALSE;
1303 //check events of LHC11a period
1304 if(run < 146858 || run > 146860) return kFALSE ;
1306 // Count number of cells with energy larger than 0.1 in SM3, cut on this number
1307 Int_t ncellsSM3 = 0;
1308 AliVCaloCells * cells = fEvent->GetEMCALCells();
1309 for(Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
1311 if(cells->GetAmplitude(icell) > 0.1 && cells->GetCellNumber(icell)/(24*48)==3) ncellsSM3++;
1314 TString triggerclasses = "";
1316 AliESDEvent *esdevent = dynamic_cast<AliESDEvent*>(fEvent);
1317 if(esdevent) triggerclasses = esdevent ->GetFiredTriggerClasses();
1318 else triggerclasses = ((AliAODEvent*)fEvent)->GetFiredTriggerClasses();
1320 Int_t ncellcut = 21;
1321 if(triggerclasses.Contains("EMC")) ncellcut = 35;
1323 if( ncellsSM3 >= ncellcut)
1325 printf("AliAnalysisTaksEMCALClusterize::IsLEDEvent() - reject event %d with ncells in SM3 %d\n",(Int_t)Entry(),ncellsSM3);
1326 if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);;
1334 //_______________________________________________________
1335 void AliAnalysisTaskEMCALClusterize::RecPoints2Clusters()
1337 // Restore clusters from recPoints
1338 // Cluster energy, global position, cells and their amplitude
1339 // fractions are restored
1342 for(Int_t i = 0; i < fClusterArr->GetEntriesFast(); i++)
1344 AliEMCALRecPoint *recPoint = (AliEMCALRecPoint*) fClusterArr->At(i);
1346 const Int_t ncells = recPoint->GetMultiplicity();
1347 Int_t ncellsTrue = 0;
1349 if(recPoint->GetEnergy() < fRecParam->GetClusteringThreshold()) continue;
1351 // cells and their amplitude fractions
1352 UShort_t absIds[ncells];
1353 Double32_t ratios[ncells];
1355 //For later check embedding
1356 AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
1358 Float_t clusterE = 0;
1359 for (Int_t c = 0; c < ncells; c++)
1361 AliEMCALDigit *digit = (AliEMCALDigit*) fDigitsArr->At(recPoint->GetDigitsList()[c]);
1363 absIds[ncellsTrue] = digit->GetId();
1364 ratios[ncellsTrue] = recPoint->GetEnergiesList()[c]/digit->GetAmplitude();
1366 // In case of unfolding, remove digits with unfolded energy too low
1369 if (recPoint->GetEnergiesList()[c] < fSelectCellMinE || ratios[ncellsTrue] < fSelectCellMinFrac)
1371 if(DebugLevel() > 1)
1373 printf("AliAnalysisTaksEMCALClusterize::RecPoints2Clusters() - Too small energy in cell of cluster: cluster cell %f, digit %f\n",
1374 recPoint->GetEnergiesList()[c],digit->GetAmplitude());
1382 //Recalculate cluster energy and number of cluster cells in case any of the cells was rejected
1384 clusterE +=recPoint->GetEnergiesList()[c];
1386 // In case of embedding, fill ratio with amount of signal,
1387 if (aodIH && aodIH->GetMergeEvents())
1389 //AliVCaloCells* inEMCALCells = InputEvent()->GetEMCALCells();
1390 AliVCaloCells* meEMCALCells = aodIH->GetEventToMerge()->GetEMCALCells();
1391 AliVCaloCells* ouEMCALCells = AODEvent()->GetEMCALCells();
1393 Float_t sigAmplitude = meEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1394 //Float_t bkgAmplitude = inEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1395 Float_t sumAmplitude = ouEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1396 //printf("\t AbsID %d, amplitude : bkg %f, sigAmplitude %f, summed %f - %f\n",absIds[ncellsTrue], bkgAmplitude, sigAmplitude, sumAmplitude, digit->GetAmplitude());
1398 if(sumAmplitude > 0) ratios[ncellsTrue] = sigAmplitude/sumAmplitude;
1399 //printf("\t \t ratio %f\n",ratios[ncellsTrue]);
1403 }// cluster cell loop
1407 if (DebugLevel() > 1)
1408 printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Skipping cluster with no cells avobe threshold E = %f, ncells %d\n",
1409 recPoint->GetEnergy(), ncells);
1413 //if(ncellsTrue != ncells) printf("Old E %f, ncells %d; New E %f, ncells %d\n",recPoint->GetEnergy(),ncells,clusterE,ncellsTrue);
1415 if(clusterE < fRecParam->GetClusteringThreshold())
1418 printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Remove cluster with energy below seed threshold %f\n",clusterE);
1425 // calculate new cluster position
1427 recPoint->EvalGlobalPosition(fRecParam->GetW0(), fDigitsArr);
1428 recPoint->GetGlobalPosition(gpos);
1431 // create a new cluster
1433 (*fCaloClusterArr)[j] = new AliAODCaloCluster() ;
1434 AliAODCaloCluster *clus = dynamic_cast<AliAODCaloCluster *>( fCaloClusterArr->At(j) ) ;
1436 clus->SetType(AliVCluster::kEMCALClusterv1);
1437 clus->SetE(clusterE);
1438 clus->SetPosition(g);
1439 clus->SetNCells(ncellsTrue);
1440 clus->SetCellsAbsId(absIds);
1441 clus->SetCellsAmplitudeFraction(ratios);
1442 clus->SetChi2(-1); //not yet implemented
1443 clus->SetTOF(recPoint->GetTime()) ; //time-of-flight
1444 clus->SetNExMax(recPoint->GetNExMax()); //number of local maxima
1445 clus->SetDistanceToBadChannel(recPoint->GetDistanceToBadTower());
1447 if(ncells == ncellsTrue)
1449 Float_t elipAxis[2];
1450 recPoint->GetElipsAxis(elipAxis);
1451 clus->SetM02(elipAxis[0]*elipAxis[0]) ;
1452 clus->SetM20(elipAxis[1]*elipAxis[1]) ;
1453 clus->SetDispersion(recPoint->GetDispersion());
1455 else if(fSelectCell)
1457 // In case some cells rejected, in unfolding case, recalculate
1458 // shower shape parameters and position
1459 if(DebugLevel() > 1)
1460 printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Cells removed from cluster (ncells %d, ncellsTrue %d), recalculate Shower Shape\n",ncells,ncellsTrue);
1462 AliVCaloCells* cells = 0x0;
1463 if (aodIH && aodIH->GetMergeEvents()) cells = AODEvent() ->GetEMCALCells();
1464 else cells = InputEvent()->GetEMCALCells();
1466 fRecoUtils->RecalculateClusterShowerShapeParameters(fGeom,cells,clus);
1467 fRecoUtils->RecalculateClusterPID(clus);
1468 fRecoUtils->RecalculateClusterPosition(fGeom,cells,clus);
1474 if ( fSetCellMCLabelFromCluster == 1 ) SetClustersMCLabelFrom2SelectedLabels(recPoint,clus) ;
1475 else if( fSetCellMCLabelFromCluster == 2 ) SetClustersMCLabelFromOriginalClusters(clus) ;
1478 // Normal case, trust what the clusterizer has found
1479 Int_t parentMult = 0;
1480 Int_t *parentList = recPoint->GetParents(parentMult);
1481 clus->SetLabel(parentList, parentMult);
1482 // printf("Label list : ");
1483 // for(Int_t ilabel = 0; ilabel < parentMult; ilabel++ ) printf(" %d ",parentList[ilabel]);
1491 //___________________________________________________________________________
1492 void AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs(Int_t & label)
1494 // MC label for Cells not remapped after ESD filtering, do it here.
1496 if(label < 0) return ;
1498 AliAODEvent * evt = dynamic_cast<AliAODEvent*> (fEvent) ;
1501 TClonesArray * arr = dynamic_cast<TClonesArray*>(evt->FindListObject("mcparticles")) ;
1504 if(label < arr->GetEntriesFast())
1506 AliAODMCParticle * particle = dynamic_cast<AliAODMCParticle *>(arr->At(label));
1507 if(!particle) return ;
1509 if(label == particle->Label()) return ; // label already OK
1510 //else printf("AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs() - Label %d - AOD stack %d \n",label, particle->Label());
1512 //else printf("AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs() - Label %d > AOD labels %d \n",label, arr->GetEntriesFast());
1514 // loop on the particles list and check if there is one with the same label
1515 for(Int_t ind = 0; ind < arr->GetEntriesFast(); ind++ )
1517 AliAODMCParticle * particle = dynamic_cast<AliAODMCParticle *>(arr->At(ind));
1518 if(!particle) continue ;
1520 if(label == particle->Label())
1523 //printf("AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs() - New Label Index %d \n",label);
1530 //printf("AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs() - Label not found set to -1 \n");
1534 //________________________________________________
1535 void AliAnalysisTaskEMCALClusterize::ResetArrays()
1537 // Reset arrays containing information for all possible cells
1538 for(Int_t j = 0; j < 12672; j++)
1540 fOrgClusterCellId[j] = -1;
1541 fCellLabels[j] = -1;
1542 fCellSecondLabels[j] = -1;
1544 fCellMatchdEta[j] = -999;
1545 fCellMatchdPhi[j] = -999;
1549 //_____________________________________________________________________________________________________
1550 void AliAnalysisTaskEMCALClusterize::SetClustersMCLabelFrom2SelectedLabels(AliEMCALRecPoint * recPoint,
1551 AliAODCaloCluster * clus)
1553 // Set the cluster MC label, the digizer was filled with most likely MC label for all cells in original cluster
1554 // Now check the second most likely MC label and add it to the new cluster
1556 Int_t parentMult = 0;
1557 Int_t *parentList = recPoint->GetParents(parentMult);
1558 clus->SetLabel(parentList, parentMult);
1560 //Write the second major contributor to each MC cluster.
1562 for ( Int_t iLoopCell = 0 ; iLoopCell < clus->GetNCells() ; iLoopCell++ )
1565 Int_t idCell = clus->GetCellAbsId(iLoopCell) ;
1568 iNewLabel = 1 ; //iNewLabel makes sure we don't write twice the same label.
1569 for ( UInt_t iLoopLabels = 0 ; iLoopLabels < clus->GetNLabels() ; iLoopLabels++ )
1571 if ( fCellSecondLabels[idCell] == -1 ) iNewLabel = 0; // -1 is never a good second label.
1572 if ( fCellSecondLabels[idCell] == clus->GetLabelAt(iLoopLabels) ) iNewLabel = 0;
1576 Int_t * newLabelArray = new Int_t[clus->GetNLabels()+1] ;
1577 for ( UInt_t iLoopNewLabels = 0 ; iLoopNewLabels < clus->GetNLabels() ; iLoopNewLabels++ )
1579 newLabelArray[iLoopNewLabels] = clus->GetLabelAt(iLoopNewLabels) ;
1582 newLabelArray[clus->GetNLabels()] = fCellSecondLabels[idCell] ;
1583 clus->SetLabel(newLabelArray,clus->GetNLabels()+1) ;
1584 delete [] newLabelArray;
1586 }//positive cell number
1590 //___________________________________________________________________________________________________
1591 void AliAnalysisTaskEMCALClusterize::SetClustersMCLabelFromOriginalClusters(AliAODCaloCluster * clus)
1593 // Get the original clusters that contribute to the new cluster, assign the labels of such clusters
1594 // to the new cluster.
1595 // Only approximatedly valid when output are V1 clusters, handle with care
1597 TArrayI clArray(100) ; //Weird if more than a few clusters are in the origin ...
1600 Int_t nLabTotOrg = 0;
1603 //Find the clusters that originally had the cells
1604 for ( Int_t iLoopCell = 0 ; iLoopCell < clus->GetNCells() ; iLoopCell++ )
1606 Int_t idCell = clus->GetCellAbsId(iLoopCell) ;
1610 Int_t idCluster = fOrgClusterCellId[idCell];
1613 for(Int_t icl =0; icl < nClu; icl++)
1615 if( idCluster == ((Int_t)clArray.GetAt(icl)) ) set = kFALSE;
1616 // printf("\t \t icell %d Cluster in array %d, IdCluster %d, in array %d, set %d\n",
1617 // iLoopCell, icl, idCluster,((Int_t)clArray.GetAt(icl)),set);
1619 if( set && idCluster >= 0)
1621 clArray.SetAt(idCluster,nClu++);
1622 //printf("******** idCluster %d \n",idCluster);
1623 nLabTotOrg+=(fEvent->GetCaloCluster(idCluster))->GetNLabels();
1625 //printf("Cluster in array %d, IdCluster %d\n",nClu-1, idCluster);
1627 //Search highest E cluster
1628 AliVCluster * clOrg = fEvent->GetCaloCluster(idCluster);
1629 //printf("\t E %f\n",clOrg->E());
1630 if(emax < clOrg->E())
1639 if(nClu==0 || nLabTotOrg == 0)
1641 if(clus->E() > 0.25) printf("AliAnalysisTaskEMCALClusterize::SetClustersMCLabelFromOriginalClusters() - Check: N org clusters %d, n tot labels %d, cluster E %f, n cells %d\n",nClu,nLabTotOrg,clus->E(), clus->GetNCells());
1642 //for(Int_t icell = 0; icell < clus->GetNCells(); icell++) printf("\t cell %d",clus->GetCellsAbsId()[icell]);
1647 // Put the first in the list the cluster with highest energy
1648 if(idMax != ((Int_t)clArray.GetAt(0))) // Max not at first position
1650 Int_t maxIndex = -1;
1651 Int_t firstCluster = ((Int_t)clArray.GetAt(0));
1652 for ( Int_t iLoopCluster = 0 ; iLoopCluster < nClu ; iLoopCluster++ )
1654 if(idMax == ((Int_t)clArray.GetAt(iLoopCluster))) maxIndex = iLoopCluster;
1657 clArray.SetAt(idMax,0);
1658 clArray.SetAt(firstCluster,maxIndex);
1662 // Get the labels list in the original clusters, assign all to the new cluster
1663 TArrayI clMCArray(nLabTotOrg) ;
1667 for ( Int_t iLoopCluster = 0 ; iLoopCluster < nClu ; iLoopCluster++ )
1669 Int_t idCluster = clArray.GetAt(iLoopCluster);
1670 //printf("New Cluster in Array %d, idCluster %d \n",iLoopCluster,idCluster);
1671 AliVCluster * clOrg = fEvent->GetCaloCluster(idCluster);
1672 Int_t nLab = clOrg->GetNLabels();
1674 for ( Int_t iLab = 0 ; iLab < nLab ; iLab++ )
1676 Int_t lab = clOrg->GetLabelAt(iLab) ;
1680 //printf("\t \t Set Label %d \n", lab);
1681 for(Int_t iLabTot =0; iLabTot < nLabTot; iLabTot++)
1683 if( lab == ((Int_t)clMCArray.GetAt(iLabTot)) ) set = kFALSE;
1684 //printf("iLoopCluster %d, Label ID in Org Cluster %d,label %d Label ID in array %d, label in array %d, set %d\n",
1685 // iLoopCluster, iLab, lab, iLabTot, ((Int_t)clMCArray.GetAt(iLabTot)),set);
1687 if( set ) clMCArray.SetAt(lab,nLabTot++);
1692 // Set the final list of labels
1694 clus->SetLabel(clMCArray.GetArray(), nLabTot);
1696 // printf("Final list of labels for new cluster : \n");
1697 // for(Int_t ice = 0; ice < clus->GetNCells() ; ice++)
1699 // printf("\t Cell %d ",clus->GetCellsAbsId()[ice]);
1700 // Int_t label = InputEvent()->GetEMCALCells()->GetCellMCLabel(clus->GetCellsAbsId()[ice]);
1701 // printf(" org %d ",label);
1702 // RemapMCLabelForAODs(label);
1703 // printf(" new %d \n",label);
1705 // for(Int_t imc = 0; imc < clus->GetNLabels(); imc++) printf("\t Label %d\n",clus->GetLabelAt(imc));
1709 //____________________________________________________________
1710 void AliAnalysisTaskEMCALClusterize::UserCreateOutputObjects()
1712 // Init geometry, create list of output clusters
1715 fOutputAODBranch = new TClonesArray("AliAODCaloCluster", 0);
1717 if(fOutputAODBranchName.Length()==0)
1719 fOutputAODBranchName = "newEMCALClustersArray";
1720 printf("Cluster branch name not set, set it to newEMCALClustersArray \n");
1723 fOutputAODBranch->SetName(fOutputAODBranchName);
1727 //fOutputAODBranch = new TClonesArray("AliAODCaloCluster", 0);
1730 //fOutputAODBranch->SetOwner(kFALSE);
1732 AddAODBranch("TClonesArray", &fOutputAODBranch);
1737 //_______________________________________________________
1738 void AliAnalysisTaskEMCALClusterize::UserExec(Option_t *)
1740 // Do clusterization event by event, execute different steps
1741 // 1) Do some checks on the kind of events (ESD, AOD) or if some filtering is needed, initializations
1742 // load things and clear arrays
1743 // 2) Clusterize a) just unfolding existing clusters (fJustUnfold)
1744 // b) recluster cells
1745 // + convert cells into digits (calibrating them)
1746 // + recluster digits into recPoints with the chosen clusterizer (V1, V1+Unfold,V2, NxN)
1747 // with methods in AliEMCALClusterizer
1748 // + transform recPoints into CaloClusters
1749 // 3) Do some final calculations in the found clusters (track-matching) and put them in an AOD branch
1755 //Remove the contents of AOD branch output list set in the previous event
1756 fOutputAODBranch->Clear("C");
1760 //Check if there is a centrality value, PbPb analysis, and if a centrality bin selection is requested
1761 //If we need a centrality bin, we select only those events in the corresponding bin.
1762 if( GetCentrality() && fCentralityBin[0] >= 0 && fCentralityBin[1] >= 0 )
1764 Float_t cen = GetEventCentrality();
1765 if(cen > fCentralityBin[1] || cen < fCentralityBin[0]) return ; //reject events out of bin.
1768 // intermediate array with new clusters : init the array only once or clear from previous event
1769 if(!fCaloClusterArr) fCaloClusterArr = new TObjArray(10000);
1770 else fCaloClusterArr->Delete();//Clear("C"); it leaks?
1772 InitGeometry(); // only once, must be done before OADB, geo OADB accessed here
1774 //Get the event, do some checks and settings
1775 CheckAndGetEvent() ;
1779 if(DebugLevel() > 0 ) printf("AliAnalysisTaksEMCALClusterize::UserExec() - Skip Event %d", (Int_t) Entry());
1783 //Init pointers, geometry, clusterizer, ocdb, aodb
1786 if(fAccessOCDB) AccessOCDB();
1787 if(fAccessOADB) AccessOADB(); // only once
1789 InitClusterization();
1795 if (fJustUnfold) ClusterUnfolding();
1796 else ClusterizeCells() ;
1801 FillCaloClusterInEvent();