Fork me on GitHub

source: svn/trunk/src/BFieldProp.cc@ 234

Last change on this file since 234 was 219, checked in by Xavier Rouby, 16 years ago

JetUtils.cc

File size: 9.9 KB
Line 
1/*
2 * ---- Delphes ----
3 * A Fast Simulator for general purpose LHC detector
4 * S. Ovyn ~~~~ severine.ovyn@uclouvain.be
5 *
6 * Center for Particle Physics and Phenomenology (CP3)
7 * Universite Catholique de Louvain (UCL)
8 * Louvain-la-Neuve, Belgium
9 * */
10
11#include "BFieldProp.h"
12#include<cmath>
13using namespace std;
14
15
16//------------------------------------------------------------------------------
17
18TrackPropagation::TrackPropagation(){
19 DET = new RESOLution();
20 init();
21}
22
23TrackPropagation::TrackPropagation(const string& DetDatacard){
24 DET = new RESOLution();
25 DET->ReadDataCard(DetDatacard);
26 init();
27}
28
29TrackPropagation::TrackPropagation(const RESOLution* DetDatacard){
30 DET= new RESOLution(*DetDatacard);
31 init();
32}
33
34TrackPropagation::TrackPropagation(const TrackPropagation & tp){
35 MAXITERATION = tp.MAXITERATION;
36 DET = new RESOLution(*(tp.DET));
37 R_max = tp.R_max; z_max = tp.z_max;
38 B_x = tp.B_x; B_y = tp.B_y; B_z = tp.B_z;
39 q = tp.q; phi_0 = tp.phi_0;
40 gammam= tp.gammam; omega = tp.omega;
41 r = tp.r; rr = tp.rr;
42 x_c = tp.x_c; y_c = tp.y_c;
43 R_c = tp.R_c; Phi_c = tp.Phi_c;
44 t = tp.t; t_z = tp.t_z; t_T = tp.t_T;
45 x_t = tp.x_t; y_t = tp.y_t; z_t = tp.z_t;
46 R_t = tp.R_t; Phi_t = tp.Phi_t;
47 Theta_t=tp.Theta_t; Eta_t = tp.Eta_t;
48 Px_t = tp.Px_t; Py_t = tp.Py_t; Pz_t = tp.Pz_t;
49 PT_t = tp.PT_t; p_t = tp.p_t; E_t = tp.E_t;
50 loop_overflow_counter = tp.loop_overflow_counter;
51}
52
53TrackPropagation& TrackPropagation::operator=(const TrackPropagation & tp) {
54 if(this==&tp) return *this;
55 MAXITERATION = tp.MAXITERATION;
56 DET = new RESOLution(*(tp.DET));
57 R_max = tp.R_max; z_max = tp.z_max;
58 B_x = tp.B_x; B_y = tp.B_y; B_z = tp.B_z;
59 q = tp.q; phi_0 = tp.phi_0;
60 gammam= tp.gammam; omega = tp.omega;
61 r = tp.r; rr = tp.rr;
62 x_c = tp.x_c; y_c = tp.y_c;
63 R_c = tp.R_c; Phi_c = tp.Phi_c;
64 t = tp.t; t_z = tp.t_z; t_T = tp.t_T;
65 x_t = tp.x_t; y_t = tp.y_t; z_t = tp.z_t;
66 R_t = tp.R_t; Phi_t = tp.Phi_t;
67 Theta_t=tp.Theta_t; Eta_t = tp.Eta_t;
68 Px_t = tp.Px_t; Py_t = tp.Py_t; Pz_t = tp.Pz_t;
69 PT_t = tp.PT_t; p_t = tp.p_t; E_t = tp.E_t;
70 loop_overflow_counter = tp.loop_overflow_counter;
71 return *this;
72}
73
74
75
76void TrackPropagation::init() {
77 MAXITERATION = 10000;
78 q= UNDEFINED; phi_0= UNDEFINED; gammam= UNDEFINED; omega=UNDEFINED; r=UNDEFINED;
79 x_c=UNDEFINED; y_c=UNDEFINED; R_c=UNDEFINED; Phi_c=UNDEFINED;
80 rr=UNDEFINED; t=UNDEFINED; t_z=UNDEFINED; t_T=UNDEFINED;
81 x_t=UNDEFINED; y_t=UNDEFINED; z_t=UNDEFINED;
82 R_t=UNDEFINED; Phi_t=UNDEFINED; Theta_t=UNDEFINED; Eta_t=UNDEFINED;
83 Px_t=UNDEFINED; Py_t=UNDEFINED; Pz_t=UNDEFINED; PT_t=UNDEFINED; p_t=UNDEFINED; E_t=UNDEFINED;
84
85 // DET has been initialised in the constructors
86 // magnetic field parameters
87 R_max = DET->TRACK_radius;
88 z_max = DET->TRACK_length/2.;
89 B_x = DET->TRACK_bfield_x;
90 B_y = DET->TRACK_bfield_y;
91 B_z = DET->TRACK_bfield_z;
92
93 loop_overflow_counter=0;
94}
95
96
97
98void TrackPropagation::Propagation(const TRootGenParticle *Part,TLorentzVector &momentum) {
99
100 q = Charge(Part->PID);
101 if(q==0) return;
102
103 if(R_max ==0) { cout << "ERROR: magnetic field has no lateral extention\n"; return;}
104 if(z_max==0) { cout << "ERROR: magnetic field has no longitudinal extention\n"; return;}
105
106 if (B_x== 0 && B_y== 0) { // faster if only B_z
107 if (B_z==0) return; // nothing to do
108
109 // initial conditions:
110 // p_X0 = Part->Px, p_Y0 = Part->Py, p_Z0 = Part->Pz, p_T0 = Part->PT;
111 // X_0 = Part->X, Y_0 = Part->Y, Z_0 = Part->Z;
112
113 // 1. initial transverse momentum p_{T0} : Part->PT
114 // initial transverse momentum direction \phi_0 = -atan(p_X0/p_Y0)
115 // relativistic gamma : gamma = E/mc² ; gammam = gamma \times m
116 // giration frequency \omega = q/(gamma m) B_z
117 // helix radius r = p_T0 / (omega gamma m)
118 phi_0 = -atan2(Part->Px,Part->Py);
119 gammam = Part->E; // here c==1
120 //cout << "gammam" << gammam << "\t gamma" << gammam/Part->M << endl;
121 omega = q * B_z /gammam;
122 r = Part->PT / (omega * gammam);
123
124 // 2. Helix parameters : center coordinates in transverse plane
125 // x_c = x_0 - r*cos(phi_0) and y_c = y_0 - r*sin(phi_0)
126 // R_c = \sqrt{x_c² + y_c²} and \Phi_c = atan{y_c/x_c}
127 x_c = Part->X - r*cos(phi_0); /// TEST !!
128 y_c = Part->Y - r*sin(phi_0);
129 R_c = sqrt(pow(x_c,2.) + pow(y_c,2.) );
130 Phi_c = atan2(y_c,x_c);
131
132 // 3. time evaluation t = min(t_T, t_z)
133 // t_T : time to exit from the sides
134 // t_T= [ Phi_c - phi_0 + atan( (R_max^2 - (R_c^2 + r^2))/(2rR_c) ) ]/omega
135 // t_z : time to exit from the front or the back
136 // t_z = gamma * m /p_z0 \times (-z_0 + z_max * sign(p_z0))
137 rr = sqrt( pow(R_c,2.) + pow(r,2.) ); // temp variable
138 t_T=0;
139 int sign_pz= (Part->Pz >0) ? 1 : -1;
140 t_z = gammam / Part->Pz * (-Part->Z + z_max*sign_pz ) ;
141 if ( fabs(R_c - r) > R_max || R_c + r < R_max ) t = t_z;
142 else {
143 t_T = (Phi_c - phi_0 + atan2( (R_max + rr)*(R_max - rr) , 2*r*R_c ) ) / omega;
144 t = min(t_T,t_z);
145 }
146
147 // 4. position in terms of x(t), y(t), z(t)
148 // x(t) = x_c + r cos (omega t + phi_0)
149 // y(t) = y_c + r sin (omega t + phi_0)
150 // z(t) = z_0 + (p_Z0/gammam) t
151 x_t = x_c + r * cos(omega * t + phi_0);
152 y_t = y_c + r * sin(omega * t + phi_0);
153 z_t = Part->Z + Part->Pz / gammam * t;
154
155 // 5. position in terms of Theta(t), Phi(t), R(t), Eta(t)
156 // R(t) = sqrt(x(t)² + y(t)²)
157 // Phi(t) = atan(y(t)/x(t))
158 // Theta(t) = atan(R(t)/z(t))
159 // Eta(t) = -ln tan (Theta(t)/2)
160 R_t = sqrt( pow(x_t,2.) + pow(y_t,2.) );
161 Phi_t = atan2( y_t, x_t);
162 if(R_t>0) {
163 Theta_t = acos( z_t / sqrt(z_t*z_t+ R_t*R_t));
164 Eta_t = - log(tan(Theta_t/2.));
165 } else{
166 Theta_t=0; Eta_t = 9999;
167 }
168
169 Px_t = - Part->PT * sin(omega*t + phi_0);
170 Py_t = Part->PT * cos(omega*t + phi_0);
171 Pz_t = Part->Pz;
172 PT_t = sqrt(Px_t*Px_t + Py_t*Py_t);
173 p_t = sqrt(PT_t*PT_t + Pz_t*Pz_t);
174 E_t=sqrt(Part->M*Part->M +p_t);
175 //if(p_t != fabs(Pz_t) ) Eta_t = log( (p_t+Pz_t)/(p_t-Pz_t) )/2.;
176 //if(p_t>0) Theta_t = acos(Pz_t/p_t);
177 momentum.SetPxPyPzE(Px_t,Py_t,Pz_t,E_t);
178
179// test zone ---
180/*
181 cout << cos(atan(R_t/z_t)) << "\t" << cos(Theta_t) << "\t" << cos(momentum.Theta()) << "\t" << Pz_t/temp_p << endl;
182 double Eta_t1 = log( (E+Pz_t)/(E-Pz_t) )/2.;
183 double Eta_t2 = log( (temp_p+Pz_t)/(temp_p-Pz_t) )/2.;
184 if(0 && fabs(Eta_t -Eta_t2)>1e-310) {
185 cout << "ERROR-BUG: Eta_t != Eta_t2\n";
186 cout << "Eta_t= " << Eta_t << "\t Eta_t1= " << Eta_t1 << "\t Eta_t2= " << Eta_t2 << endl;
187 }
188
189 double R_t2 = sqrt( pow(R_c,2.) + pow(r,2.) + 2*r*R_c*cos(phi_0 + omega*t - Phi_c) ); // cross-check
190 if(fabs(R_t - R_t2) > 1e-7)
191 cout << "ERROR-BUG: R_t != R_t2: R_t=" << R_t << " R_t2=" << R_t2 << " R_t - R_t2 =" << R_t - R_t2 << endl;
192 if( fabs(E - gammam) > 1e-3 ) {
193 cout << "ERROR-BUG: energy is not conserved in src/BFieldProp.cc\n";
194 cout << "E - momentum.E() = " << fabs(E - momentum.E()) << " gammam - E " << fabs(gammam -E) << endl; }
195 if( fabs(PT_t - Part->PT) > 1e-10 ) {
196 cout << "ERROR-BUG: PT is not conversed in src/BFieldProp.cc. ";
197 cout << "(at " << 100*(PT_t - Part->PT) << "%)\n";
198 }
199 if(momentum.Pz() != Pz_t)
200 cout << "ERROR-BUG: Pz is not conserved in src/BFieldProp.cc\n";
201
202 double temp_p0=sqrt(Part->PT*Part->PT + Part->Pz*Part->Pz);
203 if(fabs( (temp_p-temp_p0)*(temp_p+temp_p0) )>1e-10 ) {
204 cout << "ERROR-BUG: momentum |vec{p}| is not conserved in src/BFieldProp.cc\n";
205 cout << temp_p << "\t" << temp_p0 << endl;
206 }
207
208 // if x_c == y_c ==0 (set it by hand!), easy cross-check
209 //cout << "tan(phi_p)= " << momentum.Py()/momentum.Px() << "\t -1/tan(phi_x)= " << -x_t/y_t << endl;
210*/
211
212 } else { // if B_x or B_y are non zero: longer computation
213
214 float Xvertex1 = Part->X;
215 float Yvertex1 = Part->Y;
216 float Zvertex1 = Part->Z;
217
218 //out of tracking coverage?
219 if(sqrt(Xvertex1*Xvertex1+Yvertex1*Yvertex1) > R_max){return;}
220 if(fabs(Zvertex1) > z_max){return;}
221
222 double px = Part->Px / 0.003;
223 double py = Part->Py / 0.003;
224 double pz = Part->Pz / 0.003;
225 double pt = Part->PT / 0.003; // sqrt(px*px+py*py);
226 double p = sqrt(pz*pz + pt*pt); //sqrt(px*px+py*py+pz*pz);
227
228 double M = Part->M;
229 double vx = px/M;
230 double vy = py/M;
231 double vz = pz/M;
232 double qm = q/M;
233
234 double ax = qm*(B_z*vy - B_y*vz);
235 double ay = qm*(B_x*vz - B_z*vx);
236 double az = qm*(B_y*vx - B_x*vy);
237 double dt = 1/p;
238 if(pt<266 && vz < 0.0012) dt = fabs(0.001/vz); // ?????
239
240 double xold=Xvertex1; double x=xold;
241 double yold=Yvertex1; double y=yold;
242 double zold=Zvertex1; double z=zold;
243
244 double VTold = pt/M; //=sqrt(vx*vx+vy*vy);
245
246 unsigned int k = 0;
247 double VTratio=0;
248 double R_max2 = R_max*R_max;
249 double r2=0; // will be x*x+y*y
250
251 while(k < MAXITERATION){
252 k++;
253
254 vx += ax*dt;
255 vy += ay*dt;
256 vz += az*dt;
257
258 VTratio = VTold/sqrt(vx*vx+vy*vy);
259 vx *= VTratio;
260 vy *= VTratio;
261
262 ax = qm*(B_z*vy - B_y*vz);
263 ay = qm*(B_x*vz - B_z*vx);
264 az = qm*(B_y*vx - B_x*vy);
265
266 x += vx*dt;
267 y += vy*dt;
268 z += vz*dt;
269 r2 = x*x + y*y;
270
271 if( r2 > R_max2 ){
272 x /= r2/R_max2;
273 y /= r2/R_max2;
274 break;
275 }
276 if( fabs(z)>z_max)break;
277
278 xold = x;
279 yold = y;
280 zold = z;
281 } // while loop
282
283 if(k == MAXITERATION) loop_overflow_counter++;
284 //cout << "too short loop in " << loop_overflow_counter << " cases" << endl;
285
286 if(x!=0 && y!=0 && z!=0) {
287 float Theta = atan2(sqrt(r2),z);
288 double eta = -log(tan(Theta/2.));
289 double phi = atan2(y,x);
290 momentum.SetPtEtaPhiE(Part->PT,eta,phi,Part->E);
291 }
292
293 } // if b_x or b_y non zero
294}
Note: See TracBrowser for help on using the repository browser.