]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliGenInfo.h
70d9dbcdc9b190deb3518cab5f7b154814f08fcb
[u/mrichter/AliRoot.git] / PWG1 / AliGenInfo.h
1 #ifndef ALIGENINFO_H
2 #define ALIGENINFO_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6
7
8 //////////////////////////////////////////////////////////////////////////////
9 //                          Class AliGenInfo                               //
10 //   collect together MC info for comparison purposes - effieciency studies and so on//                                                                 //
11 //   marian.ivanov@cern.ch                                                  //
12 //////////////////////////////////////////////////////////////////////////////
13
14
15
16 ////////////////////////////////////////////////////////////////////////
17 //
18 // Start of implementation of the class AliTPCdigitRow
19 //
20 ////////////////////////////////////////////////////////////////////////
21
22 #include <TParticle.h>
23 #include "AliTrackReference.h"
24
25 class TFile;
26 class AliRunLoader;
27 class AliStack;
28 class AliTPCParam;
29
30 const Int_t kgRowBytes = 32;
31
32 class AliTPCdigitRow: public TObject {
33 public:
34   AliTPCdigitRow();
35   virtual ~AliTPCdigitRow(){;}
36   void SetRow(Int_t row);
37   Bool_t TestRow(Int_t row);
38   AliTPCdigitRow & operator=(const AliTPCdigitRow &digOld);
39   Int_t RowsOn(Int_t upto=8*kgRowBytes);
40   Int_t Last();
41   Int_t First();
42   void Reset();
43
44 //private:
45   UChar_t fDig[kgRowBytes];
46   ClassDef(AliTPCdigitRow,1)  // container for digit pattern
47 };
48
49
50 ////////////////////////////////////////////////////////////////////////
51 //
52 // Start of implementation of the class AliMCInfo
53 //
54 ////////////////////////////////////////////////////////////////////////
55
56 class AliMCInfo: public TObject {
57
58 public:
59   AliMCInfo();
60   ~AliMCInfo();
61   void Update();
62   AliTrackReference  fTrackRef;      // track reference saved in the output tree
63   AliTrackReference  fTrackRefOut;   // decay track reference saved in the output tree
64   AliTrackReference  fTRdecay;       // track reference at decay point
65   //
66   Int_t     fPrimPart;               // index of primary particle in TreeH
67   TParticle fParticle;               // generated particle 
68   Float_t   fMass;                   // mass of the particle
69   Float_t   fCharge;                 //
70   Int_t     fLabel;                  // track label
71   Int_t     fEventNr;                // event number
72   Int_t     fMCtracks;               // indication of how many times the track is retuturned back
73   Int_t fPdg;                        //pdg code
74   Float_t fDecayCoord[3];            // position of particle decay
75   Double_t fVDist[4];                //distance of the particle vertex from primary vertex
76   Bool_t fTPCdecay;                  //indicates decay in TPC
77   Int_t fRowsWithDigitsInn;          // number of rows with digits in the inner sectors
78   Int_t fRowsWithDigits;             // number of rows with digits in the outer sectors
79   Int_t fRowsTrackLength;            // last - first row with digit
80   Float_t fPrim;                     // theoretical dedx in tpc according particle momenta and mass
81   AliTPCdigitRow fTPCRow;                  // information about digits row pattern
82   Int_t fNTPCRef;                    // tpc references counter
83   Int_t fNITSRef;                    // ITS references counter
84   Int_t fNTRDRef;                    // TRD references counter
85   Int_t fNTOFRef;                    // TOF references counter
86   TClonesArray * fTPCReferences;     //containner with all track references -in the TPC
87   TClonesArray * fITSReferences;     //container with ITS references
88   TClonesArray * fTRDReferences;     //container with TRD references  
89   TClonesArray * fTOFReferences;     //container with TRD references  
90   //
91   ClassDef(AliMCInfo,1)  // container for 
92 };
93
94
95
96 class AliGenV0Info: public TObject {
97 public:
98   AliMCInfo fMCd;                   //info about daughter particle - second particle for V0
99   AliMCInfo fMCm;      //info about mother particle   - first particle for V0
100   TParticle fMotherP;   //particle info about mother particle
101   void Update(Float_t vertex[3]);        // put some derived info to special field 
102   Double_t    fMCDist1;    //info about closest distance according closest MC - linear DCA
103   Double_t    fMCDist2;    //info about closest distance parabolic DCA
104   //
105   Double_t     fMCPdr[3];    //momentum at vertex daughter  - according approx at DCA
106   Double_t     fMCPd[4];     //exact momentum from MC info
107   Double_t     fMCX[3];      //exact position of the vertex
108   Double_t     fMCXr[3];     //rec. position according helix
109   //
110   Double_t     fMCPm[3];    //momentum at the vertex mother
111   Double_t     fMCAngle[3]; //three angels
112   Double_t     fMCRr;       // rec position of the vertex 
113   Double_t     fMCR;        //exact r position of the vertex
114   Int_t        fPdg[2];   //pdg code of mother and daugter particles
115   Int_t        fLab[2];   //MC label of the partecle  
116   //
117   Double_t       fInvMass;  //reconstructed invariant mass -
118   Float_t        fPointAngleFi; //point angle fi
119   Float_t        fPointAngleTh; //point angle theta
120   Float_t        fPointAngle;   //point angle full
121   //
122   ClassDef(AliGenV0Info,1)  // container for  
123 };
124
125
126
127 class AliGenKinkInfo: public TObject {
128 public:
129   AliMCInfo fMCd;            //info about daughter particle - second particle for V0
130   AliMCInfo fMCm;            //info about mother particle   - first particle for V0
131   void Update();             // put some derived info to special field 
132   Float_t GetQt();           //
133   Double_t    fMCDist1;      //info about closest distance according closest MC - linear DCA
134   Double_t    fMCDist2;      //info about closest distance parabolic DCA
135   //
136   Double_t     fMCPdr[3];    //momentum at vertex daughter  - according approx at DCA
137   Double_t     fMCPd[4];     //exact momentum from MC info
138   Double_t     fMCX[3];      //exact position of the vertex
139   Double_t     fMCXr[3];     //rec. position according helix
140   //
141   Double_t     fMCPm[3];     //momentum at the vertex mother
142   Double_t     fMCAngle[3];  //three angels
143   Double_t     fMCRr;        // rec position of the vertex 
144   Double_t     fMCR;         //exact r position of the vertex
145   Int_t        fPdg[2];      //pdg code of mother and daugter particles
146   Int_t        fLab[2];      //MC label of the partecle
147   ClassDef(AliGenKinkInfo,1) // container for  
148 };
149
150
151
152 ////////////////////////////////////////////////////////////////////////
153 // 
154 // Start of implementation of the class AliGenInfoMaker
155 //
156 ////////////////////////////////////////////////////////////////////////
157
158 class AliGenInfoMaker {
159
160 public:
161   AliGenInfoMaker();
162   AliGenInfoMaker(const char * fnGalice, const char* fnRes    ="genTracks.root",
163                    Int_t nEvents=1, Int_t firstEvent=0);
164   virtual ~AliGenInfoMaker();
165   void Reset();
166   Int_t Exec();
167   Int_t Exec(Int_t nEvents, Int_t firstEventNr);
168   void CreateTreeGenTracks();
169   void CloseOutputFile();
170   Int_t TreeKLoop();
171   Int_t TreeTRLoop();
172   Int_t TreeDLoop();
173   Int_t BuildKinkInfo();  // build information about MC kinks
174   Int_t BuildV0Info();  // build information about MC kinks
175   Int_t BuildHitLines();  // build information about MC kinks
176   //void  FillInfo(Int_t iParticle);
177   void SetFirstEventNr(Int_t i) {fFirstEventNr = i;}
178   void SetNEvents(Int_t i) {fNEvents = i;}
179   void SetDebug(Int_t level) {fDebug = level;}
180   Int_t SetIO(Int_t eventNr);
181   Int_t CloseIOEvent();
182   Int_t CloseIO();
183   Int_t SetIO();
184   Float_t TR2LocalX(AliTrackReference *trackRef,
185                     AliTPCParam *paramTPC);
186   AliMCInfo * GetInfo(UInt_t i){return (i<fNParticles)? fGenInfo[i]:0;}
187   AliMCInfo * MakeInfo(UInt_t i);
188
189 public:
190   Int_t fDebug;                   //! debug flag  
191   Int_t fEventNr;                 //! current event number
192   Int_t fLabel;                   //! track label
193   Int_t fNEvents;                 //! number of events to process
194   Int_t fFirstEventNr;            //! first event to process
195   UInt_t fNParticles;              //! number of particles in TreeK
196   TTree * fTreeGenTracks;          //! output tree with generated tracks
197   TTree * fTreeKinks;             //!  output tree with Kinks
198   TTree * fTreeV0;                //!  output tree with V0
199   TTree * fTreeHitLines;          //! tree with hit lines
200   char  fFnRes[1000];             //! output file name with stored tracks
201   TFile *fFileGenTracks;          //! output file with stored fTreeGenTracks
202   //
203   AliRunLoader * fLoader;         //! pointer to the run loader
204   TTree * fTreeD;                 //! current tree with digits
205   TTree * fTreeTR;                //! current tree with TR
206   AliStack *fStack;               //! current stack
207   // 
208   AliMCInfo **   fGenInfo;    //! array with pointers to gen info
209   Int_t fNInfos;                  //! number of tracks with infos
210   //
211   AliTPCParam* fParamTPC;         //! AliTPCParam
212   Float_t fVPrim[3];             //! primary vertex position
213                                   // the fVDist[3] contains size of the 3-vector
214
215 private:
216   static const Int_t seedRow11 = 158;  // nRowUp - 1
217   static const Int_t seedRow12 = 139;  // nRowUp - 1 - (Int_t) 0.125*nRowUp
218   static const Int_t seedRow21 = 149;  // seedRow11 - shift
219   static const Int_t seedRow22 = 130;  // seedRow12 - shift
220   static const Double_t kRaddeg = 180./3.14159265358979312;
221   // 
222   static const Double_t fgTPCPtCut = 0.03; // do not store particles with generated pT less than this
223   static const Double_t fgITSPtCut = 0.2; // do not store particles with generated pT less than this
224   static const Double_t fgTRDPtCut = 0.2; // do not store particles with generated pT less than this
225   static const Double_t fgTOFPtCut = 0.15; // do not store particles with generated pT less than this
226  
227   ClassDef(AliGenInfoMaker,1)    // class which creates and fills tree with TPCGenTrack objects
228 };
229
230
231
232
233
234 AliTPCParam * GetTPCParam();
235 Float_t TPCBetheBloch(Float_t bg);
236
237 #endif