forked from vmc-project/geant3
-
Notifications
You must be signed in to change notification settings - Fork 0
/
README_geane
266 lines (207 loc) · 11.5 KB
/
README_geane
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
The Geane package allows the user to calculate the
average trajectories of particles and to calculate the
transport matrix as well as the propagated error covariance
matrix. Geane is a set of routines worked out
by the European Muon Collaboration [1] and it is integrated
to the GEANT3 system [2].
See geant3/doc/gedoc for more details.
The user should invoke the routine Ertrak and Eufill or Eufilp or
Eufilv, to carry out the tracking.
In addition to this a series of utilities are available for the user
(e.g. to transform the track representation from one
system to another or to carry out 5 X 5 matrix multiplication
in an optimal way).
Track variables, Representations
--------------------------------
The particle trajectory is characterized by 5 independent
variables as a function of one parameter (e.g. the pathlength).
Among the 5 variables 1 is related to the curvature (to the absolute
value of the momentum, p), 2 are related to the direction of the
particle and the other 2 are related to the spatial location.
The most usual representation of these 5 parameters are:
I. 1/p, lambda, phi, y_perp, z_perp
where lambda and phi are the dip and azimuthal angles related
to the momentum components in the following way:
p_x = p cos(lambda) cos(phi)
p_y = p cos(lambda) sin(phi)
p_z = p sin(lambda)
y_perp and z_perp are the coordinates of the trajectory in a
local orthonormal reference frame with the x_perp axis along the
particle direction, the y_perp being parallel to the x-y plane.
This representation is usually applied in the overall reference
frame. (In the EMC code this reference frame is labelled by 'SC'
since the overall system was identified with that of the Streamer
Chamber.)
II. 1/p, y', z', y, z
where y'=dy/dx and z'=dz/dx. This representation is particularly
useful in fixed target experiments, where the trajectory is evaluated
on successive parallel planes (which are perpendicular to the x-axis).
(In the EMC code this representation is labelled by 'SP' since a
convenient mathematical description of a trajectory being approxima-
tely parallel to the x-axis is a 'spline'.)
III. 1/p, v', w', v, w
where v'=dv/du and w'=dw/du in an orthonormal coordinate system with
axis u, v and w. This representation is paricularly useful when the
trajectory has to be evaluated on different detector planes
in a colliding beam experiment, where the planes can take a great
variety of directions.(In the EMC code this representation is
labelled by 'SD' as System of Detection.)
Of course, all the above representations of the trajectory
are equivalent and one can go from one representation to the
other by calculating the corresponding Jacobian. These Jacobians
are provided by the following EMC routines:
S/R Trscsp from I to II
S/R Trspsc from II to I
S/R Trscsd from I to III
S/R Trsdsc from III to I
User Interface
---------------------------------------------------
Ertrak(const Float_t *x1, const Float_t *p1,
const Float_t *x2, const Float_t *p2,
Int_t ipa, Option_t *chopt)
---------------------------------------------------
Performs the tracking of the track from point x1 to point x2
(Before calling this routine the user should call Eufil(l/p/v)
x1 - Starting coordinates (Cartesian)
p1 - Starting 3-momentum (Cartesian)
x2 - Final coordinates (Cartesian)
p2 - Final 3-momentum (Cartesian)
ipa - Particle code (a la GEANT) of the track
chopt
'B' 'Backward tracking' - i.e. energy loss
added to the current energy
'E' 'Exact' calculation of errors assuming
helix (i.e. pathlength not
assumed as infinitesimal)
'L' Tracking upto prescribed Lengths reached
'M' 'Mixed' prediction (not yet coded)
'O' Tracking 'Only' without calculating errors
'P' Tracking upto prescribed Planes reached
'V' Tracking upto prescribed Volumes reached
'X' Tracking upto prescribed Point approached
------------------------------------------
Eufill(Int_t n,Float_t *ein,Float_t *xlf);
------------------------------------------
User routine to set the input values for chopt = 'L'
n Number of predictions where to store results
ein Input error matrix
xlf Defines the tracklengths which if passed the
result should be stored
----------------------------------------------------------------
Eufilp(const Int_t n, Float_t *ein, Float_t *pli, Float_t *plf);
----------------------------------------------------------------
User routine to set the input values for chopt = 'P'
n Number of predictions where to store results
ein Input error matrix (in the 'Plane' system )
pli Defines the start plane
PLI(3,1) - and
PLI(3,2) - 2 unit vectors in the plane
plf Defines the end plane
PLF(3,1,I) - and
PLF(3,2,I) - 2 unit vectors in the plane
PLF(3,3,I) - point on the plane
at intermediate point I
---------------------------------------------------------------------
Eufilv(Int_t n, Float_t *ein, Char_t *namv, Int_t *numv,Int_t *iovl);
---------------------------------------------------------------------
User routine to set the input values for chopt = 'V'
n Number of predictions where to store results
ein Input error matrix
cnamv Volume name of the prediction
numv Volume number (if 0 = all volumes)
iovl = 1 prediction when entering in the volume
= 2 prediction when leaving the volume
--------------------------------------------------------------------------
void Trscsd(Float_t *pc,Float_t *rc,Float_t *pd,Float_t *rd,Float_t *h,
Float_t *ch,Int_t *ierr,Float_t *spu,Float_t *dj,Float_t *dk);
--------------------------------------------------------------------------
transforms error matrix from sc variables (1/p,lambda,phi,yt,zt)
to variables (1/p,v',w',v,w)
pc[3] 1/p,lambda,phi input
pd[3] 1/p,v',w' output
h [3] magnetic field input
rc[15] error matrix in sc variables input (triangle)
rd[15] error matrix in 1/p,v',w',v,w output (triangle)
ch charge of particle input
charge and magnetic field are needed
for correlation terms (v',yt),(v',zt),(w',yt),(w',zt)
these correlation terms appear because rc is assumed
to be the error matrix for fixed s (path length)
and rd for fixed u
dj[3] unit vector in v-direction
dk[3] unit vector in w-direction of detector system
ierr = 1 particle moves perpendicular to u-axis
( v',w' are not defined )
spu sign of u-component of particle momentum output
--------------------------------------------------------------------------
void Trsdsc(Float_t *pd,Float_t *rd,Float_t *pc,Float_t *rc,Float_t *h,
Float_t *ch,Int_t *ierr,Float_t *spu,Float_t *dj,Float_t *dk);
--------------------------------------------------------------------------
transforms error matrix from variables (1/p,v',w',v,w)
to sc variables (1/p,lambda,phi,yt,zt)
pd[3] 1/p,v',w' input
pc[3] 1/p,lambda,phi output
h[3] magnetic field input
rd[15] error matrix in 1/p,v',w',v,w input (triangle)
rc[15] error matrix in sc variables output (triangle)
ch charge of particle input
charge and magnetic field are needed
for correlation terms (lambda,v),(lambda,w),(phi,v),(phi,w)
these correlation terms appear because rc is assumed
to be the error matrix for fixed s (path length)
and rd for fixed u
dj[3] unit vector in v-direction
dk[3] unit vector in w-direction of detector system
ierr not used
spu sign of u-component of particle momentum input
--------------------------------------------------------------------------
void Trscsp(Float_t *ps,Float_t *rs,Float_t *pc,Float_t *rc,Float_t *h,
Float_t *ch,Int_t *ierr,Float_t *spx);
--------------------------------------------------------------------------
transforms error matrix from sc variables (1/p,lambda,phi,yt,zt)
to spline variables (1/p,y',z',y,z)
pc[3] 1/p,lambda,phi input
ps[3] 1/p,y',z' output
h[3] magnetic field input
rc[15] error matrix in sc variables input (triangle)
rs[15] error matrix in spline variables output (triangle)
ch charge of particle input
charge and magnetic field are needed
for correlation terms (y',yt),(y',zt),(z',yt),(z',zt)
these correlation terms appear because rc is assumed
to be the error matrix for fixed s (path length)
and rs for fixed x
ierr = 1 particle moves perpendicular to x-axis
( y',z' are not defined )
spx sign of x-component of particle momentum output
--------------------------------------------------------------------------
void Trspsc(Float_t *ps,Float_t *rs,Float_t *pc,Float_t *rc,Float_t *h,
Float_t *ch,Int_t *ierr,Float_t *spx);
--------------------------------------------------------------------------
transforms error matrix from spline variables (1/p,y',z',y,z)
to sc variables (1/p,lambda,phi,yt,zt)
ps[3] 1/p,y',z' input
pc[3] 1/p,lambda,phi output
h[3] magnetic field input
rs[15] error matrix in spline variables input (triangle)
rc[15] error matrix in sc variables output (triangle)
ch charge of particle input
charge and magnetic field are needed
for correlation terms (lambda,y),(lambda,z),(phi,y),(phi,z)
these correlation terms appear because rc is assumed
to be the error matrix for fixed s (path length)
and rs for fixed x
ierr not used
spx sign of x-component of particle momentum input
Examples of Application
=======================
1. Representing the trajectory at another point
2. Joining track elements in different parts of the detector
3. Prediction of the trajectory
4. Fitting trajectory parameters
See geant3/doc/gedoc for details.
[1] W.Wittek, EMC Internal Reports: EMC/80/15, EMCSW/80/39,
EMCSW/81/13, EMCSW/81/18
A.Haas, The EMC Utility Package: UTIL42
[2] R.Brun, F.Bruyant, M.Maire, A.C.McPherson, P.Zanarini
DD/EE/84-1, May 1986