]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSgeom.cxx
Introducing new frame 1099
[u/mrichter/AliRoot.git] / ITS / AliITSgeom.cxx
CommitLineData
4c039060 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/*
17$Log$
ad0e60d9 18Revision 1.2 1999/09/29 09:24:20 fca
19Introduction of the Copyright and cvs Log
20
4c039060 21*/
22
58005f18 23///////////////////////////////////////////////////////////////////////
24// ITS geometry manimulaiton routines. //
25// Created April 15 1999. //
26// version: 0.0.0 //
27// By: Bjorn S. Nilsen //
28// version: 0.0.1 //
29// Updated May 27 1999. //
30// Added Cylinderical random and global based changes. //
31// Added function PrintComparison. //
32///////////////////////////////////////////////////////////////////////
33#include <iostream.h>
34#include <fstream.h>
35#include <iomanip.h>
36#include <stdio.h>
37#include "AliITSgeom.h"
38#include "TRandom.h"
39
40ClassImp(AliITSgeom)
41
42//_____________________________________________________________________
43AliITSgeom::AliITSgeom(){
44////////////////////////////////////////////////////////////////////////
45// The default constructor for the AliITSgeom class. It, by default,
46// sets fNlayers to zero and zeros all pointers.
47////////////////////////////////////////////////////////////////////////
48 // Default constructor.
49 // Do not allocate anything zero everything
50 fNlayers = 0;
51 fNlad = 0;
52 fNdet = 0;
53 fg = 0;
54 fShape = 0;
55 return;
56}
57
58//_____________________________________________________________________
59AliITSgeom::~AliITSgeom(){
60////////////////////////////////////////////////////////////////////////
61// The destructor for the AliITSgeom class. If the arrays fNlad,
62// fNdet, or fg have had memory allocated to them, there pointer values
63// are non zero, then this memory space is freed and they are set
64// to zero. In addition, fNlayers is set to zero. The destruction of
65// TObjArray fShape is, by default, handled by the TObjArray destructor.
66////////////////////////////////////////////////////////////////////////
67 // Default destructor.
68 // if arrays exist delet them. Then set everything to zero.
69 if(fg!=0){
70 for(Int_t i=0;i<fNlayers;i++) delete[] fg[i];
71 delete[] fg;
72 } // end if fg!=0
73 if(fNlad!=0) delete[] fNlad;
74 if(fNdet!=0) delete[] fNdet;
75 fNlayers = 0;
76 fNlad = 0;
77 fNdet = 0;
78 fg = 0;
79 return;
80}
81
82//_____________________________________________________________________
83AliITSgeom::AliITSgeom(const char *filename){
84////////////////////////////////////////////////////////////////////////
85// The constructor for the AliITSgeom class. All of the data to fill
86// this structure is read in from the file given my the input filename.
87////////////////////////////////////////////////////////////////////////
88 FILE *pf;
89 Int_t i;
90 ITS_geom *g;
91 Int_t l,a,d;
92 Float_t x,y,z,o,p,q,r,s,t;
ad0e60d9 93 Double_t oor,pr,qr,rr,sr,tr; // Radians
58005f18 94 Double_t lr[9];
95 Double_t si; // sin(angle)
96 Double_t PI = TMath::Pi(), byPI = PI/180.;
97
98 pf = fopen(filename,"r");
99
100 fNlayers = 6; // set default number of ladders
101 fNlad = new Int_t[fNlayers];
102 fNdet = new Int_t[fNlayers];
103 // find the number of laders and detectors in this geometry.
104 for(i=0;i<fNlayers;i++){fNlad[i]=fNdet[i]=0;} // zero out arrays
105 for(;;){ // for ever loop
106 i = fscanf(pf,"%d %d %d %f %f %f %f %f %f %f %f %f",
107 &l,&a,&d,&x,&y,&z,&o,&p,&q,&r,&s,&t);
108 if(i==EOF) break;
109 if(l<1 || l>fNlayers) {
110 printf("error in file %s layer=%d min is 1 max is %d/n",
111 filename,l,fNlayers);
112 continue;
113 }// end if l
114 if(fNlad[l-1]<a) fNlad[l-1] = a;
115 if(fNdet[l-1]<d) fNdet[l-1] = d;
116 } // end for ever loop
117 // counted the number of laders and detectors now allocate space.
118 fg = new ITS_geom* [fNlayers];
119 for(i=0;i<fNlayers;i++){
120 fg[i] = 0;
121 l = fNlad[i]*fNdet[i];
122 fg[i] = new ITS_geom[l]; // allocate space for transforms
123 } // end for i
124
125 // Set up Shapes for a default configuration of 6 layers.
126 fShape = new TObjArray;
127 AddShape((TObject *) new AliITSgeomSPD()); // shape 0
128 AddShape((TObject *) new AliITSgeomSDD()); // shape 1
129 AddShape((TObject *) new AliITSgeomSPD()); // shape 2
130
131 // prepair to read in transforms
132 rewind(pf); // start over reading file
133 for(;;){ // for ever loop
134 i = fscanf(pf,"%d %d %d %f %f %f %f %f %f %f %f %f",
135 &l,&a,&d,&x,&y,&z,&o,&p,&q,&r,&s,&t);
136 if(i==EOF) break;
137 if(l<1 || l>fNlayers) {
138 printf("error in file %s layer=%d min is 1 max is %d/n",
139 filename,l,fNlayers);
140 continue;
141 }// end if l
142 l--; a--; d--; // shift layer, lader, and detector counters to zero base
143 i = d + a*fNdet[l]; // position of this detector
144 g = &(fg[l][i]);
145
ad0e60d9 146 oor = byPI*o;
58005f18 147 pr = byPI*p;
148 qr = byPI*q;
149 rr = byPI*r;
150 sr = byPI*s;
151 tr = byPI*t;
152
153 g->fx0 = x;
154 g->fy0 = y;
155 g->fz0 = z;
ad0e60d9 156 si = sin(oor);if(o== 90.0) si = +1.0;
58005f18 157 if(o==270.0) si = -1.0;
158 if(o== 0.0||o==180.) si = 0.0;
159 lr[0] = si * cos(pr);
160 lr[1] = si * sin(pr);
ad0e60d9 161 lr[2] = cos(oor);if(o== 90.0||o==270.) lr[2] = 0.0;
58005f18 162 if(o== 0.0) lr[2] = +1.0;
163 if(o==180.0) lr[2] = -1.0;
164 si = sin(qr);if(q== 90.0) si = +1.0;
165 if(q==270.0) si = -1.0;
166 if(q== 0.0||q==180.) si = 0.0;
167 lr[3] = si * cos(rr);
168 lr[4] = si * sin(rr);
169 lr[5] = cos(qr);if(q== 90.0||q==270.) lr[5] = 0.0;
170 if(q== 0.0) lr[5] = +1.0;
171 if(q==180.0) lr[5] = -1.0;
172 si = sin(sr);if(r== 90.0) si = +1.0;
173 if(r==270.0) si = -1.0;
174 if(r== 0.0||r==180.) si = 0.0;
175 lr[6] = si * cos(tr);
176 lr[7] = si * sin(tr);
177 lr[8] = cos(sr);if(r== 90.0||r==270.0) lr[8] = 0.0;
178 if(r== 0.0) lr[8] = +1.0;
179 if(r==180.0) lr[8] = -1.0;
180 // Normalize these elements
181 for(a=0;a<3;a++){// reuse float si and integers a and d.
182 si = 0.0;
183 for(d=0;d<3;d++) si += lr[3*a+d]*lr[3*a+d];
184 si = TMath::Sqrt(1./si);
185 for(d=0;d<3;d++) g->fr[3*a+d] = lr[3*a+d] = si*lr[3*a+d];
186 } // end for a
187 // get angles from matrix up to a phase of 180 degrees.
ad0e60d9 188 oor = atan2(lr[7],lr[8]);if(oor<0.0) oor += 2.0*PI;
58005f18 189 pr = asin(lr[2]); if(pr<0.0) pr += 2.0*PI;
190 qr = atan2(lr[3],lr[0]);if(qr<0.0) qr += 2.0*PI;
ad0e60d9 191 g->frx = oor;
58005f18 192 g->fry = pr;
193 g->frz = qr;
194 // l = layer-1 at this point.
195 if(l==0||l==1) g->fShapeIndex = 0; // SPD's
196 else if(l==2||l==3) g->fShapeIndex = 1; // SDD's
197 else if(l==4||l==5) g->fShapeIndex = 2; // SSD's
198 } // end for ever loop
199 fclose(pf);
200}
201
202//________________________________________________________________________
203AliITSgeom::AliITSgeom(AliITSgeom &source){
204////////////////////////////////////////////////////////////////////////
205// The copy constructor for the AliITSgeom class. It calls the
206// = operator function. See the = operator function for more details.
207////////////////////////////////////////////////////////////////////////
208 source = *this; // Just use the = operator for now.
209 return;
210}
211
212//________________________________________________________________________
213void AliITSgeom::operator=(AliITSgeom &source){
214////////////////////////////////////////////////////////////////////////
215// The = operator function for the AliITSgeom class. It makes an
216// independent copy of the class in such a way that any changes made
217// to the copied class will not affect the source class in any way.
218// This is required for many ITS alignment studies where the copied
219// class is then modified by introducing some misalignment.
220////////////////////////////////////////////////////////////////////////
221 Int_t i,j,k;
222
223 if(this == &source) return; // don't assign to ones self.
224
225 // if there is an old structure allocated delete it first.
226 if(fg != 0){
227 for(i=0;i<fNlayers;i++) delete[] fg[i];
228 delete[] fg;
229 } // end if fg != 0
230 if(fNlad != 0) delete[] fNlad;
231 if(fNdet != 0) delete[] fNdet;
232
233 fNlayers = source.fNlayers;
234 fNlad = new Int_t[fNlayers];
235 for(i=0;i<fNlayers;i++) fNlad[i] = source.fNlad[i];
236 fNdet = new Int_t[fNlayers];
237 for(i=0;i<fNlayers;i++) fNdet[i] = source.fNdet[i];
238 fShape = new TObjArray(*(source.fShape));//This does not make a proper copy.
239 fg = new ITS_geom* [fNlayers];
240 for(i=0;i<fNlayers;i++){
241 fg[i] = new ITS_geom[fNlad[i]*fNdet[i]];
242 for(j=0;j<(fNlad[i]*fNdet[i]);j++){
243 fg[i][j].fShapeIndex = source.fg[i][j].fShapeIndex;
244 fg[i][j].fx0 = source.fg[i][j].fx0;
245 fg[i][j].fy0 = source.fg[i][j].fy0;
246 fg[i][j].fz0 = source.fg[i][j].fz0;
247 fg[i][j].frx = source.fg[i][j].frx;
248 fg[i][j].fry = source.fg[i][j].fry;
249 fg[i][j].frz = source.fg[i][j].frz;
250 for(k=0;k<9;k++) fg[i][j].fr[k] = source.fg[i][j].fr[k];
251 } // end for j
252 } // end for i
253 return;
254}
255
256
257//________________________________________________________________________
258void AliITSgeom::GtoL(Int_t lay,Int_t lad,Int_t det,
259 const Float_t *g,Float_t *l){
260////////////////////////////////////////////////////////////////////////
261// The function that does the global ALICE Cartesian coordinate
262// to local active volume detector Cartesian coordinate transformation.
263// The local detector coordinate system is determined by the layer,
264// ladder, and detector numbers. The global coordinates are entered by
265// the three element Float_t array g and the local coordinate values
266// are returned by the three element Float_t array l. The order of the
267// three elements are g[0]=x, g[1]=y, and g[2]=z, similarly for l.
268////////////////////////////////////////////////////////////////////////
269 Double_t x,y,z;
270 ITS_geom *gl;
271
272 lay--; lad--; det--;
273 gl = &(fg[lay][fNdet[lay]*lad+det]);
274
275 x = g[0] - gl->fx0;
276 y = g[1] - gl->fy0;
277 z = g[2] - gl->fz0;
278 l[0] = gl->fr[0]*x + gl->fr[1]*y + gl->fr[2]*z;
279 l[1] = gl->fr[3]*x + gl->fr[4]*y + gl->fr[5]*z;
280 l[2] = gl->fr[6]*x + gl->fr[7]*y + gl->fr[8]*z;
281 return;
282}
283
284//________________________________________________________________________
285void AliITSgeom::GtoL(const Int_t *id,const Float_t *g,Float_t *l){
286////////////////////////////////////////////////////////////////////////
287// The function that does the local active volume detector Cartesian
288// coordinate to global ALICE Cartesian coordinate transformation.
289// The local detector coordinate system is determined by the layer,
290// ladder, and detector numbers. The local coordinates are entered by
291// the three element Float_t array l and the global coordinate values
292// are returned by the three element Float_t array g. The order of the
293// three elements are l[0]=x, l[1]=y, and l[2]=z, similarly for g.
294////////////////////////////////////////////////////////////////////////
295 Int_t lay,lad,det;
296 Double_t x,y,z;
297 ITS_geom *gl;
298
299 lay = id[0]; lad = id[1]; det = id[2];
300 lay--; lad--; det--;
301 gl = &(fg[lay][fNdet[lay]*lad+det]);
302
303 x = g[0] - gl->fx0;
304 y = g[1] - gl->fy0;
305 z = g[2] - gl->fz0;
306 l[0] = gl->fr[0]*x + gl->fr[1]*y + gl->fr[2]*z;
307 l[1] = gl->fr[3]*x + gl->fr[4]*y + gl->fr[5]*z;
308 l[2] = gl->fr[6]*x + gl->fr[7]*y + gl->fr[8]*z;
309 return;
310}
311//________________________________________________________________________
ad0e60d9 312void AliITSgeom::GtoL(const Int_t index,const Float_t *g,Float_t *l){
58005f18 313////////////////////////////////////////////////////////////////////////
314// The function that does the local active volume detector Cartesian
315// coordinate to global ALICE Cartesian coordinate transformation.
316// The local detector coordinate system is determined by the detector
317// index numbers (see GetModuleIndex and GetModuleID). The local
318// coordinates are entered by the three element Float_t array l and the
319// global coordinate values are returned by the three element Float_t array g.
320// The order of the three elements are l[0]=x, l[1]=y, and l[2]=z, similarly
321// for g.
322////////////////////////////////////////////////////////////////////////
323 Int_t lay,lad,det;
324 Double_t x,y,z;
325 ITS_geom *gl;
326
327 this->GetModuleId(index,lay,lad,det);
328 lay--; lad--; det--;
329 gl = &(fg[lay][fNdet[lay]*lad+det]);
330
331 x = g[0] - gl->fx0;
332 y = g[1] - gl->fy0;
333 z = g[2] - gl->fz0;
334 l[0] = gl->fr[0]*x + gl->fr[1]*y + gl->fr[2]*z;
335 l[1] = gl->fr[3]*x + gl->fr[4]*y + gl->fr[5]*z;
336 l[2] = gl->fr[6]*x + gl->fr[7]*y + gl->fr[8]*z;
337 return;
338}
339
340//________________________________________________________________________
341void AliITSgeom::LtoG(Int_t lay,Int_t lad,Int_t det,
342 const Float_t *l,Float_t *g){
343////////////////////////////////////////////////////////////////////////
344// The function that does the local active volume detector Cartesian
345// coordinate to global ALICE Cartesian coordinate transformation.
346// The local detector coordinate system is determined by the layer,
347// ladder, and detector numbers. The local coordinates are entered by
348// the three element Float_t array l and the global coordinate values
349// are returned by the three element Float_t array g. The order of the
350// three elements are l[0]=x, l[1]=y, and l[2]=z, similarly for g.
351////////////////////////////////////////////////////////////////////////
352 Double_t x,y,z;
353 ITS_geom *gl;
354
355 lay--; lad--; det--;
356 gl = &(fg[lay][fNdet[lay]*lad+det]);
357
358 x = gl->fr[0]*l[0] + gl->fr[3]*l[1] + gl->fr[6]*l[2];
359 y = gl->fr[1]*l[0] + gl->fr[4]*l[1] + gl->fr[7]*l[2];
360 z = gl->fr[2]*l[0] + gl->fr[5]*l[1] + gl->fr[8]*l[2];
361 g[0] = x + gl->fx0;
362 g[1] = y + gl->fy0;
363 g[2] = z + gl->fz0;
364 return;
365}
366
367//________________________________________________________________________
368void AliITSgeom::LtoG(const Int_t *id,const Float_t *l,Float_t *g){
369////////////////////////////////////////////////////////////////////////
370// The function that does the local active volume detector Cartesian
371// coordinate to global ALICE Cartesian coordinate transformation.
372// The local detector coordinate system is determined by the three
373// element array Id containing as it's three elements Id[0]=layer,
374// Id[1]=ladder, and Id[2]=detector numbers. The local coordinates
375// are entered by the three element Float_t array l and the global
376// coordinate values are returned by the three element Float_t array g.
377// The order of the three elements are l[0]=x, l[1]=y, and l[2]=z,
378// similarly for g.
379////////////////////////////////////////////////////////////////////////
380 Int_t lay,lad,det;
381 Double_t x,y,z;
382 ITS_geom *gl;
383
384 lay = id[0]; lad = id[1]; det = id[2];
385 lay--; lad--; det--;
386 gl = &(fg[lay][fNdet[lay]*lad+det]);
387
388 x = gl->fr[0]*l[0] + gl->fr[3]*l[1] + gl->fr[6]*l[2];
389 y = gl->fr[1]*l[0] + gl->fr[4]*l[1] + gl->fr[7]*l[2];
390 z = gl->fr[2]*l[0] + gl->fr[5]*l[1] + gl->fr[8]*l[2];
391 g[0] = x + gl->fx0;
392 g[1] = y + gl->fy0;
393 g[2] = z + gl->fz0;
394 return;
395}
396//________________________________________________________________________
ad0e60d9 397void AliITSgeom::LtoG(const Int_t index,const Float_t *l,Float_t *g){
58005f18 398////////////////////////////////////////////////////////////////////////
399// The function that does the local active volume detector Cartesian
400// coordinate to global ALICE Cartesian coordinate transformation.
401// The local detector coordinate system is determined by the detector
402// index number (see GetModuleIndex and GetModuleId). The local coordinates
403// are entered by the three element Float_t array l and the global
404// coordinate values are returned by the three element Float_t array g.
405// The order of the three elements are l[0]=x, l[1]=y, and l[2]=z,
406// similarly for g.
407////////////////////////////////////////////////////////////////////////
408 Int_t lay,lad,det;
409 Double_t x,y,z;
410 ITS_geom *gl;
411
412 this->GetModuleId(index,lay,lad,det);
413 lay--; lad--; det--;
414 gl = &(fg[lay][fNdet[lay]*lad+det]);
415
416 x = gl->fr[0]*l[0] + gl->fr[3]*l[1] + gl->fr[6]*l[2];
417 y = gl->fr[1]*l[0] + gl->fr[4]*l[1] + gl->fr[7]*l[2];
418 z = gl->fr[2]*l[0] + gl->fr[5]*l[1] + gl->fr[8]*l[2];
419 g[0] = x + gl->fx0;
420 g[1] = y + gl->fy0;
421 g[2] = z + gl->fz0;
422 return;
423}
424//________________________________________________________________________
425void AliITSgeom::GtoLMomentum(Int_t lay,Int_t lad,Int_t det,
426 const Float_t *g,Float_t *l){
427////////////////////////////////////////////////////////////////////////
428// The function that does the global ALICE Cartesian momentum
429// to local active volume detector Cartesian momentum transformation.
430// The local detector coordinate system is determined by the layer,
431// ladder, and detector numbers. The global momentums are entered by
432// the three element Float_t array g and the local momentums values
433// are returned by the three element Float_t array l. The order of the
434// three elements are g[0]=x, g[1]=y, and g[2]=z, similarly for l.
435////////////////////////////////////////////////////////////////////////
436 Double_t px,py,pz;
437 ITS_geom *gl;
438
439 lay--; lad--; det--;
440 gl = &(fg[lay][fNdet[lay]*lad+det]);
441
442 px = g[0];
443 py = g[1];
444 pz = g[2];
445 l[0] = gl->fr[0]*px + gl->fr[1]*py + gl->fr[2]*pz;
446 l[1] = gl->fr[3]*px + gl->fr[4]*py + gl->fr[5]*pz;
447 l[2] = gl->fr[6]*px + gl->fr[7]*py + gl->fr[8]*pz;
448 return;
449}
450//________________________________________________________________________
451void AliITSgeom::LtoGMomentum(Int_t lay,Int_t lad,Int_t det,
452 const Float_t *l,Float_t *g){
453////////////////////////////////////////////////////////////////////////
454// The function that does the local active volume detector Cartesian
455// momentum to global ALICE Cartesian momentum transformation.
456// The local detector momentum system is determined by the layer,
457// ladder, and detector numbers. The locall momentums are entered by
458// the three element Float_t array l and the global momentum values
459// are returned by the three element Float_t array g. The order of the
460// three elements are l[0]=x, l[1]=y, and l[2]=z, similarly for g.
461////////////////////////////////////////////////////////////////////////
462 Double_t px,py,pz;
463 ITS_geom *gl;
464
465 lay--; lad--; det--;
466 gl = &(fg[lay][fNdet[lay]*lad+det]);
467
468 px = gl->fr[0]*l[0] + gl->fr[3]*l[1] + gl->fr[6]*l[2];
469 py = gl->fr[1]*l[0] + gl->fr[4]*l[1] + gl->fr[7]*l[2];
470 pz = gl->fr[2]*l[0] + gl->fr[5]*l[1] + gl->fr[8]*l[2];
471 g[0] = px;
472 g[1] = py;
473 g[2] = pz;
474 return;
475}
476//___________________________________________________________________________
477Int_t AliITSgeom::GetModuleIndex(Int_t lay,Int_t lad,Int_t det){
478 Int_t i,j,k;
479
480 i = fNdet[lay-1] * (lad-1) + det - 1;
481 j = 0;
482 for(k=0;k<lay-1;k++) j += fNdet[k]*fNlad[k];
483 return (i+j);
484}
485//___________________________________________________________________________
486void AliITSgeom::GetModuleId(Int_t index,Int_t &lay,Int_t &lad,Int_t &det){
487 Int_t i,j,k;
488
489 j = 0;
490 for(k=0;k<fNlayers;k++){
491 j += fNdet[k]*fNlad[k];
492 if(index>j)break;
493 } // end for k
494 lay = k+1;
495 i = index -j + fNdet[k]*fNlad[k];
496 j = 0;
497 for(k=0;k<fNlad[lay-1];k++){
498 j += fNdet[k];
499 if(i>fNdet[k])break;
500 } // end for k
501 lad = k+1;
502 det = 1+i-fNdet[lay-1]*k;
503 return;
504}
505//___________________________________________________________________________
506void AliITSgeom::GlobalChange(Float_t *tran,Float_t *rot){
507////////////////////////////////////////////////////////////////////////
508// This function performs a Cartesian translation and rotation of
509// the full ITS from its default position by an amount determined by
510// the three element arrays dtranslation and drotation. If every element
511// of dtranslation and drotation are zero then there is no change made
512// the geometry. The change is global in that the exact same translation
513// and rotation is done to every detector element in the exact same way.
514// The units of the translation are those of the Monte Carlo, usually cm,
515// and those of the rotation are in radians. The elements of dtranslation
516// are dtranslation[0] = x, dtranslation[1] = y, and dtranslation[2] = z.
517// The elements of drotation are drotation[0] = rx, drotation[1] = ry, and
518// drotation[2] = rz. A change in x will move the hole ITS in the ALICE
519// global x direction, the same for a change in y. A change in z will
520// result in a translation of the ITS as a hole up or down the beam line.
521// A change in the angles will result in the inclination of the ITS with
522// respect to the beam line, except for an effective rotation about the
523// beam axis which will just rotate the ITS as a hole about the beam axis.
524////////////////////////////////////////////////////////////////////////
525 Int_t i,j,k,l;
526 Double_t rx,ry,rz;
527 Double_t sx,cx,sy,cy,sz,cz;
528 ITS_geom *gl;
529
530 for(i=0;i<fNlayers;i++){
531 for(j=0;j<fNlad[i];j++) for(k=0;k<fNdet[i];k++){
532 l = fNdet[i]*j+k; // resolved index
533 gl = &(fg[i][l]);
534 gl->fx0 += tran[0];
535 gl->fy0 += tran[1];
536 gl->fz0 += tran[2];
537 gl->frx += rot[0];
538 gl->fry += rot[1];
539 gl->frz += rot[2];
540 rx = gl->frx; ry = gl->fry; rz = gl->frz;
541 sx = sin(rx); cx = cos(rx);
542 sy = sin(ry); cy = cos(ry);
543 sz = sin(rz); cz = cos(rz);
544 gl->fr[0] = cz*cy;
545 gl->fr[1] = -cz*sy*sx - sz*cx;
546 gl->fr[2] = -cz*sy*cx + sz*sx;
547 gl->fr[3] = sz*cy;
548 gl->fr[4] = -sz*sy*sx + cz*cx;
549 gl->fr[5] = -sz*sy*cx - cz*sx;
550 gl->fr[6] = sy;
551 gl->fr[7] = cy*sx;
552 gl->fr[8] = cy*cx;
553 } // end for j,k
554 } // end for i
555 return;
556}
557
558//___________________________________________________________________________
559void AliITSgeom::GlobalCylindericalChange(Float_t *tran,Float_t *rot){
560////////////////////////////////////////////////////////////////////////
561// This function performs a cylindrical translation and rotation of
562// each ITS element by a fixed about in radius, rphi, and z from its
563// default position by an amount determined by the three element arrays
564// dtranslation and drotation. If every element of dtranslation and
565// drotation are zero then there is no change made the geometry. The
566// change is global in that the exact same distance change in translation
567// and rotation is done to every detector element in the exact same way.
568// The units of the translation are those of the Monte Carlo, usually cm,
569// and those of the rotation are in radians. The elements of dtranslation
570// are dtranslation[0] = r, dtranslation[1] = rphi, and dtranslation[2] = z.
571// The elements of drotation are drotation[0] = rx, drotation[1] = ry, and
572// drotation[2] = rz. A change in r will results in the increase of the
573// radius of each layer by the same about. A change in rphi will results in
574// the rotation of each layer by a different angle but by the same
575// circumferential distance. A change in z will result in a translation
576// of the ITS as a hole up or down the beam line. A change in the angles
577// will result in the inclination of the ITS with respect to the beam
578// line, except for an effective rotation about the beam axis which will
579// just rotate the ITS as a hole about the beam axis.
580////////////////////////////////////////////////////////////////////////
581 Int_t i,j,k,l;
582 Double_t rx,ry,rz,r,phi,rphi; // phi in radians
583 Double_t sx,cx,sy,cy,sz,cz,r0;
584 ITS_geom *gl;
585
586// printf("trans=%f %f %f rot=%f %f %f\n",tran[0],tran[1],tran[2],
587// rot[0],rot[1],rot[2]);
588 for(i=0;i<fNlayers;i++){
589 for(j=0;j<fNlad[i];j++) for(k=0;k<fNdet[i];k++){
590 l = fNdet[i]*j+k; // resolved index
591 gl = &(fg[i][l]);
592 r = r0= TMath::Hypot(gl->fy0,gl->fx0);
593 phi = atan2(gl->fy0,gl->fx0);
594 rphi = r0*phi;
595 r += tran[0];
596 rphi += tran[1];
597 phi = rphi/r0;
598 gl->fx0 = r*TMath::Cos(phi);
599 gl->fy0 = r*TMath::Sin(phi);
600 gl->fz0 += tran[2];
601 gl->frx += rot[0];
602 gl->fry += rot[1];
603 gl->frz += rot[2];
604 rx = gl->frx; ry = gl->fry; rz = gl->frz;
605 sx = sin(rx); cx = cos(rx);
606 sy = sin(ry); cy = cos(ry);
607 sz = sin(rz); cz = cos(rz);
608 gl->fr[0] = cz*cy;
609 gl->fr[1] = -cz*sy*sx - sz*cx;
610 gl->fr[2] = -cz*sy*cx + sz*sx;
611 gl->fr[3] = sz*cy;
612 gl->fr[4] = -sz*sy*sx + cz*cx;
613 gl->fr[5] = -sz*sy*cx - cz*sx;
614 gl->fr[6] = sy;
615 gl->fr[7] = cy*sx;
616 gl->fr[8] = cy*cx;
617 } // end for j,k
618 } // end for i
619 return;
620}
621
622//___________________________________________________________________________
623void AliITSgeom::RandomChange(Float_t *stran,Float_t *srot){
624////////////////////////////////////////////////////////////////////////
625// This function performs a Gaussian random displacement and/or
626// rotation about the present global position of each active
627// volume/detector of the ITS. The sigma of the random displacement
628// is determined by the three element array stranslation, for the
629// x y and z translations, and the three element array srotation,
630// for the three rotation about the axis x y and z.
631////////////////////////////////////////////////////////////////////////
632 Int_t i,j,k,l;
633 Double_t rx,ry,rz;
634 Double_t sx,cx,sy,cy,sz,cz;
635 TRandom ran;
636 ITS_geom *gl;
637
638 for(i=0;i<fNlayers;i++){
639 for(j=0;j<fNlad[i];j++) for(k=0;k<fNdet[i];k++){
640 l = fNdet[i]*j+k; // resolved index
641 gl = &(fg[i][l]);
642 gl->fx0 += ran.Gaus(0.0,stran[0]);
643 gl->fy0 += ran.Gaus(0.0,stran[1]);
644 gl->fz0 += ran.Gaus(0.0,stran[2]);
645 gl->frx += ran.Gaus(0.0, srot[0]);
646 gl->fry += ran.Gaus(0.0, srot[1]);
647 gl->frz += ran.Gaus(0.0, srot[2]);
648 rx = gl->frx; ry = gl->fry; rz = gl->frz;
649 sx = sin(rx); cx = cos(rx);
650 sy = sin(ry); cy = cos(ry);
651 sz = sin(rz); cz = cos(rz);
652 gl->fr[0] = cz*cy;
653 gl->fr[1] = -cz*sy*sx - sz*cx;
654 gl->fr[2] = -cz*sy*cx + sz*sx;
655 gl->fr[3] = sz*cy;
656 gl->fr[4] = -sz*sy*sx + cz*cx;
657 gl->fr[5] = -sz*sy*cx - cz*sx;
658 gl->fr[6] = sy;
659 gl->fr[7] = cy*sx;
660 gl->fr[8] = cy*cx;
661 } // end for j,k
662 } // end for i
663 return;
664}
665
666//___________________________________________________________________________
667void AliITSgeom::RandomCylindericalChange(Float_t *stran,Float_t *srot){
668////////////////////////////////////////////////////////////////////////
669// This function performs a Gaussian random displacement and/or
670// rotation about the present global position of each active
671// volume/detector of the ITS. The sigma of the random displacement
672// is determined by the three element array stranslation, for the
673// r rphi and z translations, and the three element array srotation,
674// for the three rotation about the axis x y and z. This random change
675// in detector position allow for the simulation of a random uncertainty
676// in the detector positions of the ITS.
677////////////////////////////////////////////////////////////////////////
678 Int_t i,j,k,l;
679 Double_t rx,ry,rz,r,phi,x,y; // phi in radians
680 Double_t sx,cx,sy,cy,sz,cz,r0;
681 TRandom ran;
682 ITS_geom *gl;
683
684// printf("trans=%f %f %f rot=%f %f %f\n",stran[0],stran[1],stran[2],
685// srot[0],srot[1],srot[2]);
686 for(i=0;i<fNlayers;i++){
687 for(j=0;j<fNlad[i];j++) for(k=0;k<fNdet[i];k++){
688 l = fNdet[i]*j+k; // resolved index
689 gl = &(fg[i][l]);
690 x = gl->fx0;
691 y = gl->fy0;
692 r = r0= TMath::Hypot(y,x);
693 phi = TMath::ATan2(y,x);
694// if(phi<0.0) phi += 2.0*TMath::Pi();
695 r += ran.Gaus(0.0,stran[0]);
696 phi += ran.Gaus(0.0,stran[1])/r0;
697// printf("fx0=%f fy0=%f rcos(phi)=%f rsin(phi)=%f\n",gl->fx0,gl->fy0,
698// r*TMath::Cos(phi),r*TMath::Sin(phi));
699 gl->fx0 = r*TMath::Cos(phi);
700 gl->fy0 = r*TMath::Sin(phi);
701// printf("r0=%f r=%f hypot=%f phi0=%f phi=%f ATan2=%f\n",
702// r0,r,TMath::Hypot(gl->fy0,gl->fx0),
703// phi0,phi,TMath::ATan2(gl->fy0,gl->fx0));
704 gl->fz0 += ran.Gaus(0.0,stran[2]);
705 gl->frx += ran.Gaus(0.0, srot[0]);
706 gl->fry += ran.Gaus(0.0, srot[1]);
707 gl->frz += ran.Gaus(0.0, srot[2]);
708 rx = gl->frx; ry = gl->fry; rz = gl->frz;
709 sx = sin(rx); cx = cos(rx);
710 sy = sin(ry); cy = cos(ry);
711 sz = sin(rz); cz = cos(rz);
712 gl->fr[0] = cz*cy;
713 gl->fr[1] = -cz*sy*sx - sz*cx;
714 gl->fr[2] = -cz*sy*cx + sz*sx;
715 gl->fr[3] = sz*cy;
716 gl->fr[4] = -sz*sy*sx + cz*cx;
717 gl->fr[5] = -sz*sy*cx - cz*sx;
718 gl->fr[6] = sy;
719 gl->fr[7] = cy*sx;
720 gl->fr[8] = cy*cx;
721 } // end for j,k
722 } // end for i
723 return;
724}
725
726//___________________________________________________________________________
727void AliITSgeom::SetByAngles(Int_t lay,Int_t lad,Int_t det,
728 Float_t rx,Float_t ry,Float_t rz){
729////////////////////////////////////////////////////////////////////////
730// This function computes a new rotation matrix based on the angles
731// rx, ry, and rz (in radians) for a give detector on the give ladder
732// in the give layer. A new
733// fg[layer-1][(fNlad[layer-1]*(ladder-1)+detector-1)].fr[] array is
734// computed.
735////////////////////////////////////////////////////////////////////////
736 ITS_geom *g;
737 Double_t sx,cx,sy,cy,sz,cz;
738
739 lay--; lad--; det--; // set to zero base now.
740 g = &(fg[lay][fNdet[lay]*lad+det]);
741
742 sx = sin(rx); cx = cos(rx);
743 sy = sin(ry); cy = cos(ry);
744 sz = sin(rz); cz = cos(rz);
745 g->frx = rx;
746 g->fry = ry;
747 g->frz = rz;
748 g->fr[0] = cz*cy;
749 g->fr[1] = -cz*sy*sx - sz*cx;
750 g->fr[2] = -cz*sy*cx + sz*sx;
751 g->fr[3] = sz*cy;
752 g->fr[4] = -sz*sy*sx + cz*cx;
753 g->fr[5] = -sz*sy*cx - cz*sx;
754 g->fr[6] = sy;
755 g->fr[7] = cy*sx;
756 g->fr[8] = cy*cx;
757 return;
758}
759
760//___________________________________________________________________________
761void AliITSgeom::GetRotMatrix(Int_t lay,Int_t lad,Int_t det,Float_t *mat){
762////////////////////////////////////////////////////////////////////////
763// Returns, in the Float_t array pointed to by mat, the full rotation
764// matrix for the give detector defined by layer, ladder, and detector.
765// It returns all nine elements of fr in the ITS_geom structure. See the
766// description of the ITS_geom structure for further details of this
767// rotation matrix.
768////////////////////////////////////////////////////////////////////////
769 Int_t i;
770 ITS_geom *g;
771
772 lay--; lad--; det--; // shift to base 0
773 g = &(fg[lay][fNdet[lay]*lad+det]);
774 for(i=0;i<9;i++) mat[i] = g->fr[i];
775 return;
776}
777
778//___________________________________________________________________________
779void AliITSgeom::PrintComparison(FILE *fp,AliITSgeom *other){
780////////////////////////////////////////////////////////////////////////
781// This function was primarily created for diagnostic reasons. It
782// print to a file pointed to by the file pointer fp the difference
783// between two AliITSgeom classes. The format of the file is basicly,
784// define d? to be the difference between the same element of the two
785// classes. For example dfrx = this->fg[i][j].frx - other->fg[i][j].frx.
786// if(at least one of dfx0, dfy0, dfz0,dfrx,dfry,dfrz are non zero) then print
787// layer ladder detector dfx0 dfy0 dfz0 dfrx dfry dfrz
788// if(at least one of the 9 elements of dfr[] are non zero) then print
789// layer ladder detector dfr[0] dfr[1] dfr[2]
790// dfr[3] dfr[4] dfr[5]
791// dfr[6] dfr[7] dfr[8]
792// Only non zero values are printed to save space. The differences are
793// typical written to a file because there are usually a lot of numbers
794// printed out and it is usually easier to read them in some nice editor
795// rather than zooming quickly past you on a screen. fprintf is used to
796// do the printing. The fShapeIndex difference is not printed at this time.
797////////////////////////////////////////////////////////////////////////
798 Int_t i,j,k,l;
799 Double_t xt,yt,zt,xo,yo,zo;
800 Double_t rxt,ryt,rzt,rxo,ryo,rzo; // phi in radians
801 ITS_geom *gt,*go;
802 Bool_t t;
803
804 for(i=0;i<this->fNlayers;i++){
805 for(j=0;j<this->fNlad[i];j++) for(k=0;k<this->fNdet[i];k++){
806 l = this->fNdet[i]*j+k; // resolved index
807 gt = &(this->fg[i][l]);
808 go = &(other->fg[i][l]);
809 xt = gt->fx0; yt = gt->fy0; zt = gt->fz0;
810 xo = go->fx0; yo = go->fy0; zo = go->fz0;
811 rxt = gt->frx; ryt = gt->fry; rzt = gt->frz;
812 rxo = go->frx; ryo = go->fry; rzo = go->frz;
813 if(!(xt==xo&&yt==yo&&zt==zo&&rxt==rxo&&ryt==ryo&&rzt==rzo))
814 fprintf(fp,"%1.1d %2.2d %2.2d dTrans=%f %f %f drot=%f %f %f\n",
815 i+1,j+1,k+1,xt-xo,yt-yo,zt-zo,rxt-rxo,ryt-ryo,rzt-rzo);
816 t = kFALSE;
817 for(i=0;i<9;i++) t = gt->fr[i] != go->fr[i];
818 if(t){
819 fprintf(fp,"%1.1d %2.2d %2.2d dfr= %e %e %e\n",i+1,j+1,k+1,
820 gt->fr[0]-go->fr[0],gt->fr[1]-go->fr[1],gt->fr[2]-go->fr[2]);
821 fprintf(fp," dfr= %e %e %e\n",
822 gt->fr[3]-go->fr[3],gt->fr[4]-go->fr[4],gt->fr[5]-go->fr[5]);
823 fprintf(fp," dfr= %e %e %e\n",
824 gt->fr[6]-go->fr[6],gt->fr[7]-go->fr[7],gt->fr[8]-go->fr[8]);
825 }
826 } // end for j,k
827 } // end for i
828 return;
829}
830
831//___________________________________________________________________________
832void AliITSgeom::PrintData(FILE *fp,Int_t lay,Int_t lad,Int_t det){
833////////////////////////////////////////////////////////////////////////
834// This function prints out the coordinate transformations for
835// the particular detector defined by layer, ladder, and detector
836// to the file pointed to by the File pointer fp. fprinf statements
837// are used to print out the numbers. The format is
838// layer ladder detector Trans= fx0 fy0 fz0 rot= frx fry frz Shape=fShapeIndex
839// dfr= fr[0] fr[1] fr[2]
840// dfr= fr[3] fr[4] fr[5]
841// dfr= fr[6] fr[7] fr[8]
842// By indicating which detector, some control over the information
843// is given to the user. The output it written to the file pointed
844// to by the file pointer fp. This can be set to stdout if you want.
845////////////////////////////////////////////////////////////////////////
846 Int_t i,j,k,l;
847 ITS_geom *gt;
848
849 i = lay-1;
850 j = lad-1;
851 k = det-1;
852 l = this->fNdet[i]*j+k; // resolved index
853 gt = &(this->fg[i][l]);
854 fprintf(fp,"%1.1d %2.2d %2.2d Trans=%f %f %f rot=%f %f %f Shape=%d\n",
855 i+1,j+1,k+1,gt->fx0,gt->fy0,gt->fz0,gt->frx,gt->fry,gt->frz,
856 gt->fShapeIndex);
857 fprintf(fp," dfr= %e %e %e\n",gt->fr[0],gt->fr[1],gt->fr[2]);
858 fprintf(fp," dfr= %e %e %e\n",gt->fr[3],gt->fr[4],gt->fr[5]);
859 fprintf(fp," dfr= %e %e %e\n",gt->fr[6],gt->fr[7],gt->fr[8]);
860 return;
861}
862//___________________________________________________________________________
863void AliITSgeom::Streamer(TBuffer &R__b){
864////////////////////////////////////////////////////////////////////////
865// The default Streamer function "written by ROOT" doesn't write out
866// the arrays referenced by pointers. Therefore, a specific Streamer function
867// has to be written. This function should not be modified but instead added
868// on to so that older versions can still be read. The proper handling of
869// the version dependent streamer function hasn't been written do to the lack
870// of finding an example at the time of writting.
871////////////////////////////////////////////////////////////////////////
872 // Stream an object of class AliITSgeom.
873 Int_t i,j,k;
874
875 if (R__b.IsReading()) {
876 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
877 TObject::Streamer(R__b);
878 R__b >> fNlayers;
879 if(fNlad!=0) delete[] fNlad;
880 if(fNdet!=0) delete[] fNdet;
881 fNlad = new Int_t[fNlayers];
882 fNdet = new Int_t[fNlayers];
883 for(i=0;i<fNlayers;i++) R__b >> fNlad[i];
884 for(i=0;i<fNlayers;i++) R__b >> fNdet[i];
885 if(fg!=0){
886 for(i=0;i<fNlayers;i++) delete[] fg[i];
887 delete[] fg;
888 } // end if fg!=0
889 fg = new ITS_geom*[fNlayers];
890 for(i=0;i<fNlayers;i++){
891 fg[i] = new ITS_geom[fNlad[i]*fNdet[i]];
892 for(j=0;j<fNlad[i]*fNdet[i];j++){
893 R__b >> fg[i][j].fShapeIndex;
894 R__b >> fg[i][j].fx0;
895 R__b >> fg[i][j].fy0;
896 R__b >> fg[i][j].fz0;
897 R__b >> fg[i][j].frx;
898 R__b >> fg[i][j].fry;
899 R__b >> fg[i][j].frz;
900 for(k=0;k<9;k++) R__b >> fg[i][j].fr[k];
901 } // end for j
902 } // end for i
903 R__b >> fShape;
904 } else {
905 R__b.WriteVersion(AliITSgeom::IsA());
906 TObject::Streamer(R__b);
907 R__b << fNlayers;
908 for(i=0;i<fNlayers;i++) R__b << fNlad[i];
909 for(i=0;i<fNlayers;i++) R__b << fNdet[i];
910 for(i=0;i<fNlayers;i++) for(j=0;j<fNlad[i]*fNdet[i];j++){
911 R__b << fg[i][j].fShapeIndex;
912 R__b << fg[i][j].fx0;
913 R__b << fg[i][j].fy0;
914 R__b << fg[i][j].fz0;
915 R__b << fg[i][j].frx;
916 R__b << fg[i][j].fry;
917 R__b << fg[i][j].frz;
918 for(k=0;k<9;k++) R__b << fg[i][j].fr[k];
919 } // end for i,j
920 R__b << fShape;
921 }
922}
923
924//___________________________________________________________________________
925ofstream & AliITSgeom::PrintGeom(ofstream &R__b){
926////////////////////////////////////////////////////////////////////////
927// The default Streamer function "written by ROOT" doesn't write out
928// the arrays referenced by pointers. Therefore, a specific Streamer function
929// has to be written. This function should not be modified but instead added
930// on to so that older versions can still be read. The proper handling of
931// the version dependent streamer function hasn't been written do to the lack
932// of finding an example at the time of writting.
933////////////////////////////////////////////////////////////////////////
934 // Stream an object of class AliITSgeom.
935 Int_t i,j,k;
936
937 R__b.setf(ios::scientific);
938 R__b << fNlayers << " ";
939 for(i=0;i<fNlayers;i++) R__b << fNlad[i] << " ";
940 for(i=0;i<fNlayers;i++) R__b << fNdet[i] << "\n";
941 for(i=0;i<fNlayers;i++) for(j=0;j<fNlad[i]*fNdet[i];j++){
942 R__b <<setprecision(16) << fg[i][j].fShapeIndex << " ";
943 R__b <<setprecision(16) << fg[i][j].fx0 << " ";
944 R__b <<setprecision(16) << fg[i][j].fy0 << " ";
945 R__b <<setprecision(16) << fg[i][j].fz0 << " ";
946 R__b <<setprecision(16) << fg[i][j].frx << " ";
947 R__b <<setprecision(16) << fg[i][j].fry << " ";
948 R__b <<setprecision(16) << fg[i][j].frz << "\n";
949 for(k=0;k<9;k++) R__b <<setprecision(16) << fg[i][j].fr[k] << " ";
950 R__b << "\n";
951 } // end for i,j
952// R__b << fShape;
953 return R__b;
954}
955
956//___________________________________________________________________________
957ifstream & AliITSgeom::ReadGeom(ifstream &R__b){
958////////////////////////////////////////////////////////////////////////
959// The default Streamer function "written by ROOT" doesn't write out
960// the arrays referenced by pointers. Therefore, a specific Streamer function
961// has to be written. This function should not be modified but instead added
962// on to so that older versions can still be read. The proper handling of
963// the version dependent streamer function hasn't been written do to the lack
964// of finding an example at the time of writting.
965////////////////////////////////////////////////////////////////////////
966 // Stream an object of class AliITSgeom.
967 Int_t i,j,k;
968
969 R__b >> fNlayers;
970 if(fNlad!=0) delete[] fNlad;
971 if(fNdet!=0) delete[] fNdet;
972 fNlad = new Int_t[fNlayers];
973 fNdet = new Int_t[fNlayers];
974 for(i=0;i<fNlayers;i++) R__b >> fNlad[i];
975 for(i=0;i<fNlayers;i++) R__b >> fNdet[i];
976 if(fg!=0){
977 for(i=0;i<fNlayers;i++) delete[] fg[i];
978 delete[] fg;
979 } // end if fg!=0
980 fg = new ITS_geom*[fNlayers];
981 for(i=0;i<fNlayers;i++){
982 fg[i] = new ITS_geom[fNlad[i]*fNdet[i]];
983 for(j=0;j<fNlad[i]*fNdet[i];j++){
984 R__b >> fg[i][j].fShapeIndex;
985 R__b >> fg[i][j].fx0;
986 R__b >> fg[i][j].fy0;
987 R__b >> fg[i][j].fz0;
988 R__b >> fg[i][j].frx;
989 R__b >> fg[i][j].fry;
990 R__b >> fg[i][j].frz;
991 for(k=0;k<9;k++) R__b >> fg[i][j].fr[k];
992 } // end for j
993 } // end for i
994// R__b >> fShape;
995 return R__b;
996}