1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
19 // Geometry class for PHOS : singleton
20 // PHOS consists of the electromagnetic calorimeter (EMCA)
21 // and a charged particle veto either in the Subatech's version (PPSD)
22 // or in the IHEP's one (CPV).
23 // The EMCA/PPSD/CPV modules are parametrized so that any configuration
24 // can be easily implemented
25 // The title is used to identify the version of CPV used.
27 // -- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC "KI" & SUBATECH)
29 // --- ROOT system ---
32 #include "TRotation.h"
33 #include "TParticle.h"
34 #include <TGeoManager.h>
36 // --- Standard library ---
38 // --- AliRoot header files ---
40 #include "AliPHOSGeometry.h"
41 #include "AliPHOSEMCAGeometry.h"
42 #include "AliPHOSRecPoint.h"
44 ClassImp(AliPHOSGeometry)
46 // these initialisations are needed for a singleton
47 AliPHOSGeometry * AliPHOSGeometry::fgGeom = 0 ;
48 Bool_t AliPHOSGeometry::fgInit = kFALSE ;
50 //____________________________________________________________________________
51 AliPHOSGeometry::AliPHOSGeometry() :
55 fIPtoUpperCPVsurface(0),
62 // must be kept public for root persistency purposes, but should never be called by the outside world
66 //____________________________________________________________________________
67 AliPHOSGeometry::AliPHOSGeometry(const AliPHOSGeometry & rhs)
69 fNModules(rhs.fNModules),
72 fIPtoUpperCPVsurface(rhs.fIPtoUpperCPVsurface),
78 Fatal("cpy ctor", "not implemented") ;
81 //____________________________________________________________________________
82 AliPHOSGeometry::AliPHOSGeometry(const Text_t* name, const Text_t* title)
83 : AliGeometry(name, title),
87 fIPtoUpperCPVsurface(0),
93 // ctor only for internal usage (singleton)
98 //____________________________________________________________________________
99 AliPHOSGeometry::~AliPHOSGeometry(void)
103 if (fRotMatrixArray) fRotMatrixArray->Delete() ;
104 if (fRotMatrixArray) delete fRotMatrixArray ;
105 if (fPHOSAngle ) delete[] fPHOSAngle ;
108 //____________________________________________________________________________
109 void AliPHOSGeometry::Init(void)
111 // Initializes the PHOS parameters :
112 // IHEP is the Protvino CPV (cathode pad chambers)
114 TString test(GetName()) ;
115 if (test != "IHEP" && test != "noCPV") {
116 AliFatal(Form("%s is not a known geometry (choose among IHEP)",
125 fGeometryEMCA = new AliPHOSEMCAGeometry();
127 fGeometryCPV = new AliPHOSCPVGeometry ();
129 fGeometrySUPP = new AliPHOSSupportGeometry();
131 fPHOSAngle = new Float_t[fNModules] ;
133 Float_t * emcParams = fGeometryEMCA->GetEMCParams() ;
135 fPHOSParams[0] = TMath::Max((Double_t)fGeometryCPV->GetCPVBoxSize(0)/2.,
136 (Double_t)(emcParams[0] - (emcParams[1]-emcParams[0])*
137 fGeometryCPV->GetCPVBoxSize(1)/2/emcParams[3]));
138 fPHOSParams[1] = emcParams[1] ;
139 fPHOSParams[2] = TMath::Max((Double_t)emcParams[2], (Double_t)fGeometryCPV->GetCPVBoxSize(2)/2.);
140 fPHOSParams[3] = emcParams[3] + fGeometryCPV->GetCPVBoxSize(1)/2. ;
142 fIPtoUpperCPVsurface = fGeometryEMCA->GetIPtoOuterCoverDistance() - fGeometryCPV->GetCPVBoxSize(1) ;
144 //calculate offset to crystal surface
145 Float_t * inthermo = fGeometryEMCA->GetInnerThermoHalfSize() ;
146 Float_t * strip = fGeometryEMCA->GetStripHalfSize() ;
147 Float_t* splate = fGeometryEMCA->GetSupportPlateHalfSize();
148 Float_t * crystal = fGeometryEMCA->GetCrystalHalfSize() ;
149 Float_t * pin = fGeometryEMCA->GetAPDHalfSize() ;
150 Float_t * preamp = fGeometryEMCA->GetPreampHalfSize() ;
151 fCrystalShift=-inthermo[1]+strip[1]+splate[1]+crystal[1]-fGeometryEMCA->GetAirGapLed()/2.+pin[1]+preamp[1] ;
152 fCryCellShift=crystal[1]-(fGeometryEMCA->GetAirGapLed()-2*pin[1]-2*preamp[1])/2;
155 for ( index = 0; index < fNModules; index++ )
156 fPHOSAngle[index] = 0.0 ; // Module position angles are set in CreateGeometry()
158 fRotMatrixArray = new TObjArray(fNModules) ;
160 // Geometry parameters are calculated
163 Double_t const kRADDEG = 180.0 / TMath::Pi() ;
164 Float_t r = GetIPtoOuterCoverDistance() + fPHOSParams[3] - GetCPVBoxSize(1) ;
165 for (Int_t iModule=0; iModule<fNModules; iModule++) {
166 fModuleCenter[iModule][0] = r * TMath::Sin(fPHOSAngle[iModule] / kRADDEG );
167 fModuleCenter[iModule][1] =-r * TMath::Cos(fPHOSAngle[iModule] / kRADDEG );
168 fModuleCenter[iModule][2] = 0.;
170 fModuleAngle[iModule][0][0] = 90;
171 fModuleAngle[iModule][0][1] = fPHOSAngle[iModule];
172 fModuleAngle[iModule][1][0] = 0;
173 fModuleAngle[iModule][1][1] = 0;
174 fModuleAngle[iModule][2][0] = 90;
175 fModuleAngle[iModule][2][1] = 270 + fPHOSAngle[iModule];
180 //____________________________________________________________________________
181 AliPHOSGeometry * AliPHOSGeometry::GetInstance()
183 // Returns the pointer of the unique instance; singleton specific
185 return static_cast<AliPHOSGeometry *>( fgGeom ) ;
188 //____________________________________________________________________________
189 AliPHOSGeometry * AliPHOSGeometry::GetInstance(const Text_t* name, const Text_t* title)
191 // Returns the pointer of the unique instance
192 // Creates it with the specified options (name, title) if it does not exist yet
194 AliPHOSGeometry * rv = 0 ;
196 if ( strcmp(name,"") == 0 )
199 fgGeom = new AliPHOSGeometry(name, title) ;
201 rv = (AliPHOSGeometry * ) fgGeom ;
210 if ( strcmp(fgGeom->GetName(), name) != 0 )
211 ::Error("GetInstance", "Current geometry is %s. You cannot call %s",
212 fgGeom->GetName(), name) ;
214 rv = (AliPHOSGeometry *) fgGeom ;
219 //____________________________________________________________________________
220 void AliPHOSGeometry::SetPHOSAngles()
222 // Calculates the position of the PHOS modules in ALICE global coordinate system
225 Double_t const kRADDEG = 180.0 / TMath::Pi() ;
226 Float_t pphi = 2 * TMath::ATan( GetOuterBoxSize(0) / ( 2.0 * GetIPtoUpperCPVsurface() ) ) ;
229 AliError(Form("PHOS modules overlap!\n pphi = %f fAngle = %f",
235 for( Int_t i = 1; i <= fNModules ; i++ ) {
236 Float_t angle = pphi * ( i - fNModules / 2.0 - 0.5 ) ;
237 fPHOSAngle[i-1] = - angle ;
241 //____________________________________________________________________________
242 Bool_t AliPHOSGeometry::AbsToRelNumbering(Int_t absId, Int_t * relid) const
244 // Converts the absolute numbering into the following array
245 // relid[0] = PHOS Module number 1:fNModules
246 // relid[1] = 0 if PbW04
248 // relid[2] = Row number inside a PHOS module
249 // relid[3] = Column number inside a PHOS module
254 Int_t phosmodulenumber = (Int_t)TMath:: Ceil( id / GetNCristalsInModule() ) ;
256 if ( phosmodulenumber > GetNModules() ) { // it is a CPV pad
258 id -= GetNPhi() * GetNZ() * GetNModules() ;
259 Float_t nCPV = GetNumberOfCPVPadsPhi() * GetNumberOfCPVPadsZ() ;
260 relid[0] = (Int_t) TMath::Ceil( id / nCPV ) ;
262 id -= ( relid[0] - 1 ) * nCPV ;
263 relid[2] = (Int_t) TMath::Ceil( id / GetNumberOfCPVPadsZ() ) ;
264 relid[3] = (Int_t) ( id - ( relid[2] - 1 ) * GetNumberOfCPVPadsZ() ) ;
266 else { // it is a PW04 crystal
268 relid[0] = phosmodulenumber ;
270 id -= ( phosmodulenumber - 1 ) * GetNPhi() * GetNZ() ;
271 relid[2] = (Int_t)TMath::Ceil( id / GetNZ() ) ;
272 relid[3] = (Int_t)( id - ( relid[2] - 1 ) * GetNZ() ) ;
276 //____________________________________________________________________________
277 void AliPHOSGeometry::GetGlobal(const AliRecPoint* recPoint, TVector3 & gpos) const
279 AliFatal(Form("Please use GetGlobalPHOS(recPoint,gpos) instead of GetGlobal!"));
282 //____________________________________________________________________________
283 void AliPHOSGeometry::GetGlobalPHOS(const AliPHOSRecPoint* recPoint, TVector3 & gpos) const
285 // Calculates the coordinates of a RecPoint and the error matrix in the ALICE global coordinate system
287 const AliPHOSRecPoint * tmpPHOS = recPoint ;
288 TVector3 localposition ;
290 tmpPHOS->GetLocalPosition(gpos) ;
293 AliFatal("Geo manager not initialized\n");
295 //construct module name
298 if(tmpPHOS->IsEmc()){
299 sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1",tmpPHOS->GetPHOSMod()) ;
300 // sprintf(path,"/ALIC_1/PHOS_%d",tmpPHOS->GetPHOSMod()) ;
304 sprintf(path,"/ALIC_1/PHOS_%d/PCPV_1",tmpPHOS->GetPHOSMod()) ;
305 dy= GetCPVBoxSize(1)/2. ; //center of CPV module
307 Double_t pos[3]={gpos.X(),gpos.Y()-dy,gpos.Z()} ;
309 pos[2]=-pos[2] ; //Opposite z directions in EMC matrix and local frame!!!
311 //now apply possible shifts and rotations
312 if (!gGeoManager->cd(path)){
313 AliFatal("Geo manager can not find path \n");
315 TGeoHMatrix *m = gGeoManager->GetCurrentMatrix();
317 m->LocalToMaster(pos,posC);
320 AliFatal("Geo matrixes are not loaded \n") ;
322 gpos.SetXYZ(posC[0],posC[1],posC[2]) ;
325 //____________________________________________________________________________
326 void AliPHOSGeometry::ImpactOnEmc(Double_t * vtx, Double_t theta, Double_t phi,
327 Int_t & moduleNumber, Double_t & z, Double_t & x) const
329 // calculates the impact coordinates on PHOS of a neutral particle
330 // emitted in the vertex vtx[3] with direction theta and phi in the ALICE global coordinate system
331 TVector3 p(TMath::Sin(theta)*TMath::Cos(phi),TMath::Sin(theta)*TMath::Sin(phi),TMath::Cos(theta)) ;
332 TVector3 v(vtx[0],vtx[1],vtx[2]) ;
335 AliFatal("Geo manager not initialized\n");
338 for(Int_t imod=1; imod<=GetNModules() ; imod++){
339 //create vector from (0,0,0) to center of crystal surface of imod module
340 Double_t tmp[3]={0.,-fCrystalShift,0.} ;
343 sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1",imod) ;
344 if (!gGeoManager->cd(path)){
345 AliFatal("Geo manager can not find path \n");
347 TGeoHMatrix *m = gGeoManager->GetCurrentMatrix();
348 Double_t posG[3]={0.,0.,0.} ;
349 if (m) m->LocalToMaster(tmp,posG);
350 TVector3 n(posG[0],posG[1],-posG[2]) ;
351 Double_t direction=n.Dot(p) ;
353 continue ; //momentum directed FROM module
354 Double_t fr = (n.Mag2()-n.Dot(v))/direction ;
355 //Calculate direction in module plain
358 Float_t * sz = fGeometryEMCA->GetInnerThermoHalfSize() ; //Wery close to the zise of the Xtl set
359 if(TMath::Abs(TMath::Abs(n.Z())<sz[2]) && n.Pt()<sz[0]){
360 moduleNumber = imod ;
362 x=TMath::Sign(n.Pt(),n.X()) ;
372 //____________________________________________________________________________
373 Bool_t AliPHOSGeometry::Impact(const TParticle * particle) const
375 // Tells if a particle enters PHOS
377 Int_t moduleNumber=0;
378 Double_t vtx[3]={particle->Vx(),particle->Vy(),particle->Vz()} ;
380 ImpactOnEmc(vtx,particle->Theta(),particle->Phi(),moduleNumber,z,x);
386 //____________________________________________________________________________
387 Bool_t AliPHOSGeometry::RelToAbsNumbering(const Int_t * relid, Int_t & absId) const
389 // Converts the relative numbering into the absolute numbering
391 // absId = from 1 to fNModules * fNPhi * fNZ
393 // absId = from N(total PHOS crystals) + 1
394 // to NCPVModules * fNumberOfCPVPadsPhi * fNumberOfCPVPadsZ
398 if ( relid[1] == 0 ) { // it is a Phos crystal
400 ( relid[0] - 1 ) * GetNPhi() * GetNZ() // the offset of PHOS modules
401 + ( relid[2] - 1 ) * GetNZ() // the offset along phi
402 + relid[3] ; // the offset along z
404 else { // it is a CPV pad
405 absId = GetNPhi() * GetNZ() * GetNModules() // the offset to separate EMCA crystals from CPV pads
406 + ( relid[0] - 1 ) * GetNumberOfCPVPadsPhi() * GetNumberOfCPVPadsZ() // the pads offset of PHOS modules
407 + ( relid[2] - 1 ) * GetNumberOfCPVPadsZ() // the pads offset of a CPV row
408 + relid[3] ; // the column number
414 //____________________________________________________________________________
415 void AliPHOSGeometry::RelPosInAlice(Int_t id, TVector3 & pos ) const
417 // Converts the absolute numbering into the global ALICE coordinate system
420 AliFatal("Geo manager not initialized\n");
425 AbsToRelNumbering(id , relid) ;
427 //construct module name
429 if(relid[1]==0){ //this is EMC
431 Double_t ps[3]= {0.0,-fCryCellShift,0.}; //Position incide the crystal
432 Double_t psC[3]={0.0,0.0,0.}; //Global position
434 //Shift and possibly apply misalignment corrections
435 Int_t nCellsXInStrip=fGeometryEMCA->GetNCellsXInStrip() ;
436 Int_t nCellsZInStrip=fGeometryEMCA->GetNCellsZInStrip() ;
437 Int_t strip=1+((Int_t) TMath::Ceil((Double_t)relid[2]/nCellsXInStrip))*fGeometryEMCA->GetNStripZ()-
438 (Int_t) TMath::Ceil((Double_t)relid[3]/nCellsZInStrip) ;
439 Int_t cellraw= relid[3]%nCellsZInStrip ;
440 if(cellraw==0)cellraw=nCellsZInStrip ;
441 Int_t cell= ((relid[2]-1)%nCellsXInStrip)*nCellsZInStrip + cellraw ;
442 sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1/PSTR_%d/PCEL_%d",
443 relid[0],strip,cell) ;
444 if (!gGeoManager->cd(path)){
445 AliFatal("Geo manager can not find path \n");
447 TGeoHMatrix *m = gGeoManager->GetCurrentMatrix();
448 if (m) m->LocalToMaster(ps,psC);
450 AliFatal("Geo matrixes are not loaded \n") ;
452 pos.SetXYZ(psC[0],psC[1],psC[2]) ;
455 //first calculate position with respect to CPV plain
456 Int_t row = relid[2] ; //offset along x axis
457 Int_t column = relid[3] ; //offset along z axis
458 Double_t ps[3]= {0.0,GetCPVBoxSize(1)/2.,0.}; //Position on top of CPV
459 Double_t psC[3]={0.0,0.0,0.}; //Global position
460 pos[0] = - ( GetNumberOfCPVPadsPhi()/2. - row - 0.5 ) * GetPadSizePhi() ; // position of pad with respect
461 pos[2] = - ( GetNumberOfCPVPadsZ() /2. - column - 0.5 ) * GetPadSizeZ() ; // of center of PHOS module
463 //now apply possible shifts and rotations
464 sprintf(path,"/ALIC_1/PHOS_%d/PCPV_1",relid[0]) ;
465 if (!gGeoManager->cd(path)){
466 AliFatal("Geo manager can not find path \n");
468 TGeoHMatrix *m = gGeoManager->GetCurrentMatrix();
469 if (m) m->LocalToMaster(ps,psC);
471 AliFatal("Geo matrixes are not loaded \n") ;
473 pos.SetXYZ(psC[0],psC[1],-psC[2]) ;
477 //____________________________________________________________________________
478 void AliPHOSGeometry::RelPosToAbsId(Int_t module, Double_t x, Double_t z, Int_t & absId) const
480 // converts local PHOS-module (x, z) coordinates to absId
482 //find Global position
484 AliFatal("Geo manager not initialized\n");
486 Double_t posL[3]={x,-fCrystalShift,-z} ; //Only for EMC!!!
489 sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1",module) ;
490 if (!gGeoManager->cd(path)){
491 AliFatal("Geo manager can not find path \n");
493 TGeoHMatrix *mPHOS = gGeoManager->GetCurrentMatrix();
495 mPHOS->LocalToMaster(posL,posG);
498 AliFatal("Geo matrixes are not loaded \n") ;
502 gGeoManager->FindNode(posG[0],posG[1],posG[2]) ;
503 //Check that path contains PSTR and extract strip number
504 TString cpath(gGeoManager->GetPath()) ;
505 Int_t indx = cpath.Index("PCEL") ;
506 if(indx==-1){ //for the few events when particle hits between srips use ideal geometry
509 relid[2] = static_cast<Int_t>(TMath::Ceil( x/ GetCellStep() + GetNPhi() / 2.) );
510 relid[3] = static_cast<Int_t>(TMath::Ceil(-z/ GetCellStep() + GetNZ() / 2.) ) ;
511 if(relid[2]<1)relid[2]=1 ;
512 if(relid[3]<1)relid[3]=1 ;
513 if(relid[2]>GetNPhi())relid[2]=GetNPhi() ;
514 if(relid[3]>GetNZ())relid[3]=GetNZ() ;
515 RelToAbsNumbering(relid,absId) ;
518 Int_t indx2 = cpath.Index("/",indx) ;
520 indx2=cpath.Length() ;
521 TString cell=cpath(indx+5,indx2-indx-5) ;
522 Int_t icell=cell.Atoi() ;
523 indx = cpath.Index("PSTR") ;
524 indx2 = cpath.Index("/",indx) ;
525 TString strip=cpath(indx+5,indx2-indx-5) ;
526 Int_t iStrip = strip.Atoi() ;
528 Int_t row = fGeometryEMCA->GetNStripZ() - (iStrip - 1) % (fGeometryEMCA->GetNStripZ()) ;
529 Int_t col = (Int_t) TMath::Ceil((Double_t) iStrip/(fGeometryEMCA->GetNStripZ())) -1 ;
531 // Absid for 8x2-strips. Looks nice :)
532 absId = (module-1)*GetNCristalsInModule() +
533 row * 2 + (col*fGeometryEMCA->GetNCellsXInStrip() + (icell - 1) / 2)*GetNZ() - (icell & 1 ? 1 : 0);
539 //____________________________________________________________________________
540 void AliPHOSGeometry::RelPosInModule(const Int_t * relid, Float_t & x, Float_t & z) const
542 // Converts the relative numbering into the local PHOS-module (x, z) coordinates
543 // Note: sign of z differs from that in the previous version (Yu.Kharlov, 12 Oct 2000)
547 AliFatal("Geo manager not initialized\n");
549 //construct module name
551 if(relid[1]==0){ //this is PHOS
553 // Calculations using ideal geometry (obsolete)
554 // x = - ( GetNPhi()/2. - relid[2] + 0.5 ) * GetCellStep() ; // position of Xtal with respect
555 // z = - ( GetNZ() /2. - relid[3] + 0.5 ) * GetCellStep() ; // of center of PHOS module
557 Double_t pos[3]= {0.0,-fCryCellShift,0.}; //Position incide the crystal
558 Double_t posC[3]={0.0,0.0,0.}; //Global position
560 //Shift and possibly apply misalignment corrections
561 Int_t nCellsXInStrip=fGeometryEMCA->GetNCellsXInStrip() ;
562 Int_t nCellsZInStrip=fGeometryEMCA->GetNCellsZInStrip() ;
563 // Int_t strip=1+(relid[3]-1)/fGeometryEMCA->GetNCellsZInStrip()+((relid[2]-1)/nCellsInStrip)*fGeometryEMCA->GetNStripZ() ;
564 Int_t strip=1+((Int_t) TMath::Ceil((Double_t)relid[2]/nCellsXInStrip))*fGeometryEMCA->GetNStripZ()-
565 (Int_t) TMath::Ceil((Double_t)relid[3]/nCellsZInStrip) ;
566 Int_t cellraw= relid[3]%nCellsZInStrip ;
567 if(cellraw==0)cellraw=nCellsZInStrip ;
568 Int_t cell= ((relid[2]-1)%nCellsXInStrip)*nCellsZInStrip + cellraw ;
569 sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1/PSTR_%d/PCEL_%d",
570 relid[0],strip,cell) ;
571 if (!gGeoManager->cd(path)){
572 AliFatal("Geo manager can not find path \n");
574 TGeoHMatrix *m = gGeoManager->GetCurrentMatrix();
575 if (m) m->LocalToMaster(pos,posC);
577 AliFatal("Geo matrixes are not loaded \n") ;
579 //Return to PHOS local system
580 Double_t posL[3]={posC[0],posC[1],posC[2]};
581 sprintf(path,"/ALIC_1/PHOS_%d",relid[0]) ;
582 if (!gGeoManager->cd(path)){
583 AliFatal("Geo manager can not find path \n");
585 TGeoHMatrix *mPHOS = gGeoManager->GetCurrentMatrix();
586 if (mPHOS) mPHOS->MasterToLocal(posC,posL);
588 AliFatal("Geo matrixes are not loaded \n") ;
590 //printf("RelPosInMod: posL=[%f,%f,%f]\n",posL[0],posL[1],posL[2]) ;
591 //printf("old: x=%f, z=%f \n",x,z);
597 //first calculate position with respect to CPV plain
598 Int_t row = relid[2] ; //offset along x axis
599 Int_t column = relid[3] ; //offset along z axis
600 Double_t pos[3]= {0.0,0.0,0.}; //Position incide the CPV printed circuit
601 Double_t posC[3]={0.0,0.0,0.}; //Global position
602 pos[0] = - ( GetNumberOfCPVPadsPhi()/2. - row - 0.5 ) * GetPadSizePhi() ; // position of pad with respect
603 pos[2] = - ( GetNumberOfCPVPadsZ() /2. - column - 0.5 ) * GetPadSizeZ() ; // of center of PHOS module
605 //now apply possible shifts and rotations
606 sprintf(path,"/ALIC_1/PHOS_%d/PCPV_1",relid[0]) ;
607 if (!gGeoManager->cd(path)){
608 AliFatal("Geo manager can not find path \n");
610 TGeoHMatrix *m = gGeoManager->GetCurrentMatrix();
611 if (m) m->LocalToMaster(pos,posC);
613 AliFatal("Geo matrixes are not loaded \n") ;
615 //Return to PHOS local system
616 Double_t posL[3]={0.,0.,0.,} ;
617 sprintf(path,"/ALIC_1/PHOS_%d",relid[0]) ;
618 if (!gGeoManager->cd(path)){
619 AliFatal("Geo manager can not find path \n");
621 TGeoHMatrix *mPHOS = gGeoManager->GetCurrentMatrix();
622 if (mPHOS) mPHOS->MasterToLocal(posC,posL);
624 AliFatal("Geo matrixes are not loaded \n") ;
634 //____________________________________________________________________________
636 void AliPHOSGeometry::GetModuleCenter(TVector3& center,
640 // Returns a position of the center of the CPV or EMC module
641 // in ideal (not misaligned) geometry
643 if (strcmp(det,"CPV") == 0) rDet = GetIPtoCPVDistance ();
644 else if (strcmp(det,"EMC") == 0) rDet = GetIPtoCrystalSurface();
646 AliFatal(Form("Wrong detector name %s",det));
648 Float_t angle = GetPHOSAngle(module); // (40,20,0,-20,-40) degrees
649 angle *= TMath::Pi()/180;
650 angle += 3*TMath::Pi()/2.;
651 center.SetXYZ(rDet*TMath::Cos(angle), rDet*TMath::Sin(angle), 0.);
654 //____________________________________________________________________________
656 void AliPHOSGeometry::Global2Local(TVector3& localPosition,
657 const TVector3& globalPosition,
660 // Transforms a global position of the rec.point to the local coordinate system
661 //Return to PHOS local system
662 Double_t posG[3]={globalPosition.X(),globalPosition.Y(),globalPosition.Z()} ;
663 Double_t posL[3]={0.,0.,0.} ;
665 sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1",module) ;
666 // sprintf(path,"/ALIC_1/PHOS_%d",module) ;
667 if (!gGeoManager->cd(path)){
668 AliFatal("Geo manager can not find path \n");
670 TGeoHMatrix *mPHOS = gGeoManager->GetCurrentMatrix();
671 if (mPHOS) mPHOS->MasterToLocal(posG,posL);
673 AliFatal("Geo matrixes are not loaded \n") ;
675 localPosition.SetXYZ(posL[0],posL[1]+fCrystalShift,-posL[2]) ;
678 Float_t angle = GetPHOSAngle(module); // (40,20,0,-20,-40) degrees
679 angle *= TMath::Pi()/180;
680 angle += 3*TMath::Pi()/2.;
681 localPosition = globalPosition;
682 localPosition.RotateZ(-angle);
685 //____________________________________________________________________________
686 void AliPHOSGeometry::Local2Global(Int_t mod, Float_t x, Float_t z,
687 TVector3& globalPosition) const
690 sprintf(path,"/ALIC_1/PHOS_%d/PEMC_1/PCOL_1/PTIO_1/PCOR_1/PAGA_1/PTII_1",mod) ;
691 // sprintf(path,"/ALIC_1/PHOS_%d",mod) ;
692 if (!gGeoManager->cd(path)){
693 AliFatal("Geo manager can not find path \n");
695 Double_t posL[3]={x,-fCrystalShift,-z} ; //Only for EMC!!!
697 TGeoHMatrix *mPHOS = gGeoManager->GetCurrentMatrix();
699 mPHOS->LocalToMaster(posL,posG);
702 AliFatal("Geo matrixes are not loaded \n") ;
704 globalPosition.SetXYZ(posG[0],posG[1],posG[2]) ;
706 //____________________________________________________________________________
707 void AliPHOSGeometry::GetIncidentVector(const TVector3 &vtx, Int_t module, Float_t x,Float_t z, TVector3 &vInc) const {
708 //Calculates vector pointing from vertex to current poisition in module local frame
709 //Note that PHOS local system and ALICE global have opposite z directions
711 Global2Local(vInc,vtx,module) ;
712 vInc.SetXYZ(vInc.X()+x,vInc.Y(),vInc.Z()+z) ;