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