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