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