]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Readers for EMCAL primary particle input.
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 8 Jan 2002 09:59:34 +0000 (09:59 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 8 Jan 2002 09:59:34 +0000 (09:59 +0000)
EVGEN/AliGenReaderEcalHijing.cxx [new file with mode: 0644]
EVGEN/AliGenReaderEcalHijing.h [new file with mode: 0644]
EVGEN/AliGenReaderEcalJets.cxx [new file with mode: 0644]
EVGEN/AliGenReaderEcalJets.h [new file with mode: 0644]

diff --git a/EVGEN/AliGenReaderEcalHijing.cxx b/EVGEN/AliGenReaderEcalHijing.cxx
new file mode 100644 (file)
index 0000000..a614b37
--- /dev/null
@@ -0,0 +1,119 @@
+/**************************************************************************
+ * 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.                  *
+ **************************************************************************/
+
+/*
+$Log$
+*/
+
+#include <TFile.h>
+#include <TTree.h>
+#include <TParticle.h>
+
+#include "AliGenReaderEcalHijing.h"
+#include "AliMC.h"
+ClassImp(AliGenReaderEcalHijing)
+
+
+AliGenReaderEcalHijing::AliGenReaderEcalHijing() 
+{
+// Default constructor
+    fNcurrent   = 0;
+    fTreeNtuple = 0;
+}
+
+void AliGenReaderEcalHijing::Init() 
+{
+//
+// reset the existing file environment and open a new root file if
+// the pointer to the Fluka tree is null
+    
+    TFile *pFile=0;
+    if (!pFile) {
+       pFile = new TFile(fFileName);
+       pFile->cd();
+       printf("\n I have opened %s file \n", fFileName);
+    }
+// get the tree address in the Fluka boundary source file
+    fTreeNtuple = (TTree*)gDirectory->Get("h2");
+    TTree *h2=fTreeNtuple;
+    h2->SetMakeClass(1);
+//Set branch addresses
+    h2->SetBranchAddress("njatt", &fNjatt);
+    h2->SetBranchAddress("nahij", &fNahij);
+    h2->SetBranchAddress("nphij", &fNphij);
+    h2->SetBranchAddress("khij",   fKhij) ;
+    h2->SetBranchAddress("pxhij",  fPxhij);
+    h2->SetBranchAddress("pyhij",  fPyhij);
+    h2->SetBranchAddress("pzhij",  fPzhij);
+    h2->SetBranchAddress("ehij",   fEhij) ;
+}
+
+Int_t AliGenReaderEcalHijing::NextEvent() 
+{
+// Read the next event  
+    Int_t nTracks, nread;
+    
+    TFile* pFile = fTreeNtuple->GetCurrentFile();
+    pFile->cd();
+
+    Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
+    if (fNcurrent < nentries) {
+       Int_t nb = (Int_t)fTreeNtuple->GetEvent(fNcurrent);
+       nread += nb;
+       fNcurrent++;
+       printf("\n Next event contains %d tracks! \n", fNphij);
+       nTracks    = fNphij;
+       fNparticle = 0;
+       return nTracks;
+    }
+    return 0;
+}
+
+TParticle* AliGenReaderEcalHijing::NextParticle() 
+{
+    Float_t p[4];
+// Read the next particle
+    Int_t ipart = fKhij[fNparticle];
+    p[0] = fPxhij[fNparticle];
+    p[1] = fPyhij[fNparticle];      
+    p[2] = fPzhij[fNparticle];
+    p[3] = fEhij[fNparticle];
+    
+    Double_t amass = TDatabasePDG::Instance()->GetParticle(ipart)->Mass();
+
+    if(p[3] <= amass) {
+       Warning("Generate","Particle %d  E = %f mass = %f \n",
+               ipart, p[3], amass);
+    } 
+    TParticle* particle = 
+       new TParticle(ipart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], 
+                     0., 0., 0., 0.);
+    fNparticle++;
+    return particle;
+}
+
+
+
+AliGenReaderEcalHijing& AliGenReaderEcalHijing::operator=(const  AliGenReaderEcalHijing& rhs)
+{
+// Assignment operator
+    return *this;
+}
+
+
+
+
+
+
diff --git a/EVGEN/AliGenReaderEcalHijing.h b/EVGEN/AliGenReaderEcalHijing.h
new file mode 100644 (file)
index 0000000..e680223
--- /dev/null
@@ -0,0 +1,46 @@
+#ifndef ALIGENREADERECALHIJING_H
+#define ALIGENREADERECALHIJING_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include "AliGenReader.h"
+
+
+class AliGenReaderEcalHijing : public AliGenReader
+{
+ public:
+    AliGenReaderEcalHijing();
+    
+    AliGenReaderEcalHijing(const AliGenReaderEcalHijing &reader){;}
+    virtual ~AliGenReaderEcalHijing(){;}
+    // Initialise 
+    virtual void Init();
+    // Read
+    virtual Int_t NextEvent();
+    virtual TParticle*  NextParticle();
+    AliGenReaderEcalHijing & operator=(const AliGenReaderEcalHijing & rhs);
+ protected:
+    Int_t             fNcurrent;      // points to the next entry
+    Int_t             fNparticle;     // 
+    
+    TTree            *fTreeNtuple;    // pointer to the TTree
+    //Declaration of leaves types
+    Int_t           fNjatt;           // Number of particles
+    Int_t           fNahij;           // Number of particles in alice accept. 
+    Int_t           fNphij;           // ?
+    Int_t           fKhij[10000];     // particle code
+    Float_t         fPxhij[10000];    // px
+    Float_t         fPyhij[10000];    // py
+    Float_t         fPzhij[10000];    // pz
+    Float_t         fEhij[10000];     // energy
+    ClassDef(AliGenReaderEcalHijing,1) // Read particles from cwn-ntuple
+};
+#endif
+
+
+
+
+
+
diff --git a/EVGEN/AliGenReaderEcalJets.cxx b/EVGEN/AliGenReaderEcalJets.cxx
new file mode 100644 (file)
index 0000000..3367267
--- /dev/null
@@ -0,0 +1,132 @@
+/**************************************************************************
+ * 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.                  *
+ **************************************************************************/
+
+/*
+$Log$
+*/
+
+#include <TFile.h>
+#include <TTree.h>
+#include <TParticle.h>
+
+#include "AliGenReaderEcalJets.h"
+#include "AliMC.h"
+ClassImp(AliGenReaderEcalJets)
+
+
+AliGenReaderEcalJets::AliGenReaderEcalJets() 
+{
+// Default constructor
+    fNcurrent   = 0;
+    fTreeNtuple = 0;
+}
+
+void AliGenReaderEcalJets::Init() 
+{
+//
+// reset the existing file environment and open a new root file if
+// the pointer to the Fluka tree is null
+    
+    TFile *pFile=0;
+    if (!pFile) {
+       pFile = new TFile(fFileName);
+       pFile->cd();
+       printf("\n I have opened %s file \n", fFileName);
+    }
+// get the tree address in the Fluka boundary source file
+    fTreeNtuple = (TTree*)gDirectory->Get("h1");
+    TTree *h1=fTreeNtuple;
+    h1->SetMakeClass(1);
+//Set branch addresses
+    h1->SetBranchAddress("nev",   &fNev);
+    h1->SetBranchAddress("x",      fX);
+    h1->SetBranchAddress("xtyp",   fXtyp);
+    h1->SetBranchAddress("npart", &fNpart);
+    h1->SetBranchAddress("xpt",    fXpt);
+    h1->SetBranchAddress("xeta",   fXeta);
+    h1->SetBranchAddress("xphi",   fXphi);
+    h1->SetBranchAddress("xid",    fXid);
+    h1->SetBranchAddress("njet",  &fNjet);
+    h1->SetBranchAddress("jet",    fJet);
+    h1->SetBranchAddress("jeta",   fJeta);
+    h1->SetBranchAddress("jphi",   fJphi);
+    h1->SetBranchAddress("nsjet", &fNsjet);
+    h1->SetBranchAddress("jset",   fJset);
+    h1->SetBranchAddress("jseta",  fJseta);
+    h1->SetBranchAddress("jsphi",  fJsphi);
+    h1->SetBranchAddress("npjet", &fNpjet);
+    h1->SetBranchAddress("jpet",   fJpet);
+    h1->SetBranchAddress("jpeta",  fJpeta);
+    h1->SetBranchAddress("jpphi",  fJpphi);
+}
+
+Int_t AliGenReaderEcalJets::NextEvent() 
+{
+// Read the next event  
+    Int_t nTracks, nread;
+    
+    TFile* pFile = fTreeNtuple->GetCurrentFile();
+    pFile->cd();
+
+    Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
+    if (fNcurrent < nentries) {
+       Int_t nb = (Int_t)fTreeNtuple->GetEvent(fNcurrent);
+       nread += nb;
+       fNcurrent++;
+       printf("\n Next event contains %d tracks! \n", fNpjet);
+       nTracks    = fNpjet;
+       fNparticle = 0;
+       return nTracks;
+    }
+    return 0;
+}
+
+TParticle* AliGenReaderEcalJets::NextParticle() 
+{
+    Float_t p[4];
+// Read the next particle
+    Int_t    ipart  = fXid[fNparticle];
+    Float_t  pt     = fXpt[fNparticle];
+    Float_t  eta    = fXeta[fNparticle];
+    Float_t  phi    = fXphi[fNparticle];
+    Float_t  theta  = 2.*TMath::ATan(TMath::Exp(-eta));
+    Double_t amass  = TDatabasePDG::Instance()->GetParticle(ipart)->Mass();
+
+    p[0] = pt*TMath::Sin(phi);
+    p[1] = pt*TMath::Cos(phi);      
+    p[2] = pt/TMath::Cos(theta);
+    p[3] = TMath::Sqrt(pt*pt+p[2]*p[2]+amass*amass);
+    
+
+    TParticle* particle = 
+       new TParticle(ipart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], 
+                     0., 0., 0., 0.);
+    fNparticle++;
+    return particle;
+}
+
+
+
+AliGenReaderEcalJets& AliGenReaderEcalJets::operator=(const  AliGenReaderEcalJets& rhs)
+{
+// Assignment operator
+    return *this;
+}
+
+
+
+
+
+
diff --git a/EVGEN/AliGenReaderEcalJets.h b/EVGEN/AliGenReaderEcalJets.h
new file mode 100644 (file)
index 0000000..fe60a76
--- /dev/null
@@ -0,0 +1,58 @@
+#ifndef ALIGENREADERECALJETS_H
+#define ALIGENREADERECALJETS_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include "AliGenReader.h"
+
+
+class AliGenReaderEcalJets : public AliGenReader
+{
+ public:
+    AliGenReaderEcalJets();
+    
+    AliGenReaderEcalJets(const AliGenReaderEcalJets &reader){;}
+    virtual ~AliGenReaderEcalJets(){;}
+    // Initialise 
+    virtual void Init();
+    // Read
+    virtual Int_t NextEvent();
+    virtual TParticle*  NextParticle();
+    AliGenReaderEcalJets & operator=(const AliGenReaderEcalJets & rhs);
+ protected:
+    Int_t           fNcurrent;      // points to the next event
+    Int_t           fNparticle;     // points to the next particle 
+    Int_t           fNev;           // event number
+    Float_t         fX[2];          // 
+    Int_t           fXtyp[2];       // parton type
+    Int_t           fNpart;         // number of particles  
+    Float_t         fXpt[200];      // pt of particle
+    Float_t         fXeta[200];     // eta of particle
+    Float_t         fXphi[200];     // phi of particle
+    Int_t           fXid[200];      // id of particle
+    Int_t           fNjet;          // number of jets 
+    Float_t         fJet[10];       // E_t of jet
+    Float_t         fJeta[10];      // eta of jet
+    Float_t         fJphi[10];      // phi of jet
+    Int_t           fNsjet;         // number of clusters
+    Float_t         fJset[10];      // E_t of cluster 
+    Float_t         fJseta[10];     // eta of cluster
+    Float_t         fJsphi[10];     // phi of cluster
+    Int_t           fNpjet;         // ?
+    Float_t         fJpet[10];      // ?
+    Float_t         fJpeta[10];     // ?
+    Float_t         fJpphi[10];     // ?
+
+    TTree            *fTreeNtuple;    // pointer to the TTree
+    //Declaration of leaves types
+    ClassDef(AliGenReaderEcalJets,1) // Read particles from cwn-ntuple
+};
+#endif
+
+
+
+
+
+