]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
updates from pdsf train
authorloizides <loizides@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 21 May 2012 10:18:53 +0000 (10:18 +0000)
committerloizides <loizides@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 21 May 2012 10:18:53 +0000 (10:18 +0000)
PWGGA/EMCALJetTasks/AliAnalysisTaskScale.cxx
PWGGA/EMCALJetTasks/AliAnalysisTaskScale.h

index 19eaeb2b671a9fe892b87f2692faa060bfa0206e..829ae50875eb15caf3ea304811442f4d44e323e4 100644 (file)
@@ -39,7 +39,13 @@ AliAnalysisTaskScale::AliAnalysisTaskScale(const char *name) :
   fHistPtEMCALvsNtrack(0), 
   fHistEtvsNtrack(0),  
   fHistScalevsNtrack(0),  
-  fHistDeltaScalevsNtrack(0)
+  fHistDeltaScalevsNtrack(0),
+  fHistTrackPtvsCent(0), 
+  fHistClusterPtvsCent(0),
+  fHistTrackEtaPhi(0), 
+  fHistClusterEtaPhi(0),
+  fMinTrackPt(0.15),
+  fMinClusterPt(0.15)
 {
   // Constructor.
 
@@ -55,20 +61,22 @@ void AliAnalysisTaskScale::UserCreateOutputObjects()
   fOutputList = new TList();
   fOutputList->SetOwner();
 
-  fHistCentrality         = new TH1F("fHistCentrality","Centrality",         101, -1, 100);
-  fHistPtTPCvsCent        = new TH2F("fHistPtTPCvsCent","rho vs cent",       101, -1, 100, 500,   0, 1000);
-  fHistPtEMCALvsCent      = new TH2F("fHistPtEMCALvsCent","rho vs cent",     101, -1, 100, 500,   0, 1000);
-  fHistEtvsCent           = new TH2F("fHistEtvsCent","rho vs cent",          101, -1, 100, 500,   0, 1000);
-  fHistScalevsCent        = new TH2F("fHistScalevsCent","rho vs cent",       101, -1, 100, 400,   0, 4);
-  fHistDeltaScalevsCent   = new TH2F("fHistDeltaScalevsCent","rho vs cent",  101, -1, 100, 400,  -2, 2);
-  fHistPtTPCvsNtrack      = new TH2F("fHistPtTPCvsNtrack","rho vs cent",     500,  0, 2500, 500,  0, 1000);
-  fHistPtEMCALvsNtrack    = new TH2F("fHistPtEMCALvsNtrack","rho vs cent",   500,  0, 2500, 500,  0, 1000);
-  fHistEtvsNtrack         = new TH2F("fHistEtvsNtrack","rho vs cent",        500,  0, 2500, 500,  0, 1000);
-  fHistScalevsNtrack      = new TH2F("fHistScalevsNtrack","rho vs cent",     500,  0, 2500, 400,  0, 4);
-  fHistDeltaScalevsNtrack = new TH2F("fHistDeltaScalevsNtrack","rho vs cent",500,  0, 2500, 400, -2, 2);
-
-  fNewScaleFunction       = new TF1("sfunc","pol2", -1, 100);
-    
+  fHistCentrality         = new TH1F("Centrality","Centrality",              101, -1, 100);
+  fHistPtTPCvsCent        = new TH2F("PtTPCvsCent","rho vs cent",            101, -1, 100, 500,   0, 1000);
+  fHistPtEMCALvsCent      = new TH2F("PtEMCALvsCent","rho vs cent",          101, -1, 100, 500,   0, 1000);
+  fHistEtvsCent           = new TH2F("EtvsCent","rho vs cent",               101, -1, 100, 500,   0, 1000);
+  fHistScalevsCent        = new TH2F("ScalevsCent","rho vs cent",            101, -1, 100, 400,   0, 4);
+  fHistDeltaScalevsCent   = new TH2F("DeltaScalevsCent","rho vs cent",       101, -1, 100, 400,  -2, 2);
+  fHistPtTPCvsNtrack      = new TH2F("PtTPCvsNtrack","rho vs cent",          500,  0, 2500, 500,  0, 1000);
+  fHistPtEMCALvsNtrack    = new TH2F("PtEMCALvsNtrack","rho vs cent",        500,  0, 2500, 500,  0, 1000);
+  fHistEtvsNtrack         = new TH2F("EtvsNtrack","rho vs cent",             500,  0, 2500, 500,  0, 1000);
+  fHistScalevsNtrack      = new TH2F("ScalevsNtrack","rho vs cent",          500,  0, 2500, 400,  0, 4);
+  fHistDeltaScalevsNtrack = new TH2F("DeltaScalevsNtrack","rho vs cent",     500,  0, 2500, 400, -2, 2);
+  fHistTrackPtvsCent      = new TH2F("TrackPtvsCent","Track pt vs cent",     101, -1, 100,  500,  0, 100);
+  fHistClusterPtvsCent    = new TH2F("ClusterPtvsCent","Cluster pt vs cent", 101, -1, 100,  500,  0, 100);
+  fHistTrackEtaPhi        = new TH2F("TrackEtaPhi","Track eta phi",          100, -1.0, 1.0, 64,  0, 6.4);
+  fHistClusterEtaPhi      = new TH2F("ClusterEtaPhi","Cluster eta phi",      100, -1.0, 1.0, 64, -3.2, 3.2);
+
   fOutputList->Add(fHistCentrality);
   fOutputList->Add(fHistPtTPCvsCent);
   fOutputList->Add(fHistPtEMCALvsCent);
@@ -80,7 +88,10 @@ void AliAnalysisTaskScale::UserCreateOutputObjects()
   fOutputList->Add(fHistEtvsNtrack);
   fOutputList->Add(fHistScalevsNtrack);
   fOutputList->Add(fHistDeltaScalevsNtrack);
-  //fOutputList->Add(fNewScaleFunction);
+  fOutputList->Add(fHistTrackPtvsCent);
+  fOutputList->Add(fHistClusterPtvsCent);
+  fOutputList->Add(fHistTrackEtaPhi);
+  fOutputList->Add(fHistClusterEtaPhi);
 
   PostData(1, fOutputList);
 }
@@ -151,6 +162,12 @@ void AliAnalysisTaskScale::UserExec(Option_t *)
     if (TMath::Abs(track->Eta()) > 0.7)   // only accept tracks in the EMCal eta range
       continue;
 
+    fHistTrackPtvsCent->Fill(cent,track->Pt());
+    fHistTrackEtaPhi->Fill(track->Eta(),track->Phi());
+
+    if (track->Pt()< fMinTrackPt)
+      continue;
+
     ptTPC += track->Pt();
     if ((track->Phi() > EmcalMaxPhi) || (track->Phi() < EmcalMinPhi))
       continue;
@@ -168,6 +185,13 @@ void AliAnalysisTaskScale::UserExec(Option_t *)
       continue;
     TLorentzVector nPart;
     c->GetMomentum(nPart, vertex);
+
+    fHistClusterPtvsCent->Fill(cent,nPart.Pt());
+    fHistClusterEtaPhi->Fill(nPart.Eta(),nPart.Phi());
+
+    if (nPart.Pt()< fMinClusterPt)
+      continue;
+
     Et += nPart.Pt();
   }
   
index 7e90411e535d3e2dbcd7acf79d5ed870dc286845..907119d8d380891a75dbb58b6aecd3c2363609a2 100644 (file)
@@ -12,10 +12,12 @@ class TF1;
 
 class AliAnalysisTaskScale : public AliAnalysisTaskSE {
  public:
-    AliAnalysisTaskScale() : AliAnalysisTaskSE(), fTracksName(), fClustersName(), fScaleFunction(0),
-      fOutputList(0), fHistCentrality(0), fHistPtTPCvsCent(0), fHistPtEMCALvsCent(0), fHistEtvsCent(0),  
-      fHistScalevsCent(0),  fHistDeltaScalevsCent(0), fHistPtTPCvsNtrack(0), fHistPtEMCALvsNtrack(0), 
-      fHistEtvsNtrack(0), fHistScalevsNtrack(0), fHistDeltaScalevsNtrack(0) {}
+  AliAnalysisTaskScale() : AliAnalysisTaskSE(), fTracksName(), fClustersName(), fScaleFunction(0),
+    fOutputList(0), fHistCentrality(0), fHistPtTPCvsCent(0), fHistPtEMCALvsCent(0), fHistEtvsCent(0),  
+    fHistScalevsCent(0),  fHistDeltaScalevsCent(0), fHistPtTPCvsNtrack(0), fHistPtEMCALvsNtrack(0), 
+    fHistEtvsNtrack(0), fHistScalevsNtrack(0), fHistDeltaScalevsNtrack(0), fHistTrackPtvsCent(0), 
+    fHistClusterPtvsCent(0), fHistTrackEtaPhi(0), fHistClusterEtaPhi(0),  fMinTrackPt(0.15), 
+    fMinClusterPt(0.15) {}
   AliAnalysisTaskScale(const char *name);
   virtual ~AliAnalysisTaskScale() {}
   
@@ -26,6 +28,8 @@ class AliAnalysisTaskScale : public AliAnalysisTaskSE {
   void                   SetTracksName(const char *n)                          { fTracksName    = n    ; }
   void                   SetClustersName(const char *n)                        { fClustersName  = n    ; }
   void                   SetScaleFunction(TF1* sf)                             { fScaleFunction = sf   ; }
+  void                   SetMinTrackPt(Double_t min)                           { fMinTrackPt    = min  ; }
+  void                   SetMinClusterPt(Double_t min)                         { fMinClusterPt  = min  ; }
   
  protected:
   virtual Double_t       GetScaleFactor(Double_t cent);
@@ -47,10 +51,16 @@ class AliAnalysisTaskScale : public AliAnalysisTaskSE {
   TH2F                  *fHistEtvsNtrack;         //!output histogram
   TH2F                  *fHistScalevsNtrack;      //!output histogram
   TH2F                  *fHistDeltaScalevsNtrack; //!output histogram
+  TH2F                  *fHistTrackPtvsCent;      //!output histogram
+  TH2F                  *fHistClusterPtvsCent;    //!output histogram
+  TH2F                  *fHistTrackEtaPhi;        //!output histogram
+  TH2F                  *fHistClusterEtaPhi;      //!output histogram
+  Double_t               fMinTrackPt;             //pt cut for scale factor calculation
+  Double_t               fMinClusterPt;           //pt cut for scale factor calculation
 
   AliAnalysisTaskScale(const AliAnalysisTaskScale&); // not implemented
   AliAnalysisTaskScale& operator=(const AliAnalysisTaskScale&); // not implemented
   
-  ClassDef(AliAnalysisTaskScale, 3); // Scale task
+  ClassDef(AliAnalysisTaskScale, 4); // Scale task
 };
 #endif