]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added vertex distribution histos
authordainese <dainese@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 14 Dec 2009 01:10:11 +0000 (01:10 +0000)
committerdainese <dainese@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 14 Dec 2009 01:10:11 +0000 (01:10 +0000)
PWG1/global/AliAnalysisTaskVertexESD.cxx
PWG1/global/AliAnalysisTaskVertexESD.h

index ffc76f0d9dbe3a872a85f0c99534ce1c3c5b146c..90d1db5dafa4ce738cc76cd090b242338746739f 100644 (file)
@@ -68,9 +68,19 @@ fRecoVtxTPC(kFALSE),
 fRecoVtxITSTPC(kFALSE),
 fOnlyITSTPCTracks(kFALSE),
 fOnlyITSSATracks(kFALSE),
+fFillNtuple(kFALSE),
 fESD(0), 
 fOutput(0), 
 fNtupleVertexESD(0),
+fhSPDVertexX(0),
+fhSPDVertexY(0),
+fhSPDVertexZ(0),
+fhTRKVertexX(0),
+fhTRKVertexY(0),
+fhTRKVertexZ(0),
+fhTPCVertexX(0),
+fhTPCVertexY(0),
+fhTPCVertexZ(0),
 fhTrackRefs(0)
 {
   // Constructor
@@ -133,6 +143,25 @@ void AliAnalysisTaskVertexESD::CreateOutputObjects()
 
   fOutput->Add(fNtupleVertexESD);
 
+  fhSPDVertexX = new TH1F("fhSPDVertexX","SPDVertex x; x vertex [cm]; events",200,-1,1);
+  fOutput->Add(fhSPDVertexX);
+  fhSPDVertexY = new TH1F("fhSPDVertexY","SPDVertex y; y vertex [cm]; events",200,-1,1);
+  fOutput->Add(fhSPDVertexY);
+  fhSPDVertexZ = new TH1F("fhSPDVertexZ","SPDVertex z; z vertex [cm]; events",200,-20,20);
+  fOutput->Add(fhSPDVertexZ);
+  fhTRKVertexX = new TH1F("fhTRKVertexX","TRKVertex x; x vertex [cm]; events",200,-1,1);
+  fOutput->Add(fhTRKVertexX);
+  fhTRKVertexY = new TH1F("fhTRKVertexY","TRKVertex y; y vertex [cm]; events",200,-1,1);
+  fOutput->Add(fhTRKVertexY);
+  fhTRKVertexZ = new TH1F("fhTRKVertexZ","TRKVertex z; z vertex [cm]; events",200,-20,20);
+  fOutput->Add(fhTRKVertexZ);
+  fhTPCVertexX = new TH1F("fhTPCVertexX","TPCVertex x; x vertex [cm]; events",200,-3,3);
+  fOutput->Add(fhTPCVertexX);
+  fhTPCVertexY = new TH1F("fhTPCVertexY","TPCVertex y; y vertex [cm]; events",200,-3,3);
+  fOutput->Add(fhTPCVertexY);
+  fhTPCVertexZ = new TH1F("fhTPCVertexZ","TPCVertex z; z vertex [cm]; events",200,-20,20);
+  fOutput->Add(fhTPCVertexZ);
+
   fhTrackRefs = new TH2F("fhTrackRefs","Track references; x; y",1000,-4,4,1000,-4,4);
   fOutput->Add(fhTrackRefs);
 
@@ -270,8 +299,38 @@ void AliAnalysisTaskVertexESD::Exec(Option_t *)
     }
     spd0cls = alimult->GetNumberOfSingleClusters()+ntrklets;
   }
+  
+  // fill histos
+  
+  if(spdv) {
+    if(spdv->GetNContributors()>0) {
+      TString title=spdv->GetTitle();
+      if(title.Contains("3D")) {
+       fhSPDVertexX->Fill(spdv->GetXv());
+       fhSPDVertexY->Fill(spdv->GetYv());
+      }
+      fhSPDVertexZ->Fill(spdv->GetZv());
+    }
+  }
+  
+  if(trkv) {
+    if(trkv->GetNContributors()>0) {
+      fhTRKVertexX->Fill(trkv->GetXv());
+      fhTRKVertexY->Fill(trkv->GetYv());
+      fhTRKVertexZ->Fill(trkv->GetZv());
+    }
+  }
+  
+  if(tpcv) {
+    if(tpcv->GetNContributors()>0) {
+      fhTPCVertexX->Fill(tpcv->GetXv());
+      fhTPCVertexY->Fill(tpcv->GetYv());
+      fhTPCVertexZ->Fill(tpcv->GetZv());
+    }
+  } 
+  
 
-
+  // fill ntuple
   Int_t isize=37;
   Float_t xnt[37];
   
@@ -333,7 +392,7 @@ void AliAnalysisTaskVertexESD::Exec(Option_t *)
 
   if(index!=isize) printf("AliAnalysisTaskVertexESD: ERROR, index!=isize\n");
 
-  fNtupleVertexESD->Fill(xnt);
+  if(fFillNtuple) fNtupleVertexESD->Fill(xnt);
   
   // Post the data already here
   PostData(0, fOutput);
index 4cfb3bb4b5a3e539b722b5315bdccd8e25ce991e..555d59d031b6cfc0a6f74d47c362808f341d506f 100644 (file)
@@ -42,7 +42,8 @@ class AliAnalysisTaskVertexESD : public AliAnalysisTask
   void           SetRerecoVertexITSTPC(Bool_t flag=kTRUE) { fRecoVtxITSTPC=flag; }
   void           SetOnlyITSTPCTracks() {fOnlyITSTPCTracks=kTRUE;}
   void           SetOnlyITSSATracks() {fOnlyITSSATracks=kTRUE;}
-  
+  void           SetFillNtuple(Bool_t fill=kTRUE) {fFillNtuple=fill;}  
+
  protected:
   Bool_t       fCheckEventType; // read only events of type 7
   Bool_t       fReadMC;         // read Monte Carlo
@@ -50,9 +51,19 @@ class AliAnalysisTaskVertexESD : public AliAnalysisTask
   Bool_t       fRecoVtxITSTPC;  // reco ITS+TPC vertex on the flight
   Bool_t       fOnlyITSTPCTracks; // only ITS-TPC tracks to redo ITSTPC vertex
   Bool_t       fOnlyITSSATracks;  // only ITS-SA tracks to redo ITSTPC vertex
+  Bool_t       fFillNtuple;      // fill ntuple 
   AliESDEvent *fESD;            // ESD object
   TList       *fOutput;         //! list send on output slot 0
   TNtuple     *fNtupleVertexESD;//! output ntuple
+  TH1F        *fhSPDVertexX; //! output histo
+  TH1F        *fhSPDVertexY; //! output histo
+  TH1F        *fhSPDVertexZ; //! output histo
+  TH1F        *fhTRKVertexX; //! output histo
+  TH1F        *fhTRKVertexY; //! output histo
+  TH1F        *fhTRKVertexZ; //! output histo
+  TH1F        *fhTPCVertexX; //! output histo
+  TH1F        *fhTPCVertexY; //! output histo
+  TH1F        *fhTPCVertexZ; //! output histo
   TH2F        *fhTrackRefs;     //! output histo
 
  private:    
@@ -62,7 +73,7 @@ class AliAnalysisTaskVertexESD : public AliAnalysisTask
   AliESDVertex* ReconstructPrimaryVertexTPC() const;
   AliESDVertex* ReconstructPrimaryVertexITSTPC() const;
   
-  ClassDef(AliAnalysisTaskVertexESD,6); // primary vertex analysis
+  ClassDef(AliAnalysisTaskVertexESD,7); // primary vertex analysis
 };
 
 #endif