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 ///////////////////////////////////////////////////////////////////////////////
18 // EMCAL tender, apply corrections to EMCAL clusters and do track matching. //
20 // Author: Deepa Thomas (Utrecht University) //
21 // Later mods/rewrite: Jiri Kral (University of Jyvaskyla) //
22 // S. Aiola, C. Loizides : Make it work for AODs //
24 ///////////////////////////////////////////////////////////////////////////////
29 #include <TInterpreter.h>
30 #include <TObjArray.h>
31 #include <TClonesArray.h>
32 #include <TGeoGlobalMagField.h>
33 #include "AliEMCALAfterBurnerUF.h"
34 #include "AliEMCALClusterizer.h"
35 #include "AliEMCALClusterizerNxN.h"
36 #include "AliEMCALClusterizerv1.h"
37 #include "AliEMCALClusterizerv2.h"
38 #include "AliEMCALDigit.h"
39 #include "AliEMCALGeometry.h"
40 #include "AliEMCALRecParam.h"
41 #include "AliEMCALRecParam.h"
42 #include "AliEMCALRecPoint.h"
43 #include "AliEMCALRecoUtils.h"
44 #include "AliEMCALTenderSupply.h"
45 #include "AliESDCaloCluster.h"
47 #include "AliOADBContainer.h"
48 #include "AliAODEvent.h"
49 #include "AliAnalysisManager.h"
50 #include "AliESDEvent.h"
52 #include "AliTender.h"
54 ClassImp(AliEMCALTenderSupply)
56 AliEMCALTenderSupply::AliEMCALTenderSupply() :
66 ,fNonLinearThreshold(-1)
67 ,fReCalibCluster(kFALSE)
69 ,fCalibrateEnergy(kFALSE)
70 ,fCalibrateTime(kFALSE)
71 ,fDoNonLinearity(kFALSE)
72 ,fBadCellRemove(kFALSE)
73 ,fRejectExoticCells(kFALSE)
74 ,fRejectExoticClusters(kFALSE)
75 ,fClusterBadChannelCheck(kFALSE)
76 ,fRecalClusPos(kFALSE)
78 ,fNCellsFromEMCALBorder(-1)
79 ,fRecalDistToBadChannels(kFALSE)
80 ,fRecalShowerShape(kFALSE)
87 ,fCutEtaPhiSeparate(kFALSE)
92 ,fReClusterize(kFALSE)
94 ,fGeomMatrixSet(kFALSE)
95 ,fLoadGeomMatrices(kFALSE)
97 ,fDoTrackMatch(kFALSE)
98 ,fDoUpdateOnly(kFALSE)
102 ,fMisalignSurvey(kdefault)
103 ,fExoticCellFraction(-1)
104 ,fExoticCellDiffTime(-1)
105 ,fExoticCellMinAmplitude(-1)
107 // Default constructor.
109 for(Int_t i = 0; i < 12; i++) fEMCALMatrix[i] = 0 ;
112 //_____________________________________________________
113 AliEMCALTenderSupply::AliEMCALTenderSupply(const char *name, const AliTender *tender) :
114 AliTenderSupply(name,tender)
123 ,fNonLinearThreshold(-1)
124 ,fReCalibCluster(kFALSE)
126 ,fCalibrateEnergy(kFALSE)
127 ,fCalibrateTime(kFALSE)
128 ,fDoNonLinearity(kFALSE)
129 ,fBadCellRemove(kFALSE)
130 ,fRejectExoticCells(kFALSE)
131 ,fRejectExoticClusters(kFALSE)
132 ,fClusterBadChannelCheck(kFALSE)
133 ,fRecalClusPos(kFALSE)
135 ,fNCellsFromEMCALBorder(-1)
136 ,fRecalDistToBadChannels(kFALSE)
137 ,fRecalShowerShape(kFALSE)
143 ,fCutEtaPhiSum(kTRUE)
144 ,fCutEtaPhiSeparate(kFALSE)
149 ,fReClusterize(kFALSE)
151 ,fGeomMatrixSet(kFALSE)
152 ,fLoadGeomMatrices(kFALSE)
154 ,fDoTrackMatch(kFALSE)
155 ,fDoUpdateOnly(kFALSE)
159 ,fMisalignSurvey(kdefault)
160 ,fExoticCellFraction(-1)
161 ,fExoticCellDiffTime(-1)
162 ,fExoticCellMinAmplitude(-1)
166 for(Int_t i = 0; i < 12; i++) fEMCALMatrix[i] = 0 ;
169 //_____________________________________________________
170 AliEMCALTenderSupply::AliEMCALTenderSupply(const char *name, AliAnalysisTaskSE *task) :
171 AliTenderSupply(name)
180 ,fNonLinearThreshold(-1)
181 ,fReCalibCluster(kFALSE)
183 ,fCalibrateEnergy(kFALSE)
184 ,fCalibrateTime(kFALSE)
185 ,fDoNonLinearity(kFALSE)
186 ,fBadCellRemove(kFALSE)
187 ,fRejectExoticCells(kFALSE)
188 ,fRejectExoticClusters(kFALSE)
189 ,fClusterBadChannelCheck(kFALSE)
190 ,fRecalClusPos(kFALSE)
192 ,fNCellsFromEMCALBorder(-1)
193 ,fRecalDistToBadChannels(kFALSE)
194 ,fRecalShowerShape(kFALSE)
200 ,fCutEtaPhiSum(kTRUE)
201 ,fCutEtaPhiSeparate(kFALSE)
206 ,fReClusterize(kFALSE)
208 ,fGeomMatrixSet(kFALSE)
209 ,fLoadGeomMatrices(kFALSE)
211 ,fDoTrackMatch(kFALSE)
212 ,fDoUpdateOnly(kFALSE)
216 ,fMisalignSurvey(kdefault)
217 ,fExoticCellFraction(-1)
218 ,fExoticCellDiffTime(-1)
219 ,fExoticCellMinAmplitude(-1)
221 // Named constructor.
223 for(Int_t i = 0; i < 12; i++) fEMCALMatrix[i] = 0 ;
226 //_____________________________________________________
227 AliEMCALTenderSupply::~AliEMCALTenderSupply()
231 if (!AliAnalysisManager::GetAnalysisManager()->IsProofMode())
233 delete fEMCALRecoUtils;
238 fDigitsArr->Clear("C");
248 //_____________________________________________________
249 void AliEMCALTenderSupply::SetDefaults()
251 // Set default settings.
253 SwitchOnBadCellRemove();
254 SwitchOnExoticCellRemove();
255 SwitchOnCalibrateEnergy();
256 SwitchOnCalibrateTime();
257 SwitchOnUpdateCell();
258 SwitchOnReclustering();
259 SwitchOnClusterBadChannelCheck();
260 SwitchOnClusterExoticChannelCheck();
261 SwitchOnTrackMatch();
264 //_____________________________________________________
265 Bool_t AliEMCALTenderSupply::RunChanged() const
269 return (fTender && fTender->RunChanged()) || (fTask && fRun != fTask->InputEvent()->GetRunNumber());
272 //_____________________________________________________
273 void AliEMCALTenderSupply::Init()
275 // Initialise EMCAL tender.
278 AliWarning("Init EMCAL Tender supply");
280 if (fConfigName.Length()>0 && gROOT->LoadMacro(fConfigName) >=0) {
281 AliDebug(1, Form("Loading settings from macro %s", fConfigName.Data()));
282 AliEMCALTenderSupply *tender = (AliEMCALTenderSupply*)gInterpreter->ProcessLine("ConfigEMCALTenderSupply()");
283 fDebugLevel = tender->fDebugLevel;
284 fEMCALGeoName = tender->fEMCALGeoName;
285 fEMCALRecoUtils = tender->fEMCALRecoUtils;
286 fConfigName = tender->fConfigName;
287 fNonLinearFunc = tender->fNonLinearFunc;
288 fNonLinearThreshold = tender->fNonLinearThreshold;
289 fReCalibCluster = tender->fReCalibCluster;
290 fUpdateCell = tender->fUpdateCell;
291 fRecalClusPos = tender->fRecalClusPos;
292 fCalibrateEnergy = tender->fCalibrateEnergy;
293 fCalibrateTime = tender->fCalibrateTime;
294 fFiducial = tender->fFiducial;
295 fNCellsFromEMCALBorder = tender->fNCellsFromEMCALBorder;
296 fRecalDistToBadChannels = tender->fRecalDistToBadChannels;
297 fRecalShowerShape = tender->fRecalShowerShape;
298 fClusterBadChannelCheck = tender->fClusterBadChannelCheck;
299 fBadCellRemove = tender->fBadCellRemove;
300 fRejectExoticCells = tender->fRejectExoticCells;
301 fRejectExoticClusters = tender->fRejectExoticClusters;
302 fMass = tender->fMass;
303 fStep = tender->fStep;
304 fCutEtaPhiSum = tender->fCutEtaPhiSum;
305 fCutEtaPhiSeparate = tender->fCutEtaPhiSeparate;
306 fRcut = tender->fRcut;
307 fEtacut = tender->fEtacut;
308 fPhicut = tender->fPhicut;
309 fReClusterize = tender->fReClusterize;
310 fLoadGeomMatrices = tender->fLoadGeomMatrices;
311 fRecParam = tender->fRecParam;
312 fDoNonLinearity = tender->fDoNonLinearity;
313 fDoTrackMatch = tender->fDoTrackMatch;
314 fDoUpdateOnly = tender->fDoUpdateOnly;
315 fMisalignSurvey = tender->fMisalignSurvey;
316 fExoticCellFraction = tender->fExoticCellFraction;
317 fExoticCellDiffTime = tender->fExoticCellDiffTime;
318 fExoticCellMinAmplitude = tender->fExoticCellMinAmplitude;
320 for(Int_t i = 0; i < 12; i++)
321 fEMCALMatrix[i] = tender->fEMCALMatrix[i] ;
325 AliInfo( "Emcal Tender settings: ======================================" );
326 AliInfo( "------------ Switches --------------------------" );
327 AliInfo( Form( "BadCellRemove : %d", fBadCellRemove ));
328 AliInfo( Form( "ExoticCellRemove : %d", fRejectExoticCells ));
329 AliInfo( Form( "CalibrateEnergy : %d", fCalibrateEnergy ));
330 AliInfo( Form( "CalibrateTime : %d", fCalibrateTime ));
331 AliInfo( Form( "UpdateCell : %d", fUpdateCell ));
332 AliInfo( Form( "DoUpdateOnly : %d", fDoUpdateOnly ));
333 AliInfo( Form( "Reclustering : %d", fReClusterize ));
334 AliInfo( Form( "ClusterBadChannelCheck : %d", fClusterBadChannelCheck ));
335 AliInfo( Form( "ClusterExoticChannelCheck : %d", fRejectExoticClusters ));
336 AliInfo( Form( "CellFiducialRegion : %d", fFiducial ));
337 AliInfo( Form( "ReCalibrateCluster : %d", fReCalibCluster ));
338 AliInfo( Form( "RecalculateClusPos : %d", fRecalClusPos ));
339 AliInfo( Form( "RecalShowerShape : %d", fRecalShowerShape ));
340 AliInfo( Form( "NonLinearityCorrection : %d", fDoNonLinearity ));
341 AliInfo( Form( "RecalDistBadChannel : %d", fRecalDistToBadChannels ));
342 AliInfo( Form( "TrackMatch : %d", fDoTrackMatch ));
343 AliInfo( "------------ Variables -------------------------" );
344 AliInfo( Form( "DebugLevel : %d", fDebugLevel ));
345 AliInfo( Form( "BasePath : %s", fBasePath.Data() ));
346 AliInfo( Form( "ConfigFileName : %s", fConfigName.Data() ));
347 AliInfo( Form( "EMCALGeometryName : %s", fEMCALGeoName.Data() ));
348 AliInfo( Form( "NonLinearityFunction : %d", fNonLinearFunc ));
349 AliInfo( Form( "NonLinearityThreshold : %d", fNonLinearThreshold ));
350 AliInfo( Form( "MisalignmentMatrixSurvey : %d", fMisalignSurvey ));
351 AliInfo( Form( "NumberOfCellsFromEMCALBorder : %d", fNCellsFromEMCALBorder ));
352 AliInfo( Form( "RCut : %f", fRcut ));
353 AliInfo( Form( "Mass : %f", fMass ));
354 AliInfo( Form( "Step : %f", fStep ));
355 AliInfo( Form( "EtaCut : %f", fEtacut ));
356 AliInfo( Form( "PhiCut : %f", fPhicut ));
357 AliInfo( Form( "ExoticCellFraction : %f", fExoticCellFraction ));
358 AliInfo( Form( "ExoticCellDiffTime : %f", fExoticCellDiffTime ));
359 AliInfo( Form( "ExoticCellMinAmplitude : %f", fExoticCellMinAmplitude ));
360 AliInfo( "=============================================================" );
366 fEMCALRecoUtils = new AliEMCALRecoUtils;
368 // init geometry if requested
369 if (fEMCALGeoName.Length()>0)
370 fEMCALGeo = AliEMCALGeometry::GetInstance(fEMCALGeoName) ;
373 fDigitsArr = new TClonesArray("AliEMCALDigit",1000);
375 // initialising non-linearity parameters
376 if( fNonLinearThreshold != -1 )
377 fEMCALRecoUtils->SetNonLinearityThreshold(fNonLinearThreshold);
378 if( fNonLinearFunc != -1 )
379 fEMCALRecoUtils->SetNonLinearityFunction(fNonLinearFunc);
381 // missalignment function
382 fEMCALRecoUtils->SetPositionAlgorithm(AliEMCALRecoUtils::kPosTowerGlobal);
385 // do not do the eta0 fiducial cut
386 if( fNCellsFromEMCALBorder != -1 )
387 fEMCALRecoUtils->SetNumberOfCellsFromEMCALBorder(fNCellsFromEMCALBorder);
388 fEMCALRecoUtils->SwitchOnNoFiducialBorderInEMCALEta0();
390 // exotic cell rejection
391 if( fExoticCellFraction != -1 )
392 fEMCALRecoUtils->SetExoticCellFractionCut( fExoticCellFraction );
393 if( fExoticCellDiffTime != -1 )
394 fEMCALRecoUtils->SetExoticCellDiffTimeCut( fExoticCellDiffTime );
395 if( fExoticCellMinAmplitude != -1 )
396 fEMCALRecoUtils->SetExoticCellMinAmplitudeCut( fExoticCellMinAmplitude );
398 // setting track matching parameters ... mass, step size and residual cut
400 fEMCALRecoUtils->SetMass(fMass);
402 fEMCALRecoUtils->SetStep(fStep);
404 // spatial cut based on separate eta/phi or common processing
406 fEMCALRecoUtils->SwitchOnCutEtaPhiSum();
408 fEMCALRecoUtils->SetCutR(fRcut);
409 } else if (fCutEtaPhiSeparate) {
410 fEMCALRecoUtils->SwitchOnCutEtaPhiSeparate();
412 fEMCALRecoUtils->SetCutEta(fEtacut);
414 fEMCALRecoUtils->SetCutPhi(fPhicut);
418 //_____________________________________________________
419 AliVEvent* AliEMCALTenderSupply::GetEvent()
421 // Return the event pointer.
424 return fTender->GetEvent();
427 return fTask->InputEvent();
433 //_____________________________________________________
434 void AliEMCALTenderSupply::ProcessEvent()
438 AliVEvent *event = GetEvent();
441 AliError("Event ptr = 0, returning");
445 // Initialising parameters once per run number
448 AliWarning( "Run changed, initializing parameters" );
449 fRun = event->GetRunNumber();
451 // init geometry if not already done
452 if (fEMCALGeoName.Length()==0) {
453 fEMCALGeoName = "EMCAL_FIRSTYEARV1";
455 fEMCALGeoName = "EMCAL_COMPLETEV1";
458 fEMCALGeoName = "EMCAL_COMPLETE12SMV1";
460 fEMCALGeo = AliEMCALGeometry::GetInstance(fEMCALGeoName);
462 AliFatal(Form("Can not create geometry: %s", fEMCALGeoName.Data()));
470 // define what recalib parameters are needed for various switches
471 // this is based on implementation in AliEMCALRecoUtils
472 Bool_t needRecoParam = fReClusterize;
473 Bool_t needBadChannels = fBadCellRemove | fClusterBadChannelCheck | fRecalDistToBadChannels | fReClusterize;
474 Bool_t needRecalib = fCalibrateEnergy | fReClusterize;
475 Bool_t needTimecalib = fCalibrateTime | fReClusterize;
476 Bool_t needMisalign = fRecalClusPos | fReClusterize;
477 Bool_t needClusterizer = fReClusterize;
479 // initiate reco params with some defaults
480 // will not overwrite, if those have been provided by user
482 Int_t initRC = InitRecParam();
485 AliInfo("Defaults reco params loaded.");
487 AliWarning("User defined reco params.");
491 if( needBadChannels ){
492 Int_t fInitBC = InitBadChannels();
494 AliError("InitBadChannels returned false, returning");
496 AliWarning("InitBadChannels OK");
498 AliWarning(Form("No external hot channel set: %d - %s", event->GetRunNumber(), fFilepass.Data()));
501 // init recalibration factors
504 Int_t fInitRecalib = InitRecalib();
506 AliError("InitRecalib returned false, returning");
508 AliWarning("InitRecalib OK");
510 AliWarning(Form("No recalibration available: %d - %s", event->GetRunNumber(), fFilepass.Data()));
512 Int_t fInitRunDepRecalib = InitRunDepRecalib();
513 if (fInitRunDepRecalib==0)
514 AliError("InitrunDepRecalib returned false, returning");
515 if (fInitRunDepRecalib==1)
516 AliWarning("InitRecalib OK");
517 if (fInitRunDepRecalib>1)
518 AliWarning(Form("No Temperature recalibration available: %d - %s", event->GetRunNumber(), fFilepass.Data()));
524 // init time calibration
526 Int_t initTC = InitTimeCalibration();
528 AliError("InitTimeCalibration returned false, returning");
530 AliWarning("InitTimeCalib OK");
532 AliWarning(Form("No external time calibration set: %d - %s", event->GetRunNumber(), fFilepass.Data()));
535 // init misalignment matrix
537 if (!InitMisalignMatrix())
538 AliError("InitMisalignmentMatrix returned false, returning");
540 AliWarning("InitMisalignMatrix OK");
544 if( needClusterizer ) {
545 if (!InitClusterization())
546 AliError("InitClusterization returned false, returning");
548 AliWarning("InitClusterization OK");
552 fEMCALRecoUtils->Print("");
555 // disable implied switches -------------------------------------------------
556 // AliEMCALRecoUtils or clusterizer functions alredy include some
557 // recalibration so based on those implied callibration te switches
558 // here are disabled to avoid duplication
560 // clusterizer does cluster energy recalibration, position recomputation
563 fReCalibCluster = kFALSE;
564 fRecalClusPos = kFALSE;
565 fRecalShowerShape = kFALSE;
568 // CONFIGURE THE RECO UTILS -------------------------------------------------
569 // configure the reco utils
570 // this option does energy recalibration
571 if( fCalibrateEnergy )
572 fEMCALRecoUtils->SwitchOnRecalibration();
574 fEMCALRecoUtils->SwitchOffRecalibration();
576 // allows time calibration
578 fEMCALRecoUtils->SwitchOnTimeRecalibration();
580 fEMCALRecoUtils->SwitchOffTimeRecalibration();
582 // allows to zero bad cells
584 fEMCALRecoUtils->SwitchOnBadChannelsRemoval();
586 fEMCALRecoUtils->SwitchOffBadChannelsRemoval();
588 // distance to bad channel recalibration
589 if( fRecalDistToBadChannels )
590 fEMCALRecoUtils->SwitchOnDistToBadChannelRecalculation();
592 fEMCALRecoUtils->SwitchOffDistToBadChannelRecalculation();
594 // exclude exotic cells
595 if( fRejectExoticCells )
596 fEMCALRecoUtils->SwitchOnRejectExoticCell();
598 fEMCALRecoUtils->SwitchOffRejectExoticCell();
600 // exclude clusters with exotic cells
601 if( fRejectExoticClusters )
602 fEMCALRecoUtils->SwitchOnRejectExoticCluster();
604 fEMCALRecoUtils->SwitchOffRejectExoticCluster();
606 // TODO: not implemented switches
607 // SwitchOnClusterEnergySmearing
608 // SwitchOnRunDepCorrection
610 // START PROCESSING ---------------------------------------------------------
611 // Test if cells present
612 AliVCaloCells *cells= event->GetEMCALCells();
613 if (cells->GetNumberOfCells()<=0)
616 AliWarning(Form("Number of EMCAL cells = %d, returning", cells->GetNumberOfCells()));
621 AliInfo(Form("Re-calibrate cluster %d\n",fReCalibCluster));
623 // mark the cells not recalibrated in case of selected
624 // time, energy recalibration or bad channel removal
625 if( fCalibrateEnergy || fCalibrateTime || fBadCellRemove )
626 fEMCALRecoUtils->ResetCellsCalibrated();
628 // CELL RECALIBRATION -------------------------------------------------------
629 // cell objects will be updated
630 // the cell calibrations are also processed locally any time those are needed
631 // in case that the cell objects are not to be updated here for later use
635 // includes exotic cell check in the UpdateCells function - is not provided
639 // switch off recalibrations so those are not done multiple times
640 // this is just for safety, the recalibrated flag of cell object
641 // should not allow for farther processing anyways
642 fEMCALRecoUtils->SwitchOffRecalibration();
643 fEMCALRecoUtils->SwitchOffTimeRecalibration();
649 // RECLUSTERIZATION ---------------------------------------------------------
657 // Store good clusters
658 TClonesArray *clusArr = dynamic_cast<TClonesArray*>(event->FindListObject("caloClusters"));
660 clusArr = dynamic_cast<TClonesArray*>(event->FindListObject("CaloClusters"));
662 AliWarning(Form("No cluster array, number of cells in event = %d, returning", cells->GetNumberOfCells()));
666 // PROCESS SINGLE CLUSTER RECALIBRATION -------------------------------------
667 // now go through clusters one by one and process separate correction
668 // as those were defined or not
669 Int_t nclusters = clusArr->GetEntriesFast();
670 for (Int_t icluster=0; icluster < nclusters; ++icluster)
672 AliVCluster *clust = static_cast<AliVCluster*>(clusArr->At(icluster));
675 if (!clust->IsEMCAL())
678 // REMOVE CLUSTERS WITH BAD CELLS -----------------------------
679 if( fClusterBadChannelCheck )
681 // careful, the the ClusterContainsBadChannel is dependent on
682 // SwitchOnBadChannelsRemoval, switching it ON automatically
683 // and returning to original value after processing
684 Bool_t badRemoval = fEMCALRecoUtils->IsBadChannelsRemovalSwitchedOn();
685 fEMCALRecoUtils->SwitchOnBadChannelsRemoval();
687 Bool_t badResult = fEMCALRecoUtils->ClusterContainsBadChannel(fEMCALGeo, clust->GetCellsAbsId(), clust->GetNCells());
689 // switch the bad channels removal back
691 fEMCALRecoUtils->SwitchOffBadChannelsRemoval();
695 delete clusArr->RemoveAt(icluster);
696 continue; //TODO is it really needed to remove it? Or should we flag it?
700 // REMOVE EXOTIC CLUSTERS -------------------------------------
701 // does process local cell recalibration energy and time without replacing
702 // the global cell values, in case of no cell recalib done yet
703 if( fRejectExoticClusters )
705 // careful, the IsExoticCluster is dependent on
706 // SwitchOnRejectExoticCell, switching it ON automatically
707 // and returning to original value after processing
708 Bool_t exRemoval = fEMCALRecoUtils->IsRejectExoticCell();
709 fEMCALRecoUtils->SwitchOnRejectExoticCell();
711 // get bunch crossing
712 Int_t bunchCrossNo = event->GetBunchCrossNumber();
714 Bool_t exResult = fEMCALRecoUtils->IsExoticCluster(clust, cells, bunchCrossNo );
716 // switch the exotic channels removal back
718 fEMCALRecoUtils->SwitchOffRejectExoticCell();
722 delete clusArr->RemoveAt(icluster);
723 continue; //TODO is it really needed to remove it? Or should we flag it?
727 // FIDUCIAL CUT -----------------------------------------------
730 // depends on SetNumberOfCellsFromEMCALBorder
731 // SwitchOnNoFiducialBorderInEMCALEta0
732 if (!fEMCALRecoUtils->CheckCellFiducialRegion(fEMCALGeo, clust, cells)){
733 delete clusArr->RemoveAt(icluster);
734 continue; //TODO it would be nice to store the distance
738 // CLUSTER ENERGY ---------------------------------------------
739 // does process local cell recalibration energy and time without replacing
740 // the global cell values, in case of no cell recalib done yet
741 if( fReCalibCluster ) {
742 fEMCALRecoUtils->RecalibrateClusterEnergy(fEMCALGeo, clust, cells);
743 if (clust->E() < 1e-9) {
744 delete clusArr->RemoveAt(icluster);
749 // CLUSTER POSITION -------------------------------------------
750 // does process local cell energy recalibration, if enabled and cells
751 // not calibrated yet
753 fEMCALRecoUtils->RecalculateClusterPosition(fEMCALGeo, cells, clust);
755 // SHOWER SHAPE -----------------------------------------------
756 if( fRecalShowerShape )
757 fEMCALRecoUtils->RecalculateClusterShowerShapeParameters(fEMCALGeo, cells, clust);
759 // NONLINEARITY -----------------------------------------------
760 if( fDoNonLinearity )
762 Float_t correctedEnergy = fEMCALRecoUtils->CorrectClusterEnergyLinearity(clust);
763 clust->SetE(correctedEnergy);
766 // DISTANCE TO BAD CHANNELS -----------------------------------
767 if( fRecalDistToBadChannels )
768 fEMCALRecoUtils->RecalculateClusterDistanceToBadChannel(fEMCALGeo, cells, clust);
776 // TRACK MATCHING -----------------------------------------------------------
777 if (!TGeoGlobalMagField::Instance()->GetField()) {
778 AliESDEvent *esd = dynamic_cast<AliESDEvent*>(event);
780 esd->InitMagneticField();
782 AliAODEvent *aod = dynamic_cast<AliAODEvent*>(event);
783 Double_t curSol = 30000*aod->GetMagneticField()/5.00668;
784 Double_t curDip = 6000 *aod->GetMuonMagFieldScale();
785 AliMagF *field = AliMagF::CreateFieldMap(curSol,curDip);
786 TGeoGlobalMagField::Instance()->SetField(field);
790 fEMCALRecoUtils->FindMatches(event,0x0,fEMCALGeo);
791 fEMCALRecoUtils->SetClusterMatchedToTrack(event);
792 fEMCALRecoUtils->SetTracksMatchedToCluster(event);
795 //_____________________________________________________
796 Bool_t AliEMCALTenderSupply::InitMisalignMatrix()
798 // Initialising misalignment matrices
800 AliVEvent *event = GetEvent();
807 AliInfo("Misalignment matrix already set");
812 AliInfo("Initialising misalignment matrix");
814 if (fLoadGeomMatrices) {
815 for(Int_t mod=0; mod < fEMCALGeo->GetNumberOfSuperModules(); ++mod)
817 if (fEMCALMatrix[mod]){
819 fEMCALMatrix[mod]->Print();
820 fEMCALGeo->SetMisalMatrix(fEMCALMatrix[mod],mod);
823 fGeomMatrixSet = kTRUE;
827 Int_t runGM = event->GetRunNumber();
830 if(fMisalignSurvey == kdefault)
831 { //take default alignment corresponding to run no
832 AliOADBContainer emcalgeoCont(Form("emcal"));
833 emcalgeoCont.InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALlocal2master.root",Form("AliEMCALgeo"));
834 mobj=(TObjArray*)emcalgeoCont.GetObject(runGM,"EmcalMatrices");
837 if(fMisalignSurvey == kSurveybyS)
838 { //take alignment at sector level
839 if (runGM <= 140000) { //2010 data
840 AliOADBContainer emcalgeoCont(Form("emcal2010"));
841 emcalgeoCont.InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALlocal2master.root",Form("AliEMCALgeo"));
842 mobj=(TObjArray*)emcalgeoCont.GetObject(100,"survey10");
844 else if (runGM>140000)
845 { // 2011 LHC11a pass1 data
846 AliOADBContainer emcalgeoCont(Form("emcal2011"));
847 emcalgeoCont.InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALlocal2master.root",Form("AliEMCALgeo"));
848 mobj=(TObjArray*)emcalgeoCont.GetObject(100,"survey11byS");
852 if(fMisalignSurvey == kSurveybyM)
853 { //take alignment at module level
854 if (runGM <= 140000) { //2010 data
855 AliOADBContainer emcalgeoCont(Form("emcal2010"));
856 emcalgeoCont.InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALlocal2master.root",Form("AliEMCALgeo"));
857 mobj=(TObjArray*)emcalgeoCont.GetObject(100,"survey10");
859 else if (runGM>140000)
860 { // 2011 LHC11a pass1 data
861 AliOADBContainer emcalgeoCont(Form("emcal2011"));
862 emcalgeoCont.InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALlocal2master.root",Form("AliEMCALgeo"));
863 mobj=(TObjArray*)emcalgeoCont.GetObject(100,"survey11byM");
869 AliFatal("Geometry matrix array not found");
873 for(Int_t mod=0; mod < (fEMCALGeo->GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
875 fEMCALMatrix[mod] = (TGeoHMatrix*) mobj->At(mod);
876 fEMCALGeo->SetMisalMatrix(fEMCALMatrix[mod],mod);
877 fEMCALMatrix[mod]->Print();
883 //_____________________________________________________
884 Int_t AliEMCALTenderSupply::InitBadChannels()
886 // Initialising bad channel maps
888 AliVEvent *event = GetEvent();
894 AliInfo("Initialising Bad channel map");
896 // init default maps first
897 if( !fEMCALRecoUtils->GetEMCALBadChannelStatusMapArray() )
898 fEMCALRecoUtils->InitEMCALBadChannelStatusMap() ;
900 Int_t runBC = event->GetRunNumber();
902 AliOADBContainer *contBC = new AliOADBContainer("");
904 { //if fBasePath specified in the ->SetBasePath()
905 if (fDebugLevel>0) AliInfo(Form("Loading Bad Channels OADB from given path %s",fBasePath.Data()));
907 TFile *fbad=new TFile(Form("%s/EMCALBadChannels.root",fBasePath.Data()),"read");
908 if (!fbad || fbad->IsZombie())
910 AliFatal(Form("EMCALBadChannels.root was not found in the path provided: %s",fBasePath.Data()));
914 if (fbad) delete fbad;
916 contBC->InitFromFile(Form("%s/EMCALBadChannels.root",fBasePath.Data()),"AliEMCALBadChannels");
919 { // Else choose the one in the $ALICE_ROOT directory
920 if (fDebugLevel>0) AliInfo("Loading Bad Channels OADB from $ALICE_ROOT/OADB/EMCAL");
922 TFile *fbad=new TFile("$ALICE_ROOT/OADB/EMCAL/EMCALBadChannels.root","read");
923 if (!fbad || fbad->IsZombie())
925 AliFatal("$ALICE_ROOT/OADB/EMCAL/EMCALBadChannels.root was not found");
929 if (fbad) delete fbad;
931 contBC->InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALBadChannels.root","AliEMCALBadChannels");
934 TObjArray *arrayBC=(TObjArray*)contBC->GetObject(runBC);
937 AliError(Form("No external hot channel set for run number: %d", runBC));
941 Int_t sms = fEMCALGeo->GetEMCGeometry()->GetNumberOfSuperModules();
942 for (Int_t i=0; i<sms; ++i)
944 TH2I *h = fEMCALRecoUtils->GetEMCALChannelStatusMap(i);
947 h=(TH2I*)arrayBC->FindObject(Form("EMCALBadChannelMap_Mod%d",i));
951 AliError(Form("Can not get EMCALBadChannelMap_Mod%d",i));
955 fEMCALRecoUtils->SetEMCALChannelStatusMap(i,h);
960 //_____________________________________________________
961 Int_t AliEMCALTenderSupply::InitRecalib()
963 // Initialising recalibration factors.
965 AliVEvent *event = GetEvent();
971 AliInfo("Initialising recalibration factors");
973 // init default maps first
974 if( !fEMCALRecoUtils->GetEMCALRecalibrationFactorsArray() )
975 fEMCALRecoUtils->InitEMCALRecalibrationFactors() ;
977 Int_t runRC = event->GetRunNumber();
979 AliOADBContainer *contRF=new AliOADBContainer("");
981 { //if fBasePath specified in the ->SetBasePath()
982 if (fDebugLevel>0) AliInfo(Form("Loading Recalib OADB from given path %s",fBasePath.Data()));
984 TFile *fRecalib= new TFile(Form("%s/EMCALRecalib.root",fBasePath.Data()),"read");
985 if (!fRecalib || fRecalib->IsZombie())
987 AliFatal(Form("EMCALRecalib.root not found in %s",fBasePath.Data()));
991 if (fRecalib) delete fRecalib;
993 contRF->InitFromFile(Form("%s/EMCALRecalib.root",fBasePath.Data()),"AliEMCALRecalib");
996 { // Else choose the one in the $ALICE_ROOT directory
997 if (fDebugLevel>0) AliInfo("Loading Recalib OADB from $ALICE_ROOT/OADB/EMCAL");
999 TFile *fRecalib= new TFile("$ALICE_ROOT/OADB/EMCAL/EMCALRecalib.root","read");
1000 if (!fRecalib || fRecalib->IsZombie())
1002 AliFatal("$ALICE_ROOT/OADB/EMCAL/EMCALRecalib.root was not found");
1006 if (fRecalib) delete fRecalib;
1008 contRF->InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALRecalib.root","AliEMCALRecalib");
1011 TObjArray *recal=(TObjArray*)contRF->GetObject(runRC);
1014 AliError(Form("No Objects for run: %d",runRC));
1018 TObjArray *recalpass=(TObjArray*)recal->FindObject(fFilepass);
1021 AliError(Form("No Objects for run: %d - %s",runRC,fFilepass.Data()));
1025 TObjArray *recalib=(TObjArray*)recalpass->FindObject("Recalib");
1028 AliError(Form("No Recalib histos found for %d - %s",runRC,fFilepass.Data()));
1032 if (fDebugLevel>0) recalib->Print();
1034 Int_t sms = fEMCALGeo->GetEMCGeometry()->GetNumberOfSuperModules();
1035 for (Int_t i=0; i<sms; ++i)
1037 TH2F *h = fEMCALRecoUtils->GetEMCALChannelRecalibrationFactors(i);
1040 h = (TH2F*)recalib->FindObject(Form("EMCALRecalFactors_SM%d",i));
1043 AliError(Form("Could not load EMCALRecalFactors_SM%d",i));
1047 fEMCALRecoUtils->SetEMCALChannelRecalibrationFactors(i,h);
1052 //_____________________________________________________
1053 Int_t AliEMCALTenderSupply::InitRunDepRecalib()
1055 // Initialising recalibration factors.
1057 AliVEvent *event = GetEvent();
1063 AliInfo("Initialising recalibration factors");
1065 // init default maps first
1066 if( !fEMCALRecoUtils->GetEMCALRecalibrationFactorsArray() )
1067 fEMCALRecoUtils->InitEMCALRecalibrationFactors() ;
1069 Int_t runRC = event->GetRunNumber();
1071 AliOADBContainer *contRF=new AliOADBContainer("");
1073 { //if fBasePath specified in the ->SetBasePath()
1074 if (fDebugLevel>0) AliInfo(Form("Loading Recalib OADB from given path %s",fBasePath.Data()));
1076 TFile *fRunDepRecalib= new TFile(Form("%s/EMCALTemperatureCorrCalib.root",fBasePath.Data()),"read");
1077 if (!fRunDepRecalib || fRunDepRecalib->IsZombie())
1079 AliFatal(Form("EMCALTemperatureCorrCalib.root not found in %s",fBasePath.Data()));
1083 if (fRunDepRecalib) delete fRunDepRecalib;
1085 contRF->InitFromFile(Form("%s/EMCALTemperatureCorrCalib.root",fBasePath.Data()),"AliEMCALRunDepTempCalibCorrections");
1088 { // Else choose the one in the $ALICE_ROOT directory
1089 if (fDebugLevel>0) AliInfo("Loading Recalib OADB from $ALICE_ROOT/OADB/EMCAL");
1091 TFile *fRunDepRecalib= new TFile("$ALICE_ROOT/OADB/EMCAL/EMCALTemperatureCorrCalib.root","read");
1092 if (!fRunDepRecalib || fRunDepRecalib->IsZombie())
1094 AliFatal("$ALICE_ROOT/OADB/EMCAL/EMCALTemperatureCorrCalib.root was not found");
1098 if (fRunDepRecalib) delete fRunDepRecalib;
1100 contRF->InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALTemperatureCorrCalib.root","AliEMCALRunDepTempCalibCorrections");
1103 TH1S *rundeprecal=(TH1S*)contRF->GetObject(runRC);
1106 AliWarning(Form("No TemperatureCorrCalib Objects for run: %d",runRC));
1107 // let's get the closest runnumber instead then..
1110 Int_t maxEntry = contRF->GetNumberOfEntries();
1112 while ( (ic < maxEntry) && (contRF->UpperLimit(ic) < runRC) ) {
1117 Int_t closest = lower;
1118 if ( (ic<maxEntry) &&
1119 (contRF->LowerLimit(ic)-runRC) < (runRC - contRF->UpperLimit(lower)) ) {
1123 AliWarning(Form("TemperatureCorrCalib Objects found closest id %d from run: %d", closest, contRF->LowerLimit(closest)));
1124 rundeprecal = (TH1S*) contRF->GetObjectByIndex(closest);
1127 if (fDebugLevel>0) rundeprecal->Print();
1129 Int_t nSM = fEMCALGeo->GetEMCGeometry()->GetNumberOfSuperModules();
1131 for (Int_t ism=0; ism<nSM; ++ism)
1133 for (Int_t icol=0; icol<48; ++icol)
1135 for (Int_t irow=0; irow<24; ++irow)
1137 Float_t factor = fEMCALRecoUtils->GetEMCALChannelRecalibrationFactor(ism,icol,irow);
1139 Int_t absID = fEMCALGeo->GetAbsCellIdFromCellIndexes(ism, irow, icol); // original calibration factor
1140 factor *= rundeprecal->GetBinContent(absID) / 10000. ; // correction dependent on T
1141 //printf("\t ism %d, icol %d, irow %d,absID %d, corrA %2.3f, corrB %2.3f, corrAB %2.3f\n",ism, icol, irow, absID,
1142 // GetEMCALChannelRecalibrationFactor(ism,icol,irow) , rundeprecal->GetBinContent(absID) / 10000., factor);
1143 fEMCALRecoUtils->SetEMCALChannelRecalibrationFactor(ism,icol,irow,factor);
1152 //_____________________________________________________
1153 Int_t AliEMCALTenderSupply::InitTimeCalibration()
1155 // Initialising bad channel maps
1156 AliVEvent *event = GetEvent();
1162 AliInfo("Initialising time calibration map");
1164 // init default maps first
1165 if ( !fEMCALRecoUtils->GetEMCALTimeRecalibrationFactorsArray() )
1166 fEMCALRecoUtils->InitEMCALTimeRecalibrationFactors() ;
1168 Int_t runBC = event->GetRunNumber();
1170 AliOADBContainer *contBC = new AliOADBContainer("");
1172 { //if fBasePath specified in the ->SetBasePath()
1173 if (fDebugLevel>0) AliInfo(Form("Loading time calibration OADB from given path %s",fBasePath.Data()));
1175 TFile *fbad=new TFile(Form("%s/EMCALTimeCalib.root",fBasePath.Data()),"read");
1176 if (!fbad || fbad->IsZombie())
1178 AliFatal(Form("EMCALTimeCalib.root was not found in the path provided: %s",fBasePath.Data()));
1182 if (fbad) delete fbad;
1184 contBC->InitFromFile(Form("%s/EMCALTimeCalib.root",fBasePath.Data()),"AliEMCALTimeCalib");
1187 { // Else choose the one in the $ALICE_ROOT directory
1188 if (fDebugLevel>0) AliInfo("Loading time calibration OADB from $ALICE_ROOT/OADB/EMCAL");
1190 TFile *fbad=new TFile("$ALICE_ROOT/OADB/EMCAL/EMCALTimeCalib.root","read");
1191 if (!fbad || fbad->IsZombie())
1193 AliFatal("$ALICE_ROOT/OADB/EMCAL/EMCALTimeCalib.root was not found");
1197 if (fbad) delete fbad;
1199 contBC->InitFromFile("$ALICE_ROOT/OADB/EMCAL/EMCALTimeCalib.root","AliEMCALTimeCalib");
1202 TObjArray *arrayBC=(TObjArray*)contBC->GetObject(runBC);
1205 AliError(Form("No external time calibration set for run number: %d", runBC));
1209 // Here, it looks for a specific pass
1210 TObjArray *arrayBCpass=(TObjArray*)arrayBC->FindObject(fFilepass);
1213 AliError(Form("No external time calibration set for: %d -%s", runBC,fFilepass.Data()));
1217 if (fDebugLevel>0) arrayBCpass->Print();
1219 for( Int_t i = 0; i < 4; i++ )
1221 TH1F *h = fEMCALRecoUtils->GetEMCALChannelTimeRecalibrationFactors( i );
1225 h = (TH1F*)arrayBCpass->FindObject(Form("hAllTimeAvBC%d",i));
1229 AliError(Form("Can not get hAllTimeAvBC%d",i));
1233 fEMCALRecoUtils->SetEMCALChannelTimeRecalibrationFactors(i,h);
1238 //_____________________________________________________
1239 void AliEMCALTenderSupply::UpdateCells()
1241 //Remove bad cells from the cell list
1242 //Recalibrate energy and time cells
1243 //This is required for later reclusterization
1245 AliVEvent *event = GetEvent();
1249 AliVCaloCells *cells = event->GetEMCALCells();
1250 Int_t bunchCrossNo = event->GetBunchCrossNumber();
1252 fEMCALRecoUtils->RecalibrateCells(cells, bunchCrossNo);
1254 // remove exotic cells - loop through cells and zero the exotic ones
1255 // just like with bad cell rejection in reco utils (inside RecalibrateCells)
1256 if( fRejectExoticCells )
1262 Short_t mclabel = -1;
1263 Bool_t isExot = kFALSE;
1265 // loop through cells
1266 Int_t nEMcell = cells->GetNumberOfCells() ;
1267 for (Int_t iCell = 0; iCell < nEMcell; iCell++)
1269 cells->GetCell( iCell, absId, ecell, tcell, mclabel, efrac );
1271 isExot = fEMCALRecoUtils->IsExoticCell( absId, cells, bunchCrossNo );
1274 cells->SetCell( iCell, absId, 0.0, -1.0, mclabel, efrac );
1276 } // reject exotic cells
1281 //_____________________________________________________
1282 TString AliEMCALTenderSupply::GetBeamType()
1284 // Get beam type : pp-AA-pA
1285 // ESDs have it directly, AODs get it from hardcoded run number ranges
1287 AliVEvent *event = GetEvent();
1290 AliError("Couldn't retrieve event!");
1296 AliESDEvent *esd = dynamic_cast<AliESDEvent*>(event);
1298 const AliESDRun *run = esd->GetESDRun();
1299 beamType = run->GetBeamType();
1303 Int_t runNumber = event->GetRunNumber();
1304 if ((runNumber >= 136851 && runNumber <= 139517) // LHC10h
1305 || (runNumber >= 166529 && runNumber <= 170593)) // LHC11h
1318 //_____________________________________________________
1319 Int_t AliEMCALTenderSupply::InitRecParam()
1321 // Init reco params if not yet exist (probably shipped by the user already)
1323 if( fRecParam != 0 )
1326 TString beamType = GetBeamType();
1328 // set some default reco params
1329 fRecParam = new AliEMCALRecParam();
1330 fRecParam->SetClusteringThreshold(0.100);
1331 fRecParam->SetMinECut(0.050);
1332 fRecParam->SetTimeCut(250);
1333 fRecParam->SetTimeMin(425);
1334 fRecParam->SetTimeMax(825);
1335 if ( beamType == "A-A") {
1336 fRecParam->SetClusterizerFlag(AliEMCALRecParam::kClusterizerv2);
1339 fRecParam->SetClusterizerFlag(AliEMCALRecParam::kClusterizerv1);
1345 //_____________________________________________________
1346 Bool_t AliEMCALTenderSupply::InitClusterization()
1348 // Initialing clusterization/unfolding algorithm and set all the needed parameters.
1350 AliVEvent *event = GetEvent();
1356 AliInfo(Form("Initialising reclustering parameters: Clusterizer type: %d",fRecParam->GetClusterizerFlag()));
1358 //---setup clusterizer
1359 delete fClusterizer;
1360 if (fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv1)
1361 fClusterizer = new AliEMCALClusterizerv1 (fEMCALGeo);
1362 else if (fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv2)
1363 fClusterizer = new AliEMCALClusterizerv2(fEMCALGeo);
1364 else if (fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerNxN)
1366 AliEMCALClusterizerNxN *clusterizer = new AliEMCALClusterizerNxN(fEMCALGeo);
1367 clusterizer->SetNRowDiff(fRecParam->GetNRowDiff());
1368 clusterizer->SetNColDiff(fRecParam->GetNColDiff());
1369 fClusterizer = clusterizer;
1373 AliFatal(Form("Clusterizer < %d > not available", fRecParam->GetClusterizerFlag()));
1377 // Set the clustering parameters
1378 fClusterizer->SetECAClusteringThreshold( fRecParam->GetClusteringThreshold() );
1379 fClusterizer->SetECALogWeight ( fRecParam->GetW0() );
1380 fClusterizer->SetMinECut ( fRecParam->GetMinECut() );
1381 fClusterizer->SetUnfolding ( fRecParam->GetUnfold() );
1382 fClusterizer->SetECALocalMaxCut ( fRecParam->GetLocMaxCut() );
1383 fClusterizer->SetTimeCut ( fRecParam->GetTimeCut() );
1384 fClusterizer->SetTimeMin ( fRecParam->GetTimeMin() );
1385 fClusterizer->SetTimeMax ( fRecParam->GetTimeMax() );
1386 fClusterizer->SetInputCalibrated ( kTRUE );
1387 fClusterizer->SetJustClusters ( kTRUE );
1389 // In case of unfolding after clusterization is requested, set the corresponding parameters
1390 if (fRecParam->GetUnfold())
1392 for (Int_t i = 0; i < 8; ++i)
1394 fClusterizer->SetSSPars(i, fRecParam->GetSSPars(i));
1396 for (Int_t i = 0; i < 3; ++i)
1398 fClusterizer->SetPar5 (i, fRecParam->GetPar5(i));
1399 fClusterizer->SetPar6 (i, fRecParam->GetPar6(i));
1401 fClusterizer->InitClusterUnfolding();
1404 fClusterizer->SetDigitsArr(fDigitsArr);
1405 fClusterizer->SetOutput(0);
1406 fClusterArr = const_cast<TObjArray *>(fClusterizer->GetRecPoints());
1410 //_____________________________________________________
1411 void AliEMCALTenderSupply::FillDigitsArray()
1413 // Fill digits from cells to a TClonesArray.
1415 AliVEvent *event = GetEvent();
1420 fDigitsArr->Clear("C");
1421 AliVCaloCells *cells = event->GetEMCALCells();
1422 Int_t ncells = cells->GetNumberOfCells();
1423 for (Int_t icell = 0, idigit = 0; icell < ncells; ++icell)
1425 Double_t cellAmplitude=0, cellTime=0, efrac = 0;
1426 Short_t cellNumber=0, mcLabel=-1;
1428 if (cells->GetCell(icell, cellNumber, cellAmplitude, cellTime, mcLabel, efrac) != kTRUE)
1431 // Do not add if energy already too low (some cells set to 0 if bad channels)
1432 if (cellAmplitude < fRecParam->GetMinECut())
1435 // If requested, do not include exotic cells
1436 if (fEMCALRecoUtils->IsExoticCell(cellNumber,cells,event->GetBunchCrossNumber()))
1439 new((*fDigitsArr)[idigit]) AliEMCALDigit(mcLabel, mcLabel, cellNumber,
1440 (Float_t)cellAmplitude, (Float_t)cellTime,
1441 AliEMCALDigit::kHG,idigit, 0, 0, 1);
1443 // AliEMCALDigit *digit = static_cast<AliEMCALDigit*>(fDigitsArr->New(idigit));
1444 // digit->SetId(cellNumber);
1445 // digit->SetTime((Float_t)cellTime);
1446 // digit->SetTimeR((Float_t)cellTime);
1447 // digit->SetIndexInList(idigit);
1448 // digit->SetType(AliEMCALDigit::kHG);
1449 // digit->SetAmplitude((Float_t)cellAmplitude);
1454 //_____________________________________________________
1455 void AliEMCALTenderSupply::Clusterize()
1459 fClusterizer->Digits2Clusters("");
1462 //_____________________________________________________
1463 void AliEMCALTenderSupply::UpdateClusters()
1465 // Update ESD cluster list.
1467 AliVEvent *event = GetEvent();
1472 TClonesArray *clus = dynamic_cast<TClonesArray*>(event->FindListObject("caloClusters"));
1474 clus = dynamic_cast<TClonesArray*>(event->FindListObject("CaloClusters"));
1477 AliError(" Null pointer to calo clusters array, returning");
1481 Int_t nents = clus->GetEntriesFast();
1482 for (Int_t i=0; i < nents; ++i)
1484 AliVCluster *c = dynamic_cast<AliVCluster*>(clus->At(i));
1488 delete clus->RemoveAt(i);
1493 RecPoints2Clusters(clus);
1496 //_____________________________________________________
1497 void AliEMCALTenderSupply::RecPoints2Clusters(TClonesArray *clus)
1499 // Convert AliEMCALRecoPoints to AliESDCaloClusters/AliAODCaloClusters.
1500 // Cluster energy, global position, cells and their amplitude fractions are restored.
1502 AliVEvent *event = GetEvent();
1507 Int_t ncls = fClusterArr->GetEntriesFast();
1508 for(Int_t i=0, nout=clus->GetEntriesFast(); i < ncls; ++i)
1510 AliEMCALRecPoint *recpoint = static_cast<AliEMCALRecPoint*>(fClusterArr->At(i));
1512 Int_t ncellsTrue = 0;
1513 const Int_t ncells = recpoint->GetMultiplicity();
1514 UShort_t absIds[ncells];
1515 Double32_t ratios[ncells];
1516 Int_t *dlist = recpoint->GetDigitsList();
1517 Float_t *elist = recpoint->GetEnergiesList();
1518 for (Int_t c = 0; c < ncells; ++c)
1520 AliEMCALDigit *digit = static_cast<AliEMCALDigit*>(fDigitsArr->At(dlist[c]));
1521 absIds[ncellsTrue] = digit->GetId();
1522 ratios[ncellsTrue] = elist[c]/digit->GetAmplitude();
1523 if (ratios[ncellsTrue] < 0.001)
1530 AliWarning("Skipping cluster with no cells");
1534 // calculate new cluster position
1536 recpoint->GetGlobalPosition(gpos);
1540 AliVCluster *c = static_cast<AliVCluster*>(clus->New(nout++));
1542 c->SetType(AliVCluster::kEMCALClusterv1);
1543 c->SetE(recpoint->GetEnergy());
1545 c->SetNCells(ncellsTrue);
1546 c->SetDispersion(recpoint->GetDispersion());
1547 c->SetEmcCpvDistance(-1); //not yet implemented
1548 c->SetChi2(-1); //not yet implemented
1549 c->SetTOF(recpoint->GetTime()) ; //time-of-flight
1550 c->SetNExMax(recpoint->GetNExMax()); //number of local maxima
1551 Float_t elipAxis[2];
1552 recpoint->GetElipsAxis(elipAxis);
1553 c->SetM02(elipAxis[0]*elipAxis[0]) ;
1554 c->SetM20(elipAxis[1]*elipAxis[1]) ;
1555 c->SetCellsAbsId(absIds);
1556 c->SetCellsAmplitudeFraction(ratios);
1560 //_____________________________________________________
1561 void AliEMCALTenderSupply::GetPass()
1563 // Get passx from filename.
1565 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1566 fInputTree = mgr->GetTree();
1570 AliError("Pointer to tree = 0, returning");
1574 fInputFile = fInputTree->GetCurrentFile();
1576 AliError("Null pointer input file, returning");
1580 TString fname(fInputFile->GetName());
1581 if (fname.Contains("pass1")) fFilepass = TString("pass1");
1582 else if (fname.Contains("pass2")) fFilepass = TString("pass2");
1583 else if (fname.Contains("pass3")) fFilepass = TString("pass3");
1584 else if (fname.Contains("pass4")) fFilepass = TString("pass4");
1587 AliError(Form("Pass number string not found: %s", fname.Data()));