// $Id$ // // Emcal jet finder task. // // Authors: C.Loizides, S.Aiola #include #include "AliEmcalJetTask.h" #include #include #include #include #include "AliAnalysisManager.h" #include "AliCentrality.h" #include "AliEMCALGeometry.h" #include "AliESDEvent.h" #include "AliEmcalJet.h" #include "AliEmcalParticle.h" #include "AliFJWrapper.h" #include "AliMCEvent.h" #include "AliVCluster.h" #include "AliVEvent.h" #include "AliVParticle.h" ClassImp(AliEmcalJetTask) //________________________________________________________________________ AliEmcalJetTask::AliEmcalJetTask() : AliAnalysisTaskSE("AliEmcalJetTask"), fTracksName("Tracks"), fCaloName("CaloClusters"), fJetsName("Jets"), fAlgo(1), fRadius(0.4), fType(0), fMinJetTrackPt(0.15), fMinJetClusPt(0.15), fPhiMin(0), fPhiMax(TMath::TwoPi()), fEtaMin(-1), fEtaMax(1), fMinJetArea(0.01), fMinJetPt(1.0), fIsInit(0), fIsMcPart(0), fIsEmcPart(0), fJets(0), fEvent(0), fTracks(0), fClus(0) { // Default constructor. } //________________________________________________________________________ AliEmcalJetTask::AliEmcalJetTask(const char *name) : AliAnalysisTaskSE(name), fTracksName("Tracks"), fCaloName("CaloClusters"), fJetsName("Jets"), fAlgo(1), fRadius(0.4), fType(0), fMinJetTrackPt(0.15), fMinJetClusPt(0.15), fPhiMin(0), fPhiMax(TMath::TwoPi()), fEtaMin(-1), fEtaMax(1), fMinJetArea(0.01), fMinJetPt(1.0), fIsInit(0), fIsMcPart(0), fIsEmcPart(0), fJets(0), fEvent(0), fTracks(0), fClus(0) { // Standard constructor. fBranchNames="ESD:AliESDRun.,AliESDHeader.,PrimaryVertex."; } //________________________________________________________________________ AliEmcalJetTask::~AliEmcalJetTask() { // Destructor } //________________________________________________________________________ void AliEmcalJetTask::UserCreateOutputObjects() { // Create user objects. fJets = new TClonesArray("AliEmcalJet"); fJets->SetName(fJetsName); } //________________________________________________________________________ void AliEmcalJetTask::UserExec(Option_t *) { // Main loop, called for each event. if (!fIsInit) { if (!DoInit()) return; fIsInit = kTRUE; } FindJets(); } //________________________________________________________________________ void AliEmcalJetTask::Terminate(Option_t *) { // Called once at the end of the analysis. } //________________________________________________________________________ void AliEmcalJetTask::FindJets() { // Find jets. if (!fTracks && !fClus) return; TString name("kt"); fastjet::JetAlgorithm jalgo(fastjet::kt_algorithm); if (fAlgo>=1) { name = "antikt"; jalgo = fastjet::antikt_algorithm; } // setup fj wrapper AliFJWrapper fjw(name, name); fjw.SetAreaType(fastjet::active_area_explicit_ghosts); fjw.SetR(fRadius); fjw.SetAlgorithm(jalgo); fjw.SetMaxRap(1); fjw.Clear(); // get primary vertex Double_t vertex[3] = {0, 0, 0}; fEvent->GetPrimaryVertex()->GetXYZ(vertex); if ((fIsMcPart || fType == 0 || fType == 1) && fTracks) { const Int_t Ntracks = fTracks->GetEntries(); for (Int_t iTracks = 0; iTracks < Ntracks; ++iTracks) { AliVParticle *t = static_cast(fTracks->At(iTracks)); if (!t) continue; if (fIsMcPart && fType == 1 && t->Charge() == 0) continue; if (fIsMcPart && fType == 2 && t->Charge() != 0) continue; if (t->Pt() < fMinJetTrackPt) continue; Double_t eta = t->Eta(); Double_t phi = t->Phi(); if ((etafEtaMax) && (phiPx(), t->Py(), t->Pz(), t->P(), iTracks + 100); } } if ((fType == 0 || fType == 2) && fClus) { const Int_t Nclus = fClus->GetEntries(); for (Int_t iClus = 0; iClus < Nclus; ++iClus) { AliVCluster *c = 0; Double_t cEta=0,cPhi=0,cPt=0; Double_t cPx=0,cPy=0,cPz=0; if (fIsEmcPart) { AliEmcalParticle *ep = static_cast(fClus->At(iClus)); if (!ep) continue; c = ep->GetCluster(); if (!c) continue; cEta = ep->Eta(); cPhi = ep->Phi(); cPt = ep->Pt(); cPx = ep->Px(); cPy = ep->Py(); cPz = ep->Pz(); } else { c = static_cast(fClus->At(iClus)); if (!c) continue; TLorentzVector nP; c->GetMomentum(nP, vertex); cEta = nP.Eta(); cPhi = nP.Phi(); cPt = nP.Pt(); cPx = nP.Px(); cPy = nP.Py(); cPz = nP.Pz(); } if (!c->IsEMCAL()) continue; if (cPt < fMinJetClusPt) continue; if ((cEtafEtaMax) && (cPhi jets_incl = fjw.GetInclusiveJets(); for (UInt_t ij=0, jetCount=0; ij constituents(fjw.GetJetConstituents(ij)); jet->SetNumberOfTracks(constituents.size()); jet->SetNumberOfClusters(constituents.size()); Int_t nt = 0; Int_t nc = 0; Double_t neutralE = 0; Double_t maxCh = 0; Double_t maxNe = 0; Int_t gall = 0; Int_t gemc = 0; Int_t cemc = 0; Int_t ncharged = 0; Int_t nneutral = 0; Double_t mcpt = 0; Double_t emcpt = 0; for(UInt_t ic = 0; ic < constituents.size(); ++ic) { Int_t uid = constituents[ic].user_index(); if ((uid == -1) /*&& (constituents[ic].kt2() < 1e-25)*/) { //ghost particle ++gall; Double_t gphi = constituents[ic].phi(); if (gphi<0) gphi += TMath::TwoPi(); gphi *= TMath::RadToDeg(); Double_t geta = constituents[ic].eta(); if ((gphi > geom->GetArm1PhiMin()) && (gphi < geom->GetArm1PhiMax()) && (geta > geom->GetArm1EtaMin()) && (geta < geom->GetArm1EtaMax())) ++gemc; } else if ((uid > 0) && fTracks) { // track constituent Int_t tid = uid - 100; AliVParticle *t = static_cast(fTracks->At(tid)); if (!t) continue; Double_t cEta = t->Eta(); Double_t cPhi = t->Phi(); Double_t cPt = t->Pt(); Double_t cP = t->P(); if (t->Charge() == 0) { neutralE += cP; ++nneutral; if (cPt > maxNe) maxNe = cPt; } else { ++ncharged; if (cPt > maxCh) maxCh = cPt; } if (fIsMcPart || t->GetLabel() == 100) // check if MC particle mcpt += cPt; if (cPhi<0) cPhi += TMath::TwoPi(); cPhi *= TMath::RadToDeg(); if ((cPhi > geom->GetArm1PhiMin()) && (cPhi < geom->GetArm1PhiMax()) && (cEta > geom->GetArm1EtaMin()) && (cEta < geom->GetArm1EtaMax())) { emcpt += cPt; ++cemc; } jet->AddTrackAt(tid, nt); ++nt; } else if (fClus) { // cluster constituent Int_t cid = -uid - 100; AliVCluster *c = 0; Double_t cEta=0,cPhi=0,cPt=0,cP=0; if (fIsEmcPart) { AliEmcalParticle *ep = static_cast(fClus->At(cid)); if (!ep) continue; c = ep->GetCluster(); if (!c) continue; cEta = ep->Eta(); cPhi = ep->Phi(); cPt = ep->Pt(); cP = ep->P(); } else { c = static_cast(fClus->At(cid)); if (!c) continue; TLorentzVector nP; c->GetMomentum(nP, vertex); cEta = nP.Eta(); cPhi = nP.Phi(); cPt = nP.Pt(); cP = nP.P(); } neutralE += cP; if (cPt > maxNe) maxNe = cPt; if (c->Chi2() == 100) // MC particle mcpt += cPt; if (cPhi<0) cPhi += TMath::TwoPi(); cPhi *= TMath::RadToDeg(); if ((cPhi > geom->GetArm1PhiMin()) && (cPhi < geom->GetArm1PhiMax()) && (cEta > geom->GetArm1EtaMin()) && (cEta < geom->GetArm1EtaMax())) { emcpt += cPt; ++cemc; } jet->AddClusterAt(cid, nc); ++nc; ++nneutral; } else { AliError(Form("%s: No logical way to end up here.", GetName())); continue; } } /* end of constituent loop */ jet->SetNumberOfTracks(nt); jet->SetNumberOfClusters(nc); jet->SortConstituents(); jet->SetMaxNeutralPt(maxNe); jet->SetMaxChargedPt(maxCh); jet->SetNEF(neutralE / jet->E()); fastjet::PseudoJet area(fjw.GetJetAreaVector(ij)); jet->SetArea(area.perp()); jet->SetAreaEta(area.eta()); jet->SetAreaPhi(area.phi()); jet->SetNumberOfCharged(ncharged); jet->SetNumberOfNeutrals(nneutral); jet->SetMCPt(mcpt); jet->SetNEmc(cemc); jet->SetPtEmc(emcpt); if (gall > 0) jet->SetAreaEmc(fjw.GetJetArea(ij) * gemc / gall); else jet->SetAreaEmc(-1); if ((jet->Phi() > geom->GetArm1PhiMin() * TMath::DegToRad()) && (jet->Phi() < geom->GetArm1PhiMax() * TMath::DegToRad()) && (jet->Eta() > geom->GetArm1EtaMin()) && (jet->Eta() < geom->GetArm1EtaMax())) jet->SetAxisInEmcal(kTRUE); jetCount++; } fJets->Sort(); } //________________________________________________________________________ Bool_t AliEmcalJetTask::DoInit() { // Init. Return true if successful. // get input collections AliAnalysisManager *am = AliAnalysisManager::GetAnalysisManager(); // get the event fEvent = InputEvent(); if (!fEvent) { AliError(Form("%s: Could not retrieve event! Returning", GetName())); return 0; } // add jets to event if not yet there fJets->Delete(); if (!(fEvent->FindListObject(fJetsName))) fEvent->AddObject(fJets); else { AliError(Form("%s: Object with name %s already in event! Returning", GetName(), fJetsName.Data())); return 0; } if (fTracksName == "Tracks") am->LoadBranch("Tracks"); if (!fTracks && !fTracksName.IsNull()) { fTracks = dynamic_cast(fEvent->FindListObject(fTracksName)); if (!fTracks) { AliError(Form("%s: Pointer to tracks %s == 0", GetName(), fTracksName.Data())); return 0; } } if (fTracks) { TClass cls(fTracks->GetClass()->GetName()); if (cls.InheritsFrom("AliMCParticle")) fIsMcPart = 1; } if (fCaloName == "CaloClusters") am->LoadBranch("CaloClusters"); if (!fClus && !fCaloName.IsNull()) { fClus = dynamic_cast(fEvent->FindListObject(fCaloName)); if (!fClus) { AliError(Form("%s: Pointer to clus %s == 0", GetName(), fCaloName.Data())); return 0; } } if (fClus) { TClass cls(fClus->GetClass()->GetName()); if (cls.InheritsFrom("AliEmcalParticle")) fIsEmcPart = 1; } return 1; }