Logo ROOT   6.08/07
Reference Guide
JetEvent.cxx
Go to the documentation of this file.
1 // A JetEvent emulates 2 detectors A and B producing each
2 // a TClonesArray of Hit objects.
3 // A TClonesArray of Track objects is built with Hits objects
4 // of detectors A and B. Eack Track object has a TRefArray of hits.
5 // A TClonesArray of Jets is made with a subset of the Track objects
6 // also stored in a TRefArray.
7 // see $ROOTSYS/tutorials/jets.C for an example creating a Tree
8 // with JetEvents.
9 
10 #include "TMath.h"
11 #include "TRandom.h"
12 #include "JetEvent.h"
13 
18 
23 
24 ////////////////////////////////////////////////////////////////////////////////
25 /// Create a JetEvent object.
26 /// When the constructor is invoked for the first time, the class static
27 /// variables fgxxx are 0 and the TClonesArray fgxxx are created.
28 
30 {
31  if (!fgTracks) fgTracks = new TClonesArray("Track", 100);
32  if (!fgJets) fgJets = new TClonesArray("Jet", 10);
33  if (!fgHitsA) fgHitsA = new TClonesArray("Hit", 10000);
34  if (!fgHitsB) fgHitsB = new TClonesArray("Hit", 1000);
35  fJets = fgJets;
36  fTracks = fgTracks;
37  fHitsA = fgHitsA;
38  fHitsB = fgHitsB;
39 }
40 
41 ////////////////////////////////////////////////////////////////////////////////
42 
44 {
45  Reset();
46 }
47 
48 ////////////////////////////////////////////////////////////////////////////////
49 ///Build one event
50 
51 void JetEvent::Build(Int_t jetm, Int_t trackm, Int_t hitam, Int_t hitbm) {
52  //Save current Object count
53  Int_t ObjectNumber = TProcessID::GetObjectCount();
54  Clear();
55 
56  Hit *hit;
57  Track *track;
58  Jet *jet;
59  fNjet = fNtrack = fNhitA = fNhitB = 0;
60 
61  fVertex.SetXYZ(gRandom->Gaus(0,0.1),
62  gRandom->Gaus(0,0.2),
63  gRandom->Gaus(0,10));
64 
65  Int_t njets = (Int_t)gRandom->Gaus(jetm,1); if (njets < 1) njets = 1;
66  for (Int_t j=0;j<njets;j++) {
67  jet = AddJet();
68  jet->fPt = gRandom->Gaus(0,10);
69  jet->fPhi = 2*TMath::Pi()*gRandom->Rndm();
70  Int_t ntracks = (Int_t)gRandom->Gaus(trackm,3); if (ntracks < 1) ntracks = 1;
71  for (Int_t t=0;t<ntracks;t++) {
72  track = AddTrack();
73  track->fPx = gRandom->Gaus(0,1);
74  track->fPy = gRandom->Gaus(0,1);
75  track->fPz = gRandom->Gaus(0,5);
76  jet->fTracks.Add(track);
77  Int_t nhitsA = (Int_t)gRandom->Gaus(hitam,5);
78  for (Int_t ha=0;ha<nhitsA;ha++) {
79  hit = AddHitA();
80  hit->fX = 10000*j + 100*t +ha;
81  hit->fY = 10000*j + 100*t +ha+0.1;
82  hit->fZ = 10000*j + 100*t +ha+0.2;
83  track->fHits.Add(hit);
84  }
85  Int_t nhitsB = (Int_t)gRandom->Gaus(hitbm,2);
86  for (Int_t hb=0;hb<nhitsB;hb++) {
87  hit = AddHitB();
88  hit->fX = 20000*j + 100*t +hb+0.3;
89  hit->fY = 20000*j + 100*t +hb+0.4;
90  hit->fZ = 20000*j + 100*t +hb+0.5;
91  track->fHits.Add(hit);
92  }
93  track->fNhit = nhitsA + nhitsB;
94  }
95  }
96  //Restore Object count
97  //To save space in the table keeping track of all referenced objects
98  //we assume that our events do not address each other. We reset the
99  //object count to what it was at the beginning of the event.
100  TProcessID::SetObjectCount(ObjectNumber);
101 }
102 
103 
104 ////////////////////////////////////////////////////////////////////////////////
105 /// Add a new Jet to the list of tracks for this event.
106 
108 {
109  TClonesArray &jets = *fJets;
110  Jet *jet = new(jets[fNjet++]) Jet();
111  return jet;
112 }
113 
114 
115 ////////////////////////////////////////////////////////////////////////////////
116 /// Add a new track to the list of tracks for this event.
117 
119 {
120  TClonesArray &tracks = *fTracks;
121  Track *track = new(tracks[fNtrack++]) Track();
122  return track;
123 }
124 
125 
126 ////////////////////////////////////////////////////////////////////////////////
127 /// Add a new hit to the list of hits in detector A
128 
130 {
131  TClonesArray &hitsA = *fHitsA;
132  Hit *hit = new(hitsA[fNhitA++]) Hit();
133  return hit;
134 }
135 
136 ////////////////////////////////////////////////////////////////////////////////
137 /// Add a new hit to the list of hits in detector B
138 
140 {
141  TClonesArray &hitsB = *fHitsB;
142  Hit *hit = new(hitsB[fNhitB++]) Hit();
143  return hit;
144 }
145 
146 ////////////////////////////////////////////////////////////////////////////////
147 
149 {
150  fJets->Clear(option);
151  fTracks->Clear(option);
152  fHitsA->Clear(option);
153  fHitsB->Clear(option);
154 }
155 
156 ////////////////////////////////////////////////////////////////////////////////
157 /// Static function to reset all static objects for this event
158 
160 {
161  delete fgJets; fgJets = 0;
162  delete fgTracks; fgTracks = 0;
163  delete fgHitsA; fgHitsA = 0;
164  delete fgHitsB; fgHitsB = 0;
165 }
166 
167 
168 
169 
170 
171 
TRefArray fTracks
Definition: JetEvent.h:54
Jet * AddJet()
Add a new Jet to the list of tracks for this event.
Definition: JetEvent.cxx:107
TRefArray fHits
Definition: JetEvent.h:37
Int_t fNtrack
Definition: JetEvent.h:69
Double_t fPt
Definition: JetEvent.h:52
Int_t fNhitA
Definition: JetEvent.h:70
static TClonesArray * fgJets
Definition: JetEvent.h:77
Float_t fPx
Definition: JetEvent.h:33
void Add(TObject *obj)
Definition: TRefArray.h:84
const char Option_t
Definition: RtypesCore.h:62
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:235
TClonesArray * fTracks
Definition: JetEvent.h:73
JetEvent()
Create a JetEvent object.
Definition: JetEvent.cxx:29
int Int_t
Definition: RtypesCore.h:41
Int_t fNhit
Definition: JetEvent.h:36
TVector3 fVertex
Definition: JetEvent.h:67
Track * AddTrack()
Add a new track to the list of tracks for this event.
Definition: JetEvent.cxx:118
void Build(Int_t jetm=3, Int_t trackm=10, Int_t hitam=100, Int_t hitbm=10)
Build one event.
Definition: JetEvent.cxx:51
void Clear(Option_t *option="")
Definition: JetEvent.cxx:148
Float_t fY
Definition: JetEvent.h:20
friend class TClonesArray
Definition: TObject.h:198
static TClonesArray * fgHitsA
Definition: JetEvent.h:79
TClonesArray * fHitsA
Definition: JetEvent.h:74
TClonesArray * fHitsB
Definition: JetEvent.h:75
TClonesArray * fJets
Definition: JetEvent.h:72
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition: TVector3.h:235
static void SetObjectCount(UInt_t number)
static function to set the current referenced object count fgNumber is incremented every time a new o...
Definition: TProcessID.cxx:392
static TClonesArray * fgTracks
Definition: JetEvent.h:78
Definition: JetEvent.h:16
Float_t fX
Definition: JetEvent.h:19
virtual void Clear(Option_t *option="")
Clear the clones array.
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:512
Int_t fNjet
Definition: JetEvent.h:68
Float_t fPz
Definition: JetEvent.h:35
virtual ~JetEvent()
Definition: JetEvent.cxx:43
Hit * AddHitB()
Add a new hit to the list of hits in detector B.
Definition: JetEvent.cxx:139
Int_t fNhitB
Definition: JetEvent.h:71
R__EXTERN TRandom * gRandom
Definition: TRandom.h:66
Double_t Pi()
Definition: TMath.h:44
#define ClassImp(name)
Definition: Rtypes.h:279
Double_t fPhi
Definition: JetEvent.h:53
Float_t fPy
Definition: JetEvent.h:34
void Reset(Option_t *option="")
Static function to reset all static objects for this event.
Definition: JetEvent.cxx:159
Hit * AddHitA()
Add a new hit to the list of hits in detector A.
Definition: JetEvent.cxx:129
An array of clone (identical) objects.
Definition: TClonesArray.h:32
static TClonesArray * fgHitsB
Definition: JetEvent.h:80
static UInt_t GetObjectCount()
Return the current referenced object count fgNumber is incremented every time a new object is referen...
Definition: TProcessID.cxx:294
Definition: JetEvent.h:30
Definition: JetEvent.h:49
Float_t fZ
Definition: JetEvent.h:21