ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/

/* $Id$ */

///////////////////////////////////////////////////////////////////////
// ITS geometry manipulation routines.                               //
// Created April 15 1999.                                            //
// version: 0.0.0                                                    //
// By: Bjorn S. Nilsen                                               //
// version: 0.0.1                                                    //
// Updated May 27 1999.                                              //
// Added Cylindrical random and global based changes.                //
//                                                                   //
// Modified and added functions Feb. 7 2006                          //
///////////////////////////////////////////////////////////////////////


////////////////////////////////////////////////////////////////////////
//     The local coordinate system by, default, is show in the following
// figures. Also shown are the ladder numbering scheme.
//Begin_Html
/*
<img src="picts/ITS/AliITSgeomMatrix_L1.gif">
</pre>
<br clear=left>
<font size=+2 color=blue>
<p>This shows the relative geometry differences between the ALICE Global
coordinate system and the local detector coordinate system.
</font>
<pre>

<pre>
<img src="picts/ITS/its1+2_convention_front_5.gif">
</pre>
<br clear=left>
<font size=+2 color=blue>
<p>This shows the front view of the SPDs and the orientation of the local
pixel coordinate system. Note that the inner pixel layer has its y coordinate
in the opposite direction from all of the other layers.
</font>
<pre>

<pre>
<img src="picts/ITS/its3+4_convention_front_5.gif">
</pre>
<br clear=left>
<font size=+2 color=blue>
<p>This shows the front view of the SDDs and the orientation of the local
pixel coordinate system.
</font>
<pre>

<pre>
<img src="picts/ITS/its5+6_convention_front_5.gif">
</pre>
<br clear=left>
<font size=+2 color=blue>
<p>This shows the front view of the SSDs and the orientation of the local
pixel coordinate system.
</font>
<pre>
*/
//End_Html
//
////////////////////////////////////////////////////////////////////////

////////////////////////////////////////////////////////////////////////
//
// version: 0
// Written by Bjorn S. Nilsen
//
// Data Members:
//
// TString    fVersion 
//     Transformation version.
// Int_t      fTrans
//     Flag to keep track of which transformation 
// Int_t      fNmodules
//      The total number of modules
// Int_t fNlayers
//     The number of ITS layers for this geometry. By default this
//  is 6, but can be modified by the creator function if there are
// more layers defined.
//
// TArrayI fNlad
//     A pointer to an array fNlayers long containing the number of 
// ladders for each layer. This array is typically created and filled 
// by the AliITSgeom creator function.
//
// TArrayI fNdet
//     A pointer to an array fNlayers long containing the number of
// active detector volumes for each ladder. This array is typically
// created and filled by the AliITSgeom creator function.
//
// TObjArray fGm containing objects of type AliITSgeomMatrix
//     A pointer to an array of AliITSgeomMatrix classes. One element 
// per module (detector) in the ITS. AliITSgeomMatrix basicly contains
// all of the necessary information about the detector and it's coordinate
// transformations.
//
////////////////////////////////////////////////////////////////////////
#include <Riostream.h>
#include <ctype.h>

#include <TRandom.h>
#include <TSystem.h>
#include <TArrayI.h>

#include "AliITSgeom.h"
#include "AliLog.h"

ClassImp(AliITSgeom)

//______________________________________________________________________
AliITSgeom::AliITSgeom():
TObject(),
fVersion("GEANT"),// Transformation version.
fTrans(0),       // Flag to keep track of which transformation 
fNmodules(0),    // The total number of modules
fNlayers(0),     // The number of layers.
fNlad(),         //[] Array of the number of ladders/layer(layer)
fNdet(),         //[] Array of the number of detector/ladder(layer)
fGm(0,0)        // Structure of translation. and rotation.
{
    //     The default constructor for the AliITSgeom class. It, by default,
    // sets fNlayers to zero and zeros all pointers.
    // Do not allocate anything zero everything.
    // Inputs:
    //    none.
    // Outputs:
    //    none.
    // Return:
    //    a zeroed AliITSgeom object.

    fGm.SetOwner(kTRUE);
    return;
}

//______________________________________________________________________
AliITSgeom::AliITSgeom(Int_t itype,Int_t nlayers,const Int_t *nlads,
                       const Int_t *ndets,Int_t mods):
TObject(),
fVersion("GEANT"),    // Transformation version.
fTrans(itype),       // Flag to keep track of which transformation 
fNmodules(mods),     // The total number of modules
fNlayers(nlayers),   // The number of layers.
fNlad(nlayers,nlads),//[] Array of the number of ladders/layer(layer)
fNdet(nlayers,ndets),//[] Array of the number of detector/ladder(layer)
fGm(mods,0)         // Structure of translation. and rotation.
{
    //     A simple constructor to set basic geometry class variables
    // Inputs:
    //      Int_t itype   the type of transformation kept.
    //                    bit 0 => Standard GEANT
    //                    bit 1 => ITS tracking
    //                    bit 2 => A change in the coordinate system 
    //                    has been made. others are still to be defined 
    //                    as needed.
    //      Int_t nlayers The number of ITS layers also set the size of 
    //                    the arrays
    //      Int_t *nlads  an array of the number of ladders for each 
    //                    layer. This array must be nlayers long.
    //      Int_t *ndets  an array of the number of detectors per ladder
    //                    for each layer. This array must be nlayers long.
    //      Int_t mods    The number of modules. Typically the sum of all the 
    //                    detectors on every layer and ladder.
    // Outputs:
    //     none
    // Return:
    //     A properly inilized AliITSgeom object.

    fGm.SetOwner(kTRUE);
    return;
}
//______________________________________________________________________
void AliITSgeom::Init(Int_t itype,Int_t nlayers,const Int_t *nlads,
                      const Int_t *ndets,Int_t mods){
    //     A simple Inilizer to set basic geometry class variables
    // Inputs:
    //      Int_t itype   the type of transformation kept.
    //                    bit 0 => Standard GEANT
    //                    bit 1 => ITS tracking
    //                    bit 2 => A change in the coordinate system 
    //                    has been made. others are still to be defined 
    //                    as needed.
    //      Int_t nlayers The number of ITS layers also set the size of 
    //                    the arrays
    //      Int_t *nlads  an array of the number of ladders for each 
    //                    layer. This array must be nlayers long.
    //      Int_t *ndets  an array of the number of detectors per ladder 
    //                    for each layer. This array must be nlayers long.
    //      Int_t mods    The number of modules. Typically the sum of all the 
    //                    detectors on every layer and ladder.
    // Outputs:
    //     none
    // Return:
    //     A properly inilized AliITSgeom object.

    fVersion  = "GEANT";     // Transformation version.
    fTrans    = itype;       // Flag to keep track of which transformation 
    fNmodules = mods;        // The total number of modules
    fNlayers  = nlayers;     // The number of layers.
    fNlad.Set(nlayers,nlads);//[] Array of the number of ladders/layer(layer)
    fNdet.Set(nlayers,ndets);//[] Array of the number of detector/ladder(layer)
    fGm.Clear();
    fGm.Expand(mods);        // Structure of translation. and rotation.
    fGm.SetOwner(kTRUE);
    return;
}
//______________________________________________________________________
void AliITSgeom::CreateMatrix(Int_t mod,Int_t lay,Int_t lad,Int_t det,
                              AliITSDetector idet,const Double_t tran[3],
                              const Double_t rot[10]){
    // Given the translation vector tran[3] and the rotation matrix rot[1],
    // this function creates and adds to the TObject Array fGm the
    // AliITSgeomMatrix object.
    // The rot[10] matrix is set up like:
    /*   / rot[0]  rot[1]  rot[2] \
    //  |  rot[3]  rot[4]  rot[5]  |
    //   \ rot[6]  rot[7]  rot[8] /  if(rot[9]!=0) then the Identity matrix
    // is used regardless of the values in rot[0]-rot[8].
    */
    // Inputs:
    //    Int_t           mod     The module number. The location in TObjArray
    //    Int_t           lay     The layer where this module is
    //    Int_t           lad     On which ladder this module is
    //    Int_t           det     Which detector on this ladder this module is
    //    AliITSDetector idet     The type of detector see AliITSgeom.h
    //    Double_t       tran[3]  The translation vector
    //    Double_t       rot[10]  The rotation matrix.
    // Outputs:
    //    none
    // Return:
    //    none.
    Int_t id[3];
    Double_t r[3][3] = {{1.0,0.0,0.0},{0.0,1.0,0.0},{0.0,0.0,1.0}};

    if(mod<0||mod>=fGm.GetSize()){ 
	Error("CreateMatrix","mod=%d is out of bounds max value=%d",mod,
	      fGm.GetSize());
	return;
    } // end if
    delete fGm.At(mod);
    id[0] = lay; id[1] = lad; id[2] = det;
    if(rot[9]!=0.0) { // null rotation
        r[0][0] = rot[0]; r[0][1] = rot[1]; r[0][2] = rot[2];
        r[1][0] = rot[3]; r[1][1] = rot[4]; r[1][2] = rot[5];
        r[2][0] = rot[6]; r[2][1] = rot[7]; r[2][2] = rot[8];
    } // end if
    fGm.AddAt(new AliITSgeomMatrix(idet,id,r,tran),mod);
}
//______________________________________________________________________
AliITSgeom::~AliITSgeom(){
    //     The destructor for the AliITSgeom class. If the arrays fNlad,
    // fNdet, or fGm have had memory allocated to them, there pointer values
    // are non zero, then this memory space is freed and they are set
    // to zero. In addition, fNlayers is set to zero. The destruction of
    // Inputs:
    //    none.
    // Outputs:
    //    none.
    // Return:
    //    none.

    return;
}
//______________________________________________________________________
AliITSgeom::AliITSgeom(const AliITSgeom &source) : 
TObject(source),
fVersion(source.fVersion), // Transformation version.
fTrans(source.fTrans),   // Flag to keep track of which transformation
fNmodules(source.fNmodules),// The total number of modules
fNlayers(source.fNlayers), // The number of layers.
fNlad(source.fNlad),    // Array of the number of ladders/layer(layer)
fNdet(source.fNdet),    // Array of the number of detector/ladder(layer)
fGm(source.fGm.GetSize(),source.fGm.LowerBound())// Structure of 
                                                  // translation and rotation.
{
    //     The copy constructor for the AliITSgeom class. It calls the
    // = operator function. See the = operator function for more details.
    // Inputs:
    //     AliITSgeom &source  The AliITSgeom class with which to make this
    //                         a copy of.
    // Outputs:
    //     none.
    // Return:
    //     none.
    Int_t i,n;

    n = source.fGm.GetLast()+1;
    for(i=source.fGm.LowerBound();i<n;i++){
        fGm.AddAt(new AliITSgeomMatrix(*((AliITSgeomMatrix*)(
                                             source.fGm.At(i)))),i);
    } // end for i
    fGm.SetOwner(kTRUE);
    return;
}
//______________________________________________________________________
AliITSgeom& AliITSgeom::operator=(const AliITSgeom &source){
    //     The = operator function for the AliITSgeom class. It makes an
    // independent copy of the class in such a way that any changes made
    // to the copied class will not affect the source class in any way.
    // This is required for many ITS alignment studies where the copied
    // class is then modified by introducing some misalignment.
    // Inputs:
    //     AliITSgeom &source  The AliITSgeom class with which to make this
    //                         a copy of.
    // Outputs:
    //     none.
    // Return:
    //     *this The a new copy of source.
    Int_t i;

    if(this == &source) return *this; // don't assign to ones self.

    // if there is an old structure allocated delete it first.
    this->fGm.Clear();

    this->fVersion  = source.fVersion;
    this->fTrans    = source.fTrans;
    this->fNmodules = source.fNmodules;
    this->fNlayers  = source.fNlayers;
    this->fNlad     = source.fNlad;
    this->fNdet     = source.fNdet;
    this->fGm.Expand(this->fNmodules);
    for(i=source.fGm.LowerBound();i<source.fGm.GetLast();i++){
        fGm.AddAt(new AliITSgeomMatrix(*((AliITSgeomMatrix*)(
                                             source.fGm.At(i)))),i);
    } // end for i
    fGm.SetOwner(kTRUE);
    return *this;
}
//______________________________________________________________________
Int_t AliITSgeom::GetModuleIndex(Int_t lay,Int_t lad,Int_t det)const{
    //      This routine computes the module index number from the layer,
    // ladder, and detector numbers. The number of ladders and detectors
    // per layer is determined when this geometry package is constructed,
    // see AliITSgeom(const char *filename) for specifics.
    // Inputs:
    //    Int_t lay  The layer number. Starting from 1.
    //    Int_t lad  The ladder number. Starting from 1.
    //    Int_t det  The detector number. Starting from 1.
    // Outputs:
    //    none.
    // Return:
    //    the module index number, starting from zero.
    Int_t i,j,k,id[3];

    i = fNdet[lay-1] * (lad-1) + det - 1;
    j = 0;
    for(k=0;k<lay-1;k++) j += fNdet[k]*fNlad[k];
    i = i+j;
    if(i>=fNmodules) return -1;
    GetGeomMatrix(i)->GetIndex(id);
    if(id[0]==lay&&id[1]==lad&&id[2]==det) return i;
    // Array of modules fGm is not in expected order. Search for this index
    for(i=0;i<fNmodules;i++){
        GetGeomMatrix(i)->GetIndex(id);
        if(id[0]==lay&&id[1]==lad&&id[2]==det) return i;
    } // end for i
    // This layer ladder and detector combination does not exist return -1.
    return -1;
}
//______________________________________________________________________
void AliITSgeom::GetModuleId(Int_t index,Int_t &lay,Int_t &lad,Int_t &det)
const{
    //      This routine computes the layer, ladder and detector number 
    // given the module index number. The number of ladders and detectors
    // per layer is determined when this geometry package is constructed,
    // see AliITSgeom(const char *filename) for specifics.
    // Inputs:
    //     Int_t index  The module index number, starting from zero.
    // Outputs:
    //     Int_t lay    The layer number. Starting from 1.
    //     Int_t lad    The ladder number. Starting from 1.
    //     Int_t det    The detector number. Starting from 1.
    // Return:
    //     none.
    Int_t id[3];
    AliITSgeomMatrix *g = GetGeomMatrix(index);

    if (g == 0x0){
        Error("GetModuleId","Can not get GeoMatrix for index = %d",index);
        lay = -1; lad = -1; det = -1;
    }else{
        g->GetIndex(id);
        lay = id[0]; lad = id[1]; det = id[2];
    }// End if
    return;
    // The old way kept for posterity.
/*
    Int_t i,j,k;
    j = 0;
    for(k=0;k<fNlayers;k++){
	j += fNdet[k]*fNlad[k];
	if(j>index)break;
    } // end for k
    lay = k+1;
    i = index -j + fNdet[k]*fNlad[k];
    j = 0;
    for(k=0;k<fNlad[lay-1];k++){
	j += fNdet[lay-1];
	if(j>i)break;
    } // end for k
    lad = k+1;
    det = 1+i-fNdet[lay-1]*k;
    return;
*/
}
//______________________________________________________________________
Int_t AliITSgeom::GetNDetTypes(Int_t &max)const{
    // Finds and returns the number of detector types used and the
    // maximum detector type value. Only counts id >=0 (no undefined
    // values. See AliITSgeom.h for list of AliITSDetecor enumerated types.
    // Inputs:
    //    none.
    // Outputs:
    //    The maximum detector type used
    // Return:
    //    The number of detector types used
    Int_t i,*n,id;

    max = -1;
    for(i=0;i<GetIndexMax();i++){
        id = GetModuleType(i);
        if(id>max) max=id;
    } // end for i
    n = new Int_t[max+1];
    for(i=0;i<max;i++) n[i] = 0;
    for(i=0;i<GetIndexMax();i++){
        id = GetModuleType(i);
        if(id>-1)n[id]++; // note id=-1 => undefined.
    } // end for i
    id = 0;
    for(i=0;i<max;i++) if(n[i]!=0) id++;
    delete[] n;
    return id+1;
}
//______________________________________________________________________
Int_t AliITSgeom::GetNDetTypes(TArrayI &maxs,AliITSDetector *types)const{
    // Finds and returns the number of detector types used and the
    // number of each detector type. Only counts id >=0 (no undefined
    // values. See AliITSgeom.h for list of AliITSDetecor enumerated types.
    // Inputs:
    //    none.
    // Outputs:
    //    The maximum detector type used
    // Return:
    //    The number of detector types used
    Int_t i,j,*n,id,max;

    max = -1;
    for(i=0;i<GetIndexMax();i++){
        id = GetModuleType(i);
        if(id>max) max=id;
    } // end for i
    n = new Int_t[max+1];
    for(i=0;i<max;i++) n[i] = 0;
    for(i=0;i<GetIndexMax();i++){
        id = GetModuleType(i);
        if(id>-1)n[id]++; // note id=-1 => undefined.
    } // end for i
    id = 0;
    for(i=0;i<=max;i++) if(n[i]!=0) id++;
    maxs.Set(id);
    j = 0;
    for(i=0;i<=max;i++) if(n[i]!=0){
        maxs[j] = n[i];
        types[j++] = (AliITSDetector) i;
    } // end for i/end if
    delete[] n;
    return id;
}
//______________________________________________________________________
Int_t AliITSgeom::GetStartDet(Int_t dtype)const{
    // returns the starting module index value for a give type of detector id.
    // This assumes that the detector types are different on different layers
    // and that they are not mixed up.
    // Inputs:
    //    Int_t dtype A detector type number. 0 for SPD, 1 for SDD, 
    //                and 2 for SSD.
    // Outputs:
    //    none.
    // Return:
    //    the module index for the first occurrence of that detector type.

    switch(dtype){
    case 0:
        return GetModuleIndex(1,1,1);
        break;
    case 1:
        return GetModuleIndex(3,1,1);
        break;
    case 2:
        return GetModuleIndex(5,1,1);
        break;
    default:
        Warning("GetStartDet","undefined detector type %d",dtype);
        return 0;
    } // end switch

    Warning("GetStartDet","undefined detector type %d",dtype);
    return 0;
}
//______________________________________________________________________
Int_t AliITSgeom::GetLastDet(Int_t dtype)const{
    // returns the last module index value for a give type of detector id.
    // This assumes that the detector types are different on different layers
    // and that they are not mixed up.
    // Inputs:
    //     Int_t dtype A detector type number. 0 for SPD, 1 for SDD, 
    //                 and 2 for SSD.
    // Outputs:
    // Return:
    //     the module index for the last occurrence of that detector type.

    switch((AliITSDetector)dtype){
    case kSPD:
        return GetModuleIndex(3,1,1)-1;
        break;
    case kSDD:
        return GetModuleIndex(5,1,1)-1;
        break;
    case kSSD:
        return GetIndexMax()-1;
        break;
    case kSSDp: case kSDDp: case kND:
    default:
        Warning("GetLastDet","undefined detector type %d",dtype);
        return 0;
    } // end switch

    Warning("GetLastDet","undefined detector type %d",dtype);
    return 0;
}

//______________________________________________________________________
void AliITSgeom::PrintData(FILE *fp,Int_t lay,Int_t lad,Int_t det)const{
    //     This function prints out the coordinate transformations for
    // the particular detector defined by layer, ladder, and detector
    // to the file pointed to by the File pointer fp. fprintf statements
    // are used to print out the numbers. The format is
    // layer ladder detector Trans= fx0 fy0 fz0 rot= frx fry frz 
    // Shape=fShapeIndex
    //                         dfr= fr[0] fr[1] fr[2]
    //                         dfr= fr[3] fr[4] fr[5]
    //                         dfr= fr[6] fr[7] fr[8]
    // By indicating which detector, some control over the information 
    // is given to the user. The output it written to the file pointed
    // to by the file pointer fp. This can be set to stdout if you want.
    // Inputs:
    //     FILE *fp           A file pointer to an opened file for 
    //                        writing in which the results of the 
    //                        comparison will be written.
    //     Int_t lay          The layer number. Starting from 1.
    //     Int_t lad          The ladder number. Starting from 1.
    //     Int_t det          The detector number. Starting from 1.
    // Outputs:
    //     none
    // Return:
    //     none.
    AliITSgeomMatrix *gt;
    Double_t t[3],r[3],m[3][3];

    gt = this->GetGeomMatrix(GetModuleIndex(lay,lad,det));
    gt->GetTranslation(t);
    gt->GetAngles(r);
    fprintf(fp,"%1.1d %2.2d %2.2d Trans=%f %f %f rot=%f %f %f Shape=%d\n",
            lay,lad,det,t[0],t[1],t[2],r[0],r[1],r[2],
            gt->GetDetectorIndex());
    gt->GetMatrix(m);
    fprintf(fp,"        dfr= %e %e %e\n",m[0][0],m[0][1],m[0][2]);
    fprintf(fp,"        dfr= %e %e %e\n",m[1][0],m[1][1],m[1][2]);
    fprintf(fp,"        dfr= %e %e %e\n",m[2][0],m[2][1],m[2][2]);
    return;
}

//______________________________________________________________________
Int_t AliITSgeom::GetNearest(const Double_t g[3],Int_t lay)const{
    //      Finds the Detector (Module) that is nearest the point g [cm] in
    // ALICE Global coordinates. If layer !=0 then the search is restricted
    // to Detectors (Modules) in that particular layer.
    // Inputs:
    //     Double_t g[3]  The ALICE Cartesian global coordinate from which the
    //                    distance is to be calculated with.
    //     Int_t lay      The layer to restrict the search to. If layer=0 then
    //                    all layers are searched. Default is lay=0.
    // Output:
    //     none.
    // Return:
    //     The module number representing the nearest module.
    Int_t    i,l,a,e,in=0;
    Double_t d,dn=1.0e10;
    Bool_t   t=lay!=0; // skip if lay = 0 default value check all layers.

    for(i=0;i<fNmodules;i++){
        if(t){GetModuleId(i,l,a,e);if(l!=lay) continue;}
        if((d=GetGeomMatrix(i)->Distance2(g))<dn){
            dn = d;
            in = i;
        } // end if
    } // end for i
    return in;
}
//______________________________________________________________________
void AliITSgeom::GetNearest27(const Double_t g[3],Int_t n[27],Int_t lay)const{
    //      Finds 27 Detectors (Modules) that are nearest the point g [cm] in
    // ALICE Global coordinates. If layer !=0 then the search is restricted
    // to Detectors (Modules) in that particular layer. The number 27 comes 
    // from including the nearest detector and all those around it (up, down,
    // left, right, forwards, backwards, and the corners).
    // Input:
    //     Double_t g[3]  The ALICE Cartesian global coordinate from which the
    //                    distance is to be calculated with.
    //     Int_t lay      The layer to restrict the search to. If layer=0 then
    //                    all layers are searched. Default is lay=0.
    // Output:
    //     Int_t n[27]    The module number representing the nearest 27 modules
    //                    in order.
    // Return:
    //     none.
    Int_t    i,l,a,e,in[27]={0,0,0,0,0,0,0,0,0,
                             0,0,0,0,0,0,0,0,0,
                             0,0,0,0,0,0,0,0,0,};
    Double_t d,dn[27]={1.0e10,1.0e10,1.0e10,1.0e10,1.0e10,1.0e10,
                       1.0e10,1.0e10,1.0e10,1.0e10,1.0e10,1.0e10,
                       1.0e10,1.0e10,1.0e10,1.0e10,1.0e10,1.0e10,
                       1.0e10,1.0e10,1.0e10,1.0e10,1.0e10,1.0e10,
                       1.0e10,1.0e10,1.0e10};
    Bool_t   t=(lay!=0); // skip if lay = 0 default value check all layers.

    for(i=0;i<fNmodules;i++){
        if(t){GetModuleId(i,l,a,e);if(l!=lay) continue;}
        for(a=0;a<27;a++){
            d = GetGeomMatrix(i)->Distance2(g);
            if(d<dn[a]){
                for(e=26;e>a;e--){dn[e] = dn[e-1];in[e] = in[e-1];}
                dn[a] = d; in[a] = i;
            } // end if d<dn[i]
        } // end for a
    } // end for i
    for(i=0;i<27;i++) n[i] = in[i];
}
//_______________________________________________________________________
void AliITSgeom::DetLToTrackingV2(Int_t md,Float_t xin,Float_t zin,
                                  Float_t &yout,Float_t &zout) const {

    //Conversion from local coordinates on detectors to local
    //coordinates used for tracking ("v2")
    // Inputs:
    //   Int_t   md      Module number
    //   Float_t xin     Standard local coordinate x
    //   Float_t zin     Standard local coordinate z
    // Output:
    //   Float_t yout    Tracking local coordinate y
    //   Float_t zout    Tracking local coordinate z
    // Return:
    //   none.
    Float_t x,y,z;
    Double_t rt[9],al;

    GetTrans(md,x,y,z);
    GetRotMatrix(md,rt);
    al = TMath::ATan2(rt[1],rt[0])+TMath::Pi();
    yout = -(-xin+(x*((Float_t)TMath::Cos(al))+y*((Float_t)TMath::Sin(al))));
    if(md<(GetModuleIndex(2,1,1))) yout *= -1; 
    zout = -zin+z; 
}
//_______________________________________________________________________
void AliITSgeom::TrackingV2ToDetL(Int_t md,Float_t yin,Float_t zin,
                                  Float_t &xout,Float_t &zout) const {
    //Conversion from local coordinates used for tracking ("v2") to
    //local detector coordinates  
    // Inputs:
    //   Int_t   md      Module number
    //   Float_t yin     Tracking local coordinate y
    //   Float_t zin     Tracking local coordinate z
    // Output:
    //   Float_t xout    Standard local coordinate x
    //   Float_t zout    Standard local coordinate z
    // Return:
    //   none.
    Float_t x,y,z;
    Double_t rt[9],al;

    GetTrans(md,x,y,z);
    GetRotMatrix(md,rt);
    al = TMath::ATan2(rt[1],rt[0])+TMath::Pi();
    xout = yin;
    if(md<(GetModuleIndex(2,1,1))) xout = -xout;
    xout += (x*((Float_t)TMath::Cos(al))+y*((Float_t)TMath::Sin(al)));
    zout  = -zin+z; 
}
//----------------------------------------------------------------------

 AliITSgeom.cxx:1
 AliITSgeom.cxx:2
 AliITSgeom.cxx:3
 AliITSgeom.cxx:4
 AliITSgeom.cxx:5
 AliITSgeom.cxx:6
 AliITSgeom.cxx:7
 AliITSgeom.cxx:8
 AliITSgeom.cxx:9
 AliITSgeom.cxx:10
 AliITSgeom.cxx:11
 AliITSgeom.cxx:12
 AliITSgeom.cxx:13
 AliITSgeom.cxx:14
 AliITSgeom.cxx:15
 AliITSgeom.cxx:16
 AliITSgeom.cxx:17
 AliITSgeom.cxx:18
 AliITSgeom.cxx:19
 AliITSgeom.cxx:20
 AliITSgeom.cxx:21
 AliITSgeom.cxx:22
 AliITSgeom.cxx:23
 AliITSgeom.cxx:24
 AliITSgeom.cxx:25
 AliITSgeom.cxx:26
 AliITSgeom.cxx:27
 AliITSgeom.cxx:28
 AliITSgeom.cxx:29
 AliITSgeom.cxx:30
 AliITSgeom.cxx:31
 AliITSgeom.cxx:32
 AliITSgeom.cxx:33
 AliITSgeom.cxx:34
 AliITSgeom.cxx:35
 AliITSgeom.cxx:36
 AliITSgeom.cxx:37
 AliITSgeom.cxx:38
 AliITSgeom.cxx:39
 AliITSgeom.cxx:40
 AliITSgeom.cxx:41
 AliITSgeom.cxx:42
 AliITSgeom.cxx:43
 AliITSgeom.cxx:44
 AliITSgeom.cxx:45
 AliITSgeom.cxx:46
 AliITSgeom.cxx:47
 AliITSgeom.cxx:48
 AliITSgeom.cxx:49
 AliITSgeom.cxx:50
 AliITSgeom.cxx:51
 AliITSgeom.cxx:52
 AliITSgeom.cxx:53
 AliITSgeom.cxx:54
 AliITSgeom.cxx:55
 AliITSgeom.cxx:56
 AliITSgeom.cxx:57
 AliITSgeom.cxx:58
 AliITSgeom.cxx:59
 AliITSgeom.cxx:60
 AliITSgeom.cxx:61
 AliITSgeom.cxx:62
 AliITSgeom.cxx:63
 AliITSgeom.cxx:64
 AliITSgeom.cxx:65
 AliITSgeom.cxx:66
 AliITSgeom.cxx:67
 AliITSgeom.cxx:68
 AliITSgeom.cxx:69
 AliITSgeom.cxx:70
 AliITSgeom.cxx:71
 AliITSgeom.cxx:72
 AliITSgeom.cxx:73
 AliITSgeom.cxx:74
 AliITSgeom.cxx:75
 AliITSgeom.cxx:76
 AliITSgeom.cxx:77
 AliITSgeom.cxx:78
 AliITSgeom.cxx:79
 AliITSgeom.cxx:80
 AliITSgeom.cxx:81
 AliITSgeom.cxx:82
 AliITSgeom.cxx:83
 AliITSgeom.cxx:84
 AliITSgeom.cxx:85
 AliITSgeom.cxx:86
 AliITSgeom.cxx:87
 AliITSgeom.cxx:88
 AliITSgeom.cxx:89
 AliITSgeom.cxx:90
 AliITSgeom.cxx:91
 AliITSgeom.cxx:92
 AliITSgeom.cxx:93
 AliITSgeom.cxx:94
 AliITSgeom.cxx:95
 AliITSgeom.cxx:96
 AliITSgeom.cxx:97
 AliITSgeom.cxx:98
 AliITSgeom.cxx:99
 AliITSgeom.cxx:100
 AliITSgeom.cxx:101
 AliITSgeom.cxx:102
 AliITSgeom.cxx:103
 AliITSgeom.cxx:104
 AliITSgeom.cxx:105
 AliITSgeom.cxx:106
 AliITSgeom.cxx:107
 AliITSgeom.cxx:108
 AliITSgeom.cxx:109
 AliITSgeom.cxx:110
 AliITSgeom.cxx:111
 AliITSgeom.cxx:112
 AliITSgeom.cxx:113
 AliITSgeom.cxx:114
 AliITSgeom.cxx:115
 AliITSgeom.cxx:116
 AliITSgeom.cxx:117
 AliITSgeom.cxx:118
 AliITSgeom.cxx:119
 AliITSgeom.cxx:120
 AliITSgeom.cxx:121
 AliITSgeom.cxx:122
 AliITSgeom.cxx:123
 AliITSgeom.cxx:124
 AliITSgeom.cxx:125
 AliITSgeom.cxx:126
 AliITSgeom.cxx:127
 AliITSgeom.cxx:128
 AliITSgeom.cxx:129
 AliITSgeom.cxx:130
 AliITSgeom.cxx:131
 AliITSgeom.cxx:132
 AliITSgeom.cxx:133
 AliITSgeom.cxx:134
 AliITSgeom.cxx:135
 AliITSgeom.cxx:136
 AliITSgeom.cxx:137
 AliITSgeom.cxx:138
 AliITSgeom.cxx:139
 AliITSgeom.cxx:140
 AliITSgeom.cxx:141
 AliITSgeom.cxx:142
 AliITSgeom.cxx:143
 AliITSgeom.cxx:144
 AliITSgeom.cxx:145
 AliITSgeom.cxx:146
 AliITSgeom.cxx:147
 AliITSgeom.cxx:148
 AliITSgeom.cxx:149
 AliITSgeom.cxx:150
 AliITSgeom.cxx:151
 AliITSgeom.cxx:152
 AliITSgeom.cxx:153
 AliITSgeom.cxx:154
 AliITSgeom.cxx:155
 AliITSgeom.cxx:156
 AliITSgeom.cxx:157
 AliITSgeom.cxx:158
 AliITSgeom.cxx:159
 AliITSgeom.cxx:160
 AliITSgeom.cxx:161
 AliITSgeom.cxx:162
 AliITSgeom.cxx:163
 AliITSgeom.cxx:164
 AliITSgeom.cxx:165
 AliITSgeom.cxx:166
 AliITSgeom.cxx:167
 AliITSgeom.cxx:168
 AliITSgeom.cxx:169
 AliITSgeom.cxx:170
 AliITSgeom.cxx:171
 AliITSgeom.cxx:172
 AliITSgeom.cxx:173
 AliITSgeom.cxx:174
 AliITSgeom.cxx:175
 AliITSgeom.cxx:176
 AliITSgeom.cxx:177
 AliITSgeom.cxx:178
 AliITSgeom.cxx:179
 AliITSgeom.cxx:180
 AliITSgeom.cxx:181
 AliITSgeom.cxx:182
 AliITSgeom.cxx:183
 AliITSgeom.cxx:184
 AliITSgeom.cxx:185
 AliITSgeom.cxx:186
 AliITSgeom.cxx:187
 AliITSgeom.cxx:188
 AliITSgeom.cxx:189
 AliITSgeom.cxx:190
 AliITSgeom.cxx:191
 AliITSgeom.cxx:192
 AliITSgeom.cxx:193
 AliITSgeom.cxx:194
 AliITSgeom.cxx:195
 AliITSgeom.cxx:196
 AliITSgeom.cxx:197
 AliITSgeom.cxx:198
 AliITSgeom.cxx:199
 AliITSgeom.cxx:200
 AliITSgeom.cxx:201
 AliITSgeom.cxx:202
 AliITSgeom.cxx:203
 AliITSgeom.cxx:204
 AliITSgeom.cxx:205
 AliITSgeom.cxx:206
 AliITSgeom.cxx:207
 AliITSgeom.cxx:208
 AliITSgeom.cxx:209
 AliITSgeom.cxx:210
 AliITSgeom.cxx:211
 AliITSgeom.cxx:212
 AliITSgeom.cxx:213
 AliITSgeom.cxx:214
 AliITSgeom.cxx:215
 AliITSgeom.cxx:216
 AliITSgeom.cxx:217
 AliITSgeom.cxx:218
 AliITSgeom.cxx:219
 AliITSgeom.cxx:220
 AliITSgeom.cxx:221
 AliITSgeom.cxx:222
 AliITSgeom.cxx:223
 AliITSgeom.cxx:224
 AliITSgeom.cxx:225
 AliITSgeom.cxx:226
 AliITSgeom.cxx:227
 AliITSgeom.cxx:228
 AliITSgeom.cxx:229
 AliITSgeom.cxx:230
 AliITSgeom.cxx:231
 AliITSgeom.cxx:232
 AliITSgeom.cxx:233
 AliITSgeom.cxx:234
 AliITSgeom.cxx:235
 AliITSgeom.cxx:236
 AliITSgeom.cxx:237
 AliITSgeom.cxx:238
 AliITSgeom.cxx:239
 AliITSgeom.cxx:240
 AliITSgeom.cxx:241
 AliITSgeom.cxx:242
 AliITSgeom.cxx:243
 AliITSgeom.cxx:244
 AliITSgeom.cxx:245
 AliITSgeom.cxx:246
 AliITSgeom.cxx:247
 AliITSgeom.cxx:248
 AliITSgeom.cxx:249
 AliITSgeom.cxx:250
 AliITSgeom.cxx:251
 AliITSgeom.cxx:252
 AliITSgeom.cxx:253
 AliITSgeom.cxx:254
 AliITSgeom.cxx:255
 AliITSgeom.cxx:256
 AliITSgeom.cxx:257
 AliITSgeom.cxx:258
 AliITSgeom.cxx:259
 AliITSgeom.cxx:260
 AliITSgeom.cxx:261
 AliITSgeom.cxx:262
 AliITSgeom.cxx:263
 AliITSgeom.cxx:264
 AliITSgeom.cxx:265
 AliITSgeom.cxx:266
 AliITSgeom.cxx:267
 AliITSgeom.cxx:268
 AliITSgeom.cxx:269
 AliITSgeom.cxx:270
 AliITSgeom.cxx:271
 AliITSgeom.cxx:272
 AliITSgeom.cxx:273
 AliITSgeom.cxx:274
 AliITSgeom.cxx:275
 AliITSgeom.cxx:276
 AliITSgeom.cxx:277
 AliITSgeom.cxx:278
 AliITSgeom.cxx:279
 AliITSgeom.cxx:280
 AliITSgeom.cxx:281
 AliITSgeom.cxx:282
 AliITSgeom.cxx:283
 AliITSgeom.cxx:284
 AliITSgeom.cxx:285
 AliITSgeom.cxx:286
 AliITSgeom.cxx:287
 AliITSgeom.cxx:288
 AliITSgeom.cxx:289
 AliITSgeom.cxx:290
 AliITSgeom.cxx:291
 AliITSgeom.cxx:292
 AliITSgeom.cxx:293
 AliITSgeom.cxx:294
 AliITSgeom.cxx:295
 AliITSgeom.cxx:296
 AliITSgeom.cxx:297
 AliITSgeom.cxx:298
 AliITSgeom.cxx:299
 AliITSgeom.cxx:300
 AliITSgeom.cxx:301
 AliITSgeom.cxx:302
 AliITSgeom.cxx:303
 AliITSgeom.cxx:304
 AliITSgeom.cxx:305
 AliITSgeom.cxx:306
 AliITSgeom.cxx:307
 AliITSgeom.cxx:308
 AliITSgeom.cxx:309
 AliITSgeom.cxx:310
 AliITSgeom.cxx:311
 AliITSgeom.cxx:312
 AliITSgeom.cxx:313
 AliITSgeom.cxx:314
 AliITSgeom.cxx:315
 AliITSgeom.cxx:316
 AliITSgeom.cxx:317
 AliITSgeom.cxx:318
 AliITSgeom.cxx:319
 AliITSgeom.cxx:320
 AliITSgeom.cxx:321
 AliITSgeom.cxx:322
 AliITSgeom.cxx:323
 AliITSgeom.cxx:324
 AliITSgeom.cxx:325
 AliITSgeom.cxx:326
 AliITSgeom.cxx:327
 AliITSgeom.cxx:328
 AliITSgeom.cxx:329
 AliITSgeom.cxx:330
 AliITSgeom.cxx:331
 AliITSgeom.cxx:332
 AliITSgeom.cxx:333
 AliITSgeom.cxx:334
 AliITSgeom.cxx:335
 AliITSgeom.cxx:336
 AliITSgeom.cxx:337
 AliITSgeom.cxx:338
 AliITSgeom.cxx:339
 AliITSgeom.cxx:340
 AliITSgeom.cxx:341
 AliITSgeom.cxx:342
 AliITSgeom.cxx:343
 AliITSgeom.cxx:344
 AliITSgeom.cxx:345
 AliITSgeom.cxx:346
 AliITSgeom.cxx:347
 AliITSgeom.cxx:348
 AliITSgeom.cxx:349
 AliITSgeom.cxx:350
 AliITSgeom.cxx:351
 AliITSgeom.cxx:352
 AliITSgeom.cxx:353
 AliITSgeom.cxx:354
 AliITSgeom.cxx:355
 AliITSgeom.cxx:356
 AliITSgeom.cxx:357
 AliITSgeom.cxx:358
 AliITSgeom.cxx:359
 AliITSgeom.cxx:360
 AliITSgeom.cxx:361
 AliITSgeom.cxx:362
 AliITSgeom.cxx:363
 AliITSgeom.cxx:364
 AliITSgeom.cxx:365
 AliITSgeom.cxx:366
 AliITSgeom.cxx:367
 AliITSgeom.cxx:368
 AliITSgeom.cxx:369
 AliITSgeom.cxx:370
 AliITSgeom.cxx:371
 AliITSgeom.cxx:372
 AliITSgeom.cxx:373
 AliITSgeom.cxx:374
 AliITSgeom.cxx:375
 AliITSgeom.cxx:376
 AliITSgeom.cxx:377
 AliITSgeom.cxx:378
 AliITSgeom.cxx:379
 AliITSgeom.cxx:380
 AliITSgeom.cxx:381
 AliITSgeom.cxx:382
 AliITSgeom.cxx:383
 AliITSgeom.cxx:384
 AliITSgeom.cxx:385
 AliITSgeom.cxx:386
 AliITSgeom.cxx:387
 AliITSgeom.cxx:388
 AliITSgeom.cxx:389
 AliITSgeom.cxx:390
 AliITSgeom.cxx:391
 AliITSgeom.cxx:392
 AliITSgeom.cxx:393
 AliITSgeom.cxx:394
 AliITSgeom.cxx:395
 AliITSgeom.cxx:396
 AliITSgeom.cxx:397
 AliITSgeom.cxx:398
 AliITSgeom.cxx:399
 AliITSgeom.cxx:400
 AliITSgeom.cxx:401
 AliITSgeom.cxx:402
 AliITSgeom.cxx:403
 AliITSgeom.cxx:404
 AliITSgeom.cxx:405
 AliITSgeom.cxx:406
 AliITSgeom.cxx:407
 AliITSgeom.cxx:408
 AliITSgeom.cxx:409
 AliITSgeom.cxx:410
 AliITSgeom.cxx:411
 AliITSgeom.cxx:412
 AliITSgeom.cxx:413
 AliITSgeom.cxx:414
 AliITSgeom.cxx:415
 AliITSgeom.cxx:416
 AliITSgeom.cxx:417
 AliITSgeom.cxx:418
 AliITSgeom.cxx:419
 AliITSgeom.cxx:420
 AliITSgeom.cxx:421
 AliITSgeom.cxx:422
 AliITSgeom.cxx:423
 AliITSgeom.cxx:424
 AliITSgeom.cxx:425
 AliITSgeom.cxx:426
 AliITSgeom.cxx:427
 AliITSgeom.cxx:428
 AliITSgeom.cxx:429
 AliITSgeom.cxx:430
 AliITSgeom.cxx:431
 AliITSgeom.cxx:432
 AliITSgeom.cxx:433
 AliITSgeom.cxx:434
 AliITSgeom.cxx:435
 AliITSgeom.cxx:436
 AliITSgeom.cxx:437
 AliITSgeom.cxx:438
 AliITSgeom.cxx:439
 AliITSgeom.cxx:440
 AliITSgeom.cxx:441
 AliITSgeom.cxx:442
 AliITSgeom.cxx:443
 AliITSgeom.cxx:444
 AliITSgeom.cxx:445
 AliITSgeom.cxx:446
 AliITSgeom.cxx:447
 AliITSgeom.cxx:448
 AliITSgeom.cxx:449
 AliITSgeom.cxx:450
 AliITSgeom.cxx:451
 AliITSgeom.cxx:452
 AliITSgeom.cxx:453
 AliITSgeom.cxx:454
 AliITSgeom.cxx:455
 AliITSgeom.cxx:456
 AliITSgeom.cxx:457
 AliITSgeom.cxx:458
 AliITSgeom.cxx:459
 AliITSgeom.cxx:460
 AliITSgeom.cxx:461
 AliITSgeom.cxx:462
 AliITSgeom.cxx:463
 AliITSgeom.cxx:464
 AliITSgeom.cxx:465
 AliITSgeom.cxx:466
 AliITSgeom.cxx:467
 AliITSgeom.cxx:468
 AliITSgeom.cxx:469
 AliITSgeom.cxx:470
 AliITSgeom.cxx:471
 AliITSgeom.cxx:472
 AliITSgeom.cxx:473
 AliITSgeom.cxx:474
 AliITSgeom.cxx:475
 AliITSgeom.cxx:476
 AliITSgeom.cxx:477
 AliITSgeom.cxx:478
 AliITSgeom.cxx:479
 AliITSgeom.cxx:480
 AliITSgeom.cxx:481
 AliITSgeom.cxx:482
 AliITSgeom.cxx:483
 AliITSgeom.cxx:484
 AliITSgeom.cxx:485
 AliITSgeom.cxx:486
 AliITSgeom.cxx:487
 AliITSgeom.cxx:488
 AliITSgeom.cxx:489
 AliITSgeom.cxx:490
 AliITSgeom.cxx:491
 AliITSgeom.cxx:492
 AliITSgeom.cxx:493
 AliITSgeom.cxx:494
 AliITSgeom.cxx:495
 AliITSgeom.cxx:496
 AliITSgeom.cxx:497
 AliITSgeom.cxx:498
 AliITSgeom.cxx:499
 AliITSgeom.cxx:500
 AliITSgeom.cxx:501
 AliITSgeom.cxx:502
 AliITSgeom.cxx:503
 AliITSgeom.cxx:504
 AliITSgeom.cxx:505
 AliITSgeom.cxx:506
 AliITSgeom.cxx:507
 AliITSgeom.cxx:508
 AliITSgeom.cxx:509
 AliITSgeom.cxx:510
 AliITSgeom.cxx:511
 AliITSgeom.cxx:512
 AliITSgeom.cxx:513
 AliITSgeom.cxx:514
 AliITSgeom.cxx:515
 AliITSgeom.cxx:516
 AliITSgeom.cxx:517
 AliITSgeom.cxx:518
 AliITSgeom.cxx:519
 AliITSgeom.cxx:520
 AliITSgeom.cxx:521
 AliITSgeom.cxx:522
 AliITSgeom.cxx:523
 AliITSgeom.cxx:524
 AliITSgeom.cxx:525
 AliITSgeom.cxx:526
 AliITSgeom.cxx:527
 AliITSgeom.cxx:528
 AliITSgeom.cxx:529
 AliITSgeom.cxx:530
 AliITSgeom.cxx:531
 AliITSgeom.cxx:532
 AliITSgeom.cxx:533
 AliITSgeom.cxx:534
 AliITSgeom.cxx:535
 AliITSgeom.cxx:536
 AliITSgeom.cxx:537
 AliITSgeom.cxx:538
 AliITSgeom.cxx:539
 AliITSgeom.cxx:540
 AliITSgeom.cxx:541
 AliITSgeom.cxx:542
 AliITSgeom.cxx:543
 AliITSgeom.cxx:544
 AliITSgeom.cxx:545
 AliITSgeom.cxx:546
 AliITSgeom.cxx:547
 AliITSgeom.cxx:548
 AliITSgeom.cxx:549
 AliITSgeom.cxx:550
 AliITSgeom.cxx:551
 AliITSgeom.cxx:552
 AliITSgeom.cxx:553
 AliITSgeom.cxx:554
 AliITSgeom.cxx:555
 AliITSgeom.cxx:556
 AliITSgeom.cxx:557
 AliITSgeom.cxx:558
 AliITSgeom.cxx:559
 AliITSgeom.cxx:560
 AliITSgeom.cxx:561
 AliITSgeom.cxx:562
 AliITSgeom.cxx:563
 AliITSgeom.cxx:564
 AliITSgeom.cxx:565
 AliITSgeom.cxx:566
 AliITSgeom.cxx:567
 AliITSgeom.cxx:568
 AliITSgeom.cxx:569
 AliITSgeom.cxx:570
 AliITSgeom.cxx:571
 AliITSgeom.cxx:572
 AliITSgeom.cxx:573
 AliITSgeom.cxx:574
 AliITSgeom.cxx:575
 AliITSgeom.cxx:576
 AliITSgeom.cxx:577
 AliITSgeom.cxx:578
 AliITSgeom.cxx:579
 AliITSgeom.cxx:580
 AliITSgeom.cxx:581
 AliITSgeom.cxx:582
 AliITSgeom.cxx:583
 AliITSgeom.cxx:584
 AliITSgeom.cxx:585
 AliITSgeom.cxx:586
 AliITSgeom.cxx:587
 AliITSgeom.cxx:588
 AliITSgeom.cxx:589
 AliITSgeom.cxx:590
 AliITSgeom.cxx:591
 AliITSgeom.cxx:592
 AliITSgeom.cxx:593
 AliITSgeom.cxx:594
 AliITSgeom.cxx:595
 AliITSgeom.cxx:596
 AliITSgeom.cxx:597
 AliITSgeom.cxx:598
 AliITSgeom.cxx:599
 AliITSgeom.cxx:600
 AliITSgeom.cxx:601
 AliITSgeom.cxx:602
 AliITSgeom.cxx:603
 AliITSgeom.cxx:604
 AliITSgeom.cxx:605
 AliITSgeom.cxx:606
 AliITSgeom.cxx:607
 AliITSgeom.cxx:608
 AliITSgeom.cxx:609
 AliITSgeom.cxx:610
 AliITSgeom.cxx:611
 AliITSgeom.cxx:612
 AliITSgeom.cxx:613
 AliITSgeom.cxx:614
 AliITSgeom.cxx:615
 AliITSgeom.cxx:616
 AliITSgeom.cxx:617
 AliITSgeom.cxx:618
 AliITSgeom.cxx:619
 AliITSgeom.cxx:620
 AliITSgeom.cxx:621
 AliITSgeom.cxx:622
 AliITSgeom.cxx:623
 AliITSgeom.cxx:624
 AliITSgeom.cxx:625
 AliITSgeom.cxx:626
 AliITSgeom.cxx:627
 AliITSgeom.cxx:628
 AliITSgeom.cxx:629
 AliITSgeom.cxx:630
 AliITSgeom.cxx:631
 AliITSgeom.cxx:632
 AliITSgeom.cxx:633
 AliITSgeom.cxx:634
 AliITSgeom.cxx:635
 AliITSgeom.cxx:636
 AliITSgeom.cxx:637
 AliITSgeom.cxx:638
 AliITSgeom.cxx:639
 AliITSgeom.cxx:640
 AliITSgeom.cxx:641
 AliITSgeom.cxx:642
 AliITSgeom.cxx:643
 AliITSgeom.cxx:644
 AliITSgeom.cxx:645
 AliITSgeom.cxx:646
 AliITSgeom.cxx:647
 AliITSgeom.cxx:648
 AliITSgeom.cxx:649
 AliITSgeom.cxx:650
 AliITSgeom.cxx:651
 AliITSgeom.cxx:652
 AliITSgeom.cxx:653
 AliITSgeom.cxx:654
 AliITSgeom.cxx:655
 AliITSgeom.cxx:656
 AliITSgeom.cxx:657
 AliITSgeom.cxx:658
 AliITSgeom.cxx:659
 AliITSgeom.cxx:660
 AliITSgeom.cxx:661
 AliITSgeom.cxx:662
 AliITSgeom.cxx:663
 AliITSgeom.cxx:664
 AliITSgeom.cxx:665
 AliITSgeom.cxx:666
 AliITSgeom.cxx:667
 AliITSgeom.cxx:668
 AliITSgeom.cxx:669
 AliITSgeom.cxx:670
 AliITSgeom.cxx:671
 AliITSgeom.cxx:672
 AliITSgeom.cxx:673
 AliITSgeom.cxx:674
 AliITSgeom.cxx:675
 AliITSgeom.cxx:676
 AliITSgeom.cxx:677
 AliITSgeom.cxx:678
 AliITSgeom.cxx:679
 AliITSgeom.cxx:680
 AliITSgeom.cxx:681
 AliITSgeom.cxx:682
 AliITSgeom.cxx:683
 AliITSgeom.cxx:684
 AliITSgeom.cxx:685
 AliITSgeom.cxx:686
 AliITSgeom.cxx:687
 AliITSgeom.cxx:688
 AliITSgeom.cxx:689
 AliITSgeom.cxx:690
 AliITSgeom.cxx:691
 AliITSgeom.cxx:692
 AliITSgeom.cxx:693
 AliITSgeom.cxx:694
 AliITSgeom.cxx:695
 AliITSgeom.cxx:696
 AliITSgeom.cxx:697
 AliITSgeom.cxx:698
 AliITSgeom.cxx:699
 AliITSgeom.cxx:700
 AliITSgeom.cxx:701
 AliITSgeom.cxx:702
 AliITSgeom.cxx:703
 AliITSgeom.cxx:704
 AliITSgeom.cxx:705
 AliITSgeom.cxx:706
 AliITSgeom.cxx:707
 AliITSgeom.cxx:708