]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
add bookkeep profile
authormverweij <marta.verweij@cern.ch>
Wed, 17 Dec 2014 12:35:07 +0000 (13:35 +0100)
committermverweij <marta.verweij@cern.ch>
Wed, 17 Dec 2014 12:35:28 +0000 (13:35 +0100)
PWGJE/EMCALJetTasks/AliEmcalJetByJetCorrection.cxx
PWGJE/EMCALJetTasks/AliEmcalJetByJetCorrection.h
PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskEmcalJetMassStructure.cxx
PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskEmcalJetMassStructure.h

index a91e3494bd0dd9894b7316d31e55969982e67bd4..d381e603f2f29ee005a89a04850fead7b17e39f1 100644 (file)
@@ -6,6 +6,7 @@
 
 #include "TH2.h"
 #include "TH3.h"
+#include "TProfile.h"
 #include "TMath.h"
 #include "TRandom3.h"
 #include "TLorentzVector.h"
@@ -26,7 +27,8 @@ TNamed(),
   fCollTemplates(),
   fInitialized(kFALSE),
   fEfficiencyFixed(1.),
-  fhEfficiency(0)
+  fhEfficiency(0),
+  fpAppliedEfficiency(0)
 {
   // Dummy constructor.
   fCollTemplates.SetOwner(kTRUE);
@@ -42,10 +44,22 @@ AliEmcalJetByJetCorrection::AliEmcalJetByJetCorrection(const char* name) :
   fCollTemplates(),
   fInitialized(kFALSE),
   fEfficiencyFixed(1.),
-  fhEfficiency(0)
+  fhEfficiency(0),
+  fpAppliedEfficiency(0)
 {
   // Default constructor.
   fCollTemplates.SetOwner(kTRUE);
+
+  const Int_t nBinPt = 118; //0-2: 20 bins; 2-100: 98 bins
+  Double_t binLimitsPt[nBinPt+1];
+  for(Int_t iPt = 0;iPt <= nBinPt;iPt++){
+    if(iPt<20){
+      binLimitsPt[iPt] = 0. + (Double_t)iPt*0.15;
+    } else {// 1.0
+      binLimitsPt[iPt] =  binLimitsPt[iPt-1] + 1.0;
+    }
+  }
+  fpAppliedEfficiency = new TProfile("fpAppliedEfficiency","fpAppliedEfficiency",nBinPt,binLimitsPt);
 }
 
 //__________________________________________________________________________
@@ -58,7 +72,8 @@ AliEmcalJetByJetCorrection::AliEmcalJetByJetCorrection(const AliEmcalJetByJetCor
   fCollTemplates(other.fCollTemplates),
   fInitialized(other.fInitialized),
   fEfficiencyFixed(other.fEfficiencyFixed),
-  fhEfficiency(other.fhEfficiency)
+  fhEfficiency(other.fhEfficiency),
+  fpAppliedEfficiency(other.fpAppliedEfficiency)
 {
   // Copy constructor.
 }
@@ -77,6 +92,7 @@ AliEmcalJetByJetCorrection& AliEmcalJetByJetCorrection::operator=(const AliEmcal
   fInitialized       = other.fInitialized;
   fEfficiencyFixed   = other.fEfficiencyFixed;
   fhEfficiency       = other.fhEfficiency;
+  fpAppliedEfficiency= other.fpAppliedEfficiency;
 
   return *this;
 }
@@ -96,9 +112,10 @@ AliEmcalJet* AliEmcalJetByJetCorrection::Eval(const AliEmcalJet *jet, TClonesArr
   
   Double_t meanPt = GetMeanPtConstituents(jet,fTracks);
   Double_t eff = GetEfficiency(meanPt);
+  fpAppliedEfficiency->Fill(meanPt,eff);
 
   Int_t np = TMath::FloorNint((double)jet->GetNumberOfTracks() * (1./eff -1.));
-
+  
   TLorentzVector corrVec; corrVec.SetPtEtaPhiM(jet->Pt(),jet->Eta(),jet->Phi(),jet->M());
 
   Double_t mass = 0.13957; //pion mass
index 6c7f9e4448cefb4ae6a8c212249a8da1d2ffc9da..374ec15ad8bef61f1018250a732224707febe124 100644 (file)
@@ -7,6 +7,7 @@
 #include <TClonesArray.h>
 
 class TH3;
+class TProfile;
 class AliEmcalJet;
 
 class AliEmcalJetByJetCorrection : public TNamed
@@ -21,7 +22,7 @@ class AliEmcalJetByJetCorrection : public TNamed
 
   void SetTemplate(TH3 *h)          { fh3JetPtDRTrackPt = h; }
   void SetJetPtBinWidth(Double_t w) { fBinWidthJetPt    = w; }
-  void SetJetPtRange(Double_t min, Double_t max) {fJetPtMin = min; fJetPtMax = max;};
+  void SetJetPtRange(Double_t min, Double_t max) {fJetPtMin = min; fJetPtMax = max;}
   void SetFixedTrackEfficiency(Double_t eff) { fEfficiencyFixed = eff; }
   void SetEfficiencyHist(TH1 *h)             { fhEfficiency     = h  ; }
 
@@ -29,6 +30,8 @@ class AliEmcalJetByJetCorrection : public TNamed
   Double_t     GetEfficiency(const Double_t pt) const;
   Double_t     GetMeanPtConstituents(const AliEmcalJet *jet, TClonesArray *fTracks) const;
 
+  TProfile    *GetAppliedEfficiency() const {return fpAppliedEfficiency;}
+
   void         Init();
   AliEmcalJet *Eval(const AliEmcalJet *jet, TClonesArray *fTracks);
   
@@ -41,8 +44,11 @@ class AliEmcalJetByJetCorrection : public TNamed
   Bool_t    fInitialized;                      // status of initialization
   Double_t  fEfficiencyFixed;                  // fixed efficiency for all pT and all types of tracks
   TH1      *fhEfficiency;                      // single particle efficiency
- private:
 
-  ClassDef(AliEmcalJetByJetCorrection, 1) // jet-by-jet correction class
+  //book-keeping object filled inside Eval()
+  TProfile *fpAppliedEfficiency;               // Control profile efficiency
+
+ private:
+  ClassDef(AliEmcalJetByJetCorrection, 2) // jet-by-jet correction class
 };
 #endif
index 9ecbea37b1a966e27ed2bb0732dc1d4d2b1d384f..7842a0ff719ace59c811d4e89d0854f4439f2277 100644 (file)
@@ -58,7 +58,8 @@ AliAnalysisTaskEmcalJetMassStructure::AliAnalysisTaskEmcalJetMassStructure() :
   fh2PtMassCorr(0),
   fhnMassResponse(0),
   fhnMassResponseCorr(0),
-  fh3JetPtDRTrackPt(0)
+  fh3JetPtDRTrackPt(0),
+  fpUsedEfficiency(0)
 {
   // Default constructor.
 
@@ -101,7 +102,8 @@ AliAnalysisTaskEmcalJetMassStructure::AliAnalysisTaskEmcalJetMassStructure(const
   fh2PtMassCorr(0),
   fhnMassResponse(0),
   fhnMassResponseCorr(0),
-  fh3JetPtDRTrackPt(0)
+  fh3JetPtDRTrackPt(0),
+  fpUsedEfficiency(0)
 {
   // Standard constructor.
 
@@ -241,6 +243,9 @@ void AliAnalysisTaskEmcalJetMassStructure::UserCreateOutputObjects()
   fhnMassResponseCorr = new THnSparseF(histName.Data(),histTitle.Data(),nBinsSparse0,nBins0,xmin0,xmax0);
   fOutput->Add(fhnMassResponseCorr);
 
+  if(fEJetByJetCorr) fpUsedEfficiency = fEJetByJetCorr->GetAppliedEfficiency();
+  fOutput->Add(fpUsedEfficiency);
+
   TH1::AddDirectory(oldStatus);
 
   PostData(1, fOutput); // Post data for ALL output slots > 0 here.
@@ -308,7 +313,7 @@ Bool_t AliAnalysisTaskEmcalJetMassStructure::FillHistograms()
           curVec.SetPtEtaPhiM(vp->Pt(),vp->Eta(),vp->Phi(),vp->M());
           sumVec+=curVec;
           corrVec+=curVec;
-          //   Printf("%d  %f",i,dr[indexes[i]]);
+
           fh3PtDRMass[fCentBin]->Fill(ptJet1,dr[indexes[i]],sumVec.M()/mJet1);
           fh3PtDRRho[fCentBin]->Fill(ptJet1,dr[indexes[i]],sumVec.Pt()/ptJet1);
           
index 79f16114d29b6882b042642461efddd93f666659..b07d0c595ba4a7c69fb8d5822dda6c75195ce468 100644 (file)
@@ -4,7 +4,7 @@
 class TH1;
 class TH2;
 class TH3;
-class TProfile2D;
+class TProfile;
 class THnSparse;
 class TClonesArray;
 class TArrayI;
@@ -75,14 +75,14 @@ class AliAnalysisTaskEmcalJetMassStructure : public AliAnalysisTaskEmcalJet {
   TH2F                              **fh2PtMassCorr;               //! jet pT vs mass corrected
   THnSparse                          *fhnMassResponse;             //! response matrix
   THnSparse                          *fhnMassResponseCorr;         //! response matrix corrected
-
   TH3F                              **fh3JetPtDRTrackPt;           //! jet pt vs dr(jet axis, constituent) vs pT,track
+  TProfile                           *fpUsedEfficiency;            //! efficiency used for correction
 
  private:
   AliAnalysisTaskEmcalJetMassStructure(const AliAnalysisTaskEmcalJetMassStructure&);            // not implemented
   AliAnalysisTaskEmcalJetMassStructure &operator=(const AliAnalysisTaskEmcalJetMassStructure&); // not implemented
 
-  ClassDef(AliAnalysisTaskEmcalJetMassStructure, 2)
+  ClassDef(AliAnalysisTaskEmcalJetMassStructure, 3)
 };
 #endif