]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
flag jets with too high p_T track
authorkleinb <kleinb@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 7 Sep 2011 12:59:16 +0000 (12:59 +0000)
committerkleinb <kleinb@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 7 Sep 2011 12:59:16 +0000 (12:59 +0000)
JETAN/AliAnalysisTaskJetCluster.cxx
JETAN/AliAnalysisTaskJetCluster.h
PWG4/JetTasks/AliAnalysisTaskJetSpectrum2.cxx
PWG4/JetTasks/AliAnalysisTaskJetSpectrum2.h
PWG4/macros/AddTaskESDFilterPWG4Train.C
PWG4/macros/ConfigTrainGrid.C
PWG4/macros/CreateTrackCutsPWG4.C
STEER/AOD/AliAODJet.h

index 8cfc8e30edb8fefc77e789c3f9e8a7fa223355b2..22c7b734bb78af5de858d2f047cc8ca1862bfb41 100644 (file)
@@ -106,6 +106,7 @@ AliAnalysisTaskJetCluster::AliAnalysisTaskJetCluster():
   fRecEtaWindow(0.5),
   fTrackPtCut(0.),                                                     
   fJetOutputMinPt(0.150),
   fRecEtaWindow(0.5),
   fTrackPtCut(0.),                                                     
   fJetOutputMinPt(0.150),
+  fMaxTrackPtInJet(100.),
   fJetTriggerPtCut(0),
   fVtxZCut(8),
   fVtxR2Cut(1),
   fJetTriggerPtCut(0),
   fVtxZCut(8),
   fVtxR2Cut(1),
@@ -218,6 +219,7 @@ AliAnalysisTaskJetCluster::AliAnalysisTaskJetCluster(const char* name):
   fRecEtaWindow(0.5),
   fTrackPtCut(0.),                                                     
   fJetOutputMinPt(0.150),
   fRecEtaWindow(0.5),
   fTrackPtCut(0.),                                                     
   fJetOutputMinPt(0.150),
+  fMaxTrackPtInJet(100.),
   fJetTriggerPtCut(0),
   fVtxZCut(8),
   fVtxR2Cut(1),
   fJetTriggerPtCut(0),
   fVtxZCut(8),
   fVtxR2Cut(1),
@@ -934,6 +936,7 @@ void AliAnalysisTaskJetCluster::UserExec(Option_t */*option*/)
        fh1PtJetConstRec->Fill(part->Pt());
        if(aodOutJet){
          aodOutJet->AddTrack(fRef->At(constituents[ic].user_index()));
        fh1PtJetConstRec->Fill(part->Pt());
        if(aodOutJet){
          aodOutJet->AddTrack(fRef->At(constituents[ic].user_index()));
+         if(part->Pt()>fMaxTrackPtInJet)aodOutJet->SetTrigger(AliAODJet::kHighTrackPtTriggered);
        }
        if(j==0)fh1PtJetConstLeadingRec->Fill(part->Pt());
       }
        }
        if(j==0)fh1PtJetConstLeadingRec->Fill(part->Pt());
       }
@@ -1020,6 +1023,7 @@ void AliAnalysisTaskJetCluster::UserExec(Option_t */*option*/)
         for(int ir = 0;ir < fNRandomCones;ir++){
           AliAODJet *jC = (AliAODJet*)fTCARandomConesOut->At(ir);  
           if(jC&&jC->DeltaR(vp)<fRparam){
         for(int ir = 0;ir < fNRandomCones;ir++){
           AliAODJet *jC = (AliAODJet*)fTCARandomConesOut->At(ir);  
           if(jC&&jC->DeltaR(vp)<fRparam){
+            if(vp->Pt()>fMaxTrackPtInJet)jC->SetTrigger(AliAODJet::kHighTrackPtTriggered);
             jC->SetBgEnergy(jC->ChargedBgEnergy()+vp->Pt(),0);
           }
         }  
             jC->SetBgEnergy(jC->ChargedBgEnergy()+vp->Pt(),0);
           }
         }  
@@ -1042,6 +1046,7 @@ void AliAnalysisTaskJetCluster::UserExec(Option_t */*option*/)
           for(int ir = 0;ir < fTCARandomConesOutRan->GetEntriesFast();ir++){
             AliAODJet *jC = (AliAODJet*)fTCARandomConesOutRan->At(ir);  
             if(jC&&jC->DeltaR(&vTmpRanR)<fRparam){
           for(int ir = 0;ir < fTCARandomConesOutRan->GetEntriesFast();ir++){
             AliAODJet *jC = (AliAODJet*)fTCARandomConesOutRan->At(ir);  
             if(jC&&jC->DeltaR(&vTmpRanR)<fRparam){
+              if(vTmpRanR.Pt()>fMaxTrackPtInJet)jC->SetTrigger(AliAODJet::kHighTrackPtTriggered);
               jC->SetBgEnergy(jC->ChargedBgEnergy()+vTmpRanR.Pt(),0);
             }
           }  
               jC->SetBgEnergy(jC->ChargedBgEnergy()+vTmpRanR.Pt(),0);
             }
           }  
index f11c0ac421d969c773fb487adacd245778de4be3..493d1d14cd60879068c30ce49f8cd8f87f3527a4 100644 (file)
@@ -74,9 +74,12 @@ class AliAnalysisTaskJetCluster : public AliAnalysisTaskSE
     virtual const char* GetJetOutputBranch(){return fNonStdBranch.Data();}
     virtual void SetJetOutputFile(const char *c){fNonStdFile = c;}
     virtual const char* GetJetOutputFile(){return fNonStdFile.Data();}
     virtual const char* GetJetOutputBranch(){return fNonStdBranch.Data();}
     virtual void SetJetOutputFile(const char *c){fNonStdFile = c;}
     virtual const char* GetJetOutputFile(){return fNonStdFile.Data();}
+    virtual void SetMaxTrackPtInJet(Float_t x){fMaxTrackPtInJet = x;}
     virtual void SetJetOutputMinPt(Float_t x){fJetOutputMinPt = x;}
     virtual void SetBackgroundCalc(Bool_t b){fUseBackgroundCalc = b;} 
 
     virtual void SetJetOutputMinPt(Float_t x){fJetOutputMinPt = x;}
     virtual void SetBackgroundCalc(Bool_t b){fUseBackgroundCalc = b;} 
 
+
+
     // for Fast Jet
     fastjet::JetAlgorithm        GetAlgorithm()         const {return fAlgorithm;}
     fastjet::Strategy            GetStrategy()          const {return fStrategy;}
     // for Fast Jet
     fastjet::JetAlgorithm        GetAlgorithm()         const {return fAlgorithm;}
     fastjet::Strategy            GetStrategy()          const {return fStrategy;}
@@ -135,6 +138,7 @@ class AliAnalysisTaskJetCluster : public AliAnalysisTaskSE
     Float_t       fRecEtaWindow;          // eta window used for corraltion plots between rec and gen 
     Float_t       fTrackPtCut;            // minimum track pt to be accepted
     Float_t       fJetOutputMinPt;        // minimum p_t for jets to be written out
     Float_t       fRecEtaWindow;          // eta window used for corraltion plots between rec and gen 
     Float_t       fTrackPtCut;            // minimum track pt to be accepted
     Float_t       fJetOutputMinPt;        // minimum p_t for jets to be written out
+    Float_t       fMaxTrackPtInJet;       // maximum track pt within a jet for flagging...
     Float_t       fJetTriggerPtCut;       // minimum jwt pT for AOD to be written
     Float_t       fVtxZCut;               // zvtx cut
     Float_t       fVtxR2Cut;              // R vtx cut (squared) 
     Float_t       fJetTriggerPtCut;       // minimum jwt pT for AOD to be written
     Float_t       fVtxZCut;               // zvtx cut
     Float_t       fVtxR2Cut;              // R vtx cut (squared) 
@@ -235,7 +239,7 @@ class AliAnalysisTaskJetCluster : public AliAnalysisTaskSE
     TList *fHistList; //!leading tracks to be skipped in the randomized event Output list
    
 
     TList *fHistList; //!leading tracks to be skipped in the randomized event Output list
    
 
-    ClassDef(AliAnalysisTaskJetCluster, 18
+    ClassDef(AliAnalysisTaskJetCluster, 19
 };
  
 #endif
 };
  
 #endif
index 49fcf9007a4ff851891b6aa1f2583ae0a5fb845b..1a162dc85fe9d865fd873479ee3afa39f7214f99 100644 (file)
@@ -93,6 +93,7 @@ AliAnalysisTaskJetSpectrum2::AliAnalysisTaskJetSpectrum2():
   fDoMatching(kFALSE),
   fNMatchJets(5),
   fNRPBins(3),
   fDoMatching(kFALSE),
   fNMatchJets(5),
   fNRPBins(3),
+  fJetTriggerExcludeMask(AliAODJet::kHighTrackPtTriggered),
   fFilterMask(0),
   fEventSelectionMask(0),
   fAnalysisType(0),
   fFilterMask(0),
   fEventSelectionMask(0),
   fAnalysisType(0),
@@ -133,6 +134,7 @@ AliAnalysisTaskJetSpectrum2::AliAnalysisTaskJetSpectrum2():
     fh1NJets[ij] = 0;
     fh1SumPtTrack[ij] = 0;
     fh1PtJetsIn[ij] = 0;
     fh1NJets[ij] = 0;
     fh1SumPtTrack[ij] = 0;
     fh1PtJetsIn[ij] = 0;
+    fh1PtJetsInRej[ij] = 0;
     fh1PtTracksIn[ij] = 0;
     fh1PtTracksInLow[ij] = 0;
     fh2NJetsPt[ij]  = 0;
     fh1PtTracksIn[ij] = 0;
     fh1PtTracksInLow[ij] = 0;
     fh2NJetsPt[ij]  = 0;
@@ -183,6 +185,7 @@ AliAnalysisTaskJetSpectrum2::AliAnalysisTaskJetSpectrum2(const char* name):
   fDoMatching(kFALSE),
   fNMatchJets(5),
   fNRPBins(3),
   fDoMatching(kFALSE),
   fNMatchJets(5),
   fNRPBins(3),
+  fJetTriggerExcludeMask(AliAODJet::kHighTrackPtTriggered),
   fFilterMask(0),
   fEventSelectionMask(0),
   fAnalysisType(0),
   fFilterMask(0),
   fEventSelectionMask(0),
   fAnalysisType(0),
@@ -223,6 +226,7 @@ AliAnalysisTaskJetSpectrum2::AliAnalysisTaskJetSpectrum2(const char* name):
     fh1NJets[ij] = 0;
     fh1SumPtTrack[ij] = 0;
     fh1PtJetsIn[ij] = 0;
     fh1NJets[ij] = 0;
     fh1SumPtTrack[ij] = 0;
     fh1PtJetsIn[ij] = 0;
+    fh1PtJetsInRej[ij] = 0;
     fh1PtTracksIn[ij] = 0;
     fh1PtTracksInLow[ij] = 0;
     fh2NJetsPt[ij]  = 0;
     fh1PtTracksIn[ij] = 0;
     fh1PtTracksInLow[ij] = 0;
     fh2NJetsPt[ij]  = 0;
@@ -444,6 +448,9 @@ void AliAnalysisTaskJetSpectrum2::UserCreateOutputObjects()
     
     fh1PtJetsIn[ij]  = new TH1F(Form("fh1PtJets%sIn",cAdd.Data()),Form("%s jets p_T;p_{T} (GeV/c)",cAdd.Data()),nBinPt,binLimitsPt);
     fHistList->Add(fh1PtJetsIn[ij]);
     
     fh1PtJetsIn[ij]  = new TH1F(Form("fh1PtJets%sIn",cAdd.Data()),Form("%s jets p_T;p_{T} (GeV/c)",cAdd.Data()),nBinPt,binLimitsPt);
     fHistList->Add(fh1PtJetsIn[ij]);
+
+    fh1PtJetsInRej[ij]  = new TH1F(Form("fh1PtJets%sInRej",cAdd.Data()),Form("%s jets p_T;p_{T} (GeV/c)",cAdd.Data()),nBinPt,binLimitsPt);
+    fHistList->Add(fh1PtJetsInRej[ij]);
     
     fh1PtTracksIn[ij] = new TH1F(Form("fh1PtTracks%sIn",cAdd.Data()),Form("%s track p_T;p_{T} (GeV/c)",cAdd.Data()),nBinPt,binLimitsPt);
     fHistList->Add(fh1PtTracksIn[ij]);
     
     fh1PtTracksIn[ij] = new TH1F(Form("fh1PtTracks%sIn",cAdd.Data()),Form("%s track p_T;p_{T} (GeV/c)",cAdd.Data()),nBinPt,binLimitsPt);
     fHistList->Add(fh1PtTracksIn[ij]);
@@ -887,6 +894,10 @@ void AliAnalysisTaskJetSpectrum2::FillJetHistos(TList &jetsList,TList &particles
   for(int ij = 0;ij < nJets;ij++){
     AliAODJet *jet = (AliAODJet*)jetsList.At(ij);
     Float_t ptJet = jet->Pt();
   for(int ij = 0;ij < nJets;ij++){
     AliAODJet *jet = (AliAODJet*)jetsList.At(ij);
     Float_t ptJet = jet->Pt();
+    if(jet->Trigger()&fJetTriggerExcludeMask){
+      fh1PtJetsInRej[iType]->Fill(ptJet);
+      continue;
+    }
     fh1PtJetsIn[iType]->Fill(ptJet);
     if(ptJet>ptOld){
       Printf("%s:%d Jets Type %d Not Sorted !! %d:%.3E %d:%.3E",(char*)__FILE__,__LINE__,iType,ij,ptJet,ij-1,ptOld);
     fh1PtJetsIn[iType]->Fill(ptJet);
     if(ptJet>ptOld){
       Printf("%s:%d Jets Type %d Not Sorted !! %d:%.3E %d:%.3E",(char*)__FILE__,__LINE__,iType,ij,ptJet,ij-1,ptOld);
index 9330fdf6cb7bedc6a1f1d693a311c57fcffa5e78..8888bd0ca5fcdd601b73b107bba2f1c29472d48a 100644 (file)
@@ -74,6 +74,7 @@ class AliAnalysisTaskJetSpectrum2 : public AliAnalysisTaskSE
     virtual void SetTrackTypeGen(Int_t i){fTrackTypeGen = i;}
     virtual void SetTrackTypeRec(Int_t i){fTrackTypeRec = i;}
     virtual void SetFilterMask(UInt_t i){fFilterMask = i;}
     virtual void SetTrackTypeGen(Int_t i){fTrackTypeGen = i;}
     virtual void SetTrackTypeRec(Int_t i){fTrackTypeRec = i;}
     virtual void SetFilterMask(UInt_t i){fFilterMask = i;}
+    virtual void SetJetTriggerExclude(Char_t i){fJetTriggerExcludeMask = i;}
     virtual void SetMatching(Bool_t b = kTRUE){fDoMatching = b;}
     virtual void SetRPMethod(Int_t i){fRPMethod = i;}
     virtual void SetEventSelectionMask(UInt_t i){fEventSelectionMask = i;}
     virtual void SetMatching(Bool_t b = kTRUE){fDoMatching = b;}
     virtual void SetRPMethod(Int_t i){fRPMethod = i;}
     virtual void SetEventSelectionMask(UInt_t i){fEventSelectionMask = i;}
@@ -153,6 +154,7 @@ class AliAnalysisTaskJetSpectrum2 : public AliAnalysisTaskSE
     Bool_t        fDoMatching;            // switch on the matching between rec and gen
     Short_t       fNMatchJets;            // number of leading jets considered from the list
     Short_t       fNRPBins;               // number of bins with respect to RP
     Bool_t        fDoMatching;            // switch on the matching between rec and gen
     Short_t       fNMatchJets;            // number of leading jets considered from the list
     Short_t       fNRPBins;               // number of bins with respect to RP
+    Char_t        fJetTriggerExcludeMask; // mask for jet triggers to exclude
     UInt_t        fFilterMask;            // filter bit for slecected tracks
     UInt_t        fEventSelectionMask;    // Selection information used to filter events
     Int_t         fAnalysisType;          // Analysis type 
     UInt_t        fFilterMask;            // filter bit for slecected tracks
     UInt_t        fEventSelectionMask;    // Selection information used to filter events
     Int_t         fAnalysisType;          // Analysis type 
@@ -199,6 +201,7 @@ class AliAnalysisTaskJetSpectrum2 : public AliAnalysisTaskSE
 
     TH1F*         fh1PtIn[kJetTypes][kMaxJets+1];  //! Jet pt  
     TH1F*         fh1PtJetsIn[kJetTypes];       //! Jet pt for all jets
 
     TH1F*         fh1PtIn[kJetTypes][kMaxJets+1];  //! Jet pt  
     TH1F*         fh1PtJetsIn[kJetTypes];       //! Jet pt for all jets
+    TH1F*         fh1PtJetsInRej[kJetTypes];    //! Jet pt for all rejected jets
     TH1F*         fh1PtTracksIn[kJetTypes];     //! track pt for all tracks
     TH1F*         fh1PtTracksInLow[kJetTypes];  //! track pt for all tracks
     
     TH1F*         fh1PtTracksIn[kJetTypes];     //! track pt for all tracks
     TH1F*         fh1PtTracksInLow[kJetTypes];  //! track pt for all tracks
     
@@ -223,7 +226,7 @@ class AliAnalysisTaskJetSpectrum2 : public AliAnalysisTaskSE
     TList *fHistList;                  //! Output list
    
 
     TList *fHistList;                  //! Output list
    
 
-    ClassDef(AliAnalysisTaskJetSpectrum2, 17) // Analysis task for standard jet analysis
+    ClassDef(AliAnalysisTaskJetSpectrum2, 18) // Analysis task for standard jet analysis
 };
  
 #endif
 };
  
 #endif
index 480a18192dad27a5a0610f0e275ad95884b430aa..ceb8e3170b27a7c9551ae0eee6fc4b86cc11c79a 100644 (file)
@@ -105,7 +105,7 @@ AliAnalysisTaskESDfilter *AddTaskESDFilter(Bool_t useKineFilter=kTRUE,
 
    // ITS cuts for new jet analysis 
    gROOT->LoadMacro("$ALICE_ROOT/PWG4/macros/CreateTrackCutsPWG4.C");
 
    // ITS cuts for new jet analysis 
    gROOT->LoadMacro("$ALICE_ROOT/PWG4/macros/CreateTrackCutsPWG4.C");
-   AliESDtrackCuts* esdTrackCutsHG0 = CreateTrackCutsPWG4(10001004);
+   AliESDtrackCuts* esdTrackCutsHG0 = CreateTrackCutsPWG4(10001005);
 
    // throw out tracks with too low number of clusters in
    // the first pass (be consistent with TPC only tracks)
 
    // throw out tracks with too low number of clusters in
    // the first pass (be consistent with TPC only tracks)
@@ -117,11 +117,11 @@ AliAnalysisTaskESDfilter *AddTaskESDFilter(Bool_t useKineFilter=kTRUE,
 
 
    // the complement to the one with SPD requirement
 
 
    // the complement to the one with SPD requirement
-   AliESDtrackCuts* esdTrackCutsHG1 = CreateTrackCutsPWG4(10011004);
+   AliESDtrackCuts* esdTrackCutsHG1 = CreateTrackCutsPWG4(10011005);
 
    // the tracks that must not be taken pass this cut and
    // non HGC1 and HG
 
    // the tracks that must not be taken pass this cut and
    // non HGC1 and HG
-   AliESDtrackCuts* esdTrackCutsHG2 = CreateTrackCutsPWG4(10021004);
+   AliESDtrackCuts* esdTrackCutsHG2 = CreateTrackCutsPWG4(10021005);
 
    
 
 
    
 
@@ -132,7 +132,7 @@ AliAnalysisTaskESDfilter *AddTaskESDFilter(Bool_t useKineFilter=kTRUE,
    esdTrackCutsH2->SetMaxChi2PerClusterITS(36.);
    esdTrackCutsH2->SetPtRange(0.15,1E10);
 
    esdTrackCutsH2->SetMaxChi2PerClusterITS(36.);
    esdTrackCutsH2->SetPtRange(0.15,1E10);
 
-   AliESDtrackCuts* esdTrackCutsGCOnly = CreateTrackCutsPWG4(10041004);
+   AliESDtrackCuts* esdTrackCutsGCOnly = CreateTrackCutsPWG4(10041005);
 
    // TPC only tracks
    AliESDtrackCuts* esdTrackCutsTPCCOnly = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
 
    // TPC only tracks
    AliESDtrackCuts* esdTrackCutsTPCCOnly = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
index 4035b0b930b910866d0e4dcf864821f888008986..9076fd350b258943c6366a7cafc8f0c8a24b8ecc 100644 (file)
@@ -56,7 +56,7 @@
   if (kPluginMode.Contains("merge")){
     // currently merging this one...
     //       cDate = "110717a";
   if (kPluginMode.Contains("merge")){
     // currently merging this one...
     //       cDate = "110717a";
-    bRun = 802; Int_t bExtra = 0; cDate = "110812a";
+    //    bRun = 802; Int_t bExtra = 0; cDate = "110812a";
   }
   kUseDebug = kFALSE;
   // this is for testing just one run...
   }
   kUseDebug = kFALSE;
   // this is for testing just one run...
index 149f4da08212490857e2c1e8303aaa97d37efba3..e3dc7fb04577ca1e7328e84bca7c76eae1d31d7b 100644 (file)
@@ -145,6 +145,34 @@ AliESDtrackCuts *CreateTrackCutsPWG4(Int_t cutMode) {
     tag = "Global tracks jet analysis with ITSrefit and NclsIter1=70, noSPD requirement";
 
   }
     tag = "Global tracks jet analysis with ITSrefit and NclsIter1=70, noSPD requirement";
 
   }
+  if(stdCutMode == 1005) {
+
+    bStdCutsDefined = kTRUE;
+
+    // TPC  
+    trackCuts->SetMinNClustersTPC(70);
+    trackCuts->SetMaxChi2PerClusterTPC(4);
+    trackCuts->SetRequireTPCStandAlone(kTRUE); //cut on NClustersTPC and chi2TPC Iter1
+    trackCuts->SetAcceptKinkDaughters(kFALSE);
+    trackCuts->SetRequireTPCRefit(kTRUE);
+    trackCuts->SetMaxFractionSharedTPCClusters(0.4);
+    // ITS
+    trackCuts->SetRequireITSRefit(kTRUE);
+    //accept secondaries
+    trackCuts->SetMaxDCAToVertexXY(2.4);
+    trackCuts->SetMaxDCAToVertexZ(3.2);
+    trackCuts->SetDCAToVertex2D(kTRUE);
+    //reject fakes
+    trackCuts->SetMaxChi2PerClusterITS(36);
+
+    trackCuts->SetRequireSigmaToVertex(kFALSE);
+
+    trackCuts->SetEtaRange(-0.9,0.9);
+    trackCuts->SetPtRange(0.15, 1E+15.);
+    tag = "Global tracks jet analysis with ITSrefit and NclsIter1=70, noSPD requirement, no upper pt cut";
+
+  }
 
 
 
 
 
 
index 9f337d5951f65c5e793ec3c6c675330042444e0a..5505f48055ff7725594085ca5b40b20d46b6bacc 100644 (file)
@@ -100,7 +100,8 @@ class AliAODJet : public AliVParticle {
     // first only one bit for EMCAL and TRD, leave space for more
     // trigger types and/or other detectors
     enum {kEMCALTriggered = 1<<0,
     // first only one bit for EMCAL and TRD, leave space for more
     // trigger types and/or other detectors
     enum {kEMCALTriggered = 1<<0,
-         kTRDTriggered = 4<<0};
+         kTRDTriggered =   1<<2,
+         kHighTrackPtTriggered = 1<<7};
 
 
  private:
 
 
  private: