]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added possibility to store tracklets and muons in MC mode
authorekryshen <evgeny.kryshen@cern.ch>
Sat, 29 Mar 2014 14:38:07 +0000 (15:38 +0100)
committerekryshen <evgeny.kryshen@cern.ch>
Sat, 29 Mar 2014 14:38:07 +0000 (15:38 +0100)
PWGCF/Correlations/Base/AliAnalysisTaskCFTree.cxx
PWGCF/Correlations/Base/AliAnalysisTaskCFTree.h
PWGCF/Correlations/Base/AliCFParticle.h

index 25eb765ec9bdf723ac83f5898b704bc30e9efa73..fc6e4ec9935d29b1453fff69aa594288603b3f47 100644 (file)
@@ -86,6 +86,7 @@ AliAnalysisTaskCFTree::AliAnalysisTaskCFTree(const char* name) :
   fSharedClusterCut(0.4),
   fCrossedRowsCut(100),
   fFoundFractionCut(0.8),
+  fStoreTracks(0),
   fStoreTracklets(0),
   fStoreMuons(0)
 {
@@ -115,7 +116,7 @@ void AliAnalysisTaskCFTree::CreateOutputObjects(){
 
   fListOfHistos->Add(fEventStatistics);
 
-  fParticles = new TClonesArray("AliCFParticle",2000);
+  if (fStoreTracks)    fParticles = new TClonesArray("AliCFParticle",2000);
   if (fStoreTracklets) fTracklets = new TClonesArray("AliCFParticle",100);
   if (fStoreMuons)     fMuons     = new TClonesArray("AliCFParticle",100);
   // create file-resident tree
@@ -131,9 +132,9 @@ void AliAnalysisTaskCFTree::CreateOutputObjects(){
   fTree->Branch("orbit",&fOrbit);
   fTree->Branch("bc",&fBc);
   fTree->Branch("mask",&fSelectMask);
-  fTree->Branch("particles",&fParticles);
-  fTree->Branch("tracklets",&fTracklets);
-  fTree->Branch("muons",&fMuons);
+  if (fStoreTracks)    fTree->Branch("particles",&fParticles);
+  if (fStoreTracklets) fTree->Branch("tracklets",&fTracklets);
+  if (fStoreMuons)     fTree->Branch("muons",&fMuons);
   PostData(0,fListOfHistos);
   PostData(1,fTree);
 }
@@ -142,7 +143,9 @@ void AliAnalysisTaskCFTree::CreateOutputObjects(){
 
 //-----------------------------------------------------------------------------
 void AliAnalysisTaskCFTree::Exec(Option_t *){
-  fParticles->Clear();
+  if (fParticles) fParticles->Clear();
+  if (fTracklets) fTracklets->Clear();
+  if (fMuons)     fMuons->Clear();
   
   fEventStatistics->Fill("before cuts",1);
   
@@ -172,20 +175,21 @@ void AliAnalysisTaskCFTree::Exec(Option_t *){
     if (TMath::Abs(fZvtx) >= fZVertexCut || vertex->GetNContributors()<fnTracksVertex)  return;
     fEventStatistics->Fill("after vertex selection",1);
     
-    for (Int_t ipart=0;ipart<event->GetNumberOfTracks();ipart++){
-      AliVTrack* track = (AliVTrack*) event->GetTrack(ipart);
-      if (!track) continue;
-      UInt_t mask = GetFilterMap(track);
-      if (!(mask & fTrackFilterBit)) continue;
+    if (fParticles) {
+      for (Int_t ipart=0;ipart<event->GetNumberOfTracks();ipart++){
+        AliVTrack* track = (AliVTrack*) event->GetTrack(ipart);
+        if (!track) continue;
+        UInt_t mask = GetFilterMap(track);
+        if (!(mask & fTrackFilterBit)) continue;
 
-      if (track->InheritsFrom("AliAODTrack")) AddTrack(track,mask,0);
-      else if (track->InheritsFrom("AliESDtrack")) {
-        if (mask)                           AddTrack(track,mask,1);
-        if (mask & fHybridConstrainedMask)  AddTrack(track,mask,2);
-        if (mask & fTPConlyConstrainedMask) AddTrack(track,mask,3);
+        if (track->InheritsFrom("AliAODTrack")) AddTrack(track,mask,0);
+        else if (track->InheritsFrom("AliESDtrack")) {
+          if (mask)                           AddTrack(track,mask,1);
+          if (mask & fHybridConstrainedMask)  AddTrack(track,mask,2);
+          if (mask & fTPConlyConstrainedMask) AddTrack(track,mask,3);
+        }
       }
     }
-
     if (fIsAOD){
       AliAODEvent* aod = (AliAODEvent*) event;
       
@@ -264,51 +268,104 @@ void AliAnalysisTaskCFTree::Exec(Option_t *){
     if (vertex->GetNContributors()<fnTracksVertex)  fZvtx=-1000;
     fEventStatistics->Fill("after check for vertex contributors",1);
     
-    UInt_t* masks = new UInt_t[nProduced];
-    for (Int_t i=0;i<nProduced;i++) masks[i]=0;
-    
-    // Loop over reconstructed tracks to set masks
-    for (Int_t ipart=0;ipart<event->GetNumberOfTracks();ipart++){
-      AliVTrack* part = (AliVTrack*) event->GetTrack(ipart);
-      if (!part) continue;
-      Int_t label = TMath::Abs(part->GetLabel());
-      if (label>=nProduced) continue;
-      masks[label] |= GetFilterMap(part);
+    if (fParticles) {
+      UInt_t* masks = new UInt_t[nProduced];
+      for (Int_t i=0;i<nProduced;i++) masks[i]=0;
+      
+      // Loop over reconstructed tracks to set masks
+      for (Int_t ipart=0;ipart<event->GetNumberOfTracks();ipart++){
+        AliVTrack* part = (AliVTrack*) event->GetTrack(ipart);
+        if (!part) continue;
+        Int_t label = TMath::Abs(part->GetLabel());
+        if (label>=nProduced) continue;
+        masks[label] |= GetFilterMap(part);
+      }
+      
+      // Loop over mc tracks to be stored
+      for (Int_t ipart=0;ipart<nProduced;ipart++){
+        AliVParticle* part = 0;
+        Bool_t isPrimary = 0;
+        if (mcTracks) { // AOD analysis
+          AliAODMCParticle* particle = (AliAODMCParticle*) mcTracks->At(ipart);
+          if (!particle) continue;
+          // skip injected signals
+          AliAODMCParticle* mother = particle;
+          while (!mother->IsPrimary()) mother = (AliAODMCParticle*) mcTracks->At(mother->GetMother());
+          if (mother->GetLabel()>=nPrimGen) continue;
+          part = particle;
+          // check for primary
+          isPrimary = particle->IsPhysicalPrimary();
+        } else { // ESD analysis
+          AliMCParticle* particle = (AliMCParticle*) mcEvent->GetTrack(ipart);
+          if (!particle) continue;
+          // skip injected signals
+          AliMCParticle* mother = particle;
+          while (mother->GetMother()>=0) mother = (AliMCParticle*) mcEvent->GetTrack(mother->GetMother());
+          if (mother->GetLabel()>=nPrimGen) continue;
+          part = particle;
+          // check for primary
+          isPrimary = mcEvent->IsPhysicalPrimary(ipart);
+        }
+        
+        // store only primaries and all reconstructed (non-zero mask)
+        Int_t mask = masks[ipart];
+        if (isPrimary) mask |= (1 << 30);
+        if (!mask) continue; 
+        AddTrack(part,mask);
+      }
+      delete[] masks;
     }
     
-    // Loop over mc tracks to be stored
-    for (Int_t ipart=0;ipart<nProduced;ipart++){
-      AliVParticle* part = 0;
-      Bool_t isPrimary = 0;
-      if (mcTracks) { // AOD analysis
-        AliAODMCParticle* particle = (AliAODMCParticle*) mcTracks->At(ipart);
-        if (!particle) continue;
-        // skip injected signals
-        AliAODMCParticle* mother = particle;
-        while (!mother->IsPrimary()) mother = (AliAODMCParticle*) mcTracks->At(mother->GetMother());
-        if (mother->GetLabel()>=nPrimGen) continue;
-        part = particle;
-        // check for primary
-        isPrimary = particle->IsPhysicalPrimary();
-      } else { // ESD analysis
-        AliMCParticle* particle = (AliMCParticle*) mcEvent->GetTrack(ipart);
-        if (!particle) continue;
-        // skip injected signals
-        AliMCParticle* mother = particle;
-        while (mother->GetMother()>=0) mother = (AliMCParticle*) mcEvent->GetTrack(mother->GetMother());
-        if (mother->GetLabel()>=nPrimGen) continue;
-        part = particle;
-        // check for primary
-        isPrimary = mcEvent->IsPhysicalPrimary(ipart);
+    if (fIsAOD){
+      AliAODEvent* aod = (AliAODEvent*) event;
+    
+      if (fStoreTracklets) {
+        AliAODTracklets* tracklets = aod->GetTracklets();
+        Int_t nTracklets = tracklets->GetNumberOfTracklets();
+        for (Int_t i=0;i<nTracklets;i++){
+          Float_t phi   = tracklets->GetPhi(i);
+          Float_t theta = tracklets->GetTheta(i);
+          Float_t dphi  = tracklets->GetDeltaPhi(i);
+          AliCFParticle* tracklet = new ((*fTracklets)[fTracklets->GetEntriesFast()]) AliCFParticle(dphi,-TMath::Log(TMath::Tan(theta/2)),phi,0,0,4);
+          Int_t label1 = tracklets->GetLabel(i,0);
+          Int_t label2 = tracklets->GetLabel(i,1);
+          if (label1!=label2) continue;
+          AliAODMCParticle* particle = (AliAODMCParticle*) mcTracks->At(label1);
+          if (!particle) continue;
+          Short_t charge = particle->Charge();
+          Float_t ptMC   = particle->Pt();
+          Float_t etaMC  = particle->Eta();
+          Float_t phiMC  = particle->Phi();
+          Float_t pdg    = particle->PdgCode();
+          tracklet->SetCharge(charge);
+          tracklet->SetAt(ptMC,0);
+          tracklet->SetAt(etaMC,1);
+          tracklet->SetAt(phiMC,2);
+          tracklet->SetAt(pdg,3);
+        }
+      }
+
+      if (fStoreMuons){
+        for (Int_t iTrack = 0; iTrack < aod->GetNTracks(); iTrack++) {
+          AliAODTrack* track = aod->GetTrack(iTrack);
+          if (!track->IsMuonTrack()) continue;
+          Float_t pt     = track->Pt();
+          Float_t eta    = track->Eta();
+          Float_t phi    = track->Phi();
+          Short_t charge = track->Charge();
+          Float_t dca    = track->DCA();
+          Float_t chi2   = track->Chi2perNDF();
+          Float_t rabs   = track->GetRAtAbsorberEnd();
+          Int_t   mask   = track->GetMatchTrigger();
+          if (rabs < 17.6 || rabs > 89.5) continue;
+          if (eta < -4 || eta > -2.5) continue;
+          AliCFParticle* part = new ((*fMuons)[fMuons->GetEntriesFast()]) AliCFParticle(pt,eta,phi,charge,mask,3);
+          part->SetAt(dca,0);
+          part->SetAt(chi2,1);
+          part->SetAt(rabs,2);
+        }
       }
-      
-      // store only primaries and all reconstructed (non-zero mask)
-      Int_t mask = masks[ipart];
-      if (isPrimary) mask |= (1 << 30);
-      if (!mask) continue; 
-      AddTrack(part,mask);
     }
-    delete[] masks;
   }
   fTree->Fill();
   PostData(0,fListOfHistos);
index 154847932df6dd47a8ef64ea2713f297b777124c..7bc8ecffc18dd7504b4312dc7407017ade6b0ac4 100644 (file)
@@ -45,6 +45,7 @@ class AliAnalysisTaskCFTree : public AliAnalysisTask {
   void SetCrossedRowsCut(Int_t val)     { fCrossedRowsCut = val; }
   void SetFoundFractionCut(Float_t val) { fFoundFractionCut = val;  }
   // Switchers for additional data to be stored
+  void SetStoreTracks(Bool_t val=kTRUE)    { fStoreTracks    = val; }
   void SetStoreTracklets(Bool_t val=kTRUE) { fStoreTracklets = val; }
   void SetStoreMuons(Bool_t val=kTRUE)     { fStoreMuons     = val; }
 
@@ -92,6 +93,7 @@ class AliAnalysisTaskCFTree : public AliAnalysisTask {
   Int_t   fCrossedRowsCut;    // cut on crossed rows
   Float_t fFoundFractionCut;  // cut on crossed rows/findable clusters
   //
+  Bool_t fStoreTracks;        // if kTRUE - Barrel tracks will be stored as AliCFParticles
   Bool_t fStoreTracklets;     // if kTRUE - SPD tracklets will be stored as AliCFParticles
   Bool_t fStoreMuons;         // if kTRUE - muon tracks will be stored as AliCFParticles
   ClassDef(AliAnalysisTaskCFTree,2);
index 3defe7ce883b29cb55502a42093250a18b578d9f..10423ffc0183ef7ba3cbc7f24883c322d31f36c7 100644 (file)
@@ -52,7 +52,11 @@ class AliCFParticle : public AliVParticle, public TArrayF {
 
   virtual Bool_t IsEqual(const TObject* obj) const { return (obj->GetUniqueID() == GetUniqueID()); }
 
-  virtual void SetPhi(Double_t phi) { fPhi = phi; }
+  virtual void SetPt(Double_t pt)        { fPt     = pt;     }
+  virtual void SetEta(Double_t eta)      { fEta    = eta;    }
+  virtual void SetPhi(Double_t phi)      { fPhi    = phi;    }
+  virtual void SetCharge(Short_t charge) { fCharge = charge; }
+  virtual void SetMask(Short_t mask)     { fMask   = mask;   }
  protected:
   Float_t fPt;
   Float_t fEta;