]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenReaderTreeK.cxx
13d594442c4cf26dc62ea57e16b9b696d585859b
[u/mrichter/AliRoot.git] / EVGEN / AliGenReaderTreeK.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.1  2001/11/09 09:11:24  morsch
19 Realisation of AliGenReader that reads the kine tree (TreeK).
20
21 */
22 #include <TFile.h>
23 #include <TTree.h>
24 #include <TParticle.h>
25
26 #include "AliGenReaderTreeK.h"
27 #include "AliStack.h"
28 #include "AliHeader.h"
29 #include "AliRun.h"
30
31 ClassImp(AliGenReaderTreeK);
32
33
34 AliGenReaderTreeK::AliGenReaderTreeK():AliGenReader() 
35 {
36 //  Default constructor
37     fFileName       = NULL;
38     fStack          = 0;
39     fHeader         = 0;
40     fNcurrent       = 0;
41     fNparticle      = 0;
42     fFile           = 0;
43     fBaseFile       = 0;
44     fTreeE          = 0;
45 }
46
47 AliGenReaderTreeK::~AliGenReaderTreeK() 
48 {
49 // Destructor
50     delete fTreeE;
51 }
52
53 void AliGenReaderTreeK::Init() 
54 {
55 // Initialization
56 // Connect base file and file to read from
57
58     TTree *ali = gAlice->TreeE();
59     if (ali) {
60         fBaseFile = ali->GetCurrentFile();
61     } else {
62         printf("\n Warning: Basefile cannot be found !\n");
63     }
64     fFile = new TFile(fFileName);
65 }
66
67 Int_t AliGenReaderTreeK::NextEvent() 
68 {
69 // Read the next event  
70 //  cd to file with old kine tree    
71     if (!fBaseFile) Init();
72     if (fStack) delete fStack;
73     fStack = new AliStack(1000);
74     fFile->cd();
75 //  Connect treeE
76     if (fTreeE) {
77         delete fTreeE;
78     } else {
79         fTreeE = (TTree*)gDirectory->Get("TE");
80     }
81
82     if (fHeader) delete fHeader;
83     fHeader = 0;
84     fTreeE->SetBranchAddress("Header", &fHeader);
85 //  Get next event
86     fTreeE->GetEntry(fNcurrent);
87     fStack = fHeader->Stack();
88     fStack->GetEvent(fNcurrent);
89     
90 //  cd back to base file
91     fBaseFile->cd();
92 //
93     fNcurrent++;
94     fNparticle = 0;
95     Int_t ntrack =  fStack->GetNtrack();
96     printf("\n Next event contains %d particles", ntrack);
97     
98     return  ntrack;
99 }
100
101 TParticle* AliGenReaderTreeK::NextParticle() 
102 {
103 //  Return next particle
104     TParticle* part = fStack->Particle(fNparticle);
105     fNparticle++;
106     return part;
107 }
108
109
110
111 AliGenReaderTreeK& AliGenReaderTreeK::operator=(const  AliGenReaderTreeK& rhs)
112 {
113 // Assignment operator
114     return *this;
115 }
116
117
118