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 // : AliCaloTrackReader: Fills AOD data in 3 TObjArrays (PHOS, EMCAL, CTS)
25 // This part is commented: Mixing analysis can be done, input AOD with events
26 // is opened in the AliCaloTrackReader::Init()
28 //-- Author: Gustavo Conesa (LNF-INFN)
29 //////////////////////////////////////////////////////////////////////////////
32 // --- ROOT system ---
36 //---- ANALYSIS system ----
37 #include "AliMCEvent.h"
38 #include "AliAODMCHeader.h"
39 #include "AliGenPythiaEventHeader.h"
40 #include "AliESDEvent.h"
41 #include "AliAODEvent.h"
42 #include "AliVTrack.h"
43 #include "AliVParticle.h"
44 #include "AliMixedEvent.h"
45 #include "AliESDtrack.h"
46 #include "AliESDtrackCuts.h"
47 #include "AliTriggerAnalysis.h"
48 #include "AliESDVZERO.h";
50 //---- PartCorr/EMCAL ---
51 #include "AliEMCALRecoUtils.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), fComparePtHardAndJetPt(kFALSE), fPtHardAndJetPtFactor(7),
62 fCTSPtMin(0), fEMCALPtMin(0),fPHOSPtMin(0), fAODBranchList(new TList ),
63 fCTSTracks(new TObjArray()), fEMCALClusters(new TObjArray()), fPHOSClusters(new TObjArray()),
64 fEMCALCells(0x0), fPHOSCells(0x0),
65 fInputEvent(0x0), fOutputEvent(0x0),fMC(0x0),
66 fFillCTS(0),fFillEMCAL(0),fFillPHOS(0),
67 fFillEMCALCells(0),fFillPHOSCells(0),
68 fRemoveSuspiciousClusters(kFALSE), fSmearClusterEnergy(kFALSE), fRandom(),
69 // fSecondInputAODTree(0x0), fSecondInputAODEvent(0x0),
70 // fSecondInputFileName(""),fSecondInputFirstEvent(0),
71 // fCTSTracksNormalInputEntries(0), fEMCALClustersNormalInputEntries(0),
72 // fPHOSClustersNormalInputEntries(0),
73 fTrackStatus(0), fESDtrackCuts(0), fTrackMult(0), fTrackMultEtaCut(0.8),
74 fReadStack(kFALSE), fReadAODMCParticles(kFALSE),
75 fDeltaAODFileName("deltaAODPartCorr.root"),fFiredTriggerClassName(""),
76 fAnaLED(kFALSE),fTaskName(""),fCaloUtils(0x0),
77 fMixedEvent(NULL), fNMixedEvent(1), fVertex(NULL),
78 fWriteOutputDeltaAOD(kFALSE),fOldAOD(kFALSE),fCaloFilterPatch(kFALSE),
79 fEMCALClustersListName(""),fZvtxCut(0.),
80 fDoEventSelection(kFALSE), fDoV0ANDEventSelection(kFALSE),
81 fTriggerAnalysis (new AliTriggerAnalysis), fCentralityClass("V0M"),fCentralityOpt(10)
86 //Initialize parameters
90 //_________________________________
91 AliCaloTrackReader::~AliCaloTrackReader() {
94 if(fFiducialCut) delete fFiducialCut ;
97 fAODBranchList->Delete();
98 delete fAODBranchList ;
102 if(fDataType!=kMC)fCTSTracks->Clear() ;
103 else fCTSTracks->Delete() ;
108 if(fDataType!=kMC)fEMCALClusters->Clear("C") ;
109 else fEMCALClusters->Delete() ;
110 delete fEMCALClusters ;
114 if(fDataType!=kMC)fPHOSClusters->Clear("C") ;
115 else fPHOSClusters->Delete() ;
116 delete fPHOSClusters ;
120 // delete fEMCALCells ;
124 // delete fPHOSCells ;
128 for (Int_t i = 0; i < fNMixedEvent; i++) {
129 delete [] fVertex[i] ;
135 if(fESDtrackCuts) delete fESDtrackCuts;
136 if(fTriggerAnalysis) delete fTriggerAnalysis;
138 // Pointers not owned, done by the analysis frame
139 // if(fInputEvent) delete fInputEvent ;
140 // if(fOutputEvent) delete fOutputEvent ;
141 // if(fMC) delete fMC ;
143 // if(fSecondInputAODTree){
144 // fSecondInputAODTree->Clear();
145 // delete fSecondInputAODTree;
148 // if(fSecondInputAODEvent) delete fSecondInputAODEvent ;
150 // Pointer not owned, deleted by maker
151 //if (fCaloUtils) delete fCaloUtils ;
155 //_________________________________________________________________________
156 Bool_t AliCaloTrackReader::ComparePtHardAndJetPt(){
157 // Check the event, if the requested ptHard is much larger than the jet pT, then there is a problem.
159 if(!fReadStack) return kTRUE; //Information not filtered to AOD
161 if(!strcmp(GetGenEventHeader()->ClassName(), "AliGenPythiaEventHeader")){
163 AliGenPythiaEventHeader* pygeh= (AliGenPythiaEventHeader*) GetGenEventHeader();
164 Int_t nTriggerJets = pygeh->NTriggerJets();
165 Float_t ptHard = pygeh->GetPtHard();
167 //if(fDebug > 1) printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %f\n",nTriggerJets, ptHard);
168 Float_t tmpjet[]={0,0,0,0};
169 for(Int_t ijet = 0; ijet< nTriggerJets; ijet++){
170 pygeh->TriggerJet(ijet, tmpjet);
171 jet = new TParticle(94, 21, -1, -1, -1, -1, tmpjet[0],tmpjet[1],tmpjet[2],tmpjet[3], 0,0,0,0);
172 //Compare jet pT and pt Hard
173 //if(fDebug > 1) printf("AliMCAnalysisUtils:: %d pycell jet pT %f\n",ijet, jet->Pt());
174 if(jet->Pt() > fPtHardAndJetPtFactor * ptHard) {
175 printf("AliMCAnalysisUtils::PythiaEventHeader: Njets: %d, pT Hard %2.2f, pycell jet pT %2.2f, rejection factor %1.1f\n",
176 nTriggerJets, ptHard, jet->Pt(), fPtHardAndJetPtFactor);
187 //____________________________________________________________________________
188 AliStack* AliCaloTrackReader::GetStack() const {
189 //Return pointer to stack
193 if(fDebug > 1) printf("AliCaloTrackReader::GetStack() - Stack is not available\n");
198 //____________________________________________________________________________
199 AliHeader* AliCaloTrackReader::GetHeader() const {
200 //Return pointer to header
202 return fMC->Header();
204 printf("AliCaloTrackReader::Header is not available\n");
208 //____________________________________________________________________________
209 AliGenEventHeader* AliCaloTrackReader::GetGenEventHeader() const {
210 //Return pointer to Generated event header
212 return fMC->GenEventHeader();
214 printf("AliCaloTrackReader::GenEventHeader is not available\n");
219 //____________________________________________________________________________
220 TClonesArray* AliCaloTrackReader::GetAODMCParticles(Int_t input) const {
221 //Return list of particles in AOD. Do it for the corresponding input event.
223 TClonesArray * rv = NULL ;
224 if(fDataType == kAOD){
228 AliAODEvent * evt = dynamic_cast<AliAODEvent*> (fInputEvent) ;
230 rv = (TClonesArray*)evt->FindListObject("mcparticles");
232 printf("AliCaloTrackReader::GetAODMCParticles() - wrong AOD input index? %d, or non existing tree? \n",input);
234 } //else if(input == 1 && fSecondInputAODEvent){ //Second input AOD
235 // rv = (TClonesArray*) fSecondInputAODEvent->FindListObject("mcparticles");
239 printf("AliCaloTrackReader::GetAODMCParticles() - Input are not AODs\n");
245 //____________________________________________________________________________
246 AliAODMCHeader* AliCaloTrackReader::GetAODMCHeader(Int_t input) const {
247 //Return MC header in AOD. Do it for the corresponding input event.
248 AliAODMCHeader *mch = NULL;
249 if(fDataType == kAOD){
252 mch = (AliAODMCHeader*)((AliAODEvent*)fInputEvent)->FindListObject("mcheader");
254 // //Second input AOD
255 // else if(input == 1){
256 // mch = (AliAODMCHeader*) fSecondInputAODEvent->FindListObject("mcheader");
259 printf("AliCaloTrackReader::GetAODMCHeader() - wrong AOD input index, %d\n",input);
263 printf("AliCaloTrackReader::GetAODMCHeader() - Input are not AODs\n");
269 //_______________________________________________________________
270 void AliCaloTrackReader::Init()
272 //Init reader. Method to be called in AliAnaPartCorrMaker
274 //Get the file with second input events if the filename is given
275 //Get the tree and connect the AODEvent. Only with AODs
277 if(fReadStack && fReadAODMCParticles){
278 printf("AliCaloTrackReader::Init() - Cannot access stack and mcparticles at the same time, change them \n");
280 fReadAODMCParticles = kFALSE;
283 // if(fSecondInputFileName!=""){
284 // if(fDataType == kAOD){
285 // TFile * input2 = new TFile(fSecondInputFileName,"read");
286 // printf("AliCaloTrackReader::Init() - Second input file opened: %s, size %d \n", input2->GetName(), (Int_t) input2->GetSize());
287 // fSecondInputAODTree = (TTree*) input2->Get("aodTree");
288 // if(fSecondInputAODTree) printf("AliCaloTrackReader::Init() - Second input tree opened: %s, entries %d \n",
289 // fSecondInputAODTree->GetName(), (Int_t) fSecondInputAODTree->GetEntries());
291 // printf("AliCaloTrackReader::Init() - Second input tree not available, STOP \n");
294 // fSecondInputAODEvent = new AliAODEvent;
295 // fSecondInputAODEvent->ReadFromTree(fSecondInputAODTree);
296 // if(fSecondInputFirstEvent >= fSecondInputAODTree->GetEntriesFast()){
297 // printf("AliCaloTrackReader::Init() - Requested first event of second input %d, is larger than number of events %d, STOP\n",
298 // fSecondInputFirstEvent, (Int_t) fSecondInputAODTree->GetEntriesFast());
302 // else printf("AliCaloTrackReader::Init() - Second input not added, reader is not AOD\n");
306 //_______________________________________________________________
307 void AliCaloTrackReader::InitParameters()
309 //Initialize the parameters of the analysis.
315 //Do not filter the detectors input by default.
319 fFillEMCALCells = kFALSE;
320 fFillPHOSCells = kFALSE;
322 //fSecondInputFileName = "" ;
323 //fSecondInputFirstEvent = 0 ;
324 fReadStack = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
325 fReadAODMCParticles = kFALSE; // Check in the constructor of the other readers if it was set or in the configuration file
326 fDeltaAODFileName = "deltaAODPartCorr.root";
327 fFiredTriggerClassName = "";
331 //We want tracks fitted in the detectors:
332 //fTrackStatus=AliESDtrack::kTPCrefit;
333 //fTrackStatus|=AliESDtrack::kITSrefit;
335 fESDtrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010();
337 fV0ADC[0] = 0; fV0ADC[1] = 0;
338 fV0Mul[0] = 0; fV0Mul[1] = 0;
343 fCentralityBin[0]=fCentralityBin[1]=-1;
346 fSmearClusterEnergy = kFALSE;
347 fSmearClusterParam[0] = 0.07; // * sqrt E term
348 fSmearClusterParam[1] = 0.02; // * E term
349 fSmearClusterParam[2] = 0.00; // constant
353 //________________________________________________________________
354 void AliCaloTrackReader::Print(const Option_t * opt) const
357 //Print some relevant parameters set for the analysis
361 printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
362 printf("Task name : %s\n", fTaskName.Data()) ;
363 printf("Data type : %d\n", fDataType) ;
364 printf("CTS Min pT : %2.1f GeV/c\n", fCTSPtMin) ;
365 printf("EMCAL Min pT : %2.1f GeV/c\n", fEMCALPtMin) ;
366 printf("PHOS Min pT : %2.1f GeV/c\n", fPHOSPtMin) ;
367 printf("Use CTS = %d\n", fFillCTS) ;
368 printf("Use EMCAL = %d\n", fFillEMCAL) ;
369 printf("Use PHOS = %d\n", fFillPHOS) ;
370 printf("Use EMCAL Cells = %d\n", fFillEMCALCells) ;
371 printf("Use PHOS Cells = %d\n", fFillPHOSCells) ;
372 printf("Track status = %d\n", (Int_t) fTrackStatus) ;
373 printf("Track Mult Eta Cut = %d\n", (Int_t) fTrackMultEtaCut) ;
374 printf("Write delta AOD = %d\n", fWriteOutputDeltaAOD) ;
376 if(fComparePtHardAndJetPt)
377 printf("Compare jet pt and pt hard to accept event, factor = %2.2f",fPtHardAndJetPtFactor);
379 // if(fSecondInputFileName!="") {
380 // printf("Second Input File Name = %s\n", fSecondInputFileName.Data()) ;
381 // printf("Second Input First Event = %d\n", fSecondInputFirstEvent) ;
384 printf("Read Kine from, stack? %d, AOD ? %d \n", fReadStack, fReadAODMCParticles) ;
385 printf("Delta AOD File Name = %s\n", fDeltaAODFileName.Data()) ;
387 printf("Centrality: Class %s, Option %d, Bin [%d,%d] \n", fCentralityClass.Data(),fCentralityOpt,fCentralityBin[0], fCentralityBin[1]) ;
393 //___________________________________________________
394 Bool_t AliCaloTrackReader::FillInputEvent(const Int_t iEntry, const char * /*currentFileName*/) {
395 //Fill the event counter and input lists that are needed, called by the analysis maker.
397 fEventNumber = iEntry;
398 //fCurrentFileName = TString(currentFileName);
400 if(fDebug >= 0) printf("AliCaloTrackReader::FillInputEvent() - Input event not available, skip event analysis\n");
403 //Select events only fired by a certain trigger configuration if it is provided
405 if(fInputEvent->GetHeader())
406 eventType = ((AliVHeader*)fInputEvent->GetHeader())->GetEventType();
407 if( fFiredTriggerClassName !="" && !fAnaLED){
409 return kFALSE; //Only physics event, do not use for simulated events!!!
411 printf("AliCaloTrackReader::FillInputEvent() - FiredTriggerClass <%s>, selected class <%s>, compare name %d\n",
412 GetFiredTriggerClasses().Data(),fFiredTriggerClassName.Data(), GetFiredTriggerClasses().Contains(fFiredTriggerClassName));
413 if( !GetFiredTriggerClasses().Contains(fFiredTriggerClassName) ) return kFALSE;
416 // kStartOfRun = 1, // START_OF_RUN
417 // kEndOfRun = 2, // END_OF_RUN
418 // kStartOfRunFiles = 3, // START_OF_RUN_FILES
419 // kEndOfRunFiles = 4, // END_OF_RUN_FILES
420 // kStartOfBurst = 5, // START_OF_BURST
421 // kEndOfBurst = 6, // END_OF_BURST
422 // kPhysicsEvent = 7, // PHYSICS_EVENT
423 // kCalibrationEvent = 8, // CALIBRATION_EVENT
424 // kFormatError = 9, // EVENT_FORMAT_ERROR
425 // kStartOfData = 10, // START_OF_DATA
426 // kEndOfData = 11, // END_OF_DATA
427 // kSystemSoftwareTriggerEvent = 12, // SYSTEM_SOFTWARE_TRIGGER_EVENT
428 // kDetectorSoftwareTriggerEvent = 13 // DETECTOR_SOFTWARE_TRIGGER_EVENT
430 if(eventType!=7 && fDebug > 1 )printf("AliCaloTrackReader::FillInputEvent() - DO LED, Event Type <%d>, 8 Calibration \n", eventType);
431 if(eventType!=8)return kFALSE;
434 //In case of analysis of events with jets, skip those with jet pt > 5 pt hard
435 if(fComparePtHardAndJetPt && GetStack()) {
436 if(!ComparePtHardAndJetPt()) return kFALSE ;
439 //In case of mixing events with other AOD file
440 // if(fDataType == kAOD && fSecondInputAODTree){
443 // printf("AliCaloTrackReader::FillInputEvent() - Get event %d from second input AOD file \n", iEntry+fSecondInputFirstEvent);
444 // if(fSecondInputAODTree->GetEntriesFast() <= iEntry+fSecondInputFirstEvent) {
445 // if(fSecondInputAODTree->GetEntriesFast() == iEntry+fSecondInputFirstEvent)
446 // printf("AliCaloTrackReader::FillInputEvent() - Skip events from event %d, no more events in second AOD file \n", iEntry);
451 // Int_t nbytes = fSecondInputAODTree->GetEvent(iEntry+fSecondInputFirstEvent);
452 // if ( nbytes == 0 ) {//If nothing in AOD
453 // printf("AliCaloTrackReader::FillInputEvent() - Nothing in Second AOD input, STOP\n");
461 //Reject events with Z vertex too large, only for SE analysis, if not, cut on the analysis code
462 if(!GetMixedEvent() && TMath::Abs(fVertex[0][2]) > fZvtxCut) return kFALSE;
464 //------------------------------------------------------
465 //Event rejection depending on vertex, pileup, v0and
466 //------------------------------------------------------
467 if(fDoEventSelection){
468 if(!fCaloFilterPatch){
469 //Do not analyze events with pileup
470 Bool_t bPileup = fInputEvent->IsPileupFromSPD(3, 0.8, 3., 2., 5.); //Default values, if not it does not compile
471 //Bool_t bPileup = event->IsPileupFromSPD();
472 if(bPileup) return kFALSE;
474 if(fDoV0ANDEventSelection){
475 Bool_t bV0AND = kTRUE;
476 AliESDEvent* esd = dynamic_cast<AliESDEvent*> (fInputEvent);
478 bV0AND = fTriggerAnalysis->IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0AND);
479 //else bV0AND = //FIXME FOR AODs
480 if(!bV0AND) return kFALSE;
483 if(!CheckForPrimaryVertex()) return kFALSE;
486 if(fInputEvent->GetNumberOfCaloClusters() > 0) {
487 AliVCluster * calo = fInputEvent->GetCaloCluster(0);
488 if(calo->GetNLabels() == 4){
489 Int_t * selection = calo->GetLabels();
490 Bool_t bPileup = selection[0];
491 if(bPileup) return kFALSE;
493 Bool_t bGoodV = selection[1];
494 if(!bGoodV) return kFALSE;
496 if(fDoV0ANDEventSelection){
497 Bool_t bV0AND = selection[2];
498 if(!bV0AND) return kFALSE;
501 fTrackMult = selection[3];
502 if(fTrackMult == 0) return kFALSE;
504 //First filtered AODs, track multiplicity stored there.
505 fTrackMult = (Int_t) ((AliAODHeader*)fInputEvent->GetHeader())->GetCentrality();
506 if(fTrackMult == 0) return kFALSE;
508 }//at least one cluster
510 //printf("AliCaloTrackReader::FillInputEvent() - No clusters in event\n");
511 //Remove events with vertex (0,0,0), bad vertex reconstruction
512 if(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;
514 //First filtered AODs, track multiplicity stored there.
515 fTrackMult = (Int_t) ((AliAODHeader*)fInputEvent->GetHeader())->GetCentrality();
516 if(fTrackMult == 0) return kFALSE;
518 }// CaloFileter patch
520 //------------------------------------------------------
522 //Check if there is a centrality value, PbPb analysis, and if a centrality bin selection is requested
523 //If we need a centrality bin, we select only those events in the corresponding bin.
524 if(GetCentrality() && fCentralityBin[0]>=0 && fCentralityBin[1]>=0 && fCentralityOpt==100){
525 Int_t cen = GetEventCentrality();
526 if(cen > fCentralityBin[1] || cen < fCentralityBin[0]) return kFALSE; //reject events out of bin.
529 //Fill the arrays with cluster/tracks/cells data
531 FillInputEMCALCells();
533 FillInputPHOSCells();
537 //Accept events with at least one track
538 if(fTrackMult == 0 && fDoEventSelection) return kFALSE;
551 //__________________________________________________
552 void AliCaloTrackReader::ResetLists() {
553 // Reset lists, called by the analysis maker
555 if(fCTSTracks) fCTSTracks -> Clear();
556 if(fEMCALClusters) fEMCALClusters -> Clear("C");
557 if(fPHOSClusters) fPHOSClusters -> Clear("C");
558 // if(fEMCALCells) fEMCALCells -> Clear("");
559 // if(fPHOSCells) fPHOSCells -> Clear("");
561 fV0ADC[0] = 0; fV0ADC[1] = 0;
562 fV0Mul[0] = 0; fV0Mul[1] = 0;
566 //____________________________________________________________________________
567 void AliCaloTrackReader::SetInputEvent(AliVEvent* const input)
570 fMixedEvent = dynamic_cast<AliMixedEvent*>(GetInputEvent()) ;
572 fNMixedEvent = fMixedEvent->GetNumberOfEvents() ;
575 //Delete previous vertex
577 for (Int_t i = 0; i < fNMixedEvent; i++) {
578 delete [] fVertex[i] ;
583 fVertex = new Double_t*[fNMixedEvent] ;
584 for (Int_t i = 0; i < fNMixedEvent; i++) {
585 fVertex[i] = new Double_t[3] ;
586 fVertex[i][0] = 0.0 ;
587 fVertex[i][1] = 0.0 ;
588 fVertex[i][2] = 0.0 ;
592 //__________________________________________________
593 Int_t AliCaloTrackReader::GetEventCentrality() const {
594 //Return current event centrality
597 if(fCentralityOpt==100) return (Int_t) GetCentrality()->GetCentralityPercentile(fCentralityClass);
598 else if(fCentralityOpt==10) return GetCentrality()->GetCentralityClass10(fCentralityClass);
599 else if(fCentralityOpt==5) return GetCentrality()->GetCentralityClass5(fCentralityClass);
601 printf("AliAnaPartCorrBaseClass::Unknown centrality option %d, use 5, 10 or 100\n",fCentralityOpt);
609 //____________________________________________________________________________
610 void AliCaloTrackReader::GetVertex(Double_t vertex[3]) const {
611 //Return vertex position to be used for single event analysis
612 vertex[0]=fVertex[0][0];
613 vertex[1]=fVertex[0][1];
614 vertex[2]=fVertex[0][2];
617 //____________________________________________________________________________
618 void AliCaloTrackReader::GetVertex(Double_t vertex[3], const Int_t evtIndex) const {
619 //Return vertex position for mixed event, recover the vertex in a particular event.
621 //Int_t evtIndex = 0; // for single events only one vertex stored in position 0, default value
622 //if (fMixedEvent && clusterID >=0) {
623 // evtIndex=GetMixedEvent()->EventIndexForCaloCluster(clusterID) ;
626 vertex[0]=fVertex[evtIndex][0]; vertex[1]=fVertex[evtIndex][1]; vertex[2]=fVertex[evtIndex][2];
632 //____________________________________________________________________________
633 void AliCaloTrackReader::FillVertexArray() {
635 //Fill data member with vertex
636 //In case of Mixed event, multiple vertices
638 //Delete previous vertex
640 for (Int_t i = 0; i < fNMixedEvent; i++) {
641 delete [] fVertex[i] ;
646 fVertex = new Double_t*[fNMixedEvent] ;
647 for (Int_t i = 0; i < fNMixedEvent; i++) {
648 fVertex[i] = new Double_t[3] ;
649 fVertex[i][0] = 0.0 ;
650 fVertex[i][1] = 0.0 ;
651 fVertex[i][2] = 0.0 ;
654 if (!fMixedEvent) { //Single event analysis
657 if(fInputEvent->GetPrimaryVertex()){
658 fInputEvent->GetPrimaryVertex()->GetXYZ(fVertex[0]);
661 printf("AliCaloTrackReader::FillVertexArray() - NULL primary vertex\n");
662 fVertex[0][0]=0.; fVertex[0][1]=0.; fVertex[0][2]=0.;
663 }//Primary vertex pointer do not exist
665 } else {//MC read event
666 fVertex[0][0]=0.; fVertex[0][1]=0.; fVertex[0][2]=0.;
670 printf("AliCaloTrackReader::FillVertexArray() - Single Event Vertex : %f,%f,%f\n",fVertex[0][0],fVertex[0][1],fVertex[0][2]);
672 } else { // MultiEvent analysis
673 for (Int_t iev = 0; iev < fNMixedEvent; iev++) {
674 if (fMixedEvent->GetVertexOfEvent(iev))
675 fMixedEvent->GetVertexOfEvent(iev)->GetXYZ(fVertex[iev]);
676 else { // no vertex found !!!!
677 AliWarning("No vertex found");
681 printf("AliCaloTrackReader::FillVertexArray() - Multi Event %d Vertex : %f,%f,%f\n",iev,fVertex[iev][0],fVertex[iev][1],fVertex[iev][2]);
688 //____________________________________________________________________________
689 void AliCaloTrackReader::FillInputCTS() {
690 //Return array with Central Tracking System (CTS) tracks
692 if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputCTS()\n");
694 Int_t nTracks = fInputEvent->GetNumberOfTracks() ;
698 for (Int_t itrack = 0; itrack < nTracks; itrack++) {////////////// track loop
699 AliVTrack * track = (AliVTrack*)fInputEvent->GetTrack(itrack) ; // retrieve track from esd
701 //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
702 if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus))
707 if(fDataType==kESD && !fESDtrackCuts->AcceptTrack((AliESDtrack*)track)) continue;
709 // Track filter selection
710 //if (fTrackFilter) {
711 // selectInfo = fTrackFilter->IsSelected(esdTrack);
712 // if (!selectInfo && !(esd->GetPrimaryVertex())->UsesTrack(esdTrack->GetID())) continue;
715 //Count the tracks in eta < 0.9
716 //printf("Eta %f cut %f\n",TMath::Abs(track->Eta()),fTrackMultEtaCut);
717 if(TMath::Abs(track->Eta())< fTrackMultEtaCut) fTrackMult++;
719 track->GetPxPyPz(p) ;
720 TLorentzVector momentum(p[0],p[1],p[2],0);
722 if(fCTSPtMin < momentum.Pt()){
724 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS"))
727 if(fDebug > 2 && momentum.Pt() > 0.1)
728 printf("AliCaloTrackReader::FillInputCTS() - Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
729 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
732 track->SetID(itrack);
735 fCTSTracks->Add(track);
737 }//Pt and Fiducial cut passed.
740 //fCTSTracksNormalInputEntries = fCTSTracks->GetEntriesFast();
742 printf("AliCaloTrackReader::FillInputCTS() - aod entries %d, input tracks %d, pass status %d, multipliticy %d\n", fCTSTracks->GetEntriesFast(), nTracks, nstatus, fTrackMult);//fCTSTracksNormalInputEntries);
744 // //If second input event available, add the clusters.
745 // if(fSecondInputAODTree && fSecondInputAODEvent){
746 // nTracks = fSecondInputAODEvent->GetNumberOfTracks() ;
747 // if(fDebug > 1) printf("AliCaloTrackReader::FillInputCTS() - Add second input tracks, entries %d\n", nTracks) ;
748 // for (Int_t itrack = 0; itrack < nTracks; itrack++) {////////////// track loop
749 // AliAODTrack * track = ((AliAODEvent*)fSecondInputAODEvent)->GetTrack(itrack) ; // retrieve track from esd
751 // //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
752 // if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus)) continue ;
754 // track->GetPxPyPz(p) ;
755 // TLorentzVector momentum(p[0],p[1],p[2],0);
757 // if(fCTSPtMin < momentum.Pt()){
759 // if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) continue;
761 // if(fDebug > 2 && momentum.Pt() > 0.1) printf("AliCaloTrackReader::FillInputCTS() - Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
762 // momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
764 // fCTSTracks->Add(track);
766 // }//Pt and Fiducial cut passed.
769 // if(fDebug > 1) printf("AliCaloTrackReader::FillInputCTS() - aod normal entries %d, after second input %d\n", fCTSTracksNormalInputEntries, fCTSTracks->GetEntriesFast());
770 // } //second input loop
774 //____________________________________________________________________________
775 void AliCaloTrackReader::FillInputEMCALAlgorithm(AliVCluster * clus, const Int_t iclus) {
776 //Fill the EMCAL data in the array, do it
780 vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
782 //Check if the cluster contains any bad channel and if close to calorimeter borders
783 if(GetCaloUtils()->ClusterContainsBadChannel("EMCAL",clus->GetCellsAbsId(), clus->GetNCells()))
785 if(!GetCaloUtils()->CheckCellFiducialRegion(clus, (AliVCaloCells*)fInputEvent->GetEMCALCells(), fInputEvent, vindex))
788 //Remove suspicious clusters
789 if(fRemoveSuspiciousClusters){
790 Int_t ncells = clus->GetNCells();
791 Float_t energy = clus->E();
792 Float_t minNCells = 1+energy/3;//-x*x*0.0033
793 if(ncells < minNCells) {
794 //if(energy > 2)printf("AliCaloTrackReader::FillInputEMCALAlgorithm() - Remove cluster: e %2.2f, Ncells %d, min Ncells %2.1f\n",energy,ncells,minNCells);
798 // if(energy > 2)printf("AliCaloTrackReader::FillInputEMCALAlgorithm() - Keep cluster: e %2.2f, Ncells %d, min Ncells %2.1f\n",energy,ncells,minNCells);
803 TLorentzVector momentum ;
805 clus->GetMomentum(momentum, fVertex[vindex]);
807 if(fEMCALPtMin < momentum.Pt()){
809 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL"))
812 if(fDebug > 2 && momentum.E() > 0.1)
813 printf("AliCaloTrackReader::FillInputEMCAL() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
814 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
817 //clus->GetPosition(pos);
818 //printf("Before Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
820 //Recalibrate the cluster energy
821 if(GetCaloUtils()->IsRecalibrationOn()) {
822 Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, GetEMCALCells());
824 //printf("Recalibrated Energy %f\n",clus->E());
825 GetCaloUtils()->RecalculateClusterShowerShapeParameters(GetEMCALCells(),clus);
826 GetCaloUtils()->RecalculateClusterPID(clus);
830 //Recalculate distance to bad channels, if new list of bad channels provided
831 GetCaloUtils()->RecalculateClusterDistanceToBadChannel(GetEMCALCells(),clus);
833 //Recalculate cluster position
834 if(GetCaloUtils()->IsRecalculationOfClusterPositionOn()){
835 GetCaloUtils()->RecalculateClusterPosition(GetEMCALCells(),clus);
836 //clus->GetPosition(pos);
837 //printf("After Corrections: e %f, x %f, y %f, z %f\n",clus->E(),pos[0],pos[1],pos[2]);
840 //Correct non linearity
841 if(GetCaloUtils()->IsCorrectionOfClusterEnergyOn()){
842 GetCaloUtils()->CorrectClusterEnergy(clus) ;
843 //printf("Linearity Corrected Energy %f\n",clus->E());
846 //In case of MC analysis, to match resolution/calibration in real data
847 if(fSmearClusterEnergy){
848 Float_t energy = clus->E();
849 Float_t rdmEnergy = fRandom.Gaus(energy,fSmearClusterParam[0]*TMath::Sqrt(energy)+
850 fSmearClusterParam[1]*energy+fSmearClusterParam[2]);
851 clus->SetE(rdmEnergy);
852 if(fDebug > 2) printf("\t Energy %f, smeared %f\n", energy, clus->E());
858 fEMCALClusters->Add(clus);
862 //____________________________________________________________________________
863 void AliCaloTrackReader::FillInputEMCAL() {
864 //Return array with EMCAL clusters in aod format
866 if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputEMCAL()\n");
868 //Loop to select clusters in fiducial cut and fill container with aodClusters
869 if(fEMCALClustersListName==""){
870 Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
871 for (Int_t iclus = 0; iclus < nclusters; iclus++) {
872 AliVCluster * clus = 0;
873 if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
874 if (IsEMCALCluster(clus)){
875 FillInputEMCALAlgorithm(clus, iclus);
880 //Recalculate track matching
881 if(fDataType==kESD)GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent);
883 }//Get the clusters from the input event
885 TClonesArray * clusterList = dynamic_cast<TClonesArray*> (fOutputEvent->FindListObject(fEMCALClustersListName));
887 printf("AliCaloTrackReader::FillInputEMCAL() - Wrong name of list with clusters? <%s>\n",fEMCALClustersListName.Data());
890 Int_t nclusters = clusterList->GetEntriesFast();
891 for (Int_t iclus = 0; iclus < nclusters; iclus++) {
892 AliVCluster * clus = dynamic_cast<AliVCluster*> (clusterList->At(iclus));
893 //printf("E %f\n",clus->E());
894 if (clus) FillInputEMCALAlgorithm(clus, iclus);
895 else printf("AliCaloTrackReader::FillInputEMCAL() - Null cluster in list!\n");
899 //Recalculate track matching, not necessary, already done in the reclusterization task
900 //GetCaloUtils()->RecalculateClusterTrackMatching(fInputEvent,clusterList);
904 //fEMCALClustersNormalInputEntries = fEMCALClusters->GetEntriesFast();
905 if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - aod entries %d\n", fEMCALClusters->GetEntriesFast());//fEMCALClustersNormalInputEntries);
907 //If second input event available, add the clusters.
908 // if(fSecondInputAODTree && fSecondInputAODEvent){
909 // GetSecondInputAODVertex(v);
910 // nclusters = ((AliAODEvent*)fSecondInputAODEvent)->GetNumberOfCaloClusters();
911 // if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - Add second input clusters, entries %d\n", nclusters) ;
912 // for (Int_t iclus = 0; iclus < nclusters; iclus++) {
913 // AliAODCaloCluster * clus = 0;
914 // if ( (clus = ((AliAODEvent*)fSecondInputAODEvent)->GetCaloCluster(iclus)) ) {
915 // if (clus->IsEMCAL()){
916 // TLorentzVector momentum ;
917 // clus->GetMomentum(momentum, v);
919 // if(fEMCALPtMin < momentum.Pt()){
921 // if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) continue;
923 // if(fDebug > 2 && momentum.E() > 0.1) printf("AliCaloTrackReader::FillInputEMCAL() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
924 // momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
925 // fEMCALClusters->Add(clus);
926 // }//Pt and Fiducial cut passed.
928 // }// cluster exists
931 // if(fDebug > 1) printf("AliCaloTrackReader::FillInputEMCAL() - aod normal entries %d, after second input %d\n", fEMCALClustersNormalInputEntries, fEMCALClusters->GetEntriesFast());
933 // } //second input loop
936 //____________________________________________________________________________
937 void AliCaloTrackReader::FillInputPHOS() {
938 //Return array with PHOS clusters in aod format
940 if(fDebug > 2 ) printf("AliCaloTrackReader::FillInputPHOS()\n");
942 //Loop to select clusters in fiducial cut and fill container with aodClusters
943 Int_t nclusters = fInputEvent->GetNumberOfCaloClusters();
944 for (Int_t iclus = 0; iclus < nclusters; iclus++) {
945 AliVCluster * clus = 0;
946 if ( (clus = fInputEvent->GetCaloCluster(iclus)) ) {
947 if (IsPHOSCluster(clus)){
948 //Check if the cluster contains any bad channel and if close to calorimeter borders
951 vindex = fMixedEvent->EventIndexForCaloCluster(iclus);
952 if( GetCaloUtils()->ClusterContainsBadChannel("PHOS",clus->GetCellsAbsId(), clus->GetNCells()))
954 if(!GetCaloUtils()->CheckCellFiducialRegion(clus, fInputEvent->GetPHOSCells(), fInputEvent, vindex))
957 TLorentzVector momentum ;
959 clus->GetMomentum(momentum, fVertex[vindex]);
961 if(fPHOSPtMin < momentum.Pt()){
963 if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS"))
966 if(fDebug > 2 && momentum.E() > 0.1)
967 printf("AliCaloTrackReader::FillInputPHOS() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
968 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
970 //Recalibrate the cluster energy
971 if(GetCaloUtils()->IsRecalibrationOn()) {
972 Float_t energy = GetCaloUtils()->RecalibrateClusterEnergy(clus, (AliAODCaloCells*)GetPHOSCells());
980 fPHOSClusters->Add(clus);
982 }//Pt and Fiducial cut passed.
987 //fPHOSClustersNormalInputEntries = fPHOSClusters->GetEntriesFast() ;
988 if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS() - aod entries %d\n", fPHOSClusters->GetEntriesFast());//fPHOSClustersNormalInputEntries);
990 //If second input event available, add the clusters.
991 // if(fSecondInputAODTree && fSecondInputAODEvent){
992 // GetSecondInputAODVertex(v);
993 // nclusters = ((AliAODEvent*)fSecondInputAODEvent)->GetNumberOfCaloClusters();
994 // if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS() - Add second input clusters, entries %d\n", nclusters);
995 // for (Int_t iclus = 0; iclus < nclusters; iclus++) {
996 // AliAODCaloCluster * clus = 0;
997 // if ( (clus = ((AliAODEvent*)fSecondInputAODEvent)->GetCaloCluster(iclus)) ) {
998 // if (clus->IsPHOS()){
999 // TLorentzVector momentum ;
1000 // clus->GetMomentum(momentum, v);
1002 // if(fPHOSPtMin < momentum.Pt()){
1004 // if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"PHOS")) continue;
1006 // if(fDebug > 2 && momentum.E() > 0.1) printf("AliCaloTrackReader::FillInputPHOS() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
1007 // momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
1008 // fPHOSClusters->Add(clus);
1009 // }//Pt and Fiducial cut passed.
1011 // }// cluster exists
1013 // if(fDebug > 1) printf("AliCaloTrackReader::FillInputPHOS() - aod normal entries %d, after second input %d\n", fPHOSClustersNormalInputEntries, fPHOSClusters->GetEntriesFast());
1014 // } //second input loop
1018 //____________________________________________________________________________
1019 void AliCaloTrackReader::FillInputEMCALCells() {
1020 //Return array with EMCAL cells in aod format
1022 fEMCALCells = fInputEvent->GetEMCALCells();
1026 //____________________________________________________________________________
1027 void AliCaloTrackReader::FillInputPHOSCells() {
1028 //Return array with PHOS cells in aod format
1030 fPHOSCells = fInputEvent->GetPHOSCells();
1034 //____________________________________________________________________________
1035 void AliCaloTrackReader::FillInputVZERO(){
1036 //Fill VZERO information in data member, add all the channels information.
1037 AliVVZERO* v0 = fInputEvent->GetVZEROData();
1038 //printf("Init V0: ADC (%d,%d), Multiplicity (%d,%d) \n",fV0ADC[0],fV0ADC[1],fV0Mul[0],fV0Mul[1]);
1042 AliESDVZERO* esdV0 = dynamic_cast<AliESDVZERO*> (v0);
1043 for (Int_t i = 0; i < 32; i++)
1045 if(esdV0){//Only available in ESDs
1046 fV0ADC[0] += (Int_t)esdV0->GetAdcV0C(i);
1047 fV0ADC[1] += (Int_t)esdV0->GetAdcV0A(i);
1049 fV0Mul[0] += (Int_t)v0->GetMultiplicityV0C(i);
1050 fV0Mul[1] += (Int_t)v0->GetMultiplicityV0A(i);
1053 printf("V0: ADC (%d,%d), Multiplicity (%d,%d) \n",fV0ADC[0],fV0ADC[1],fV0Mul[0],fV0Mul[1]);
1057 printf("Cannot retrieve V0 ESD! Run w/ null V0 charges");
1062 //____________________________________________________________________________
1063 Bool_t AliCaloTrackReader::IsEMCALCluster(AliVCluster* cluster) const {
1064 // Check if it is a cluster from EMCAL. For old AODs cluster type has
1065 // different number and need to patch here
1067 if(fDataType==kAOD && fOldAOD)
1069 if (cluster->GetType() == 2) return kTRUE;
1074 return cluster->IsEMCAL();
1079 //____________________________________________________________________________
1080 Bool_t AliCaloTrackReader::IsPHOSCluster(AliVCluster * cluster) const {
1081 //Check if it is a cluster from PHOS.For old AODs cluster type has
1082 // different number and need to patch here
1084 if(fDataType==kAOD && fOldAOD)
1086 Int_t type = cluster->GetType();
1087 if (type == 0 || type == 1) return kTRUE;
1092 return cluster->IsPHOS();
1097 //____________________________________________________________________________
1098 Bool_t AliCaloTrackReader::CheckForPrimaryVertex(){
1099 //Check if the vertex was well reconstructed, copy from V0Reader of conversion group
1102 AliESDEvent * event = dynamic_cast<AliESDEvent*> (fInputEvent);
1103 if(!event) return kFALSE;
1105 if(event->GetPrimaryVertexTracks()->GetNContributors() > 0) {
1109 if(event->GetPrimaryVertexTracks()->GetNContributors() < 1) {
1111 if(event->GetPrimaryVertexSPD()->GetNContributors() > 0) {
1112 //cout<<"spd vertex type::"<< fESDEvent->GetPrimaryVertex()->GetName() << endl;
1116 if(event->GetPrimaryVertexSPD()->GetNContributors() < 1) {
1117 // cout<<"bad vertex type::"<< event->GetPrimaryVertex()->GetName() << endl;