]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSGetterLight.h
Hits primary definition changes so that exclude primaries from depth of crystal
[u/mrichter/AliRoot.git] / PHOS / AliPHOSGetterLight.h
CommitLineData
f2bde07c 1#ifndef ALIPHOSGETTERLIGHT_H
2#define ALIPHOSGETTERLIGHT_H
3/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
5
6/* $Id$ */
7
8//_________________________________________________________________________
9// Class to mask AliPHOSGetter in "on flight" reconstruction (i.e. without
10// writing to file and creation full ALICE data structure)
11//
12//*-- Author: D.Peressounko (RRC KI)
13
14
15// --- ROOT system ---
16class TClonesArray ;
17class TObjArray ;
18
19// --- Standard library ---
20
21// --- AliRoot header files ---
22#include "AliPHOSGetter.h"
23#include "AliPHOSGeometry.h"
24class AliPHOSClusterizer ;
25class AliPHOSTrackSegmentMaker ;
26class AliPHOSPID ;
27
28class AliPHOSGetterLight : public AliPHOSGetter {
29
30public:
31 AliPHOSGetterLight() ; // ctor
32
33 virtual ~AliPHOSGetterLight() ; // dtor
34
35 static AliPHOSGetterLight * Instance(const char* /*headerFile*/,
36 const char* version = AliConfig::GetDefaultEventFolderName(),
37 Option_t * openingOption = "READ" ) ;
38 static AliPHOSGetterLight * Instance() ;
39
40 //=========== General information about run ==============
41 virtual Bool_t IsLoaded(TString /*tree*/) const { Error("IsLoaded","NotDefined"); return kFALSE ; }
42 virtual void SetLoaded(TString /*tree*/) { Error("SetLoaded","NotDefined"); }
43
44 virtual Int_t MaxEvent() const {return 1 ;} //always "read" event 1
45 virtual Int_t EventNumber() const {return 0; } //always the same event
46 virtual Bool_t VersionExists(TString & /*opt*/) const {return kFALSE;}
47 virtual UShort_t EventPattern(void) const {return 0;}
48 virtual Float_t BeamEnergy(void) const {return 10.;}
49
50 //========== PHOSGeometry and PHOS =============
51 virtual AliPHOS * PHOS() const { Error("PHOS()","NotDefined"); return 0 ;}
52 virtual AliPHOSGeometry * PHOSGeometry() const {AliPHOSGeometry * geom = AliPHOSGeometry::GetInstance("IHEP","IHEP") ; return geom ;}
53
54 //========== Methods to read something from file ==========
55 virtual void Event(Int_t /*event*/, const char * /*opt = "HSDRTP"*/){} //Use data already in memory
56 virtual void Track(Int_t /*itrack*/) { Error("Track()","NotDefined");}
57
58
59 //-----------------now getter's data--------------------------------------
60 virtual AliPHOSCalibrationDB * CalibrationDB(){return fcdb; }
61 virtual void SetCalibrationDB(AliPHOSCalibrationDB * cdb) {fcdb = cdb ;}
62
63 //=========== Primaries ============
64 virtual TClonesArray * Primaries(void) {Error("Primaries()","NotDefined"); return 0;}
65 virtual TParticle * Primary(Int_t /*index*/) const {Error("Primary()","NotDefined"); return 0;}
66 virtual Int_t NPrimaries()const {Error("NPrimaries()","NotDefined"); return 0; }
67 virtual TParticle * Secondary(const TParticle * /*p*/, Int_t /*index*/) const {Error("Secondary()","NotDefined"); return 0; }
68
69 //=========== Hits =================
70 virtual TClonesArray * Hits(void) {Error("Hits()","Method not defined") ; return 0;}
71 virtual AliPHOSHit * Hit(Int_t /*index*/) {Error("Hit()","Method not defined"); return 0 ;}
72 virtual TTree * TreeH() const {Error("TreeH","Method not defined") ; return 0;}
73
74 //=========== SDigits ==============
75 virtual TClonesArray * SDigits(){Error("Sdigits","Method not defined") ; return 0;}
76 virtual AliPHOSDigit * SDigit(Int_t /*index*/) { Error("SDigit","Method not defined") ;return 0 ;}
77 virtual TTree * TreeS() const {Error("TresS","Method not defined") ; return 0 ;}
78 virtual AliPHOSSDigitizer * SDigitizer() {Error("SDigitizer","Method not defined") ; return 0 ;}
79
80 virtual TString GetSDigitsFileName() const { Error("OnFlight","Method not defined") ; return "" ; }
81 virtual Int_t LoadSDigits(Option_t* /*opt=""*/) const { Error("LoadSDigits","Method not defined") ;return 0 ; }
82 virtual Int_t LoadSDigitizer(Option_t* /*opt=""*/) const { Error("LoadSdigitizer","Method not defined") ; return 0 ; }
83 virtual Int_t WriteSDigits(Option_t* /*opt=""*/) const { Error("WriteSDigits","Method not defined") ; return 0 ; }
84 virtual Int_t WriteSDigitizer(Option_t* /*opt=""*/) const {Error("WriteSDigitizer","Method not defined") ; return 0 ; }
85
86 //========== Digits ================
87 virtual TClonesArray * Digits(){return fDigits ; }
88 virtual AliPHOSDigit * Digit(Int_t index) { return static_cast<AliPHOSDigit *>(fDigits->At(index)) ;}
89 virtual TTree * TreeD() const {Error("TreeD","Method not defined") ; return 0;}
90 virtual AliPHOSDigitizer * Digitizer(){Error("Digitizer","Method not defined") ; return 0;}
91 virtual TString GetDigitsFileName() const { Error("GetDigitsFileName","Method not defined") ; return "" ; }
92 virtual Int_t LoadDigits(Option_t* /*opt=""*/) const {Error("LoadDigits","Method not defined") ; return 0 ; }
93 virtual Int_t LoadDigitizer(Option_t* /*opt=""*/) const {Error("LoadDigitizer","Method not defined") ; return 0;}
94 virtual Int_t WriteDigits(Option_t* /*opt=""*/) const { Error("WriteDigits","Method not defined") ; return 0; }
95 virtual Int_t WriteDigitizer(Option_t* /*opt=""*/) const {Error("WriteDigitizer","Method not defined"); return 0 ;}
96
97 //Methods to distinguish raw and simulated digits
98 virtual Bool_t IsRawDigits(void) const {return fRawDigits;}
99 virtual void SetRawDigits(Bool_t isRaw = kTRUE){fRawDigits = isRaw;}
100
101 //========== RecPoints =============
102 virtual TObjArray * EmcRecPoints(){return fEmcRecPoints ;}
103 virtual AliPHOSEmcRecPoint * EmcRecPoint(Int_t index) { return static_cast<AliPHOSEmcRecPoint *>(fEmcRecPoints->At(index)) ;}
104 virtual TObjArray * CpvRecPoints(){return fCpvRecPoints ;}
105 virtual AliPHOSCpvRecPoint * CpvRecPoint(Int_t index) { return static_cast<AliPHOSCpvRecPoint *>(fCpvRecPoints->At(index)) ;}
106 virtual TTree * TreeR() const {Error("TreeR","Method not defined") ; return 0;}
107 virtual AliPHOSClusterizer * Clusterizer() { return fClusterizer;}
108 virtual TString GetRecPointsFileName() const { Error("RecPointsFileName","Method not defined") ;return "" ; }
109 virtual Int_t LoadRecPoints(Option_t* /*opt=""*/) const {Error("LoadRecPoints","Method not defined") ; return 0; }
110 virtual Int_t LoadClusterizer(Option_t* /*opt=""*/) const {Error("LoadClusterizer","Method not defined") ; return 0 ;}
111 virtual Int_t WriteRecPoints(Option_t* /*opt=""*/) const {Error("WriteRecPoints","Method not defined") ; return 0; }
112 virtual Int_t WriteClusterizer(Option_t* /*opt=""*/) const {Error("WriteClusterizer","Method not defined"); return 0 ;}
113
114 //========== TrackSegments TClonesArray * TrackSegments(const char * name = 0) {
115 virtual TClonesArray * TrackSegments(){return fTS ;} ;
116 virtual AliPHOSTrackSegment * TrackSegment(Int_t index) { return static_cast<AliPHOSTrackSegment *>(fTS->At(index)) ;}
117 virtual TTree * TreeT() const {Error("TreeT","Method not defined") ; return 0 ;}
118 virtual AliPHOSTrackSegmentMaker * TrackSegmentMaker(){ return fTSM ;}
119 virtual TString GetTracksFileName() const { Error("GetTSFileName","Method not defiled") ; return "" ; }
120 virtual Int_t LoadTracks(Option_t* /*opt=""*/) const { Error("LoadTracks","Method not defined") ;return 0; }
121 virtual Int_t LoadTrackSegementMaker(Option_t* /*opt=""*/) const {Error("LoadTSMaker","Method noe defined") ; return 0 ;}
122 virtual Int_t WriteTracks(Option_t* /*opt=""*/) const { Error("WriteTracks","Method not defined") ; return 0; }
123 virtual Int_t WriteTrackSegmentMaker(Option_t* /*opt=""*/) const {Error("WriteTSM","Method not defined") ; return 0 ;}
124
125 //========== RecParticles ===========
126 virtual TClonesArray * RecParticles(){ return fRP;}
127 virtual AliPHOSRecParticle * RecParticle(Int_t index) { return static_cast<AliPHOSRecParticle *>(fRP->At(index)) ;}
128 virtual TTree * TreeP() const {Error("TreeP","Method net defined"); return 0 ;}
129 virtual AliPHOSPID * PID(){return fPID ;}
130 virtual TString GetRecParticlesFileName() const { Error("GetRPFileName","Method not defined") ; return "" ; }
131 virtual Int_t LoadRecParticles(Option_t* /*opt=""*/) const { Error("LoadRP","Method not defined") ;return 0 ; }
132 virtual Int_t LoadPID(Option_t* /*opt=""*/) const {Error("LoadPID","Method not defined") ;return 0 ; }
133 virtual Int_t WriteRecParticles(Option_t* /*opt=""*/) const { Error("WriteRP","Method not defined"); return 0 ; }
134 virtual Int_t WritePID(Option_t* /*opt=""*/) const {Error("WritePID","Method not defined"); return 0 ; }
135
136 //========== Raw ===========
137 // virtual Int_t ReadRaw(Int_t event) ;
138
139 // virtual void SetDebug(Int_t level) {fgDebug = level;} // Set debug level
140 virtual void PostClusterizer(AliPHOSClusterizer * clu)const{((AliPHOSGetterLight*)fgObjGetter)->fClusterizer = clu;}
141 virtual void PostPID(AliPHOSPID * pid)const{((AliPHOSGetterLight*)fgObjGetter)->fPID = pid ;}
142 virtual void PostTrackSegmentMaker(AliPHOSTrackSegmentMaker * tr)const{((AliPHOSGetterLight*)fgObjGetter)->fTSM = tr ;}
143 //virtual void PostSDigitizer (AliPHOSSDigitizer * sdigitizer)
144 // const {PhosLoader()->PostSDigitizer(sdigitizer);}
145 //virtual void PostDigitizer (AliPHOSDigitizer * digitizer)
146 // const {PhosLoader()->PostDigitizer(dynamic_cast<AliDigitizer *>(digitizer));}
147
148 virtual TString Version() const { return "OnFlight" ; }
149 virtual AliPHOSLoader * PhosLoader() const { Error("PhosLoader","Method not defined") ;return 0 ; }
150 virtual void Reset(){}
151
152 virtual AliESD * ESD() const { return 0 ; }
153
154private:
155
156 AliPHOSGetterLight(const char* headerFile,
157 const char* version = AliConfig::GetDefaultEventFolderName(),
158 Option_t * openingOption = "READ") ;
159
160private :
161
162 TClonesArray * fDigits ;
163 TObjArray * fEmcRecPoints ;
164 TObjArray * fCpvRecPoints ;
165 TClonesArray * fTS ;
166 TClonesArray * fRP ;
167
168 AliPHOSCalibrationDB * fcdb ;
169
170 AliPHOSClusterizer * fClusterizer ;
171 AliPHOSTrackSegmentMaker * fTSM ;
172 AliPHOSPID * fPID ;
173
174 Bool_t fRawDigits ;
175
176 ClassDef(AliPHOSGetterLight,1) // Getter for \"on flyght\" reconstruction
177
178};
179
180#endif // AliPHOSGETTERLIGHT_H