]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALEMCGeometry.cxx
cout replaced by AliDebug (Plamen)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALEMCGeometry.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: AliEMCALEMCGeometry.cxx 29514 2008-10-26 10:24:38Z hristov $*/
17
18 //_________________________________________________________________________
19 // Geometry class  for EMCAL : singleton  
20 // EMCAL consists of layers of scintillator and lead
21 // with scintillator fiber arranged as "shish-kebab" skewers 
22 // Places the the Barrel Geometry of The EMCAL at Midrapidity
23 // between 80 and 180(or 190) degrees of Phi and
24 // -0.7 to 0.7 in eta 
25 //
26 //     EMCAL geometry tree:
27 //     EMCAL -> superModule -> module -> tower(cell)
28 //     Indexes
29 //     absId -> nSupMod     -> nModule -> (nIphi,nIeta)
30 //
31 //   Name choices: 
32 //   EMCAL_PDC06 (geometry used for PDC06 simulations, kept for backward compatibility)
33 //      = equivalent to SHISH_77_TRD1_2X2_FINAL_110DEG in old notation
34 //   EMCAL_COMPLETE (geometry for expected complete detector)
35 //      = equivalent to SHISH_77_TRD1_2X2_FINAL_110DEG scTh=0.176 pbTh=0.144
36 //          in old notation
37 //   EMCAL_WSUC (Wayne State test stand)
38 //      = no definite equivalent in old notation, was only used by
39 //          Aleksei, but kept for testing purposes
40 //
41 //   etc.
42 //
43 //
44 //
45 //*-- Author: Sahal Yacoob (LBL / UCT)
46 //     and  : Yves Schutz (SUBATECH)
47 //     and  : Jennifer Klay (LBL)
48 //     and  : Aleksei Pavlinov (WSU) 
49 //     and  : Magali Estienne (SUBATECH)
50
51 // --- Root header files ---
52 #include <TObjArray.h>
53 #include <TObjString.h>
54 #include <TRegexp.h>
55
56 // -- ALICE Headers.
57 #include "AliLog.h"
58
59 // --- EMCAL headers
60 #include "AliEMCALEMCGeometry.h"
61 #include <cassert>
62
63 ClassImp(AliEMCALEMCGeometry)
64
65 // these initialisations are needed for a singleton
66 Bool_t    AliEMCALEMCGeometry::fgInit      = kFALSE;
67 const Char_t*   AliEMCALEMCGeometry::fgkDefaultGeometryName = "EMCAL_COMPLETE";
68
69
70 AliEMCALEMCGeometry::AliEMCALEMCGeometry() 
71   : TNamed(),
72     fGeoName(0),fArrayOpts(0),fNAdditionalOpts(0),fECPbRadThickness(0.),fECScintThick(0.),
73     fNECLayers(0),fArm1PhiMin(0.),fArm1PhiMax(0.),fArm1EtaMin(0.),fArm1EtaMax(0.),fIPDistance(0.),
74     fShellThickness(0.),fZLength(0.),fNZ(0),fNPhi(0),fSampling(0.),fNumberOfSuperModules(0),
75     fFrontSteelStrip(0.),fLateralSteelStrip(0.),fPassiveScintThick(0.),fPhiModuleSize(0.),
76     fEtaModuleSize(0.),fPhiTileSize(0.),fEtaTileSize(0.),fLongModuleSize(0.),fNPhiSuperModule(0),
77     fNPHIdiv(0),fNETAdiv(0), fNCells(0),fNCellsInSupMod(0),fNCellsInModule(0),
78     // Trigger staff
79     fNTRUEta(0), fNTRUPhi(0), fNModulesInTRUEta(0), fNModulesInTRUPhi(0), fNEtaSubOfTRU(0),
80     // 
81     fTrd1Angle(0.),f2Trd1Dx2(0.),
82     fPhiGapForSM(0.),fKey110DEG(0),fPhiBoundariesOfSM(0), fPhiCentersOfSM(0),fEtaMaxOfTRD1(0),
83     fCentersOfCellsEtaDir(0), fCentersOfCellsXDir(0),fCentersOfCellsPhiDir(0),
84     fEtaCentersOfCells(0),fPhiCentersOfCells(0),fShishKebabTrd1Modules(0),
85     fILOSS(-1), fIHADR(-1),
86     //obsolete member data
87     fAlFrontThick(0.), fGap2Active(0.), fSteelFrontThick(0.), fTrd2AngleY(0.),
88     f2Trd2Dy2(0.), fEmptySpace(0.), fTubsR(0.), fTubsTurnAngle(0.)
89
90   // default ctor only for internal usage (singleton)
91   // must be kept public for root persistency purposes, 
92   // but should never be called by the outside world    
93
94   AliDebug(2, "AliEMCALEMCGeometry : default ctor ");
95 }
96 //______________________________________________________________________
97 AliEMCALEMCGeometry::AliEMCALEMCGeometry(const Text_t* name, const Text_t* title) :
98   TNamed(name,title),
99     fGeoName(0),fArrayOpts(0),fNAdditionalOpts(0),fECPbRadThickness(0.),fECScintThick(0.),
100     fNECLayers(0),fArm1PhiMin(0.),fArm1PhiMax(0.),fArm1EtaMin(0.),fArm1EtaMax(0.),fIPDistance(0.),
101     fShellThickness(0.),fZLength(0.),fNZ(0),fNPhi(0),fSampling(0.),fNumberOfSuperModules(0),
102     fFrontSteelStrip(0.),fLateralSteelStrip(0.),fPassiveScintThick(0.),fPhiModuleSize(0.),
103     fEtaModuleSize(0.),fPhiTileSize(0.),fEtaTileSize(0.),fLongModuleSize(0.),fNPhiSuperModule(0),
104     fNPHIdiv(0),fNETAdiv(0), fNCells(0),fNCellsInSupMod(0),fNCellsInModule(0),
105     // Trigger staff
106     fNTRUEta(0), fNTRUPhi(0), fNModulesInTRUEta(0), fNModulesInTRUPhi(0), fNEtaSubOfTRU(0),
107     // 
108     fTrd1Angle(0.),f2Trd1Dx2(0.),
109     fPhiGapForSM(0.),fKey110DEG(0),fPhiBoundariesOfSM(0), fPhiCentersOfSM(0), fEtaMaxOfTRD1(0),
110     fCentersOfCellsEtaDir(0),fCentersOfCellsXDir(0),fCentersOfCellsPhiDir(0),
111     fEtaCentersOfCells(0),fPhiCentersOfCells(0),fShishKebabTrd1Modules(0),
112     fILOSS(-1), fIHADR(-1), 
113     //obsolete member data
114     fAlFrontThick(0.), fGap2Active(0.), fSteelFrontThick(0.), fTrd2AngleY(0.),
115     f2Trd2Dy2(0.), fEmptySpace(0.), fTubsR(0.), fTubsTurnAngle(0.)
116 {
117   // ctor only for internal usage (singleton)
118   AliDebug(2, Form("AliEMCALEMCGeometry(%s,%s) ", name,title));
119
120   Init();
121
122   //  CreateListOfTrd1Modules();
123
124   if (AliDebugLevel()>=2) {
125     PrintGeometry();
126   }
127
128 }
129 //______________________________________________________________________
130 AliEMCALEMCGeometry::AliEMCALEMCGeometry(const AliEMCALEMCGeometry& geom)
131   : TNamed(geom),
132     fGeoName(geom.fGeoName),
133     fArrayOpts(geom.fArrayOpts),
134     fNAdditionalOpts(geom.fNAdditionalOpts),
135     fECPbRadThickness(geom.fECPbRadThickness),
136     fECScintThick(geom.fECScintThick),
137     fNECLayers(geom.fNECLayers),
138     fArm1PhiMin(geom.fArm1PhiMin),
139     fArm1PhiMax(geom.fArm1PhiMax),
140     fArm1EtaMin(geom.fArm1EtaMin),
141     fArm1EtaMax(geom.fArm1EtaMax),
142     fIPDistance(geom.fIPDistance),
143     fShellThickness(geom.fShellThickness),
144     fZLength(geom.fZLength),
145     fNZ(geom.fNZ),
146     fNPhi(geom.fNPhi),
147     fSampling(geom.fSampling),
148     fNumberOfSuperModules(geom.fNumberOfSuperModules),
149     fFrontSteelStrip(geom.fFrontSteelStrip),
150     fLateralSteelStrip(geom.fLateralSteelStrip),
151     fPassiveScintThick(geom.fPassiveScintThick),
152     fPhiModuleSize(geom.fPhiModuleSize),
153     fEtaModuleSize(geom.fEtaModuleSize),
154     fPhiTileSize(geom.fPhiTileSize),
155     fEtaTileSize(geom.fEtaTileSize),
156     fLongModuleSize(geom.fLongModuleSize),
157     fNPhiSuperModule(geom.fNPhiSuperModule),
158     fNPHIdiv(geom.fNPHIdiv),
159     fNETAdiv(geom.fNETAdiv),
160     fNCells(geom.fNCells),
161     fNCellsInSupMod(geom.fNCellsInSupMod),
162     fNCellsInModule(geom.fNCellsInModule),
163     // Trigger staff
164     fNTRUEta(geom.fNTRUEta),
165     fNTRUPhi(geom.fNTRUPhi),
166     fNModulesInTRUEta(geom.fNModulesInTRUEta),
167     fNModulesInTRUPhi(geom.fNModulesInTRUPhi),
168     fNEtaSubOfTRU(geom.fNEtaSubOfTRU),
169     //
170     fTrd1Angle(geom.fTrd1Angle),
171     f2Trd1Dx2(geom.f2Trd1Dx2),
172     fPhiGapForSM(geom.fPhiGapForSM),
173     fKey110DEG(geom.fKey110DEG),
174     fPhiBoundariesOfSM(geom.fPhiBoundariesOfSM),
175     fPhiCentersOfSM(geom.fPhiCentersOfSM),
176     fEtaMaxOfTRD1(geom.fEtaMaxOfTRD1),
177     fCentersOfCellsEtaDir(geom.fCentersOfCellsEtaDir),
178     fCentersOfCellsXDir(geom.fCentersOfCellsXDir),
179     fCentersOfCellsPhiDir(geom.fCentersOfCellsPhiDir),
180     fEtaCentersOfCells(geom.fEtaCentersOfCells),
181     fPhiCentersOfCells(geom.fPhiCentersOfCells),
182     fShishKebabTrd1Modules(geom.fShishKebabTrd1Modules),
183     fILOSS(geom.fILOSS), fIHADR(geom.fIHADR),
184     //obsolete member data
185     fAlFrontThick(geom.fAlFrontThick),
186     fGap2Active(geom.fGap2Active),
187     fSteelFrontThick(geom.fSteelFrontThick),
188     fTrd2AngleY(geom.fTrd2AngleY),
189     f2Trd2Dy2(geom.f2Trd2Dy2),
190     fEmptySpace(geom.fEmptySpace),
191     fTubsR(geom.fTubsR),
192     fTubsTurnAngle(geom.fTubsTurnAngle)
193 {
194   //copy ctor
195 }
196
197 //______________________________________________________________________
198 AliEMCALEMCGeometry::~AliEMCALEMCGeometry(void){
199     // dtor
200 }
201
202 //______________________________________________________________________
203 void AliEMCALEMCGeometry::Init(void){
204   //
205   // Initializes the EMCAL parameters based on the name
206   // Only Shashlyk geometry is available, but various combinations of
207   // layers and number of supermodules can be selected with additional
208   // options or geometry name
209   //
210
211   fkAdditionalOpts[0] = "nl=";       // number of sampling layers (fNECLayers)
212   fkAdditionalOpts[1] = "pbTh=";     // cm, Thickness of the Pb   (fECPbRadThick)
213   fkAdditionalOpts[2] = "scTh=";     // cm, Thickness of the Sc    (fECScintThick)
214   fkAdditionalOpts[3] = "latSS=";    // cm, Thickness of lateral steel strip (fLateralSteelStrip)
215   fkAdditionalOpts[4] = "allILOSS="; // = 0,1,2,3,4 (4 - energy loss without fluctuation)
216   fkAdditionalOpts[5] = "allIHADR="; // = 0,1,2 (0 - no hadronic interaction)
217
218   fNAdditionalOpts = sizeof(fkAdditionalOpts) / sizeof(char*);
219
220   // geometry
221   fgInit = kFALSE; // Assume failed until proven otherwise.
222   fGeoName   = GetName();
223   fGeoName.ToUpper();
224
225   //Convert old geometry names to new ones
226   if(fGeoName.Contains("SHISH_77_TRD1_2X2_FINAL_110DEG")) {
227     if(fGeoName.Contains("PBTH=0.144") && fGeoName.Contains("SCTH=0.176")) {
228       fGeoName = "EMCAL_COMPLETE";
229     } else {
230       fGeoName = "EMCAL_PDC06";
231     }
232   }
233   if(fGeoName.Contains("WSUC")) fGeoName = "EMCAL_WSUC";
234
235   //check that we have a valid geometry name
236   if(!(fGeoName.Contains("EMCAL_PDC06") || fGeoName.Contains("EMCAL_COMPLETE") || fGeoName.Contains("EMCAL_WSUC") || fGeoName.Contains("EMCAL_1stYear"))) {
237     Fatal("Init", "%s is an undefined geometry!", fGeoName.Data()) ; 
238   }
239
240   // Option to know whether we have the "half" supermodule(s) or not
241   fKey110DEG = 0;
242   if(fGeoName.Contains("COMPLETE") || fGeoName.Contains("PDC06")) fKey110DEG = 1; // for GetAbsCellId
243   fShishKebabTrd1Modules = 0;
244
245   // JLK 13-Apr-2008
246   //default parameters are those of EMCAL_COMPLETE geometry
247   //all others render variations from these at the end of
248   //geometry-name specific options
249
250   fNumberOfSuperModules = 12;       // 12 = 6 * 2 (6 in phi, 2 in Z)
251   fNPhi                 = 12;       // module granularity in phi within smod (azimuth)
252   fNZ                   = 24;       // module granularity along Z within smod (eta)
253   fNPHIdiv = fNETAdiv   = 2;        // tower granularity within module
254   fArm1PhiMin           = 80.0;     // degrees, Starting EMCAL Phi position
255   fArm1PhiMax           = 200.0;    // degrees, Ending EMCAL Phi position
256   fArm1EtaMin           = -0.7;     // pseudorapidity, Starting EMCAL Eta position
257   fArm1EtaMax           = +0.7;     // pseudorapidity, Ending EMCAL Eta position
258   fIPDistance           = 428.0;    // cm, radial distance to front face from nominal vertex point
259   fPhiGapForSM          = 2.;       // cm, only for final TRD1 geometry
260   fFrontSteelStrip      = 0.025;    // 0.025cm = 0.25mm  (13-may-05 from V.Petrov)
261   fPassiveScintThick    = 0.8;      // 0.8cm   = 8mm     (13-may-05 from V.Petrov)
262   fLateralSteelStrip    = 0.01;     // 0.01cm  = 0.1mm   (13-may-05 from V.Petrov) - was 0.025
263   fTrd1Angle            = 1.5;      // in degrees       
264
265   fSampling             = 1.;       // should be calculated with call to DefineSamplingFraction()
266   fNECLayers            = 77;       // (13-may-05 from V.Petrov) - can be changed with additional options
267   fECScintThick         = 0.176;    // scintillator layer thickness
268   fECPbRadThickness     = 0.144;    // lead layer thickness
269
270   fPhiModuleSize = 12.26 - fPhiGapForSM / Float_t(fNPhi); // first assumption
271   fEtaModuleSize = fPhiModuleSize;
272
273   fZLength              = 700.;     // Z coverage (cm)
274
275
276   //needs to be called for each geometry and before setting geometry
277   //parameters which can depend on the outcome
278   CheckAdditionalOptions();
279
280   //modifications to the above for PDC06 geometry
281   if(fGeoName.Contains("PDC06")){ // 18-may-05 - about common structure
282     fECScintThick  = fECPbRadThickness = 0.16;// (13-may-05 from V.Petrov)    
283     CheckAdditionalOptions();
284   }
285
286   //modifications to the above for WSUC geometry
287   if(fGeoName.Contains("WSUC")){ // 18-may-05 - about common structure
288     fPhiModuleSize = 12.5;     // 20-may-05 - rectangular shape
289     fEtaModuleSize = 11.9;
290     fECScintThick  = fECPbRadThickness = 0.16;// (13-may-05 from V.Petrov)
291     fNumberOfSuperModules = 1; // 27-may-05
292     fShellThickness = 30.;       // should be change 
293     fNPhi = fNZ = 4; 
294     CheckAdditionalOptions();
295   }
296
297   if(fGeoName.Contains("1stYear")){     
298         fNumberOfSuperModules = 2;      
299          
300         if(fGeoName.Contains("LowerEta")) {
301                 fNPhiSuperModule = 1;           
302         }
303         else if(fGeoName.Contains("LowerPhi_SideA")){
304         fNPhiSuperModule = 2;   
305         fArm1EtaMax=0;          
306         }
307         else if(fGeoName.Contains("LowerPhi_SideC")){
308         fNPhiSuperModule = 2;           
309         fArm1EtaMin=0;  
310         }
311                 
312       CheckAdditionalOptions(); 
313   }
314
315   // constant for transition absid <--> indexes
316   fNCellsInModule  = fNPHIdiv*fNETAdiv;
317   fNCellsInSupMod = fNCellsInModule*fNPhi*fNZ;
318   fNCells         = fNCellsInSupMod*fNumberOfSuperModules;
319   if(GetKey110DEG()) fNCells -= fNCellsInSupMod;
320
321   fNPhiSuperModule = fNumberOfSuperModules/2;
322   if(fNPhiSuperModule < 1) fNPhiSuperModule = 1;
323     
324   fPhiTileSize = fPhiModuleSize/double(fNPHIdiv) - fLateralSteelStrip; // 13-may-05 
325   fEtaTileSize = fEtaModuleSize/double(fNETAdiv) - fLateralSteelStrip; // 13-may-05 
326
327   fLongModuleSize = fNECLayers*(fECScintThick + fECPbRadThickness);  
328   f2Trd1Dx2 = fEtaModuleSize + 2.*fLongModuleSize*TMath::Tan(fTrd1Angle*TMath::DegToRad()/2.);
329   if(!fGeoName.Contains("WSUC")) fShellThickness  = TMath::Sqrt(fLongModuleSize*fLongModuleSize + f2Trd1Dx2*f2Trd1Dx2);
330
331   //These parameters are used to create the mother volume to hold the supermodules
332   //2cm padding added to allow for misalignments - JLK 30-May-2008
333   fEnvelop[0]     = fIPDistance - 1.; // mother volume inner radius
334   fEnvelop[1]     = fIPDistance + fShellThickness + 1.; // mother volume outer r.
335   fEnvelop[2]     = fZLength + 2.; //mother volume length 
336
337   // Local coordinates
338   fParSM[0] = GetShellThickness()/2.;        
339   fParSM[1] = GetPhiModuleSize() * GetNPhi()/2.;
340   fParSM[2] = fZLength/4.;  //divide by 4 to get half-length of SM
341
342   // SM phi boundaries - (0,1),(2,3) .. (10,11) - has the same boundaries; Nov 7, 2006 
343   fPhiBoundariesOfSM.Set(fNumberOfSuperModules);
344   fPhiCentersOfSM.Set(fNumberOfSuperModules/2);
345   fPhiBoundariesOfSM[0] = TMath::PiOver2() - TMath::ATan2(fParSM[1] , fIPDistance); // 1th and 2th modules)
346   fPhiCentersOfSM[0]     = TMath::PiOver2();
347   if(fNumberOfSuperModules > 1) 
348     fPhiBoundariesOfSM[1] = TMath::PiOver2() + TMath::ATan2(fParSM[1] , fIPDistance);
349   if(fNumberOfSuperModules > 2) {
350     for(int i=1; i<=4; i++) { // from 2th ro 9th
351       fPhiBoundariesOfSM[2*i]   = fPhiBoundariesOfSM[0] + 20.*TMath::DegToRad()*i;
352       fPhiBoundariesOfSM[2*i+1] = fPhiBoundariesOfSM[1] + 20.*TMath::DegToRad()*i;
353       fPhiCentersOfSM[i]         = fPhiCentersOfSM[0]     + 20.*TMath::DegToRad()*i;
354     }
355   }
356   if(fNumberOfSuperModules > 10) {
357     fPhiBoundariesOfSM[11] = 190.*TMath::DegToRad();
358     fPhiBoundariesOfSM[10] = fPhiBoundariesOfSM[11] - TMath::ATan2((fParSM[1]) , fIPDistance);
359     fPhiCentersOfSM[5]      = (fPhiBoundariesOfSM[10]+fPhiBoundariesOfSM[11])/2.; 
360   }
361
362   //called after setting of scintillator and lead layer parameters
363   DefineSamplingFraction();
364
365   
366   // TRU parameters - Apr 29,08 by PAI. 
367   // These parameters values was updated at Nov 05, 2007
368   // As is on Olivier  BOURRION (LPSC) ppt preasentation 
369   // at ALICE trigger meeting at 13th-14th March
370   fNTRUEta = 1;           // was 3
371   fNTRUPhi = 3;           // was 1
372   fNModulesInTRUEta = 24; // was 8
373   fNModulesInTRUPhi = 4;  // was 12
374   // Jet trigger 
375   // 3*6*10 + 2*6*2 = 204 -> matrix (nphi(17), neta(12))
376   fNEtaSubOfTRU     = 6;  
377   
378
379   fgInit = kTRUE; 
380 }
381
382 //___________________________________________________________________
383 void AliEMCALEMCGeometry::PrintGeometry()
384 {
385   // Separate routine is callable from broswer; Nov 7,2006
386   printf("\nInit: geometry of EMCAL named %s :\n", fGeoName.Data());
387   if(fArrayOpts) {
388     for(Int_t i=0; i<fArrayOpts->GetEntries(); i++){
389       TObjString *o = (TObjString*)fArrayOpts->At(i);
390       printf(" %i : %s \n", i, o->String().Data());
391     }
392   }
393   printf("Granularity: %d in eta and %d in phi\n", GetNZ(), GetNPhi()) ;
394   printf("Layout: phi = (%7.1f, %7.1f), eta = (%5.2f, %5.2f), IP = %7.2f -> for EMCAL envelope only\n",  
395            GetArm1PhiMin(), GetArm1PhiMax(),GetArm1EtaMin(), GetArm1EtaMax(), GetIPDistance() );
396
397   printf( "               ECAL      : %d x (%f cm Pb, %f cm Sc) \n", 
398   GetNECLayers(), GetECPbRadThick(), GetECScintThick() ) ; 
399   printf("                fSampling %5.2f \n",  fSampling );
400   printf(" fIPDistance       %6.3f cm \n", fIPDistance);
401   printf(" fNPhi %i   |  fNZ %i \n", fNPhi, fNZ);
402   printf(" fNCellsInModule %i : fNCellsInSupMod %i : fNCells %i\n",fNCellsInModule, fNCellsInSupMod, fNCells);
403   printf(" X:Y module size     %6.3f , %6.3f cm \n", fPhiModuleSize, fEtaModuleSize);
404   printf(" X:Y   tile size     %6.3f , %6.3f cm \n", fPhiTileSize, fEtaTileSize);
405   printf(" #of sampling layers %i(fNECLayers) \n", fNECLayers);
406   printf(" fLongModuleSize     %6.3f cm \n", fLongModuleSize);
407   printf(" #supermodule in phi direction %i \n", fNPhiSuperModule );
408   printf(" fILOSS %i : fIHADR %i \n", fILOSS, fIHADR);
409   printf(" fTrd1Angle %7.4f\n", fTrd1Angle);
410   printf(" f2Trd1Dx2  %7.4f\n",  f2Trd1Dx2);
411   printf("SM dimensions(TRD1) : dx %7.2f dy %7.2f dz %7.2f (SMOD, BOX)\n", 
412          fParSM[0],fParSM[1],fParSM[2]);
413   printf(" fPhiGapForSM  %7.4f cm (%7.4f <- phi size in degree)\n",  
414          fPhiGapForSM, TMath::ATan2(fPhiGapForSM,fIPDistance)*TMath::RadToDeg());
415   if(GetKey110DEG()) printf(" Last two modules have size 10 degree in  phi (180<phi<190)\n");
416   printf(" phi SM boundaries \n"); 
417   for(int i=0; i<fPhiBoundariesOfSM.GetSize()/2.; i++) {
418     printf(" %i : %7.5f(%7.2f) -> %7.5f(%7.2f) : center %7.5f(%7.2f) \n", i, 
419            fPhiBoundariesOfSM[2*i], fPhiBoundariesOfSM[2*i]*TMath::RadToDeg(),
420            fPhiBoundariesOfSM[2*i+1], fPhiBoundariesOfSM[2*i+1]*TMath::RadToDeg(),
421            fPhiCentersOfSM[i], fPhiCentersOfSM[i]*TMath::RadToDeg());
422   }
423
424 }
425
426 //______________________________________________________________________
427 void AliEMCALEMCGeometry::CheckAdditionalOptions()
428 {
429   // Feb 06,2006
430   // Additional options that
431   // can be used to select
432   // the specific geometry of 
433   // EMCAL to run
434   // Dec 27,2006
435   // adeed allILOSS= and allIHADR= for MIP investigation
436   fArrayOpts = new TObjArray;
437   Int_t nopt = ParseString(fGeoName, *fArrayOpts);
438   if(nopt==1) { // no aditional option(s)
439     fArrayOpts->Delete();
440     delete fArrayOpts;
441     fArrayOpts = 0; 
442     return;
443   }              
444   for(Int_t i=1; i<nopt; i++){
445     TObjString *o = (TObjString*)fArrayOpts->At(i); 
446
447     TString addOpt = o->String();
448     Int_t indj=-1;
449     for(Int_t j=0; j<fNAdditionalOpts; j++) {
450       TString opt = fkAdditionalOpts[j];
451       if(addOpt.Contains(opt,TString::kIgnoreCase)) {
452           indj = j;
453         break;
454       }
455     }
456     if(indj<0) {
457       AliDebug(2,Form("<E> option |%s| unavailable : ** look to the file AliEMCALGeometry.h **\n", 
458                       addOpt.Data()));
459       assert(0);
460     } else {
461       AliDebug(2,Form("<I> option |%s| is valid : number %i : |%s|\n", 
462                       addOpt.Data(), indj, fkAdditionalOpts[indj]));
463       if       (addOpt.Contains("NL=",TString::kIgnoreCase))   {// number of sampling layers
464         sscanf(addOpt.Data(),"NL=%i", &fNECLayers);
465         AliDebug(2,Form(" fNECLayers %i (new) \n", fNECLayers));
466       } else if(addOpt.Contains("PBTH=",TString::kIgnoreCase)) {//Thickness of the Pb(fECPbRadThicknes)
467         sscanf(addOpt.Data(),"PBTH=%f", &fECPbRadThickness);
468       } else if(addOpt.Contains("SCTH=",TString::kIgnoreCase)) {//Thickness of the Sc(fECScintThick)
469         sscanf(addOpt.Data(),"SCTH=%f", &fECScintThick);
470       } else if(addOpt.Contains("LATSS=",TString::kIgnoreCase)) {// Thickness of lateral steel strip (fLateralSteelStrip)
471         sscanf(addOpt.Data(),"LATSS=%f", &fLateralSteelStrip);
472         AliDebug(2,Form(" fLateralSteelStrip %f (new) \n", fLateralSteelStrip));
473       } else if(addOpt.Contains("ILOSS=",TString::kIgnoreCase)) {// As in Geant
474         sscanf(addOpt.Data(),"ALLILOSS=%i", &fILOSS);
475         AliDebug(2,Form(" fILOSS %i \n", fILOSS));
476       } else if(addOpt.Contains("IHADR=",TString::kIgnoreCase)) {// As in Geant
477         sscanf(addOpt.Data(),"ALLIHADR=%i", &fIHADR);
478         AliDebug(2,Form(" fIHADR %i \n", fIHADR));
479       }
480     }
481   }
482 }
483
484 //__________________________________________________________________
485 void AliEMCALEMCGeometry::DefineSamplingFraction()
486 {
487   // Jun 05,2006
488   // Look http://rhic.physics.wayne.edu/~pavlinov/ALICE/SHISHKEBAB/RES/linearityAndResolutionForTRD1.html
489   // Keep for compatibilty
490   //
491   if(fNECLayers == 69) {        // 10% layer reduction
492     fSampling = 12.55;
493   } else if(fNECLayers == 61) { // 20% layer reduction
494     fSampling = 12.80;
495   } else if(fNECLayers == 77) {
496     if       (fECScintThick>0.159 && fECScintThick<0.161) { // original sampling fraction, equal layers
497       fSampling = 12.327; // fECScintThick = fECPbRadThickness = 0.160;
498     } else if (fECScintThick>0.175 && fECScintThick<0.177) { // 10% Pb thicknes reduction
499       fSampling = 10.5; // fECScintThick = 0.176, fECPbRadThickness=0.144;
500     } else if(fECScintThick>0.191 && fECScintThick<0.193) { // 20% Pb thicknes reduction
501       fSampling = 8.93; // fECScintThick = 0.192, fECPbRadThickness=0.128;
502     }
503
504   }
505 }
506
507 //________________________________________________________________________________________________
508 Double_t AliEMCALEMCGeometry::GetPhiCenterOfSM(Int_t nsupmod) const
509 {
510   //returns center of supermodule in phi
511   int i = nsupmod/2;
512   return fPhiCentersOfSM[i];
513
514 }
515
516 //________________________________________________________________________________________________
517 Bool_t AliEMCALEMCGeometry::GetPhiBoundariesOfSM(Int_t nSupMod, Double_t &phiMin, Double_t &phiMax) const
518 {
519   // 0<= nSupMod <=11; phi in rad
520     static int i;
521   if(nSupMod<0 || nSupMod >11) return kFALSE;
522   i = nSupMod/2;
523   phiMin = (Double_t)fPhiBoundariesOfSM[2*i];
524   phiMax = (Double_t)fPhiBoundariesOfSM[2*i+1];
525   return kTRUE;
526 }
527
528 //________________________________________________________________________________________________
529 Bool_t AliEMCALEMCGeometry::GetPhiBoundariesOfSMGap(Int_t nPhiSec, Double_t &phiMin, Double_t &phiMax) const
530 {
531   // 0<= nPhiSec <=4; phi in rad
532   // 0;  gap boundaries between  0th&2th  | 1th&3th SM
533   // 1;  gap boundaries between  2th&4th  | 3th&5th SM
534   // 2;  gap boundaries between  4th&6th  | 5th&7th SM
535   // 3;  gap boundaries between  6th&8th  | 7th&9th SM
536   // 4;  gap boundaries between  8th&10th | 9th&11th SM
537   if(nPhiSec<0 || nPhiSec >4) return kFALSE;
538   phiMin = fPhiBoundariesOfSM[2*nPhiSec+1];
539   phiMax = fPhiBoundariesOfSM[2*nPhiSec+2];
540   return kTRUE;
541 }
542
543 //________________________________________________________________________________________________
544 int AliEMCALEMCGeometry::ParseString(const TString &topt, TObjArray &Opt)
545
546         //Parse string, does what? GCB 08/09
547         Ssiz_t begin, index, end, end2;
548         begin = index = end = end2 = 0;
549         TRegexp separator("[^ ;,\\t\\s/]+");
550         while ( (begin < topt.Length()) && (index != kNPOS) ) {
551                 // loop over given options
552                 index = topt.Index(separator,&end,begin);
553                 if (index >= 0 && end >= 1) {
554                         TString substring(topt(index,end));
555                         Opt.Add(new TObjString(substring.Data()));
556                 }
557                 begin += end+1;
558         }
559         return Opt.GetEntries();
560 }
561