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