]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSGeoUtils.cxx
style modifications (Markus)
[u/mrichter/AliRoot.git] / PHOS / AliPHOSGeoUtils.cxx
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 /* $Id: AliPHOSGeometry.cxx 25590 2008-05-06 07:09:11Z prsnko $ */
17
18 //_________________________________________________________________________
19 // Geometry class  for PHOS 
20 // PHOS consists of the electromagnetic calorimeter (EMCA)
21 // and a charged particle veto (CPV)
22 // The EMCA/CPV modules are parametrized so that any configuration
23 // can be easily implemented 
24 // The title is used to identify the version of CPV used.
25 //                  
26 // -- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC "KI" & SUBATECH)
27
28 // --- ROOT system ---
29
30 #include "TClonesArray.h"
31 #include "TVector3.h"
32 #include "TParticle.h"
33 #include <TGeoManager.h>
34 #include <TGeoMatrix.h>
35
36 // --- Standard library ---
37
38 // --- AliRoot header files ---
39 #include "AliLog.h"
40 #include "AliPHOSEMCAGeometry.h"
41 #include "AliPHOSCPVGeometry.h"
42 #include "AliPHOSSupportGeometry.h"
43 #include "AliPHOSGeoUtils.h"
44
45 ClassImp(AliPHOSGeoUtils)
46
47 //____________________________________________________________________________
48 AliPHOSGeoUtils::AliPHOSGeoUtils():
49   fGeometryEMCA(0x0),fGeometryCPV(0x0),fGeometrySUPP(0x0),
50   fNModules(0),fNCristalsInModule(0),fNPhi(0),fNZ(0),
51   fNumberOfCPVPadsPhi(0),fNumberOfCPVPadsZ(0),
52   fNCellsXInStrip(0),fNCellsZInStrip(0),fNStripZ(0),
53   fCrystalShift(0.),fCryCellShift(0.),fCryStripShift(0.),fCellStep(0.),
54   fPadSizePhi(0.),fPadSizeZ(0.),fCPVBoxSizeY(0.),fMisalArray(0x0)
55  
56 {
57     // default ctor 
58     // must be kept public for root persistency purposes, but should never be called by the outside world
59 }  
60
61 //____________________________________________________________________________
62 AliPHOSGeoUtils::AliPHOSGeoUtils(const AliPHOSGeoUtils & rhs)
63   : TNamed(rhs),
64   fGeometryEMCA(0x0),fGeometryCPV(0x0),fGeometrySUPP(0x0),
65   fNModules(0),fNCristalsInModule(0),fNPhi(0),fNZ(0),
66   fNumberOfCPVPadsPhi(0),fNumberOfCPVPadsZ(0),
67   fNCellsXInStrip(0),fNCellsZInStrip(0),fNStripZ(0),
68   fCrystalShift(0.),fCryCellShift(0.),fCryStripShift(0.),fCellStep(0.),
69   fPadSizePhi(0.),fPadSizeZ(0.),fCPVBoxSizeY(0.),fMisalArray(0x0)
70 {
71   Fatal("cpy ctor", "not implemented") ; 
72 }
73
74 //____________________________________________________________________________
75 AliPHOSGeoUtils::AliPHOSGeoUtils(const Text_t* name, const Text_t* title) 
76     : TNamed(name, title),
77   fGeometryEMCA(0x0),fGeometryCPV(0x0),fGeometrySUPP(0x0),
78   fNModules(0),fNCristalsInModule(0),fNPhi(0),fNZ(0),
79   fNumberOfCPVPadsPhi(0),fNumberOfCPVPadsZ(0),
80   fNCellsXInStrip(0),fNCellsZInStrip(0),fNStripZ(0),
81   fCrystalShift(0.),fCryCellShift(0.),fCryStripShift(0.),fCellStep(0.),
82   fPadSizePhi(0.),fPadSizeZ(0.),fCPVBoxSizeY(0.),fMisalArray(0x0)
83
84   // ctor only for normal usage 
85
86   fGeometryEMCA = new AliPHOSEMCAGeometry() ;
87   fGeometryCPV  = new AliPHOSCPVGeometry() ;
88   fGeometrySUPP = new AliPHOSSupportGeometry() ;
89
90   fNModules     = 5;
91   fNPhi  = fGeometryEMCA->GetNPhi() ;
92   fNZ    = fGeometryEMCA->GetNZ() ;
93   fNCristalsInModule = fNPhi*fNZ ;
94   fNCellsXInStrip= fGeometryEMCA->GetNCellsXInStrip() ;
95   fNCellsZInStrip= fGeometryEMCA->GetNCellsZInStrip() ;
96   fNStripZ = fGeometryEMCA->GetNStripZ() ;
97   fXtlArrSize[0]=fGeometryEMCA->GetInnerThermoHalfSize()[0] ; //Wery close to the zise of the Xtl set
98   fXtlArrSize[1]=fGeometryEMCA->GetInnerThermoHalfSize()[1] ; //Wery close to the zise of the Xtl set
99   fXtlArrSize[2]=fGeometryEMCA->GetInnerThermoHalfSize()[2] ; //Wery close to the zise of the Xtl set
100
101   //calculate offset to crystal surface
102   const Float_t * inthermo = fGeometryEMCA->GetInnerThermoHalfSize() ;
103   const Float_t * strip    = fGeometryEMCA->GetStripHalfSize() ;
104   const Float_t * splate   = fGeometryEMCA->GetSupportPlateHalfSize();
105   const Float_t * crystal  = fGeometryEMCA->GetCrystalHalfSize() ;
106   const Float_t * pin      = fGeometryEMCA->GetAPDHalfSize() ;
107   const Float_t * preamp   = fGeometryEMCA->GetPreampHalfSize() ;
108   fCrystalShift=-inthermo[1]+strip[1]+splate[1]+crystal[1]-fGeometryEMCA->GetAirGapLed()/2.+pin[1]+preamp[1] ;
109   fCryCellShift=crystal[1]-(fGeometryEMCA->GetAirGapLed()-2*pin[1]-2*preamp[1])/2;
110   fCryStripShift=fCryCellShift+splate[1] ;
111   fCellStep = 2.*fGeometryEMCA->GetAirCellHalfSize()[0] ;
112
113   fNumberOfCPVPadsPhi = fGeometryCPV->GetNumberOfCPVPadsPhi() ;
114   fNumberOfCPVPadsZ   = fGeometryCPV->GetNumberOfCPVPadsZ() ;
115   fPadSizePhi = fGeometryCPV->GetCPVPadSizePhi() ;
116   fPadSizeZ   = fGeometryCPV->GetCPVPadSizeZ() ; 
117   fCPVBoxSizeY= fGeometryCPV->GetCPVBoxSize(1) ;
118
119   for(Int_t mod=0; mod<5; mod++){
120     fEMCMatrix[mod]=0 ;
121     for(Int_t istrip=0; istrip<224; istrip++)
122       fStripMatrix[mod][istrip]=0 ;
123     fCPVMatrix[mod]=0;
124     fPHOSMatrix[mod]=0 ;
125   }
126  
127 }
128
129 //____________________________________________________________________________
130 AliPHOSGeoUtils & AliPHOSGeoUtils::operator = (const AliPHOSGeoUtils  & /*rvalue*/) { 
131
132   Fatal("assignment operator", "not implemented") ; 
133     return *this ;
134 }
135
136 //____________________________________________________________________________
137 AliPHOSGeoUtils::~AliPHOSGeoUtils(void)
138 {
139   // dtor
140   if(fGeometryEMCA){
141     delete fGeometryEMCA; fGeometryEMCA = 0 ;
142   }
143   if(fGeometryCPV){
144     delete fGeometryCPV; fGeometryCPV=0 ;
145   }
146   if(fGeometrySUPP){
147     delete fGeometrySUPP ; fGeometrySUPP=0 ;
148   }
149   if(fMisalArray){
150     delete fMisalArray; fMisalArray=0 ;
151   }
152
153   for(Int_t mod=0; mod<5; mod++){
154     delete fEMCMatrix[mod] ;
155     for(Int_t istrip=0; istrip<224; istrip++)
156       delete fStripMatrix[mod][istrip];
157     delete fCPVMatrix[mod];
158     //    delete fPHOSMatrix[mod];
159   }
160 }
161 //____________________________________________________________________________
162 Bool_t AliPHOSGeoUtils::AbsToRelNumbering(Int_t absId, Int_t * relid) const
163 {
164   // Converts the absolute numbering into the following array
165   //  relid[0] = PHOS Module number 1:fNModules 
166   //  relid[1] = 0 if PbW04
167   //           = -1 if CPV
168   //  relid[2] = Row number inside a PHOS module
169   //  relid[3] = Column number inside a PHOS module
170
171   Float_t id = absId ;
172
173   Int_t phosmodulenumber = (Int_t)TMath:: Ceil( id / fNCristalsInModule ) ; 
174   
175   if ( phosmodulenumber >  fNModules ) { // it is a CPV pad
176     
177     id -=  fNPhi * fNZ *  fNModules ; 
178     Float_t nCPV  = fNumberOfCPVPadsPhi * fNumberOfCPVPadsZ ;
179     relid[0] = (Int_t) TMath::Ceil( id / nCPV ) ;
180     relid[1] = -1 ;
181     id -= ( relid[0] - 1 ) * nCPV ; 
182     relid[2] = (Int_t) TMath::Ceil( id / fNumberOfCPVPadsZ ) ;
183     relid[3] = (Int_t) ( id - ( relid[2] - 1 ) * fNumberOfCPVPadsZ ) ; 
184   } 
185   else { // it is a PW04 crystal
186
187     relid[0] = phosmodulenumber ;
188     relid[1] = 0 ;
189     id -= ( phosmodulenumber - 1 ) *  fNPhi * fNZ ; 
190     relid[2] = (Int_t)TMath::Ceil( id / fNZ )  ;
191     relid[3] = (Int_t)( id - ( relid[2] - 1 ) * fNZ ) ; 
192   } 
193   return kTRUE ; 
194 }
195 //____________________________________________________________________________
196 Bool_t AliPHOSGeoUtils::RelToAbsNumbering(const Int_t * relid, Int_t &  absId) const
197 {
198   // Converts the relative numbering into the absolute numbering
199   // EMCA crystals:
200   //  absId = from 1 to fNModules * fNPhi * fNZ
201   // CPV pad:
202   //  absId = from N(total PHOS crystals) + 1
203   //          to NCPVModules * fNumberOfCPVPadsPhi * fNumberOfCPVPadsZ
204
205   if ( relid[1] ==  0 ) {                            // it is a Phos crystal
206     absId =
207       ( relid[0] - 1 ) * fNPhi * fNZ         // the offset of PHOS modules
208       + ( relid[2] - 1 ) * fNZ                   // the offset along phi
209       +   relid[3] ;                                 // the offset along z
210   }
211   else { // it is a CPV pad
212     absId =    fNPhi * fNZ *  fNModules         // the offset to separate EMCA crystals from CPV pads
213       + ( relid[0] - 1 ) * fNumberOfCPVPadsPhi * fNumberOfCPVPadsZ   // the pads offset of PHOS modules 
214       + ( relid[2] - 1 ) * fNumberOfCPVPadsZ                            // the pads offset of a CPV row
215       +   relid[3] ;                                                         // the column number
216   }
217   
218   return kTRUE ; 
219 }
220
221 //____________________________________________________________________________
222 void AliPHOSGeoUtils::RelPosInModule(const Int_t * relid, Float_t & x, Float_t & z) const 
223 {
224   // Converts the relative numbering into the local PHOS-module (x, z) coordinates
225
226   if(relid[1]==0){ //this is PHOS
227
228     Double_t pos[3]= {0.0,-fCryCellShift,0.}; //Position incide the crystal 
229     Double_t posC[3]={0.0,0.0,0.}; //Global position
230
231     //Shift and possibly apply misalignment corrections
232     Int_t strip=1+((Int_t) TMath::Ceil((Double_t)relid[2]/fNCellsXInStrip))*fNStripZ-
233                 (Int_t) TMath::Ceil((Double_t)relid[3]/fNCellsZInStrip) ;
234     pos[0]=((relid[2]-1)%fNCellsXInStrip-fNCellsXInStrip/2+0.5)*fCellStep ;
235     pos[2]=(-(relid[3]-1)%fNCellsZInStrip+fNCellsZInStrip/2-0.5)*fCellStep ;
236
237     Int_t mod = relid[0] ;
238     const TGeoHMatrix * m2 = GetMatrixForStrip(mod, strip) ;
239     m2->LocalToMaster(pos,posC);
240
241     //Return to PHOS local system  
242     Double_t posL2[3]={posC[0],posC[1],posC[2]};
243     const TGeoHMatrix *mPHOS2 = GetMatrixForModule(mod) ;
244     mPHOS2->MasterToLocal(posC,posL2);
245     x=posL2[0] ;
246     z=-posL2[2];
247     return ;
248   }
249   else{//CPV
250     //first calculate position with respect to CPV plain 
251     Int_t row        = relid[2] ; //offset along x axis
252     Int_t column     = relid[3] ; //offset along z axis
253     Double_t pos[3]= {0.0,0.0,0.}; //Position incide the CPV printed circuit
254     Double_t posC[3]={0.0,0.0,0.}; //Global position
255     pos[0] = - ( fNumberOfCPVPadsPhi/2. - row    - 0.5 ) * fPadSizePhi  ; // position of pad  with respect
256     pos[2] = - ( fNumberOfCPVPadsZ  /2. - column - 0.5 ) * fPadSizeZ  ; // of center of PHOS module
257
258     //now apply possible shifts and rotations
259     const TGeoHMatrix *m = GetMatrixForCPV(relid[0]) ;
260     m->LocalToMaster(pos,posC);
261     //Return to PHOS local system
262     Double_t posL[3]={0.,0.,0.,} ;
263     const TGeoHMatrix *mPHOS = GetMatrixForPHOS(relid[0]) ;
264     mPHOS->MasterToLocal(posC,posL);
265     x=posL[0] ;
266     z=posL[1];
267     return ;
268  
269   }
270   
271 }
272 //____________________________________________________________________________
273 void AliPHOSGeoUtils::RelPosToAbsId(Int_t module, Double_t x, Double_t z, Int_t & absId) const
274 {
275   // converts local PHOS-module (x, z) coordinates to absId 
276
277   //Calculate AbsId using ideal geometry. Should be sufficient for primary particles calculation
278   //(the only place where this method used currently)
279   Int_t relid[4]={module,0,1,1} ;
280   relid[2] = static_cast<Int_t>(TMath::Ceil( x/ fCellStep + fNPhi / 2.) );
281   relid[3] = fNZ+1-static_cast<Int_t>(TMath::Ceil(-z/ fCellStep + fNZ   / 2.) ) ;
282   if(relid[2]<1)relid[2]=1 ;
283   if(relid[3]<1)relid[3]=1 ;
284   if(relid[2]>fNPhi)relid[2]=fNPhi ;
285   if(relid[3]>fNZ)relid[3]=fNZ ;
286   RelToAbsNumbering(relid,absId) ;
287
288 /*
289   //find Global position
290   if (!gGeoManager){
291     printf("Geo manager not initialized\n");
292     abort() ;
293   }
294   Double_t posL[3]={x,-fCrystalShift,-z} ; //Only for EMC!!!
295   Double_t posG[3] ;
296   char path[100] ;
297   sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1",module) ;
298   if (!gGeoManager->cd(path)){
299     printf("Geo manager can not find path \n");
300     abort() ;
301   }
302   TGeoHMatrix *mPHOS = gGeoManager->GetCurrentMatrix();
303   if (mPHOS){
304      mPHOS->LocalToMaster(posL,posG);
305   }
306   else{
307     printf("Geo matrixes are not loaded \n") ;
308     abort() ;
309   }
310
311   Int_t relid[4] ;
312   gGeoManager->FindNode(posG[0],posG[1],posG[2]) ;
313   //Check that path contains PSTR and extract strip number
314   TString cpath(gGeoManager->GetPath()) ;
315   Int_t indx = cpath.Index("PCEL") ;
316   if(indx==-1){ //for the few events when particle hits between srips use ideal geometry
317     relid[0] = module ;
318     relid[1] = 0 ;
319     relid[2] = static_cast<Int_t>(TMath::Ceil( x/ fCellStep + fNPhi / 2.) );
320     relid[3] = static_cast<Int_t>(TMath::Ceil(-z/ fCellStep + fNZ   / 2.) ) ;
321     if(relid[2]<1)relid[2]=1 ;
322     if(relid[3]<1)relid[3]=1 ;
323     if(relid[2]>fNPhi)relid[2]=fNPhi ;
324     if(relid[3]>fNZ)relid[3]=fNZ ;
325     RelToAbsNumbering(relid,absId) ;
326   }
327   else{
328     Int_t indx2 = cpath.Index("/",indx) ;
329     if(indx2==-1)
330       indx2=cpath.Length() ;
331     TString cell=cpath(indx+5,indx2-indx-5) ;
332     Int_t icell=cell.Atoi() ;
333     indx = cpath.Index("PSTR") ;
334     indx2 = cpath.Index("/",indx) ;
335     TString strip=cpath(indx+5,indx2-indx-5) ;
336     Int_t iStrip = strip.Atoi() ; 
337
338     Int_t row = fNStripZ - (iStrip - 1) % (fNStripZ) ;
339     Int_t col = (Int_t) TMath::Ceil((Double_t) iStrip/(fNStripZ)) -1 ;
340  
341     // Absid for 8x2-strips. Looks nice :)
342     absId = (module-1)*fNCristalsInModule +
343                   row * 2 + (col*fNCellsXInStrip + (icell - 1) / 2)*fNZ - (icell & 1 ? 1 : 0);
344  
345   }
346 */
347  
348 }
349
350 //____________________________________________________________________________
351 void AliPHOSGeoUtils::RelPosToRelId(Int_t module, Double_t x, Double_t z, Int_t * relId) const
352 {
353   //Evaluates RelId of the crystall with given coordinates
354
355   Int_t absId ;
356   RelPosToAbsId(module, x,z,absId) ;
357   AbsToRelNumbering(absId,relId) ;
358 }
359
360 //____________________________________________________________________________
361 void AliPHOSGeoUtils::RelPosInAlice(Int_t id, TVector3 & pos ) const
362 {
363   // Converts the absolute numbering into the global ALICE coordinate system
364   
365   if (!gGeoManager){
366     printf("Geo manager not initialized\n");
367     abort();
368   }
369     
370   Int_t relid[4] ;
371     
372   AbsToRelNumbering(id , relid) ;
373     
374   //construct module name
375   if(relid[1]==0){ //this is EMC
376  
377     Double_t ps[3]= {0.0,-fCryStripShift,0.}; //Position incide the crystal
378     Double_t psC[3]={0.0,0.0,0.}; //Global position
379  
380     //Shift and possibly apply misalignment corrections
381     Int_t strip=1+((Int_t) TMath::Ceil((Double_t)relid[2]/fNCellsXInStrip))*fNStripZ-
382                 (Int_t) TMath::Ceil((Double_t)relid[3]/fNCellsZInStrip) ;
383     ps[0]=((relid[2]-1)%fNCellsXInStrip-fNCellsXInStrip/2+0.5)*fCellStep ;
384     ps[2]=(-(relid[3]-1)%fNCellsZInStrip+fNCellsZInStrip/2-0.5)*fCellStep ;
385  
386     Int_t mod = relid[0] ;
387     const TGeoHMatrix * m2 = GetMatrixForStrip(mod, strip) ;
388     m2->LocalToMaster(ps,psC);
389     pos.SetXYZ(psC[0],psC[1],psC[2]) ; 
390  
391   }
392   else{
393     //first calculate position with respect to CPV plain
394     Int_t row        = relid[2] ; //offset along x axis
395     Int_t column     = relid[3] ; //offset along z axis
396     Double_t ps[3]= {0.0,fCPVBoxSizeY/2.,0.}; //Position on top of CPV
397     Double_t psC[3]={0.0,0.0,0.}; //Global position
398     pos[0] = - ( fNumberOfCPVPadsPhi/2. - row    - 0.5 ) * fPadSizePhi  ; // position of pad  with respect
399     pos[2] = - ( fNumberOfCPVPadsZ  /2. - column - 0.5 ) * fPadSizeZ  ; // of center of PHOS module
400  
401     //now apply possible shifts and rotations
402     const TGeoHMatrix *m = GetMatrixForCPV(relid[0]) ;
403     m->LocalToMaster(ps,psC);
404     pos.SetXYZ(psC[0],psC[1],-psC[2]) ; 
405   }
406
407
408 //____________________________________________________________________________
409 void AliPHOSGeoUtils::Local2Global(Int_t mod, Float_t x, Float_t z,
410                                    TVector3& globalPosition) const 
411 {
412   Double_t posL[3]={x,-fCrystalShift,-z} ; //Only for EMC!!!
413   Double_t posG[3] ;
414   const TGeoHMatrix *mPHOS = GetMatrixForModule(mod) ;
415   mPHOS->LocalToMaster(posL,posG);
416   globalPosition.SetXYZ(posG[0],posG[1],posG[2]) ;
417 }
418 //____________________________________________________________________________
419 void AliPHOSGeoUtils::Global2Local(TVector3& localPosition,
420                                    const TVector3& globalPosition,
421                                    Int_t module) const
422 {
423   // Transforms a global position to the local coordinate system
424   // of the module 
425   //Return to PHOS local system
426   Double_t posG[3]={globalPosition.X(),globalPosition.Y(),globalPosition.Z()} ;
427   Double_t posL[3]={0.,0.,0.} ;
428   const TGeoHMatrix *mPHOS = GetMatrixForModule(module) ;
429   if(mPHOS){
430     mPHOS->MasterToLocal(posG,posL);
431     localPosition.SetXYZ(posL[0],posL[1]+fCrystalShift,-posL[2]) ;  
432   }
433   else{
434     localPosition.SetXYZ(999.,999.,999.) ; //module does not exist in given configuration
435   }
436  
437 }
438 //____________________________________________________________________________
439 Bool_t AliPHOSGeoUtils::GlobalPos2RelId(TVector3 & global, Int_t * relId){
440   //Converts position in global ALICE coordinates to relId 
441   //returns false if x,z coordinates are beyond PHOS
442   //distande to PHOS surface is NOT calculated 
443   TVector3 loc ;
444   for(Int_t mod=1; mod<=fNModules; mod++){
445     Global2Local(loc,global,mod) ;
446     //If in Acceptance
447     if((TMath::Abs(loc.Z())<fXtlArrSize[2]) && (TMath::Abs(loc.X())<fXtlArrSize[0])){
448        RelPosToRelId(mod,loc.X(),loc.Z(),relId);
449        return kTRUE ;
450     }
451   }
452   return kFALSE ; 
453
454 }
455 //____________________________________________________________________________
456 Bool_t AliPHOSGeoUtils::ImpactOnEmc(const TParticle * particle,
457        Int_t & moduleNumber, Double_t & z, Double_t & x) const
458 {
459   // Tells if a particle enters PHOS and evaluates hit position
460   Double_t vtx[3]={particle->Vx(),particle->Vy(),particle->Vz()} ;
461   return ImpactOnEmc(vtx,particle->Theta(),particle->Phi(),moduleNumber,z,x);
462 }
463  
464 //____________________________________________________________________________
465 Bool_t AliPHOSGeoUtils::ImpactOnEmc(const Double_t * vtx, Double_t theta, Double_t phi, 
466                                   Int_t & moduleNumber, Double_t & z, Double_t & x) const
467 {
468   // calculates the impact coordinates on PHOS of a neutral particle  
469   // emitted in the vertex vtx[3] with direction vec(p) in the ALICE global coordinate system
470   TVector3 p(TMath::Sin(theta)*TMath::Cos(phi),TMath::Sin(theta)*TMath::Sin(phi),TMath::Cos(theta)) ;
471   return ImpactOnEmc(vtx,p,moduleNumber,z,x) ;
472
473 }
474 //____________________________________________________________________________
475 Bool_t AliPHOSGeoUtils::ImpactOnEmc(const Double_t * vtx, const TVector3 &p,
476                                   Int_t & moduleNumber, Double_t & z, Double_t & x) const
477 {
478   // calculates the impact coordinates on PHOS of a neutral particle  
479   // emitted in the vertex vtx[3] with direction theta and phi in the ALICE global coordinate system
480   TVector3 v(vtx[0],vtx[1],vtx[2]) ;
481
482   for(Int_t imod=1; imod<=fNModules ; imod++){
483     //create vector from (0,0,0) to center of crystal surface of imod module
484     Double_t tmp[3]={0.,-fCrystalShift,0.} ;
485  
486     const TGeoHMatrix *m = GetMatrixForModule(imod) ;
487     if(!m) //module does not exist in given configuration
488       continue ; 
489     Double_t posG[3]={0.,0.,0.} ;
490     m->LocalToMaster(tmp,posG);
491     TVector3 n(posG[0],posG[1],posG[2]) ; 
492     Double_t direction=n.Dot(p) ;
493     if(direction<=0.)
494       continue ; //momentum directed FROM module
495     Double_t fr = (n.Mag2()-n.Dot(v))/direction ;  
496     //Calculate direction in module plain
497     n-=v+fr*p ;
498     n*=-1. ;
499     if(TMath::Abs(TMath::Abs(n.Z())<fXtlArrSize[2]) && n.Pt()<fXtlArrSize[0]){
500       moduleNumber = imod ;
501       z=n.Z() ;
502       x=TMath::Sign(n.Pt(),n.X()) ;
503       //no need to return to local system since we calcilated distance from module center
504       //and tilts can not be significant.
505       return kTRUE ;
506     }
507   }
508   //Not in acceptance
509   x=0; z=0 ;
510   moduleNumber=0 ;
511   return kFALSE ;
512
513 }
514 //____________________________________________________________________________
515 void AliPHOSGeoUtils::GetIncidentVector(const TVector3 &vtx, Int_t module, Float_t x,Float_t z, TVector3 &vInc) const {
516   //Calculates vector pointing from vertex to current poisition in module local frame
517   //Note that PHOS local system and ALICE global have opposite z directions
518
519   Global2Local(vInc,vtx,module) ; 
520   vInc.SetXYZ(vInc.X()+x,vInc.Y(),vInc.Z()+z) ;
521 }
522 //____________________________________________________________________________
523 const TGeoHMatrix * AliPHOSGeoUtils::GetMatrixForModule(Int_t mod)const {
524   //Provides shift-rotation matrix for module mod
525
526   //If GeoManager exists, take matrixes from it
527   if(gGeoManager){
528     char path[255] ;
529     sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1",mod) ;
530     //    sprintf(path,"/ALIC_1/PHOS_%d",relid[0]) ;
531     if (!gGeoManager->cd(path)){
532       AliWarning(Form("Geo manager can not find path %s \n",path));
533       return 0;
534     }
535     return gGeoManager->GetCurrentMatrix();
536   }
537   if(fEMCMatrix[mod-1]){
538     return fEMCMatrix[mod-1] ;
539   }
540   else{
541     AliWarning("Can not find PHOS misalignment matrixes\n") ;
542     AliWarning("Either import TGeoManager from geometry.root or \n");
543     AliWarning("read stored matrixes from AliESD Header: \n") ;
544     AliWarning("AliPHOSGeoUtils::SetMisalMatrixes(header->GetPHOSMisalMatrix()) \n") ; 
545     return 0 ;
546   }
547   return 0 ;
548 }
549 //____________________________________________________________________________
550 const TGeoHMatrix * AliPHOSGeoUtils::GetMatrixForStrip(Int_t mod, Int_t strip)const {
551   //Provides shift-rotation matrix for strip unit of the module mod
552
553   //If GeoManager exists, take matrixes from it
554   if(gGeoManager){
555     char path[255] ;
556     sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1/PSTR_%d",mod,strip) ;
557     if (!gGeoManager->cd(path)){
558       AliWarning(Form("Geo manager can not find path %s \n",path));
559       return 0 ;
560     }
561     return gGeoManager->GetCurrentMatrix();
562   } 
563   if(fStripMatrix[mod-1][strip-1]){
564     return fStripMatrix[mod-1][strip-1] ;
565   }
566   else{
567     AliWarning("Can not find PHOS misalignment matrixes\n") ;
568     AliWarning("Either import TGeoManager from geometry.root or \n");
569     AliWarning("read stored matrixes from AliESD Header: \n") ; 
570     AliWarning("AliPHOSGeoUtils::SetMisalMatrixes(header->GetPHOSMisalMatrix()) \n") ;
571     return 0 ;
572   } 
573   return 0 ;
574 }
575 //____________________________________________________________________________
576 const TGeoHMatrix * AliPHOSGeoUtils::GetMatrixForCPV(Int_t mod)const {
577   //Provides shift-rotation matrix for CPV of the module mod
578
579   //If GeoManager exists, take matrixes from it
580   if(gGeoManager){ 
581     char path[255] ;
582     //now apply possible shifts and rotations
583     sprintf(path,"/ALIC_1/PHOS_%d/PCPV_1",mod) ;
584     if (!gGeoManager->cd(path)){
585       AliWarning(Form("Geo manager can not find path %s \n",path));
586       return 0 ;
587     }
588     return gGeoManager->GetCurrentMatrix();
589   }
590   if(fCPVMatrix[mod-1]){
591     return fCPVMatrix[mod-1] ;
592   }
593   else{
594     AliWarning("Can not find PHOS misalignment matrixes\n") ;
595     AliWarning("Either import TGeoManager from geometry.root or \n");
596     AliWarning("read stored matrixes from AliESD Header: \n") ;  
597     AliWarning("AliPHOSGeoUtils::SetMisalMatrixes(header->GetPHOSMisalMatrix()) \n") ;
598     return 0 ;
599   }
600   return 0 ;
601
602 //____________________________________________________________________________
603 const TGeoHMatrix * AliPHOSGeoUtils::GetMatrixForPHOS(Int_t mod)const {
604   //Provides shift-rotation matrix for PHOS (EMC+CPV) 
605
606   //If GeoManager exists, take matrixes from it
607   if(gGeoManager){
608     char path[255] ;
609     sprintf(path,"/ALIC_1/PHOS_%d",mod) ;
610     if (!gGeoManager->cd(path)){
611       AliWarning(Form("Geo manager can not find path %s \n",path));
612       return 0 ;
613     }
614     return gGeoManager->GetCurrentMatrix();
615   }
616   if(fPHOSMatrix[mod-1]){
617     return fPHOSMatrix[mod-1] ;
618   }
619   else{
620     AliWarning("Can not find PHOS misalignment matrixes\n") ;
621     AliWarning("Either import TGeoManager from geometry.root or \n");
622     AliWarning("read stored matrixes from AliESD Header:  \n") ;   
623     AliWarning("AliPHOSGeoUtils::SetMisalMatrixes(header->GetPHOSMisalMatrix()) \n") ;
624     return 0 ;
625   }
626   return 0 ;
627 }
628 //____________________________________________________________________________
629 void AliPHOSGeoUtils::SetMisalMatrix(const TGeoHMatrix * m, Int_t mod){
630   //Fills pointers to geo matrixes
631  
632   fPHOSMatrix[mod]=m ;
633
634   //If module does not exist, make sure all its matrices are zero
635   if(m==NULL){
636     fEMCMatrix[mod]=NULL ;
637     Int_t istrip=0 ;
638     for(Int_t irow = 0; irow < fGeometryEMCA->GetNStripX(); irow ++){
639       for(Int_t icol = 0; icol < fGeometryEMCA->GetNStripZ(); icol ++){
640         fStripMatrix[mod][istrip]=NULL ;
641       }
642     } 
643     fCPVMatrix[mod]=NULL ;
644     return ;
645   }
646
647   //Calculate maxtrices for PTII
648   if(!fMisalArray)
649     fMisalArray = new TClonesArray("TGeoHMatrix",1120+10) ;
650   Int_t nr = fMisalArray->GetEntriesFast() ;
651   Double_t rotEMC[9]={1.,0.,0.,0.,0.,-1.,0.,1.,0.} ;
652   const Float_t * inthermo = fGeometryEMCA->GetInnerThermoHalfSize() ;
653   const Float_t * strip    = fGeometryEMCA->GetStripHalfSize() ;
654   const Float_t * covparams = fGeometryEMCA->GetAlCoverParams() ;
655   const Float_t * warmcov = fGeometryEMCA->GetWarmAlCoverHalfSize() ;
656   Float_t z = fGeometryCPV->GetCPVBoxSize(1) / 2. - warmcov[2] + covparams[3]-inthermo[1] ;
657   Double_t locTII[3]={0.,0.,z} ; 
658   Double_t globTII[3] ;
659
660   if (fEMCMatrix[mod] == NULL)
661     fEMCMatrix[mod] = new((*fMisalArray)[nr])TGeoHMatrix() ;
662   nr++ ;
663   fEMCMatrix[mod]->SetRotation(rotEMC) ;
664   fEMCMatrix[mod]->MultiplyLeft(fPHOSMatrix[mod]) ;
665   fPHOSMatrix[mod]->LocalToMaster(locTII,globTII) ;
666   fEMCMatrix[mod]->SetTranslation(globTII) ;
667  
668   //Now calculate ideal matrixes for strip misalignment.
669   //For the moment we can not store them in ESDHeader
670
671   Double_t loc[3]={0.,inthermo[1] - strip[1],0.} ; 
672   Double_t glob[3] ;
673
674   Int_t istrip=0 ;
675   for(Int_t irow = 0; irow < fGeometryEMCA->GetNStripX(); irow ++){
676     loc[0] = (2*irow + 1 - fGeometryEMCA->GetNStripX())* strip[0] ;
677     for(Int_t icol = 0; icol < fGeometryEMCA->GetNStripZ(); icol ++){
678       loc[2] = (2*icol + 1 - fGeometryEMCA->GetNStripZ()) * strip[2] ;
679       fEMCMatrix[mod]->LocalToMaster(loc,glob) ;
680       if (fStripMatrix[mod][istrip] == NULL)
681         fStripMatrix[mod][istrip] = new((*fMisalArray)[nr])TGeoHMatrix(*(fEMCMatrix[mod])) ; //Use same rotation as PHOS module
682       nr++ ;
683       fStripMatrix[mod][istrip]->SetTranslation(glob) ;
684       istrip++;
685     }
686   }
687  
688   //Now calculate CPV matrixes
689   const Float_t * emcParams = fGeometryEMCA->GetEMCParams() ;
690   Double_t globCPV[3] ;
691   Double_t locCPV[3]={0.,0.,- emcParams[3]} ;
692   Double_t rot[9]={1.,0.,0.,0.,0.,1.,0.,-1.,0.} ;
693
694   if (fCPVMatrix[mod] == NULL)
695     fCPVMatrix[mod] = new((*fMisalArray)[nr])TGeoHMatrix() ;
696   nr++ ;
697   fCPVMatrix[mod]->SetRotation(rot) ;
698   fCPVMatrix[mod]->MultiplyLeft(fPHOSMatrix[mod]) ;
699   fCPVMatrix[mod]->ReflectY(kFALSE) ;
700   fPHOSMatrix[mod]->LocalToMaster(locCPV,globCPV) ;
701   fCPVMatrix[mod]->SetTranslation(globCPV) ;
702
703 }
704