StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StPicoDst.h
1 
8 #ifndef StPicoDst_h
9 #define StPicoDst_h
10 
11 // ROOT headers
12 #include "TClonesArray.h"
13 
14 // PicoDst headers
15 #include "StPicoArrays.h"
16 
17 // Forward declarations
18 class StPicoEvent;
19 class StPicoTrack;
20 class StPicoEmcTrigger;
21 class StPicoMtdTrigger;
22 class StPicoBTowHit;
23 class StPicoBTofHit;
24 class StPicoMtdHit;
25 class StPicoEpdHit;
26 class StPicoBbcHit;
27 class StPicoFmsHit;
30 class StPicoMtdPidTraits;
32 class StPicoBEmcSmdEHit;
33 class StPicoBEmcSmdPHit;
34 class StPicoETofHit;
36 class StPicoMcVertex;
37 class StPicoMcTrack;
38 
39 //_________________
40 class StPicoDst {
41 
42  public:
43 
44 #if defined (__TFG__VERSION__)
45  StPicoDst() { fgPicoDst = this;}
46  virtual ~StPicoDst() {fgPicoDst = 0;}
47  virtual Bool_t IsGoodTrigger() const;
48 #else /* ! __TFG__VERSION__ */
49  StPicoDst() { /* emtpy */}
52  ~StPicoDst() { /* empty*/ }
53 #endif
54 
56  static void set(TClonesArray**);
58  static void unset();
60  static TClonesArray* picoArray(Int_t type) { return picoArrays[type]; }
61 
63  static StPicoEvent* event() { return (StPicoEvent*)picoArrays[StPicoArrays::Event]->UncheckedAt(0); }
65  static StPicoTrack* track(Int_t i) { return (StPicoTrack*)picoArrays[StPicoArrays::Track]->UncheckedAt(i); }
67  static StPicoEmcTrigger* emcTrigger(Int_t i) { return (StPicoEmcTrigger*)picoArrays[StPicoArrays::EmcTrigger]->UncheckedAt(i); }
69  static StPicoMtdTrigger* mtdTrigger(Int_t i) { return (StPicoMtdTrigger*)picoArrays[StPicoArrays::MtdTrigger]->UncheckedAt(i); }
71  static StPicoBTowHit* btowHit(Int_t i) { return (StPicoBTowHit*)picoArrays[StPicoArrays::BTowHit]->UncheckedAt(i); }
73  static StPicoBTofHit* btofHit(Int_t i) { return (StPicoBTofHit*)picoArrays[StPicoArrays::BTofHit]->UncheckedAt(i); }
75  static StPicoMtdHit* mtdHit(Int_t i) { return (StPicoMtdHit*)picoArrays[StPicoArrays::MtdHit]->UncheckedAt(i); }
77  static StPicoBbcHit* bbcHit(Int_t i) { return (StPicoBbcHit*)picoArrays[StPicoArrays::BbcHit]->UncheckedAt(i); }
79  static StPicoEpdHit* epdHit(Int_t i) { return (StPicoEpdHit*)picoArrays[StPicoArrays::EpdHit]->UncheckedAt(i); }
81  static StPicoFmsHit* fmsHit(Int_t i) { return (StPicoFmsHit*)picoArrays[StPicoArrays::FmsHit]->UncheckedAt(i); }
83  static StPicoBEmcPidTraits* bemcPidTraits(Int_t i) { return (StPicoBEmcPidTraits*)picoArrays[StPicoArrays::BEmcPidTraits]->UncheckedAt(i); }
85  static StPicoBTofPidTraits* btofPidTraits(Int_t i) { return (StPicoBTofPidTraits*)picoArrays[StPicoArrays::BTofPidTraits]->UncheckedAt(i); }
87  static StPicoMtdPidTraits* mtdPidTraits(Int_t i) { return (StPicoMtdPidTraits*)picoArrays[StPicoArrays::MtdPidTraits]->UncheckedAt(i); }
89  static StPicoTrackCovMatrix* trackCovMatrix(Int_t i) { return (StPicoTrackCovMatrix*)picoArrays[StPicoArrays::TrackCovMatrix]->UncheckedAt(i); }
91  static StPicoBEmcSmdEHit* bemcSmdEHit(Int_t i) { return (StPicoBEmcSmdEHit*)picoArrays[StPicoArrays::BEmcSmdEHit]->UncheckedAt(i); }
93  static StPicoBEmcSmdPHit* bemcSmdPHit(Int_t i) { return (StPicoBEmcSmdPHit*)picoArrays[StPicoArrays::BEmcSmdPHit]->UncheckedAt(i); }
95  static StPicoETofHit* etofHit(Int_t i) { return (StPicoETofHit*)picoArrays[StPicoArrays::ETofHit]->UncheckedAt(i); }
97  static StPicoETofPidTraits* etofPidTraits(Int_t i) { return (StPicoETofPidTraits*)picoArrays[StPicoArrays::ETofPidTraits]->UncheckedAt(i); }
99  static StPicoMcVertex* mcVertex(Int_t i) { return (StPicoMcVertex*)picoArrays[StPicoArrays::McVertex]->UncheckedAt(i); }
101  static StPicoMcTrack* mcTrack(Int_t i) { return (StPicoMcTrack*)picoArrays[StPicoArrays::McTrack]->UncheckedAt(i); }
102 
104  static UInt_t numberOfTracks() { return picoArrays[StPicoArrays::Track]->GetEntriesFast(); }
106  static UInt_t numberOfEmcTriggers() { return picoArrays[StPicoArrays::EmcTrigger]->GetEntriesFast(); }
108  static UInt_t numberOfMtdTriggers() { return picoArrays[StPicoArrays::MtdTrigger]->GetEntriesFast(); }
110  static UInt_t numberOfBTowHits() { return picoArrays[StPicoArrays::BTowHit]->GetEntriesFast(); }
112  static UInt_t numberOfBTofHits() { return picoArrays[StPicoArrays::BTofHit]->GetEntriesFast(); }
114  static UInt_t numberOfMtdHits() { return picoArrays[StPicoArrays::MtdHit]->GetEntriesFast(); }
116  static UInt_t numberOfBbcHits() { return picoArrays[StPicoArrays::BbcHit]->GetEntriesFast(); }
118  static UInt_t numberOfEpdHits() { return picoArrays[StPicoArrays::EpdHit]->GetEntriesFast(); }
120  static UInt_t numberOfFmsHits() { return picoArrays[StPicoArrays::FmsHit]->GetEntriesFast(); }
122  static UInt_t numberOfBEmcPidTraits() { return picoArrays[StPicoArrays::BEmcPidTraits]->GetEntriesFast(); }
124  static UInt_t numberOfBTofPidTraits() { return picoArrays[StPicoArrays::BTofPidTraits]->GetEntriesFast(); }
126  static UInt_t numberOfMtdPidTraits() { return picoArrays[StPicoArrays::MtdPidTraits]->GetEntriesFast(); }
128  static UInt_t numberOfTrackCovMatrices() { return picoArrays[StPicoArrays::TrackCovMatrix]->GetEntriesFast(); }
130  static UInt_t numberOfBEmcSmdEHits() { return picoArrays[StPicoArrays::BEmcSmdEHit]->GetEntriesFast(); }
132  static UInt_t numberOfBEmcSmdPHits() { return picoArrays[StPicoArrays::BEmcSmdPHit]->GetEntriesFast(); }
134  static UInt_t numberOfETofHits() { return picoArrays[StPicoArrays::ETofHit]->GetEntriesFast(); }
136  static UInt_t numberOfETofPidTraits() { return picoArrays[StPicoArrays::ETofPidTraits]->GetEntriesFast(); }
138  static UInt_t numberOfMcVertices() { return picoArrays[StPicoArrays::McVertex]->GetEntriesFast(); }
140  static UInt_t numberOfMcTracks() { return picoArrays[StPicoArrays::McTrack]->GetEntriesFast(); }
141 
143  void print() const;
145  static void printTracks();
147  static void printTriggers();
149  static void printBTowHits();
151  static void printBTofHits();
153  static void printMtdHits();
155  static void printFmsHits();
157  static void printBEmcPidTraits();
159  static void printBTofPidTraits();
161  static void printMtdPidTraits();
163  static void printTrackCovMatrices();
165  static void printBEmcSmdEHits();
167  static void printBEmcSmdPHits();
169  static void printETofHits();
171  static void printETofPidTraits();
173  static void printMcVertices();
175  static void printMcTracks();
176 
177 #if defined (__TFG__VERSION__)
178  static StPicoDst *instance() {return fgPicoDst;}
179 #endif /* __TFG__VERSION__ */
180 
181  private:
182 
184  static TClonesArray** picoArrays;
185 
186 #if defined (__TFG__VERSION__)
187  static StPicoDst *fgPicoDst;
188 #endif /* __TFG__VERSION__ */
189 };
190 
191 #endif
static void printTrackCovMatrices()
Print track covariance matrix info.
Definition: StPicoDst.cxx:243
static StPicoMcTrack * mcTrack(Int_t i)
Return pointer to i-th MC track.
Definition: StPicoDst.h:101
static StPicoMtdHit * mtdHit(Int_t i)
Return pointer to i-th mtd hit.
Definition: StPicoDst.h:75
static StPicoETofHit * etofHit(Int_t i)
Return pointer to i-th etof hit.
Definition: StPicoDst.h:95
static void printBEmcSmdEHits()
Print BEMC SMD eta info.
Definition: StPicoDst.cxx:261
static UInt_t numberOfFmsHits()
Return number of FMS hits.
Definition: StPicoDst.h:120
static void unset()
Reset the pointers to the TClonesArrays to 0.
Definition: StPicoDst.cxx:29
Holds information about MTD-matched track.
static UInt_t numberOfBTowHits()
Return number of BTow hits.
Definition: StPicoDst.h:110
static UInt_t numberOfEmcTriggers()
Return number of Emc triggers.
Definition: StPicoDst.h:106
static void printMtdPidTraits()
Print MTD PID trait info.
Definition: StPicoDst.cxx:225
static UInt_t numberOfMtdHits()
Return number of MTD hits.
Definition: StPicoDst.h:114
Holds information about Monte Carlo vertex.
Holds information about BEMC tower.
Definition: StPicoBTowHit.h:19
static StPicoBbcHit * bbcHit(Int_t i)
Return pointer to i-th bbc hit.
Definition: StPicoDst.h:77
static TClonesArray * picoArray(Int_t type)
Return pointer to the n-th TClonesArray.
Definition: StPicoDst.h:60
static void printETofHits()
Print ETOF hit info.
Definition: StPicoDst.cxx:297
static UInt_t numberOfBTofPidTraits()
Return number of BTof PID traits.
Definition: StPicoDst.h:124
static StPicoBTowHit * btowHit(Int_t i)
Return pointer to i-th btow hit.
Definition: StPicoDst.h:71
static StPicoBEmcSmdPHit * bemcSmdPHit(Int_t i)
Return pointer to i-th BEMC SMD phi hit.
Definition: StPicoDst.h:93
static void printMcTracks()
Print MC track info.
Definition: StPicoDst.cxx:80
static void printETofPidTraits()
Print ETOF PID trait info.
Definition: StPicoDst.cxx:315
static void printFmsHits()
Print FMS hit info.
Definition: StPicoDst.cxx:171
Class storing MTD trigger information including VPD, QT, MT101, TF201.
Holds BEMC SmdEta hit information.
static UInt_t numberOfETofPidTraits()
Return number of ETOF PID traits.
Definition: StPicoDst.h:136
static StPicoEpdHit * epdHit(Int_t i)
Return pointer to i-th epd hit.
Definition: StPicoDst.h:79
Holds information about FMS hit.
Definition: StPicoFmsHit.h:17
Holds EMC trigger information.
Hold information about eTOF-matched tracks.
static StPicoEvent * event()
Return pointer to current StPicoEvent (class holding the event wise information)
Definition: StPicoDst.h:63
static StPicoETofPidTraits * etofPidTraits(Int_t i)
Return pointer to i-th etof pidTraits.
Definition: StPicoDst.h:97
Holds information about track parameters.
Definition: StPicoTrack.h:35
~StPicoDst()
Destructor.
Definition: StPicoDst.h:52
static StPicoBEmcPidTraits * bemcPidTraits(Int_t i)
Return pointer to i-th emc pidTraits.
Definition: StPicoDst.h:83
static UInt_t numberOfEpdHits()
Return number of EPD hits.
Definition: StPicoDst.h:118
static StPicoBEmcSmdEHit * bemcSmdEHit(Int_t i)
Return pointer to i-th BEMC SMD eta hit.
Definition: StPicoDst.h:91
static void printBTowHits()
Print BTOW hit info.
Definition: StPicoDst.cxx:117
static UInt_t numberOfETofHits()
Return number of ETof hits.
Definition: StPicoDst.h:134
static StPicoMtdPidTraits * mtdPidTraits(Int_t i)
Return pointer to i-th mtd pidTraits.
Definition: StPicoDst.h:87
static UInt_t numberOfTracks()
Return number of tracks.
Definition: StPicoDst.h:104
static void set(TClonesArray **)
Set the pointers to the TClonesArrays.
Definition: StPicoDst.cxx:34
static void printBTofPidTraits()
Print BTOF PID trait info.
Definition: StPicoDst.cxx:207
static StPicoTrackCovMatrix * trackCovMatrix(Int_t i)
Return pointer to i-th track covariance matrix.
Definition: StPicoDst.h:89
static void printBEmcSmdPHits()
Print BEMC SMD phi info.
Definition: StPicoDst.cxx:279
Stores BTOF hit information.
Definition: StPicoBTofHit.h:18
Main class that keeps TClonesArrays with main classes.
Definition: StPicoDst.h:40
static UInt_t numberOfBbcHits()
Return number of BBC hits.
Definition: StPicoDst.h:116
Keep information about Barrel ElectroMagnetic Calorimeter (BEMC) matched tracks.
static StPicoEmcTrigger * emcTrigger(Int_t i)
Return pointer to i-th trigger data.
Definition: StPicoDst.h:67
static StPicoBTofHit * btofHit(Int_t i)
Return pointer to i-th btof hit.
Definition: StPicoDst.h:73
static StPicoMtdTrigger * mtdTrigger(Int_t i)
Return pointer to i-th MTD trigger data.
Definition: StPicoDst.h:69
static UInt_t numberOfMtdPidTraits()
Return number of MTD traits.
Definition: StPicoDst.h:126
Holds information about Monte Carlo track parameters.
Definition: StPicoMcTrack.h:32
static StPicoTrack * track(Int_t i)
Return pointer to i-th track.
Definition: StPicoDst.h:65
Holds information about the hit from MTD.
Definition: StPicoMtdHit.h:18
static void printBEmcPidTraits()
Print BEMC PID trait info.
Definition: StPicoDst.cxx:189
static StPicoMcVertex * mcVertex(Int_t i)
Return pointer to i-th MC vertex.
Definition: StPicoDst.h:99
static void printTracks()
Print track info.
Definition: StPicoDst.cxx:46
static void printBTofHits()
Print BTOF hit info.
Definition: StPicoDst.cxx:135
static UInt_t numberOfBEmcSmdEHits()
Return number of BEMC SMD eta hits.
Definition: StPicoDst.h:130
Stores eTOF hit information.
Definition: StPicoETofHit.h:19
Holds BEMC SmdPhi hit information.
static UInt_t numberOfBEmcSmdPHits()
Return number of BEMC SMD phi hits.
Definition: StPicoDst.h:132
Stores global information about the event.
Definition: StPicoEvent.h:24
StPicoDst()
Default constructor.
Definition: StPicoDst.h:50
static UInt_t numberOfTrackCovMatrices()
Return number of track covariance matrices.
Definition: StPicoDst.h:128
static void printTriggers()
Print trigger.
Definition: StPicoDst.cxx:97
static StPicoFmsHit * fmsHit(Int_t i)
Return pointer to i-th fms hit.
Definition: StPicoDst.h:81
static void printMcVertices()
Print MC vertex info.
Definition: StPicoDst.cxx:63
static UInt_t numberOfBTofHits()
Return number of BTof hits.
Definition: StPicoDst.h:112
static UInt_t numberOfMcVertices()
Return number of MC vertices.
Definition: StPicoDst.h:138
void print() const
Print information.
Definition: StPicoDst.cxx:39
static UInt_t numberOfMtdTriggers()
Return number of MTD triggers.
Definition: StPicoDst.h:108
static void printMtdHits()
Print MTD hit info.
Definition: StPicoDst.cxx:153
static UInt_t numberOfBEmcPidTraits()
Return number of BEMC PID traits.
Definition: StPicoDst.h:122
static StPicoBTofPidTraits * btofPidTraits(Int_t i)
Return pointer to i-th btof pidTraits.
Definition: StPicoDst.h:85
static UInt_t numberOfMcTracks()
Return number of MC tracks.
Definition: StPicoDst.h:140