]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenReaderTreeK.cxx
order of filling in 2D hitsos and GetValues mthod correceted
[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.2  2001/11/12 14:31:00  morsch
19 Memory leaks fixed. (M. Bondila)
20
21 Revision 1.1  2001/11/09 09:11:24  morsch
22 Realisation of AliGenReader that reads the kine tree (TreeK).
23
24 */
25 #include <TFile.h>
26 #include <TTree.h>
27 #include <TParticle.h>
28
29 #include "AliGenReaderTreeK.h"
30 #include "AliStack.h"
31 #include "AliHeader.h"
32 #include "AliRun.h"
33
34 ClassImp(AliGenReaderTreeK);
35
36
37 AliGenReaderTreeK::AliGenReaderTreeK():AliGenReader() 
38 {
39 //  Default constructor
40     fFileName       = NULL;
41     fStack          = 0;
42     fHeader         = 0;
43     fNcurrent       = 0;
44     fNparticle      = 0;
45     fFile           = 0;
46     fBaseFile       = 0;
47     fTreeE          = 0;
48 }
49
50 AliGenReaderTreeK::AliGenReaderTreeK(const AliGenReaderTreeK &reader)
51 {
52     ;
53 }
54
55
56 AliGenReaderTreeK::~AliGenReaderTreeK() 
57 {
58 // Destructor
59     delete fTreeE;
60 }
61
62 void AliGenReaderTreeK::Init() 
63 {
64 // Initialization
65 // Connect base file and file to read from
66
67     TTree *ali = gAlice->TreeE();
68     if (ali) {
69         fBaseFile = ali->GetCurrentFile();
70     } else {
71         printf("\n Warning: Basefile cannot be found !\n");
72     }
73     fFile = new TFile(fFileName);
74 }
75
76 Int_t AliGenReaderTreeK::NextEvent() 
77 {
78 // Read the next event  
79 //  cd to file with old kine tree    
80     if (!fBaseFile) Init();
81     if (fStack) delete fStack;
82     fStack = new AliStack(1000);
83     fFile->cd();
84 //  Connect treeE
85     if (fTreeE) {
86         delete fTreeE;
87     } else {
88         fTreeE = (TTree*)gDirectory->Get("TE");
89     }
90
91     if (fHeader) delete fHeader;
92     fHeader = 0;
93     fTreeE->SetBranchAddress("Header", &fHeader);
94 //  Get next event
95     fTreeE->GetEntry(fNcurrent);
96     fStack = fHeader->Stack();
97     fStack->GetEvent(fNcurrent);
98     
99 //  cd back to base file
100     fBaseFile->cd();
101 //
102     fNcurrent++;
103     fNparticle = 0;
104     Int_t ntrack =  fStack->GetNtrack();
105     printf("\n Next event contains %d particles", ntrack);
106     
107     return  ntrack;
108 }
109
110 TParticle* AliGenReaderTreeK::NextParticle() 
111 {
112 //  Return next particle
113     TParticle* part = fStack->Particle(fNparticle);
114     fNparticle++;
115     return part;
116 }
117
118
119
120 AliGenReaderTreeK& AliGenReaderTreeK::operator=(const  AliGenReaderTreeK& rhs)
121 {
122 // Assignment operator
123     return *this;
124 }
125
126
127