]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
New files belonging to the PPD version of the JETAN code
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 19 Jan 2006 11:20:00 +0000 (11:20 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 19 Jan 2006 11:20:00 +0000 (11:20 +0000)
JETAN/AliJetDistributions.cxx [new file with mode: 0644]
JETAN/AliJetDistributions.h [new file with mode: 0644]
JETAN/centua1.C [new file with mode: 0644]
JETAN/run.sh [new file with mode: 0644]

diff --git a/JETAN/AliJetDistributions.cxx b/JETAN/AliJetDistributions.cxx
new file mode 100644 (file)
index 0000000..da2b2f8
--- /dev/null
@@ -0,0 +1,213 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+//---------------------------------------------------------------------
+// JetDistributions class 
+// Get different basic distributions
+// Authors: mercedes.lopez.noriega@cern.ch
+//---------------------------------------------------------------------
+#include "AliJetDistributions.h"
+ClassImp(AliJetDistributions)
+////////////////////////////////////////////////////////////////////////
+// includes
+#include <Riostream.h>
+#include <TH1F.h>
+#include <TH2F.h>
+#include <TFile.h>
+#include <TTree.h>
+#include <TMath.h>
+  
+#include "AliJetProductionDataPDC2004.h"
+#include "AliJet.h"
+#include "AliJetKineReaderHeader.h"
+#include "AliJetESDReaderHeader.h"
+#include "AliLeading.h"
+  
+////////////////////////////////////////////////////////////////////////
+// constructor/destructor
+  
+AliJetDistributions::AliJetDistributions()
+{
+  // Constructor
+  fDirectory    = 0x0;   
+  fFile         = "jets.root";   
+  fEventMin     =  0;
+  fEventMax     = -1;
+  fRunMin       =  0;
+  fRunMax       = 11;
+  fPercentage   = -1.0;
+  SetReaderHeader();
+  
+  fPart = 0;
+  fGenJ = 0;
+  fRecJ = 0;
+
+  fRetaH = 0;
+  fRphiH = 0;
+  fRptH = 0;
+  fRetaphiH = 0;
+
+  fMultH = 0;
+
+  // for Analize()
+  fPythia = kFALSE;
+  fDoPart = kTRUE;
+  fDoGenJ = kTRUE;
+  fDoRecJ = kTRUE;
+
+}
+
+////////////////////////////////////////////////////////////////////////
+// define histogrames 
+
+void AliJetDistributions::DefineHistograms()
+{
+  // Define histograms to be used
+  fRetaH = new TH1F("fRetaH","Reconstructed eta",140,-1.2,1.2);
+  SetProperties(fRetaH,"#eta","entries");
+  fRphiH = new TH1F("fRphiH","Reconstructed phi",18,0,2.0*TMath::Pi());
+  SetProperties(fRphiH,"#phi","entries");
+  fRptH = new TH1F("fRptH","Reconstructed pt",150,0,150);
+  SetProperties(fRptH,"p_{T} (GeV/c)","entries");
+
+  fRetaphiH = new TH2F("fRetaphiH","Reconstructed eta vs. phi",140,-1.2,1.2,18,0.,2.0*TMath::Pi());
+  SetProperties(fRetaphiH,"#eta","#phi");
+
+  fMultH = new TH1F("fMultH","Reconstructed Multiplicity",1000,0,10000);
+  SetProperties(fMultH,"Multiplicity","entries");
+}
+
+void AliJetDistributions::SetProperties(TH1* h,const char* x, const char* y) const
+{
+  // Properties of histograms (x title, y title and error propagation)
+  h->SetXTitle(x);
+  h->SetYTitle(y);
+  h->Sumw2();
+}
+////////////////////////////////////////////////////////////////////////
+// fill histogrames 
+
+void AliJetDistributions::FillHistograms()
+{
+  // Run data 
+  AliJetProductionDataPDC2004* runData = new AliJetProductionDataPDC2004();
+  
+  // Loop over runs
+  TFile* jFile = 0x0;
+  for (Int_t iRun = fRunMin; iRun <= fRunMax; iRun++) {
+    char fn[20];
+    sprintf(fn,"%s/%s.root",fDirectory,(runData->GetRunTitle(iRun)).Data());
+    jFile = new TFile(fn);
+    printf("  Analyzing run: %d %s\n", iRun,fn);       
+    
+    // Get reader header and events to be looped over
+    AliJetReaderHeader *jReaderH = (AliJetReaderHeader*)(jFile->Get(fReaderHeader));
+    if (fEventMin == -1) fEventMin =  jReaderH->GetFirstEvent();
+    if (fEventMax == -1) {
+      fEventMax =  jReaderH->GetLastEvent();
+    } else {
+      fEventMax = TMath::Min(fEventMax, jReaderH->GetLastEvent());
+    }
+
+    //AliUA1JetHeader *jH = (AliUA1JetHeader *) (jFile->Get("AliUA1JetHeader"));
+    //jH->PrintParameters();
+
+    
+    // Loop over events
+    for (Int_t i = fEventMin; i < fEventMax; i++) {
+      //printf("  Analyzing run: %d  Event %d / %d \n",iRun, i, fEventMax);
+
+      // Get next tree with AliJet
+      char nameT[100];
+      sprintf(nameT, "TreeJ%d",i);
+      TTree *jetT =(TTree *)(jFile->Get(nameT));
+      if (fDoRecJ) jetT->SetBranchAddress("FoundJet",    &fRecJ);
+      if (fDoGenJ) jetT->SetBranchAddress("GenJet",      &fGenJ);
+      if (fDoPart) jetT->SetBranchAddress("LeadingPart", &fPart);
+      
+      jetT->GetEntry(0);
+
+      FillDistributions(fRecJ);
+
+      delete jetT;
+    } // end loop over events in one file
+    if (jFile) jFile->Close();
+    delete jFile;
+  } // end loop over files
+}
+
+void AliJetDistributions::FillDistributions(AliJet *j)
+{
+  // Fill histrograms
+  TArrayI inJet = j->GetInJet();
+  TArrayF etain = j->GetEtaIn();
+  TArrayF ptin = j->GetPtIn();
+  TArrayF phiin = j->GetPhiIn();
+
+  fMultH->Fill(inJet.GetSize(),1);
+  for(Int_t i=0;i<inJet.GetSize();i++) {
+    fRetaH->Fill(etain[i],1);
+    fRphiH->Fill(phiin[i],1);
+    fRptH->Fill(ptin[i],1);
+    fRetaphiH->Fill(etain[i],phiin[i],1);
+  }
+}
+
+////////////////////////////////////////////////////////////////////////
+// Plot histogrames 
+
+void AliJetDistributions::PlotHistograms()
+{
+  // improved!
+  fMultH->Draw();
+  fRetaH->Draw();
+  fRphiH->Draw();
+  fRetaphiH->Draw();
+  fRptH->Draw();
+
+}
+
+
+////////////////////////////////////////////////////////////////////////
+// Save histogrames 
+
+void AliJetDistributions::SaveHistograms()
+{
+  // Save histogramas
+  TFile *fOut = new TFile(fFile,"recreate");
+  fOut->cd();
+  fMultH->Write();
+  fRetaphiH->Write();
+  fRetaH->Write();
+  fRphiH->Write();
+  fRptH->Write();
+  fOut->Close();
+}
+
+// main Analysis function
+
+void AliJetDistributions::Analyze()
+{
+  // Do the analysis
+  DefineHistograms();
+  FillHistograms();
+  PlotHistograms();
+  SaveHistograms();
+}
+
+
diff --git a/JETAN/AliJetDistributions.h b/JETAN/AliJetDistributions.h
new file mode 100644 (file)
index 0000000..5981d63
--- /dev/null
@@ -0,0 +1,83 @@
+#ifndef ALIJETDISTRIBUTIONS_H
+#define ALIJETDISTRIBUTIONS_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+//---------------------------------------------------------------------
+// JetDistributions class 
+// Get different basic distributions
+// Author: mercedes.lopez.noriega@cern.ch
+//---------------------------------------------------------------------
+
+#include <TObject.h> 
+class AliLeading;
+class AliJet;
+class TH1;
+class TH1F;
+class TH2F;
+class TProfile;
+class TLorentzVector;
+
+class AliJetDistributions : public TObject
+{
+ public:
+    AliJetDistributions();
+    ~AliJetDistributions(){;}
+
+    void Analyze();
+    void DefineHistograms();
+    void FillHistograms();
+    void FillDistributions(AliJet *j);
+    void PlotHistograms();
+    void SaveHistograms();
+    
+    // Setter
+    void SetDirectory(char* directory) {fDirectory = directory;}
+    void SetOutputFile(char* file) {fFile = file;}
+    void SetPercentage(Float_t p) { fPercentage = p;}
+    void SetEventRange(Int_t imin, Int_t imax) {fEventMin = imin; fEventMax = imax;}
+    void SetRunRange(Int_t imin, Int_t imax) {fRunMin = imin; fRunMax = imax;}
+    void SetPythia(Bool_t f = kFALSE){fPythia = f;}    
+    void SetProperties(TH1* h,const char* x, const char* y) const;
+    void SetReaderHeader(char *s="AliJetKineReaderHeader") {fReaderHeader = s;}
+    void SetPartPtCut(Float_t c) { fPartPtCut = c; }
+
+    void SetDoLeadPart(Bool_t f = kTRUE) {fDoPart = f;}
+    void SetDoGenJet(Bool_t f = kTRUE) {fDoGenJ = f;}
+    void SetDoRecJet(Bool_t f = kTRUE) {fDoRecJ = f;}
+    
+ private:
+    char*  fReaderHeader;// Reader header
+    char*  fDirectory;   // Directory
+    char*  fFile     ;   // Output file name
+    Int_t  fEventMin;    // Minimum event number
+    Int_t  fEventMax;    // Maximum event number
+    Int_t  fRunMin;      // Minimum run number 
+    Int_t  fRunMax;      // Maximum run number
+    Float_t fPercentage; // percentage of pt from signal particles to accept a jet
+    Float_t fPartPtCut;  // cut in the pt of particles in histos
+
+    // user options
+    Bool_t fPythia;      // if pythia events
+    Bool_t fDoPart;      // do analysis of leading particle
+    Bool_t fDoGenJ;      // do analysis of leading generated jet
+    Bool_t fDoRecJ;      // do analysis of leading rec jet
+   
+    // leading hets and particles
+    AliLeading* fPart;   // pointer to leading particle
+    AliJet* fGenJ;       // pointer to leading generated jet
+    AliJet* fRecJ;       // pointer to leading reconstructed jet
+
+    // histos for reconstructed particles
+    TH1F* fRetaH;   // Eta of reconstructed particle
+    TH1F* fRphiH;   // Phi of reconstructed particle
+    TH1F* fRptH;    // Pt of reconstructed particle
+    TH2F* fRetaphiH;// Eta vs Phi of reconstructed particles
+
+    TH1F* fMultH;   // Multiplicity
+
+    ClassDef(AliJetDistributions,1)
+};
+#endif
diff --git a/JETAN/centua1.C b/JETAN/centua1.C
new file mode 100644 (file)
index 0000000..b7e8c9f
--- /dev/null
@@ -0,0 +1,46 @@
+void centua1(char* bin)
+{
+  gSystem->Load("$(ALICE_ROOT)/lib/tgt_$(ALICE_TARGET)/libJETAN");
+  
+  char dir[100];
+  //sprintf(dir, "/data1/jgcn/partonic_events/unquenched/%s", bin);
+  sprintf(dir, "/home/guest/alice/data/cent1_nq/%s", bin);
+
+  AliJetESDReaderHeader *krh = new AliJetESDReaderHeader(); // hijing
+  //AliJetKineReaderHeader *krh = new AliJetKineReaderHeader(); // pythia
+  krh->SetComment("");
+  krh->SetDirectory(dir);
+  krh->SetPattern("miniesd"); // hijing
+  krh->SetFirstEvent(0);
+  krh->SetLastEvent(100000);
+  krh->SetPtCut(1.0);
+  //krh->SetFastSimTPC(kFALSE); // pythia
+  //krh->SetFastSimEMCAL(kFALSE); // pythia
+  // define reader and set its header
+  AliJetESDReader *kr = new AliJetESDReader(); // hijing
+  //AliJetKineReader *kr = new AliJetKineReader(); // pythia
+  kr->SetReaderHeader(krh);
+  // define jet header
+  AliUA1JetHeader *jh=new AliUA1JetHeader();
+  jh->SetComment("UA1 jet code with radius 1");
+  jh->SetMode(1);
+  jh->SetRadius(0.4);
+  jh->SetMinCellEt(0.);
+  jh->SetEtSeed(4.);
+  jh->SetLegoNbinPhi(420.);
+  jh->SetLegoNbinEta(120.);
+  jh->SetLegoEtaMin(-0.9);
+  jh->SetLegoEtaMax(+0.9);
+    
+  // define jet finder. Set its header and reader
+  AliUA1JetFinder *jf = new AliUA1JetFinder();
+  jf->SetJetHeader(jh);
+  jf->SetJetReader(kr);
+  jf->SetPlotMode(kTRUE);
+  jf->SetOutputFile("jets.root");
+  // do the job
+  jf->Run();
+}
diff --git a/JETAN/run.sh b/JETAN/run.sh
new file mode 100644 (file)
index 0000000..d3e8f02
--- /dev/null
@@ -0,0 +1,79 @@
+aliroot -b << EOF 
+.x testua1.C("et14")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/150-180GeV.root
+
+
+aliroot -b << EOF 
+.x testua1.C("et13")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/125-150GeV.root
+
+aliroot -b << EOF 
+.x testua1.C("et12")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/104-125GeV.root
+
+aliroot -b << EOF 
+.x testua1.C("et11")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/86-104GeV.root
+
+
+aliroot -b << EOF 
+.x testua1.C("et10")
+.q
+EOF
+mv jets.root  jets_pt_0_r_1.0/72-86GeV.root
+
+aliroot -b << EOF 
+.x testua1.C("et9")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/60-72GeV.root
+
+aliroot -b << EOF 
+.x testua1.C("et8")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/50-60GeV.root
+
+
+aliroot -b << EOF 
+.x testua1.C("et7")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/42-50GeV.root
+
+aliroot -b << EOF 
+.x testua1.C("et6")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/35-42GeV.root
+
+
+aliroot -b << EOF 
+.x testua1.C("et5")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/29-35GeV.root
+
+
+aliroot -b << EOF 
+.x testua1.C("et4")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/24-29GeV.root
+
+
+aliroot -b << EOF 
+.x testua1.C("et3")
+.q
+EOF
+mv jets.root jets_pt_0_r_1.0/20-24GeV.root
+
+