20 #include <TClonesArray.h>
22 #include <TObjArray.h>
26 #include "AliAnalysisUtils.h"
27 #include "AliAODEvent.h"
28 #include "AliAODMCParticle.h"
29 #include "AliAODTrack.h"
30 #include "AliCentrality.h"
31 #include "AliESDEvent.h"
32 #include "AliESDtrack.h"
33 #include "AliInputEventHandler.h"
35 #include "AliMCEvent.h"
36 #include "AliPicoTrack.h"
38 #include "AliVCluster.h"
39 #include "AliVEvent.h"
40 #include "AliVTrack.h"
41 #include "AliVVertex.h"
42 #include "AliEmcalTriggerPatchInfo.h"
55 ClassImp(HighPtTracks::AliReducedTrackSelectionContainer)
56 ClassImp(HighPtTracks::AliReducedHighPtEventCreator)
59 namespace HighPtTracks {
64 AliReducedHighPtEventCreator::AliReducedHighPtEventCreator():
65 AliAnalysisTaskEmcal(),
68 fTrackSelections(NULL),
69 fSwapTriggerThresholds(kFALSE),
76 fCentralityMethod(
"V0A")
85 AliReducedHighPtEventCreator::AliReducedHighPtEventCreator(
const char* name):
86 AliAnalysisTaskEmcal(name, kTRUE),
89 fTrackSelections(NULL),
90 fSwapTriggerThresholds(kFALSE),
97 fCentralityMethod(
"V0A")
99 DefineOutput(2, TTree::Class());
101 fTrackSelections =
new TObjArray;
102 fTrackSelections->SetOwner(kTRUE);
108 AliReducedHighPtEventCreator::~AliReducedHighPtEventCreator() {
109 if(fTrackSelections)
delete fTrackSelections;
115 void AliReducedHighPtEventCreator::UserCreateOutputObjects() {
116 AliAnalysisTaskEmcal::UserCreateOutputObjects();
119 fOutputTree =
new TTree(
"ReducedEvent",
"A simple reduced event");
120 fOutputEvent =
new AliReducedHighPtEvent();
121 fOutputTree->Branch(
"ReducedEvent",
"AliReducedHighPtEvent", fOutputEvent, 128000, 0);
123 PostData(1, fOutput);
124 PostData(2, fOutputTree);
137 Bool_t AliReducedHighPtEventCreator::Run() {
139 AliError(
"Cluster container missing");
143 AliError(
"Track container missing");
146 if(!SelectEvent(fInputEvent))
return kFALSE;
147 new(fOutputEvent) AliReducedHighPtEvent(kTRUE);
150 fOutputEvent->SetVertexZ(fInputEvent->GetPrimaryVertex()->GetZ());
151 AliCentrality *centralityHandler = fInputEvent->GetCentrality();
152 if(centralityHandler) fOutputEvent->SetCentralityPercentile(centralityHandler->GetCentralityPercentile(fCentralityMethod.Data()));
153 ConvertTriggerPatches(fTriggerPatchInfo, fOutputEvent->GetPatchContainer());
154 TString triggerString(fInputEvent->GetFiredTriggerClasses());
155 fOutputEvent->SetDecisionFromTriggerString(triggerString.Contains(
"EG1"), triggerString.Contains(
"EG2"), triggerString.Contains(
"EJ1"), triggerString.Contains(
"EJ2"));
156 fOutputEvent->SetMinBiasEvent(fInputHandler->IsEventSelected() & AliVEvent::kINT7);
158 std::map<int, int> mcindexmap;
161 AliReducedMCHeader *mcheader =
new AliReducedMCHeader();
163 mcheader->SetCrossSection(fXsection);
164 mcheader->SetNumberOfTrials(fNTrials);
165 mcheader->SetPtHard(fPtHard);
167 fOutputEvent->SetMonteCarloHeader(mcheader);
170 AliMCEvent *mcev = MCEvent();
172 for(Int_t ipart = 0; ipart < mcev->GetNumberOfTracks(); ipart++){
173 AliVParticle *part = mcev->GetTrack(ipart);
174 Double_t pt(TMath::Abs(part->Pt())), eta(part->Eta());
175 if(pt < fMinPt || pt > fMaxPt)
continue;
176 if(eta < fMinEta || eta > fMaxEta)
continue;
177 if(!part->Charge())
continue;
178 if(part->IsA() == AliAODMCParticle::Class()){
179 AliAODMCParticle *aodpart =
static_cast<AliAODMCParticle *
>(part);
180 if(!aodpart->IsPhysicalPrimary())
continue;
183 if(!mcev->Stack()->IsPhysicalPrimary(ipart))
continue;
186 AliReducedGeneratedParticle *reducedgen =
new AliReducedGeneratedParticle(npart, part->PdgCode(), part->Px(), part->Py(), part->Pz(), part->E());
187 fOutputEvent->AddReducedGeneratedParticle(reducedgen);
188 mcindexmap.insert(std::pair<int,int>(ipart, npart));
194 std::map<int, int> clusterindexmap;
196 Double_t vtxpos[3];fInputEvent->GetPrimaryVertex()->GetXYZ(vtxpos);
197 for(TIter cliter = TIter(fCaloClusters).Begin(); cliter != TIter::End(); ++cliter){
198 AliVCluster *incluster =
static_cast<AliVCluster *
>(*cliter);
199 if(!SelectCluster(incluster))
continue;
200 TLorentzVector clustervec;
201 incluster->GetMomentum(clustervec, vtxpos);
202 AliReducedEmcalCluster *redcluster =
new AliReducedEmcalCluster(ncluster, incluster->E(), clustervec.Eta(), clustervec.Phi(), incluster->GetM02(), incluster->GetM20());
204 TArrayD cellEnergies;
205 GetCellEnergies(incluster, cellEnergies);
206 TArrayI indices(cellEnergies.GetSize());
207 TMath::Sort(cellEnergies.GetSize(), cellEnergies.GetArray(), indices.GetArray(), kTRUE);
208 redcluster->SetLeadingCellEnergies(
209 cellEnergies[indices[0]],
210 cellEnergies.GetSize() > 1 ? cellEnergies[indices[1]] : 0,
211 cellEnergies.GetSize() > 2 ? cellEnergies[indices[2]] : 0
215 for(Int_t ilab = 0; ilab < incluster->GetNLabels(); ilab++){
216 AliVParticle *assigned = MCEvent()->GetTrack(TMath::Abs(incluster->GetLabels()[ilab]));
217 if(!assigned)
continue;
218 redcluster->AddTrueContributor(assigned->PdgCode(), assigned->Px(), assigned->Py(), assigned->Pz(), assigned->E());
221 fOutputEvent->AddReducedCluster(redcluster);
222 clusterindexmap.insert(std::pair<int,int>(incluster->GetID(), ncluster));
227 for(TIter trkiter = TIter(fTracks).Begin(); trkiter != TIter::End(); ++trkiter){
228 AliVTrack *trackRaw =
static_cast<AliVTrack *
>(*trkiter), *trackToCheck(NULL);
230 if(trackRaw->IsA() == AliPicoTrack::Class()){
231 AliPicoTrack *picotrack =
static_cast<AliPicoTrack *
>(trackRaw);
232 trackToCheck = picotrack->GetTrack();
234 trackToCheck = trackRaw;
236 FixTrackInputEvent(trackToCheck);
238 if(!SelectTrack(trackToCheck,cutIndices))
continue;
239 AliReducedReconstructedTrack *rectrack =
new AliReducedReconstructedTrack;
241 trackToCheck->GetPxPyPz(pvec);
242 rectrack->SetMomentumVector(pvec[0], pvec[1], pvec[2]);
243 for(Int_t icut = 0; icut < cutIndices.GetSize(); icut++) rectrack->SetTrackCuts(cutIndices[icut]);
246 Int_t label = TMath::Abs(trackToCheck->GetLabel());
247 if(label < fMCEvent->GetNumberOfTracks()){
248 std::map<int,int>::iterator found = mcindexmap.find(label);
249 if(found != mcindexmap.end()){
250 rectrack->SetMatchedParticleIndex(found->second);
252 if(label >= 0) rectrack->SetGoodTrackLabel(kTRUE);
253 rectrack->SetTPCClusters(trackToCheck->GetTPCNcls());
254 rectrack->SetTPCCrossedRows(GetTPCCrossedRows(trackToCheck));
255 rectrack->SetTPCSharedClusters(trackToCheck->GetTPCSharedMapPtr()->CountBits());
256 rectrack->SetTPCFindableClusters(trackToCheck->GetTPCNclsF());
260 Int_t clusterIndex = trackToCheck->GetEMCALcluster();
261 if(clusterIndex >= 0 && clusterIndex < fCaloClusters->GetEntries()){
262 std::map<int,int>::iterator found = clusterindexmap.find(clusterIndex);
263 if(found != clusterindexmap.end()){
264 rectrack->SetMatchedClusterIndex(found->second);
270 PostData(2, fOutputTree);
279 void AliReducedHighPtEventCreator::AddVirtualTrackSelection(
281 fTrackSelections->Add(
new AliReducedTrackSelectionContainer(cutindex, sel));
289 Bool_t AliReducedHighPtEventCreator::SelectEvent(AliVEvent* event)
const {
290 const AliVVertex *primvtx =
event->GetPrimaryVertex();
291 if(!primvtx || !primvtx->GetNContributors())
return kFALSE;
292 if(TMath::Abs(primvtx->GetZ()) > 10.)
return kFALSE;
293 if(event->IsPileupFromSPD(3, 0.8, 3., 2., 5.))
return kFALSE;
294 AliAnalysisUtils eventSelUtil;
295 if(!eventSelUtil.IsVertexSelected2013pA(event))
return kFALSE;
304 Bool_t AliReducedHighPtEventCreator::SelectCluster(
const AliVCluster* clust)
const {
305 if(!clust->IsEMCAL())
return kFALSE;
306 if(clust->E() < fMinClusterE || clust->E() > fMaxClusterE)
return kFALSE;
316 Int_t AliReducedHighPtEventCreator::SelectTrack(AliVTrack* track, TArrayI& cutindices)
const {
317 std::vector<int> selected;
319 double pt = TMath::Abs(track->Pt()), eta = track->Eta();
320 if(pt < fMinPt || pt > fMaxPt)
return 0;
321 if(eta < fMinEta || eta > fMaxEta)
return 0;
323 for(TIter cutiter = TIter(fTrackSelections).Begin(); cutiter != TIter::End(); ++cutiter){
324 AliReducedTrackSelectionContainer *mycut =
static_cast<AliReducedTrackSelectionContainer *
>(*cutiter);
325 if(!mycut->GetTrackSelection()->IsTrackAccepted(track))
continue;
326 selected.push_back(mycut->GetIndex());
328 if(!selected.size())
return 0;
329 cutindices.Set(selected.size());
331 for(std::vector<int>::iterator inditer = selected.begin(); inditer != selected.end(); ++inditer){
332 cutindices[counter++] = *inditer;
334 return selected.size();
342 Int_t AliReducedHighPtEventCreator::GetTPCCrossedRows(
const AliVTrack* trk)
const {
343 if(trk->IsA() == AliESDtrack::Class()){
344 return (static_cast<const AliESDtrack *>(trk))->GetTPCCrossedRows();
345 }
else if(trk->IsA() == AliAODTrack::Class()){
346 return (static_cast<const AliAODTrack *>(trk))->GetTPCNCrossedRows();
356 void AliReducedHighPtEventCreator::GetCellEnergies(AliVCluster* emccluster, TArrayD& energies)
const {
357 if(!fInputEvent->GetEMCALCells()) {
358 AliError(
"No EMCAL cells");
361 AliDebug(2, Form(
"Number of cells: %d, array: %p", emccluster->GetNCells(), emccluster->GetCellsAbsId()));
362 energies.Set(emccluster->GetNCells());
363 for(
int icell = 0; icell < emccluster->GetNCells(); icell++){
365 energies[icell] = fInputEvent->GetEMCALCells()->GetCellAmplitude(emccluster->GetCellsAbsId()[icell]);
374 void AliReducedHighPtEventCreator::ConvertTriggerPatches(TClonesArray* patches,
375 AliReducedPatchContainer* cont) {
377 AliError(
"Trigger patch container not found\n");
380 for(TIter patchIter = TIter(patches).Begin(); patchIter != TIter::End(); ++patchIter){
381 AliEmcalTriggerPatchInfo *mypatch =
static_cast<AliEmcalTriggerPatchInfo *
>(*patchIter);
382 if(!mypatch->IsOfflineSimple() && mypatch->IsLevel0())
continue;
384 Bool_t isDefined(kFALSE);
385 if(mypatch->IsOfflineSimple()){
391 AliDebug(2,
"Unknown offline patch type");
394 AliDebug(2, Form(
"Adding offline patch of type %d",
int(triggertype)));
395 cont->AddTriggerPatch(kTRUE, triggertype, mypatch->GetPatchE(), mypatch->GetADCAmp(), mypatch->GetEtaGeo(), mypatch->GetPhiGeo());
402 AliDebug(2,
"Unknown online patch type");
405 AliDebug(2, Form(
"Adding online patch of type %d",
int(triggertype)));
406 cont->AddTriggerPatch(kFALSE, triggertype, mypatch->GetPatchE(), mypatch->GetADCAmp(), mypatch->GetEtaGeo(), mypatch->GetPhiGeo());
416 void AliReducedHighPtEventCreator::FixTrackInputEvent(AliVTrack* trk) {
417 if(!trk->GetEvent()){
418 if(trk->IsA() == AliESDtrack::Class())
419 (static_cast<AliESDtrack *>(trk))->SetESDEvent(static_cast<AliESDEvent *>(fInputEvent));
420 else if(trk->IsA() == AliAODTrack::Class())
421 (static_cast<AliAODTrack *>(trk))->SetAODEvent(static_cast<AliAODEvent *>(fInputEvent));
422 else if(trk->IsA() == AliPicoTrack::Class()){
423 AliPicoTrack *mytrk =
static_cast<AliPicoTrack *
>(trk);
424 if(!mytrk->GetEvent()){
425 if(mytrk->GetTrack()->IsA() == AliESDtrack::Class())
426 (static_cast<AliESDtrack *>(mytrk->GetTrack()))->SetESDEvent(static_cast<AliESDEvent *>(fInputEvent));
428 (static_cast<AliAODTrack *>(mytrk->GetTrack()))->SetAODEvent(static_cast<AliAODEvent *>(fInputEvent));
437 AliReducedTrackSelectionContainer::AliReducedTrackSelectionContainer():
439 fTrackSelection(NULL)
448 AliReducedTrackSelectionContainer::AliReducedTrackSelectionContainer(
458 AliReducedTrackSelectionContainer::~AliReducedTrackSelectionContainer() {
459 if(fTrackSelection)
delete fTrackSelection;
Reduced event structure for high- analysis.
ClassImp(AliAnalysisTaskTriggerRates) AliAnalysisTaskTriggerRates
Jet trigger, low threshold.
Interface for virtual track selection.
Declaration of class AliReducedHighPtEventCreator and AliReducedTrackSelectionContainer.
Declaration of class AlliReducedPatchInfo.
Declaration of class AliReducedPatchContainer, a container for reduced trigger patches.
Declaration of class AliReducedGeneratedParticle.
Gamma trigger, low threshold.
Reduced information about reconstructed EMCAL clusters.
Jet trigger, high threshold.
Declartion of class AliEMCalPtTaskVTrackSelection.
Gamma trigger, high threshold.