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 **************************************************************************/
17 //_________________________________________________________________________
18 // Base class for reading data: MonteCarlo, ESD or AOD, of PHOS EMCAL and
19 // Central Barrel Tracking detectors (CTS).
20 // Not all MC particles/tracks/clusters are kept, some kinematical/fiducial restrictions are done.
21 // Mother class of : AliCaloTrackESDReader: Fills ESD data in 3 TObjArrays (PHOS, EMCAL, CTS)
22 // : AliCaloTrackMCReader : Fills Kinematics data in 3 TObjArrays (PHOS, EMCAL, CTS)
23 // : AliCaloTrackAODReader: Fills AOD data in 3 TObjArrays (PHOS, EMCAL, CTS)
24 //-- Author: Gustavo Conesa (LNF-INFN)
25 //////////////////////////////////////////////////////////////////////////////
28 // --- ROOT system ---
31 // ---- ANALYSIS system ----
32 #include "AliMCEvent.h"
33 #include "AliAODMCHeader.h"
34 #include "AliGenPythiaEventHeader.h"
35 #include "AliESDEvent.h"
36 #include "AliAODEvent.h"
37 #include "AliVTrack.h"
38 #include "AliVParticle.h"
39 #include "AliMixedEvent.h"
40 #include "AliESDtrack.h"
41 #include "AliESDtrackCuts.h"
42 #include "AliTriggerAnalysis.h"
43 #include "AliESDVZERO.h"
44 #include "AliVCaloCells.h"
46 // ---- Detectors ----
47 #include "AliPHOSGeoUtils.h"
48 #include "AliEMCALGeometry.h"
51 #include "AliCalorimeterUtils.h"
52 #include "AliCaloTrackReader.h"
54 ClassImp(AliCaloTrackReader)
57 //________________________________________
58 AliCaloTrackReader::AliCaloTrackReader() :
59 TObject(), fEventNumber(-1), //fCurrentFileName(""),
60 fDataType(0), fDebug(0),
61 fFiducialCut(0x0), fCheckFidCut(kFALSE),
62 fComparePtHardAndJetPt(0), fPtHardAndJetPtFactor(7),
63 fCTSPtMin(0), fEMCALPtMin(0), fPHOSPtMin(0),
64 fCTSPtMax(1000), fEMCALPtMax(1000), fPHOSPtMax(1000),
65 fAODBranchList(new TList ),
66 fCTSTracks(new TObjArray()), fEMCALClusters(new TObjArray()), fPHOSClusters(new TObjArray()),
67 fEMCALCells(0x0), fPHOSCells(0x0),
68 fInputEvent(0x0), fOutputEvent(0x0),fMC(0x0),
69 fFillCTS(0), fFillEMCAL(0), fFillPHOS(0),
70 fFillEMCALCells(0), fFillPHOSCells(0),
71 fRecalculateClusters(kFALSE),fSelectEmbeddedClusters(kFALSE),
72 fTrackStatus(0), fTrackFilterMask(0), fESDtrackCuts(0),
73 fTrackMult(0), fTrackMultEtaCut(0.8),
74 fReadStack(kFALSE), fReadAODMCParticles(kFALSE),
75 fDeltaAODFileName("deltaAODPartCorr.root"),
76 fFiredTriggerClassName(""), fAnaLED(kFALSE),
77 fTaskName(""), fCaloUtils(0x0),
78 fMixedEvent(NULL), fNMixedEvent(1), fVertex(NULL),
79 fWriteOutputDeltaAOD(kFALSE),fOldAOD(kFALSE), fCaloFilterPatch(kFALSE),
80 fEMCALClustersListName(""), fZvtxCut(0.),
81 fAcceptFastCluster(kTRUE), fRemoveLEDEvents(kFALSE),
82 fDoEventSelection(kFALSE), fDoV0ANDEventSelection(kFALSE), fUseEventsWithPrimaryVertex(kFALSE),
83 fTriggerAnalysis (new AliTriggerAnalysis),
84 fCentralityClass("V0M"), fCentralityOpt(10),
85 fEventPlaneMethod("Q")
90 //Initialize parameters
94 //_______________________________________
95 AliCaloTrackReader::~AliCaloTrackReader()
102 fAODBranchList->Delete();
103 delete fAODBranchList ;
107 if(fDataType!=kMC)fCTSTracks->Clear() ;
108 else fCTSTracks->Delete() ;
113 if(fDataType!=kMC)fEMCALClusters->Clear("C") ;
114 else fEMCALClusters->Delete() ;
115 delete fEMCALClusters ;
119 if(fDataType!=kMC)fPHOSClusters->Clear("C") ;
120 else fPHOSClusters->Delete() ;
121 delete fPHOSClusters ;
125 for (Int_t i = 0; i < fNMixedEvent; i++) {
126 delete [] fVertex[i] ;
132 delete fESDtrackCuts;
133 delete fTriggerAnalysis;
135 // Pointers not owned, done by the analysis frame
136 // if(fInputEvent) delete fInputEvent ;
137 // if(fOutputEvent) delete fOutputEvent ;
138 // if(fMC) delete fMC ;
139 // Pointer not owned, deleted by maker
140 // if (fCaloUtils) delete fCaloUtils ;
144 //________________________________________________
145 Bool_t AliCaloTrackReader::ComparePtHardAndJetPt()
147 // Check the event, if the requested ptHard is much larger than the jet pT, then there is a problem.
149 if(!fReadStack) return kTRUE; //Information not filtered to AOD
151 if(!strcmp(GetGenEventHeader()->ClassName(), "AliGenPythiaEventHeader")){
153 AliGenPythiaEventHeader* pygeh= (AliGenPythiaEventHeader*) GetGenEventHeader();
154 Int_t nTriggerJets = pygeh->NTriggerJets();
155 Float_t ptHard = pygeh->GetPtHard();
157 //if(fDebug > 1) printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %f\n",nTriggerJets, ptHard);
158 Float_t tmpjet[]={0,0,0,0};
159 for(Int_t ijet = 0; ijet< nTriggerJets; ijet++){
160 pygeh->TriggerJet(ijet, tmpjet);
161 jet = new TParticle(94, 21, -1, -1, -1, -1, tmpjet[0],tmpjet[1],tmpjet[2],tmpjet[3], 0,0,0,0);
162 //Compare jet pT and pt Hard
163 //if(fDebug > 1) printf("AliMCAnalysisUtils:: %d pycell jet pT %f\n",ijet, jet->Pt());
164 if(jet->Pt() > fPtHardAndJetPtFactor * ptHard) {
165 printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %2.2f, pycell jet pT %2.2f, rejection factor %1.1f\n",
166 nTriggerJets, ptHard, jet->Pt(), fPtHardAndJetPtFactor);
177 //____________________________________________
178 AliStack* AliCaloTrackReader::GetStack() const
180 //Return pointer to stack
184 if(fDebug > 1) printf("AliCaloTrackReader::GetStack() - Stack is not available\n");
189 //______________________________________________
190 AliHeader* AliCaloTrackReader::GetHeader() const
192 //Return pointer to header
194 return fMC->Header();
196 printf("AliCaloTrackReader::Header is not available\n");
201 //______________________________________________________________
202 AliGenEventHeader* AliCaloTrackReader::GetGenEventHeader() const
204 //Return pointer to Generated event header
206 return fMC->GenEventHeader();
208 printf("AliCaloTrackReader::GenEventHeader is not available\n");
213 //____________________________________________________________________
214 TClonesArray* AliCaloTrackReader::GetAODMCParticles(Int_t input) const
216 //Return list of particles in AOD. Do it for the corresponding input event.
218 TClonesArray * rv = NULL ;
219 if(fDataType == kAOD){
223 AliAODEvent * evt = dynamic_cast<AliAODEvent*> (fInputEvent) ;
225 rv = (TClonesArray*)evt->FindListObject("mcparticles");
227 printf("AliCaloTrackReader::GetAODMCParticles() - wrong AOD input index? %d, or non existing tree? \n",input);
232 printf("AliCaloTrackReader::GetAODMCParticles() - Input are not AODs\n");
238 //___________________________________________________________________
239 AliAODMCHeader* AliCaloTrackReader::GetAODMCHeader(Int_t input) const
241 //Return MC header in AOD. Do it for the corresponding input event.
242 AliAODMCHeader *mch = NULL;
243 if(fDataType == kAOD){
246 mch = (AliAODMCHeader*)((AliAODEvent*)fInputEvent)->FindListObject("mcheader");
248 // //Second input AOD
249 // else if(input == 1){
250 // mch = (AliAODMCHeader*) fSecondInputAODEvent->FindListObject("mcheader");
253 printf("AliCaloTrackReader::GetAODMCHeader() - wrong AOD input index, %d\n",input);
257 printf("AliCaloTrackReader::GetAODMCHeader() - Input are not AODs\n");
263 //_____________________________
264 void AliCaloTrackReader::Init()
266 //Init reader. Method to be called in AliAnaPartCorrMaker
268 if(fReadStack && fReadAODMCParticles){
269 printf("AliCaloTrackReader::Init() - Cannot access stack and mcparticles at the same time, change them \n");
271 fReadAODMCParticles = kFALSE;
276 //_______________________________________
277 void AliCaloTrackReader::InitParameters()
279 //Initialize the parameters of the analysis.
285 fEMCALPtMax = 1000. ;
288 //Do not filter the detectors input by default.
292 fFillEMCALCells = kFALSE;
293 fFillPHOSCells = kFALSE;
295 fReadStack = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
296 fReadAODMCParticles = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
297 fDeltaAODFileName = "deltaAODPartCorr.root";
298 fFiredTriggerClassName = "";
302 //We want tracks fitted in the detectors:
303 //fTrackStatus=AliESDtrack::kTPCrefit;
304 //fTrackStatus|=AliESDtrack::kITSrefit;
305 fTrackFilterMask = 128; //For AODs, but what is the difference between fTrackStatus and fTrackFilterMask?
307 fESDtrackCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts(); //initialize with TPC only tracks
309 fV0ADC[0] = 0; fV0ADC[1] = 0;
310 fV0Mul[0] = 0; fV0Mul[1] = 0;
315 fCentralityBin[0]=fCentralityBin[1]=-1;
319 //________________________________________________________
320 void AliCaloTrackReader::Print(const Option_t * opt) const
323 //Print some relevant parameters set for the analysis
327 printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
328 printf("Task name : %s\n", fTaskName.Data()) ;
329 printf("Data type : %d\n", fDataType) ;
330 printf("CTS Min pT : %2.1f GeV/c\n", fCTSPtMin) ;
331 printf("EMCAL Min pT : %2.1f GeV/c\n", fEMCALPtMin) ;
332 printf("PHOS Min pT : %2.1f GeV/c\n", fPHOSPtMin) ;
333 printf("CTS Max pT : %2.1f GeV/c\n", fCTSPtMax) ;
334 printf("EMCAL Max pT : %2.1f GeV/c\n", fEMCALPtMax) ;
335 printf("PHOS Max pT : %2.1f GeV/c\n", fPHOSPtMax) ;
336 printf("Use CTS = %d\n", fFillCTS) ;
337 printf("Use EMCAL = %d\n", fFillEMCAL) ;
338 printf("Use PHOS = %d\n", fFillPHOS) ;
339 printf("Use EMCAL Cells = %d\n", fFillEMCALCells) ;
340 printf("Use PHOS Cells = %d\n", fFillPHOSCells) ;
341 printf("Track status = %d\n", (Int_t) fTrackStatus) ;
342 printf("Track filter mask (AODs) = %d\n", (Int_t) fTrackFilterMask) ;
343 printf("Track Mult Eta Cut = %d\n", (Int_t) fTrackMultEtaCut) ;
344 printf("Write delta AOD = %d\n", fWriteOutputDeltaAOD) ;
345 printf("Recalculate Clusters = %d\n", fRecalculateClusters) ;
347 if(fComparePtHardAndJetPt)
348 printf("Compare jet pt and pt hard to accept event, factor = %2.2f",fPtHardAndJetPtFactor);
350 printf("Read Kine from, stack? %d, AOD ? %d \n", fReadStack, fReadAODMCParticles) ;
351 printf("Delta AOD File Name = %s\n", fDeltaAODFileName.Data()) ;
352 printf("Centrality: Class %s, Option %d, Bin [%d,%d] \n", fCentralityClass.Data(),fCentralityOpt,fCentralityBin[0], fCentralityBin[1]) ;
358 //_________________________________________________________________________
359 Bool_t AliCaloTrackReader::FillInputEvent(const Int_t iEntry,
360 const char * /*currentFileName*/)
362 //Fill the event counter and input lists that are needed, called by the analysis maker.
364 fEventNumber = iEntry;
365 //fCurrentFileName = TString(currentFileName);
367 if(fDebug >= 0) printf("AliCaloTrackReader::FillInputEvent() - Input event not available, skip event analysis\n");
370 //Select events only fired by a certain trigger configuration if it is provided
372 if(fInputEvent->GetHeader())
373 eventType = ((AliVHeader*)fInputEvent->GetHeader())->GetEventType();
375 if (GetFiredTriggerClasses().Contains("FAST") && !GetFiredTriggerClasses().Contains("ALL") && !fAcceptFastCluster) {
376 if(fDebug > 0) printf("AliCaloTrackReader::FillInputEvent - Do not count events from fast cluster, trigger name %s\n",fFiredTriggerClassName.Data());
380 //-------------------------------------------------------------------------------------
381 // Reject event if large clusters with large energy
382 // Use only for LHC11a data for the moment, and if input is clusterizer V1 or V1+unfolding
383 // If clusterzer NxN or V2 it does not help
384 //-------------------------------------------------------------------------------------
385 if(fRemoveLEDEvents){
388 //printf("Event %d\n",GetEventNumber());
389 for (Int_t i = 0; i < fInputEvent->GetNumberOfCaloClusters(); i++)
391 AliVCluster *clus = fInputEvent->GetCaloCluster(i);
393 if ((clus->E() > 500 && clus->GetNCells() > 200 ) || clus->GetNCells() > 200) {
394 Int_t absID = clus->GetCellsAbsId()[0];
395 Int_t sm = GetCaloUtils()->GetEMCALGeometry()->GetSuperModuleNumber(absID);
396 if(fDebug > 0) printf("AliCaloTrackReader::FillInputEvent - reject event %d with cluster : E %f, ncells %d, absId(0) %d, SM %d\n",GetEventNumber(),clus->E(), clus->GetNCells(),absID, sm);
402 // Count number of cells with energy larger than 0.1 in SM3, cut on this number
405 for(Int_t icell = 0; icell < fInputEvent->GetEMCALCells()->GetNumberOfCells(); icell++){
406 Int_t absID = fInputEvent->GetEMCALCells()->GetCellNumber(icell);
407 Int_t sm = GetCaloUtils()->GetEMCALGeometry()->GetSuperModuleNumber(absID);
408 if(fInputEvent->GetEMCALCells()->GetAmplitude(icell) > 0.1 && sm==3) ncellsSM3++;
409 if(fInputEvent->GetEMCALCells()->GetAmplitude(icell) > 0.1 && sm==4) ncellsSM4++;
413 if(fFiredTriggerClassName.Contains("EMC")) ncellcut = 35;
415 if(ncellsSM3 >= ncellcut || ncellsSM4 >= 100) {
416 if(fDebug > 0) printf(" AliCaloTrackReader::FillInputEvent() - reject event with ncells in SM3 %d and SM4 %d\n",ncellsSM3, ncellsSM4);
419 }// Remove LED events
421 // Reject pure LED events?
422 if( fFiredTriggerClassName !="" && !fAnaLED){
424 return kFALSE; //Only physics event, do not use for simulated events!!!
426 printf("AliCaloTrackReader::FillInputEvent() - FiredTriggerClass <%s>, selected class <%s>, compare name %d\n",
427 GetFiredTriggerClasses().Data(),fFiredTriggerClassName.Data(), GetFiredTriggerClasses().Contains(fFiredTriggerClassName));
428 if( !GetFiredTriggerClasses().Contains(fFiredTriggerClassName) ) return kFALSE;
429 else if(fDebug > 0) printf("AliCaloTrackReader::FillInputEvent() - Accepted triggered event\n");
432 // kStartOfRun = 1, // START_OF_RUN
433 // kEndOfRun = 2, // END_OF_RUN
434 // kStartOfRunFiles = 3, // START_OF_RUN_FILES
435 // kEndOfRunFiles = 4, // END_OF_RUN_FILES
436 // kStartOfBurst = 5, // START_OF_BURST
437 // kEndOfBurst = 6, // END_OF_BURST
438 // kPhysicsEvent = 7, // PHYSICS_EVENT
439 // kCalibrationEvent = 8, // CALIBRATION_EVENT
440 // kFormatError = 9, // EVENT_FORMAT_ERROR
441 // kStartOfData = 10, // START_OF_DATA
442 // kEndOfData = 11, // END_OF_DATA
443 // kSystemSoftwareTriggerEvent = 12, // SYSTEM_SOFTWARE_TRIGGER_EVENT
444 // kDetectorSoftwareTriggerEvent = 13 // DETECTOR_SOFTWARE_TRIGGER_EVENT
446 if(eventType!=7 && fDebug > 1 )printf("AliCaloTrackReader::FillInputEvent() - DO LED, Event Type <%d>, 8 Calibration \n", eventType);
447 if(eventType!=8)return kFALSE;
450 //In case of analysis of events with jets, skip those with jet pt > 5 pt hard
451 if(fComparePtHardAndJetPt && GetStack()) {
452 if(!ComparePtHardAndJetPt()) return kFALSE ;
457 //Reject events with Z vertex too large, only for SE analysis, if not, cut on the analysis code
458 if(!GetMixedEvent() && TMath::Abs(fVertex[0][2]) > fZvtxCut) return kFALSE;
460 //------------------------------------------------------
461 //Event rejection depending on vertex, pileup, v0and
462 //------------------------------------------------------
463 if(fDoEventSelection){
464 if(!fCaloFilterPatch){
465 //Do not analyze events with pileup
466 Bool_t bPileup = fInputEvent->IsPileupFromSPD(3, 0.8, 3., 2., 5.); //Default values, if not it does not compile
467 //Bool_t bPileup = event->IsPileupFromSPD();
468 if(bPileup) return kFALSE;
470 if(fDoV0ANDEventSelection){
471 Bool_t bV0AND = kTRUE;
472 AliESDEvent* esd = dynamic_cast<AliESDEvent*> (fInputEvent);
474 bV0AND = fTriggerAnalysis->IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0AND);
475 //else bV0AND = //FIXME FOR AODs
476 if(!bV0AND) return kFALSE;
479 if(fUseEventsWithPrimaryVertex && !CheckForPrimaryVertex()) return kFALSE;
483 if(fInputEvent->GetNumberOfCaloClusters() > 0) {
484 AliVCluster * calo = fInputEvent->GetCaloCluster(0);
485 if(calo->GetNLabels() == 4){
486 Int_t * selection = calo->GetLabels();
487 Bool_t bPileup = selection[0];
488 if(bPileup) return kFALSE;
490 Bool_t bGoodV = selection[1];
491 if(fUseEventsWithPrimaryVertex && !bGoodV) return kFALSE;
493 if(fDoV0ANDEventSelection){
494 Bool_t bV0AND = selection[2];
495 if(!bV0AND) return kFALSE;
498 fTrackMult = selection[3];
499 if(fTrackMult == 0) return kFALSE;
501 //First filtered AODs, track multiplicity stored there.
502 fTrackMult = (Int_t) ((AliAODHeader*)fInputEvent->GetHeader())->GetCentrality();
503 if(fTrackMult == 0) return kFALSE;
505 }//at least one cluster
507 //printf("AliCaloTrackReader::FillInputEvent() - No clusters in event\n");
508 //Remove events with vertex (0,0,0), bad vertex reconstruction
509 if(fUseEventsWithPrimaryVertex && TMath::Abs(fVertex[0][0]) < 1.e-6 && TMath::Abs(fVertex[0][1]) < 1.e-6 && TMath::Abs(fVertex[0][2]) < 1.e-6) return kFALSE;
511 //First filtered AODs, track multiplicity stored there.
512 fTrackMult = (Int_t) ((AliAODHeader*)fInputEvent->GetHeader())->GetCentrality();
513 if(fTrackMult == 0) return kFALSE;
515 }// CaloFileter patch
517 //------------------------------------------------------
519 //Check if there is a centrality value, PbPb analysis, and if a centrality bin selection is requested
520 //If we need a centrality bin, we select only those events in the corresponding bin.
521 if(GetCentrality() && fCentralityBin[0]>=0 && fCentralityBin[1]>=0 && fCentralityOpt==100){
522 Int_t cen = GetEventCentrality();
523 if(cen > fCentralityBin[1] || cen < fCentralityBin[0]) return kFALSE; //reject events out of bin.
526 //Fill the arrays with cluster/tracks/cells data
528 FillInputEMCALCells();
530 FillInputPHOSCells();
534 //Accept events with at least one track
535 if(fTrackMult == 0 && fDoEventSelection) return kFALSE;
548 //___________________________________
549 void AliCaloTrackReader::ResetLists()
551 // Reset lists, called by the analysis maker
553 if(fCTSTracks) fCTSTracks -> Clear();
554 if(fEMCALClusters) fEMCALClusters -> Clear("C");
555 if(fPHOSClusters) fPHOSClusters -> Clear("C");
557 fV0ADC[0] = 0; fV0ADC[1] = 0;
558 fV0Mul[0] = 0; fV0Mul[1] = 0;
562 //____________________________________________________________
563 void AliCaloTrackReader::SetInputEvent(AliVEvent* const input)
566 fMixedEvent = dynamic_cast<AliMixedEvent*>(GetInputEvent()) ;
568 fNMixedEvent = fMixedEvent->GetNumberOfEvents() ;
571 //Delete previous vertex
573 for (Int_t i = 0; i < fNMixedEvent; i++) {
574 delete [] fVertex[i] ;
579 fVertex = new Double_t*[fNMixedEvent] ;
580 for (Int_t i = 0; i < fNMixedEvent; i++) {
581 fVertex[i] = new Double_t[3] ;
582 fVertex[i][0] = 0.0 ;
583 fVertex[i][1] = 0.0 ;
584 fVertex[i][2] = 0.0 ;
588 //__________________________________________________
589 Int_t AliCaloTrackReader::GetEventCentrality() const
591 //Return current event centrality
594 if(fCentralityOpt==100) return (Int_t) GetCentrality()->GetCentralityPercentile(fCentralityClass); // 100 bins max
595 else if(fCentralityOpt==10) return GetCentrality()->GetCentralityClass10(fCentralityClass);// 10 bins max
596 else if(fCentralityOpt==20) return GetCentrality()->GetCentralityClass5(fCentralityClass); // 20 bins max
598 printf("AliAnaPartCorrBaseClass::Unknown centrality option %d, use 10, 20 or 100\n",fCentralityOpt);
606 //__________________________________________________________
607 void AliCaloTrackReader::GetVertex(Double_t vertex[3]) const
609 //Return vertex position to be used for single event analysis
610 vertex[0]=fVertex[0][0];
611 vertex[1]=fVertex[0][1];
612 vertex[2]=fVertex[0][2];
615 //____________________________________________________________
616 void AliCaloTrackReader::GetVertex(Double_t vertex[3],
617 const Int_t evtIndex) const
619 //Return vertex position for mixed event, recover the vertex in a particular event.
621 vertex[0]=fVertex[evtIndex][0]; vertex[1]=fVertex[evtIndex][1]; vertex[2]=fVertex[evtIndex][2];
625 //________________________________________
626 void AliCaloTrackReader::FillVertexArray()
629 //Fill data member with vertex
630 //In case of Mixed event, multiple vertices
632 //Delete previous vertex
634 for (Int_t i = 0; i < fNMixedEvent; i++) {
635 delete [] fVertex[i] ;
640 fVertex = new Double_t*[fNMixedEvent] ;
641 for (Int_t i = 0; i < fNMixedEvent; i++) {
642 fVertex[i] = new Double_t[3] ;
643 fVertex[i][0] = 0.0 ;
644 fVertex[i][1] = 0.0 ;
645 fVertex[i][2] = 0.0 ;
648 if (!fMixedEvent) { //Single event analysis
651 if(fInputEvent->GetPrimaryVertex()){
652 fInputEvent->GetPrimaryVertex()->GetXYZ(fVertex[0]);
655 printf("AliCaloTrackReader::FillVertexArray() - NULL primary vertex\n");
656 fVertex[0][0]=0.; fVertex[0][1]=0.; fVertex[0][2]=0.;
657 }//Primary vertex pointer do not exist
659 } else {//MC read event
660 fVertex[0][0]=0.; fVertex[0][1]=0.; fVertex[0][2]=0.;
664 printf("AliCaloTrackReader::FillVertexArray() - Single Event Vertex : %f,%f,%f\n",fVertex[0][0],fVertex[0][1],fVertex[0][2]);
666 } else { // MultiEvent analysis
667 for (Int_t iev = 0; iev < fNMixedEvent; iev++) {
668 if (fMixedEvent->GetVertexOfEvent(iev))
669 fMixedEvent->GetVertexOfEvent(iev)->GetXYZ(fVertex[iev]);
670 else { // no vertex found !!!!
671 AliWarning("No vertex found");
675 printf("AliCaloTrackReader::FillVertexArray() - Multi Event %d Vertex : %f,%f,%f\n",iev,fVertex[iev][0],fVertex[iev][1],fVertex[iev][2]);
682 //_____________________________________
683 void AliCaloTrackReader::FillInputCTS()
685 //Return array with Central Tracking System (CTS) tracks
687 if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputCTS()\n");
689 Int_t nTracks = fInputEvent->GetNumberOfTracks() ;
693 for (Int_t itrack = 0; itrack < nTracks; itrack++) {////////////// track loop
694 AliVTrack * track = (AliVTrack*)fInputEvent->GetTrack(itrack) ; // retrieve track from esd
696 //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
697 if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus))
702 if (fDataType==kESD && !fESDtrackCuts->AcceptTrack((AliESDtrack*)track))
706 else if(fDataType==kAOD)
708 AliAODTrack *aodtrack = dynamic_cast <AliAODTrack*>(track);
711 printf("AliCaloTrackReader::FillInputCTS():AOD track type: %c \n", aodtrack->GetType());
712 if (fDataType!=kMC && aodtrack->TestFilterBit(fTrackFilterMask)==kFALSE) continue;
713 if(aodtrack->GetType()!=AliAODTrack::kPrimary) continue;
717 //Count the tracks in eta < 0.9
718 //printf("Eta %f cut %f\n",TMath::Abs(track->Eta()),fTrackMultEtaCut);
719 if(TMath::Abs(track->Eta())< fTrackMultEtaCut) fTrackMult++;
721 track->GetPxPyPz(p) ;
722 TLorentzVector momentum(p[0],p[1],p[2],0);
724 if(fCTSPtMin < momentum.Pt() && fCTSPtMax > momentum.Pt()){
726 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS"))
729 if(fDebug > 2 && momentum.Pt() > 0.1)
730 printf("AliCaloTrackReader::FillInputCTS() - Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
731 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
734 track->SetID(itrack);
737 fCTSTracks->Add(track);
739 }//Pt and Fiducial cut passed.
742 //fCTSTracksNormalInputEntries = fCTSTracks->GetEntriesFast();
744 printf("AliCaloTrackReader::FillInputCTS() - aod entries %d, input tracks %d, pass status %d, multipliticy %d\n", fCTSTracks->GetEntriesFast(), nTracks, nstatus, fTrackMult);//fCTSTracksNormalInputEntries);
748 //__________________________________________________________________
749 void AliCaloTrackReader::FillInputEMCALAlgorithm(AliVCluster * clus,
752 //Fill the EMCAL data in the array, do it
756 vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
758 //Reject clusters with bad channels, close to borders and exotic;
759 if(!GetCaloUtils()->GetEMCALRecoUtils()->IsGoodCluster(clus,GetCaloUtils()->GetEMCALGeometry(),GetEMCALCells(),fInputEvent->GetBunchCrossNumber())) return;
761 //Mask all cells in collumns facing ALICE thick material if requested
762 if(GetCaloUtils()->GetNMaskCellColumns()){
767 Bool_t shared = kFALSE;
768 GetCaloUtils()->GetEMCALRecoUtils()->GetMaxEnergyCell(GetCaloUtils()->GetEMCALGeometry(), GetEMCALCells(),clus,absId,iSupMod,ieta,iphi,shared);
769 if(GetCaloUtils()->MaskFrameCluster(iSupMod, ieta)) return;
772 if(fSelectEmbeddedClusters){
773 if(clus->GetNLabels()==0 || clus->GetLabel() < 0) return;
774 //else printf("Embedded cluster, %d, n label %d label %d \n",iclus,clus->GetNLabels(),clus->GetLabel());
778 //clus->GetPosition(pos);
779 //printf("Before Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
781 if(fRecalculateClusters){
782 //Recalibrate the cluster energy
783 if(GetCaloUtils()->IsRecalibrationOn()) {
785 Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, GetEMCALCells());
788 //printf("Recalibrated Energy %f\n",clus->E());
790 GetCaloUtils()->RecalculateClusterShowerShapeParameters(GetEMCALCells(),clus);
791 GetCaloUtils()->RecalculateClusterPID(clus);
795 //Recalculate distance to bad channels, if new list of bad channels provided
796 GetCaloUtils()->RecalculateClusterDistanceToBadChannel(GetEMCALCells(),clus);
798 //Recalculate cluster position
799 if(GetCaloUtils()->IsRecalculationOfClusterPositionOn()){
800 GetCaloUtils()->RecalculateClusterPosition(GetEMCALCells(),clus);
801 //clus->GetPosition(pos);
802 //printf("After Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
807 if(GetCaloUtils()->GetEMCALRecoUtils()->IsTimeRecalibrationOn()) {
809 Double_t tof = clus->GetTOF();
811 Int_t absIdMax = GetCaloUtils()->GetMaxEnergyCell(fEMCALCells, clus,frac);
813 if(fDataType==AliCaloTrackReader::kESD){
814 tof = fEMCALCells->GetCellTime(absIdMax);
817 GetCaloUtils()->GetEMCALRecoUtils()->RecalibrateCellTime(absIdMax,fInputEvent->GetBunchCrossNumber(),tof);
821 }// Time recalibration
823 //Correct non linearity
824 if(GetCaloUtils()->IsCorrectionOfClusterEnergyOn()){
825 GetCaloUtils()->CorrectClusterEnergy(clus) ;
826 //printf("Linearity Corrected Energy %f\n",clus->E());
828 //In case of MC analysis, to match resolution/calibration in real data
829 Float_t rdmEnergy = GetCaloUtils()->GetEMCALRecoUtils()->SmearClusterEnergy(clus);
830 // printf("\t Energy %f, smeared %f\n", clus->E(),rdmEnergy);
831 clus->SetE(rdmEnergy);
834 TLorentzVector momentum ;
836 clus->GetMomentum(momentum, fVertex[vindex]);
838 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) return;
840 if(fEMCALPtMin > momentum.E() || fEMCALPtMax < momentum.E()) return;
842 if(fDebug > 2 && momentum.E() > 0.1)
843 printf("AliCaloTrackReader::FillInputEMCAL() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
844 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
849 fEMCALClusters->Add(clus);
853 //_______________________________________
854 void AliCaloTrackReader::FillInputEMCAL()
856 //Return array with EMCAL clusters in aod format
858 if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputEMCAL()\n");
860 // First recalibrate cells, time or energy
861 // if(GetCaloUtils()->IsRecalibrationOn())
862 // GetCaloUtils()->GetEMCALRecoUtils()->RecalibrateCells(GetCaloUtils()->GetEMCALGeometry(),
864 // fInputEvent->GetBunchCrossNumber());
866 //Loop to select clusters in fiducial cut and fill container with aodClusters
867 if(fEMCALClustersListName==""){
868 Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
869 for (Int_t iclus = 0; iclus < nclusters; iclus++) {
870 AliVCluster * clus = 0;
871 if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
872 if (IsEMCALCluster(clus)){
873 FillInputEMCALAlgorithm(clus, iclus);
878 //Recalculate track matching
879 GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent);
881 }//Get the clusters from the input event
883 TClonesArray * clusterList = 0x0;
884 if(fOutputEvent) clusterList = dynamic_cast<TClonesArray*> (fOutputEvent->FindListObject(fEMCALClustersListName));
886 //printf("AliCaloTrackReader::FillInputEMCAL() - Wrong name of list with clusters? Try input event <%s>\n",fEMCALClustersListName.Data());
887 //List not in output event, try input event
888 clusterList = dynamic_cast<TClonesArray*> (fInputEvent->FindListObject(fEMCALClustersListName));
890 printf("AliCaloTrackReader::FillInputEMCAL() - Wrong name of list with clusters? <%s>\n",fEMCALClustersListName.Data());
895 Int_t nclusters = clusterList->GetEntriesFast();
896 for (Int_t iclus = 0; iclus < nclusters; iclus++) {
897 AliVCluster * clus = dynamic_cast<AliVCluster*> (clusterList->At(iclus));
898 //printf("E %f\n",clus->E());
899 if (clus) FillInputEMCALAlgorithm(clus, iclus);
900 else printf("AliCaloTrackReader::FillInputEMCAL() - Null cluster in list!\n");
904 // Recalculate track matching, not necessary if already done in the reclusterization task.
905 // in case it was not done ...
906 GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent,clusterList);
910 if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - aod entries %d\n", fEMCALClusters->GetEntriesFast());
914 //______________________________________
915 void AliCaloTrackReader::FillInputPHOS()
917 //Return array with PHOS clusters in aod format
919 if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputPHOS()\n");
921 //Loop to select clusters in fiducial cut and fill container with aodClusters
922 Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
923 for (Int_t iclus = 0; iclus < nclusters; iclus++) {
924 AliVCluster * clus = 0;
925 if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
926 if (IsPHOSCluster(clus)){
927 //Check if the cluster contains any bad channel and if close to calorimeter borders
930 vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
931 if( GetCaloUtils()->ClusterContainsBadChannel("PHOS",clus->GetCellsAbsId(), clus->GetNCells()))
933 if(!GetCaloUtils()->CheckCellFiducialRegion(clus, fInputEvent->GetPHOSCells(), fInputEvent, vindex))
936 if(fRecalculateClusters){
938 //Recalibrate the cluster energy
939 if(GetCaloUtils()->IsRecalibrationOn()) {
940 Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, (AliAODCaloCells*)GetPHOSCells());
946 TLorentzVector momentum ;
948 clus->GetMomentum(momentum, fVertex[vindex]);
950 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS")) continue;
952 if(fPHOSPtMin > momentum.E() || fPHOSPtMax < momentum.E()) continue;
954 if(fDebug > 2 && momentum.E() > 0.1)
955 printf("AliCaloTrackReader::FillInputPHOS() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
956 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
963 fPHOSClusters->Add(clus);
969 if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS() - aod entries %d\n", fPHOSClusters->GetEntriesFast());
973 //____________________________________________
974 void AliCaloTrackReader::FillInputEMCALCells()
976 //Return array with EMCAL cells in aod format
978 fEMCALCells = fInputEvent->GetEMCALCells();
982 //___________________________________________
983 void AliCaloTrackReader::FillInputPHOSCells()
985 //Return array with PHOS cells in aod format
987 fPHOSCells = fInputEvent->GetPHOSCells();
991 //_______________________________________
992 void AliCaloTrackReader::FillInputVZERO()
994 //Fill VZERO information in data member, add all the channels information.
995 AliVVZERO* v0 = fInputEvent->GetVZEROData();
996 //printf("Init V0: ADC (%d,%d), Multiplicity (%d,%d) \n",fV0ADC[0],fV0ADC[1],fV0Mul[0],fV0Mul[1]);
1000 AliESDVZERO* esdV0 = dynamic_cast<AliESDVZERO*> (v0);
1001 for (Int_t i = 0; i < 32; i++)
1003 if(esdV0){//Only available in ESDs
1004 fV0ADC[0] += (Int_t)esdV0->GetAdcV0C(i);
1005 fV0ADC[1] += (Int_t)esdV0->GetAdcV0A(i);
1007 fV0Mul[0] += (Int_t)v0->GetMultiplicityV0C(i);
1008 fV0Mul[1] += (Int_t)v0->GetMultiplicityV0A(i);
1011 printf("V0: ADC (%d,%d), Multiplicity (%d,%d) \n",fV0ADC[0],fV0ADC[1],fV0Mul[0],fV0Mul[1]);
1016 printf("Cannot retrieve V0 ESD! Run w/ null V0 charges\n ");
1021 //___________________________________________________________________
1022 Bool_t AliCaloTrackReader::IsEMCALCluster(AliVCluster* cluster) const
1024 // Check if it is a cluster from EMCAL. For old AODs cluster type has
1025 // different number and need to patch here
1027 if(fDataType==kAOD && fOldAOD)
1029 if (cluster->GetType() == 2) return kTRUE;
1034 return cluster->IsEMCAL();
1039 //___________________________________________________________________
1040 Bool_t AliCaloTrackReader::IsPHOSCluster(AliVCluster * cluster) const
1042 //Check if it is a cluster from PHOS.For old AODs cluster type has
1043 // different number and need to patch here
1045 if(fDataType==kAOD && fOldAOD)
1047 Int_t type = cluster->GetType();
1048 if (type == 0 || type == 1) return kTRUE;
1053 return cluster->IsPHOS();
1058 //________________________________________________
1059 Bool_t AliCaloTrackReader::CheckForPrimaryVertex()
1061 //Check if the vertex was well reconstructed, copy from V0Reader of conversion group
1064 AliESDEvent * event = dynamic_cast<AliESDEvent*> (fInputEvent);
1065 if(!event) return kTRUE;
1067 if(event->GetPrimaryVertexTracks()->GetNContributors() > 0) {
1071 if(event->GetPrimaryVertexTracks()->GetNContributors() < 1) {
1073 if(event->GetPrimaryVertexSPD()->GetNContributors() > 0) {
1074 //cout<<"spd vertex type::"<< fESDEvent->GetPrimaryVertex()->GetName() << endl;
1078 if(event->GetPrimaryVertexSPD()->GetNContributors() < 1) {
1079 // cout<<"bad vertex type::"<< event->GetPrimaryVertex()->GetName() << endl;
1088 //____________________________________________________________
1089 void AliCaloTrackReader::SetTrackCuts(AliESDtrackCuts * cuts)
1093 if(fESDtrackCuts) delete fESDtrackCuts ;
1095 fESDtrackCuts = cuts ;