Project
Loading...
Searching...
No Matches
Geometry.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11#include "EMCALBase/Geometry.h"
12
13#include <RtypesCore.h>
14#include <TMath.h>
15#include <TVector3.h>
16#include <TMathBase.h>
17#include <TVector2.h>
18#include <TParticle.h>
19#include <TString.h>
20#include <TGeoNode.h>
21#include <TJAlienCredentials.h>
22#include <TObjArray.h>
23#include <fairlogger/Logger.h>
24
25#include <cstring>
26#include <cctype>
27#include <cmath>
28#include <iomanip>
29#include <ostream>
30#include <string>
31#include <algorithm>
32#include <cstdio>
33#include <string_view>
34#include <tuple>
35
36#include <TGeoBBox.h>
37#include <TGeoManager.h>
38#include <TGeoMatrix.h>
39
42#include "CCDB/CcdbApi.h"
44#include "GPUROOTCartesianFwd.h"
45
46#include <boost/algorithm/string/predicate.hpp>
47
48using namespace o2::emcal;
49
50// these initialisations are needed for a singleton
51Geometry* Geometry::sGeom = nullptr;
52
54 : mGeoName(geo.mGeoName),
55 mKey110DEG(geo.mKey110DEG),
56 mnSupModInDCAL(geo.mnSupModInDCAL),
57 mNCellsInSupMod(geo.mNCellsInSupMod),
58 mNETAdiv(geo.mNETAdiv),
59 mNPHIdiv(geo.mNPHIdiv),
60 mNCellsInModule(geo.mNCellsInModule),
61 mPhiBoundariesOfSM(geo.mPhiBoundariesOfSM),
62 mPhiCentersOfSM(geo.mPhiCentersOfSM),
63 mPhiCentersOfSMSec(geo.mPhiCentersOfSMSec),
64 mPhiCentersOfCells(geo.mPhiCentersOfCells),
65 mCentersOfCellsEtaDir(geo.mCentersOfCellsEtaDir),
66 mCentersOfCellsPhiDir(geo.mCentersOfCellsPhiDir),
67 mEtaCentersOfCells(geo.mEtaCentersOfCells),
68 mNCells(geo.mNCells),
69 mNPhi(geo.mNPhi),
70 mCentersOfCellsXDir(geo.mCentersOfCellsXDir),
71 mArm1EtaMin(geo.mArm1EtaMin),
72 mArm1EtaMax(geo.mArm1EtaMax),
73 mArm1PhiMin(geo.mArm1PhiMin),
74 mArm1PhiMax(geo.mArm1PhiMax),
75 mEtaMaxOfTRD1(geo.mEtaMaxOfTRD1),
76 mDCALPhiMin(geo.mDCALPhiMin),
77 mDCALPhiMax(geo.mDCALPhiMax),
78 mEMCALPhiMax(geo.mEMCALPhiMax),
79 mDCALStandardPhiMax(geo.mDCALStandardPhiMax),
80 mDCALInnerExtandedEta(geo.mDCALInnerExtandedEta),
81 mDCALInnerEdge(geo.mDCALInnerEdge),
82 mShishKebabTrd1Modules(geo.mShishKebabTrd1Modules),
83 mPhiModuleSize(geo.mPhiModuleSize),
84 mEtaModuleSize(geo.mEtaModuleSize),
85 mPhiTileSize(geo.mPhiTileSize),
86 mEtaTileSize(geo.mEtaTileSize),
87 mNZ(geo.mNZ),
88 mIPDistance(geo.mIPDistance),
89 mLongModuleSize(geo.mLongModuleSize),
90 mShellThickness(geo.mShellThickness),
91 mZLength(geo.mZLength),
92 mSampling(geo.mSampling),
93 mECPbRadThickness(geo.mECPbRadThickness),
94 mECScintThick(geo.mECScintThick),
95 mNECLayers(geo.mNECLayers),
96 mNumberOfSuperModules(geo.mNumberOfSuperModules),
97 mEMCSMSystem(geo.mEMCSMSystem),
98 mFrontSteelStrip(geo.mFrontSteelStrip),
99 mLateralSteelStrip(geo.mLateralSteelStrip),
100 mPassiveScintThick(geo.mPassiveScintThick),
101 mPhiSuperModule(geo.mPhiSuperModule),
102 mNPhiSuperModule(geo.mNPhiSuperModule),
103 mTrd1Angle(geo.mTrd1Angle),
104 m2Trd1Dx2(geo.m2Trd1Dx2),
105 mPhiGapForSM(geo.mPhiGapForSM),
106 mTrd1AlFrontThick(geo.mTrd1AlFrontThick),
107 mTrd1BondPaperThick(geo.mTrd1BondPaperThick),
108 mILOSS(geo.mILOSS),
109 mIHADR(geo.mIHADR),
110 mSteelFrontThick(geo.mSteelFrontThick), // obsolete data member?
111 mCellIndexLookup(geo.mCellIndexLookup)
112{
113 memcpy(mEnvelop, geo.mEnvelop, sizeof(Float_t) * 3);
114 memcpy(mParSM, geo.mParSM, sizeof(Float_t) * 3);
115
116 memset(SMODULEMATRIX, 0, sizeof(TGeoHMatrix*) * EMCAL_MODULES);
117}
118
119Geometry::Geometry(const std::string_view name, const std::string_view mcname, const std::string_view mctitle)
120 : mGeoName(name),
121 mKey110DEG(0),
122 mnSupModInDCAL(0),
123 mNCellsInSupMod(0),
124 mNETAdiv(0),
125 mNPHIdiv(0),
126 mNCellsInModule(0),
127 mPhiBoundariesOfSM(),
128 mPhiCentersOfSM(),
129 mPhiCentersOfSMSec(),
130 mPhiCentersOfCells(),
131 mCentersOfCellsEtaDir(),
132 mCentersOfCellsPhiDir(),
133 mEtaCentersOfCells(),
134 mNCells(0),
135 mNPhi(0),
136 mCentersOfCellsXDir(),
137 mArm1EtaMin(0),
138 mArm1EtaMax(0),
139 mArm1PhiMin(0),
140 mArm1PhiMax(0),
141 mEtaMaxOfTRD1(0),
142 mDCALPhiMin(0),
143 mDCALPhiMax(0),
144 mEMCALPhiMax(0),
145 mDCALStandardPhiMax(0),
146 mDCALInnerExtandedEta(0),
147 mDCALInnerEdge(0.),
148 mShishKebabTrd1Modules(),
149 mPhiModuleSize(0.),
150 mEtaModuleSize(0.),
151 mPhiTileSize(0.),
152 mEtaTileSize(0.),
153 mNZ(0),
154 mIPDistance(0.),
155 mLongModuleSize(0.),
156 mShellThickness(0.),
157 mZLength(0.),
158 mSampling(0.),
159 mECPbRadThickness(0.),
160 mECScintThick(0.),
161 mNECLayers(0),
162 mNumberOfSuperModules(0),
163 mEMCSMSystem(),
164 mFrontSteelStrip(0.),
165 mLateralSteelStrip(0.),
166 mPassiveScintThick(0.),
167 mPhiSuperModule(0),
168 mNPhiSuperModule(0),
169 mTrd1Angle(0.),
170 m2Trd1Dx2(0.),
171 mPhiGapForSM(0.),
172 mTrd1AlFrontThick(0.0),
173 mTrd1BondPaperThick(0.),
174 mILOSS(-1),
175 mIHADR(-1),
176 mSteelFrontThick(0.) // obsolete data member?
177{
178 DefineEMC(mcname, mctitle);
180
182
184 for (auto icell = 0; icell < mNCells; icell++) {
185 mCellIndexLookup[icell] = CalculateCellIndex(icell);
186 }
187
188 memset(SMODULEMATRIX, 0, sizeof(TGeoHMatrix*) * EMCAL_MODULES);
189
190 LOG(debug) << "Name <<" << name << ">>";
191}
192
194{
195 LOG(fatal) << "assignment operator, not implemented";
196 return *this;
197}
198
200{
201 if (this == sGeom) {
202 LOG(error) << "Do not call delete on me";
203 return;
204 }
205
206 for (Int_t smod = 0; smod < mNumberOfSuperModules; smod++) {
207 if (SMODULEMATRIX[smod]) {
208 delete SMODULEMATRIX[smod];
209 }
210 }
211}
212
214{
215 Geometry* rv = static_cast<Geometry*>(sGeom);
216 if (!rv) {
218 }
219 return rv;
220}
221
222Geometry* Geometry::GetInstance(const std::string_view name, const std::string_view mcname,
223 const std::string_view mctitle)
224{
225 if (!sGeom) {
226 if (!name.length()) { // get default geometry
227 sGeom = new Geometry(DEFAULT_GEOMETRY, mcname, mctitle);
228 } else {
229 sGeom = new Geometry(name, mcname, mctitle);
230 } // end if strcmp(name,"")
231 return sGeom;
232 } else {
233 if (sGeom->GetName() != name) {
234 LOG(info) << "\n current geometry is " << sGeom->GetName() << " : you should not call " << name;
235 } // end
236 return sGeom;
237 } // end if sGeom
238
239 return nullptr;
240}
241
242Geometry* Geometry::GetInstanceFromRunNumber(Int_t runNumber, const std::string_view geoName,
243 const std::string_view mcname, const std::string_view mctitle)
244{
245 using boost::algorithm::contains;
246
247 // printf("AliEMCALGeometry::GetInstanceFromRunNumber() - run %d, geoName <<%s>> \n",runNumber,geoName.Data());
248
249 if (runNumber >= 104064 && runNumber < 140000) {
250 // 2009-2010 runs
251 // First year geometry, 4 SM.
252
253 if (contains(geoName, "FIRSTYEARV1") && geoName != std::string("")) {
254 LOG(info) << "o2::emcal::Geometry::GetInstanceFromRunNumber() *** ATTENTION *** \n"
255 << "\t Specified geometry name <<" << geoName << ">> for run " << runNumber
256 << " is not considered! \n"
257 << "\t In use <<EMCAL_FIRSTYEARV1>>, check run number and year";
258 } else {
259 LOG(info)
260 << "o2::emcal::Geometry::GetInstanceFromRunNumber() - Initialized geometry with name <<EMCAL_FIRSTYEARV1>>";
261 }
262
263 return Geometry::GetInstance("EMCAL_FIRSTYEARV1", mcname, mctitle);
264 } else if (runNumber >= 140000 && runNumber <= 170593) {
265 // Almost complete EMCAL geometry, 10 SM. Year 2011 configuration
266
267 if (contains(geoName, "COMPLETEV1") && geoName != std::string("")) {
268 LOG(info) << "o2::emcal::Geometry::GetInstanceFromRunNumber() *** ATTENTION *** \n"
269 << "\t Specified geometry name <<" << geoName << ">> for run " << runNumber
270 << " is not considered! \n"
271 << "\t In use <<EMCAL_COMPLETEV1>>, check run number and year";
272 } else {
273 LOG(info)
274 << "o2::emcal::Geometry::GetInstanceFromRunNumber() - Initialized geometry with name <<EMCAL_COMPLETEV1>>";
275 }
276 return Geometry::GetInstance("EMCAL_COMPLETEV1", mcname, mctitle);
277 } else if (runNumber > 176000 && runNumber <= 197692) {
278 // Complete EMCAL geometry, 12 SM. Year 2012 and on
279 // The last 2 SM were not active, anyway they were there.
280
281 if (contains(geoName, "COMPLETE12SMV1") && geoName != std::string("")) {
282 LOG(info) << "o2::emcal::Geometry::GetInstanceFromRunNumber() *** ATTENTION *** \n"
283 << "\t Specified geometry name <<" << geoName << " >> for run " << runNumber
284 << " is not considered! \n"
285 << "\t In use <<EMCAL_COMPLETE12SMV1>>, check run number and year";
286 } else {
287 LOG(info) << "o2::emcal::Geometry::GetInstanceFromRunNumber() - Initialized geometry with name "
288 "<<EMCAL_COMPLETE12SMV1>>";
289 }
290 return Geometry::GetInstance("EMCAL_COMPLETE12SMV1", mcname, mctitle);
291 } else // Run 2
292 {
293 // EMCAL + DCAL geometry, 20 SM. Year 2015 and on
294
295 if (contains(geoName, "DCAL_8SM") && geoName != std::string("")) {
296 LOG(info) << "o2::emcal::Geometry::GetInstanceFromRunNumber() *** ATTENTION *** \n"
297 << "\t Specified geometry name <<" << geoName << ">> for run " << runNumber
298 << " is not considered! \n"
299 << "\t In use <<EMCAL_COMPLETE12SMV1_DCAL_8SM>>, check run number and year";
300 } else {
301 LOG(info) << "o2::emcal::Geometry::GetInstanceFromRunNumber() - Initialized geometry with name "
302 "<<EMCAL_COMPLETE12SMV1_DCAL_8SM>>";
303 }
304 return Geometry::GetInstance("EMCAL_COMPLETE12SMV1_DCAL_8SM", mcname, mctitle);
305 }
306}
307
308void Geometry::DefineSamplingFraction(const std::string_view mcname, const std::string_view mctitle)
309{
310 // Jun 05,2006
311 // Look http://rhic.physics.wayne.edu/~pavlinov/ALICE/SHISHKEBAB/RES/linearityAndResolutionForTRD1.html
312 // Keep for compatibility
313 //
314 using boost::algorithm::contains;
315
316 // Sampling factor for G3
317 mSampling = 10.87; // Default value - Nov 25,2010
318 if (mNECLayers == 69) { // 10% layer reduction
319 mSampling = 12.55;
320 } else if (mNECLayers == 61) { // 20% layer reduction
321 mSampling = 12.80;
322 } else if (mNECLayers == 77) {
323 if (contains(mGeoName, "V1")) {
324 mSampling = 10.87; // Adding paper sheets and cover plate; Nov 25,2010
325 } else if (mECScintThick > 0.159 && mECScintThick < 0.161) { // original sampling fraction, equal layers
326 mSampling = 12.327; // fECScintThick = fECPbRadThickness = 0.160;
327 } else if (mECScintThick > 0.175 && mECScintThick < 0.177) { // 10% Pb thicknes reduction
328 mSampling = 10.5; // fECScintThick = 0.176, fECPbRadThickness=0.144;
329 } else if (mECScintThick > 0.191 && mECScintThick < 0.193) { // 20% Pb thicknes reduction
330 mSampling = 8.93; // fECScintThick = 0.192, fECPbRadThickness=0.128;
331 }
332 }
333
334 Float_t samplingFactorTranportModel = 1.;
335 // Note: The sampling factors are chosen so that results from the simulation
336 // engines correspond well with testbeam data
337
338 if (contains(mcname, "Geant3")) {
339 samplingFactorTranportModel = 1.; // 0.988 // Do nothing
340 } else if (contains(mcname, "Fluka")) {
341 samplingFactorTranportModel = 1.; // To be set
342 } else if (contains(mcname, "Geant4")) {
343 std::string physicslist = mctitle.substr(mctitle.find(":") + 2).data();
344 LOG(info) << "Selected physics list: " << physicslist;
345 // sampling factors for different Geant4 physics list
346 // GEANT4 10.7 -> EMCAL-784
347
348 // set a default (there may be many physics list strings)
349 samplingFactorTranportModel = 0.81;
350 if (physicslist == "FTFP_BERT_EMV+optical") {
351 samplingFactorTranportModel = 0.821;
352 } else if (physicslist == "FTFP_BERT_EMV+optical+biasing") {
353 samplingFactorTranportModel = 0.81;
354 } else if (physicslist == "FTFP_INCLXX_EMV+optical") {
355 samplingFactorTranportModel = 0.81;
356 }
357 }
358
359 LOG(info) << "MC modeler <" << mcname << ">, Title <" << mctitle << ">: Sampling " << std::setw(2)
360 << std::setprecision(3) << mSampling << ", model fraction with respect to G3 "
361 << samplingFactorTranportModel << ", final sampling " << mSampling * samplingFactorTranportModel;
362
363 mSampling *= samplingFactorTranportModel;
364}
365
366void Geometry::DefineEMC(std::string_view mcname, std::string_view mctitle)
367{
368 using boost::algorithm::contains;
369
370 // geometry
371 std::transform(mGeoName.begin(), mGeoName.end(), mGeoName.begin(), ::toupper);
372
373 // Convert old geometry names to new ones
374 if (contains(mGeoName, "SHISH_77_TRD1_2X2_FINAL_110DEG")) {
375 if (contains(mGeoName, "PBTH=0.144") && contains(mGeoName, "SCTH=0.176")) {
376 mGeoName = "EMCAL_COMPLETE";
377 } else {
378 mGeoName = "EMCAL_PDC06";
379 }
380 }
381
382 if (contains(mGeoName, "WSUC")) {
383 mGeoName = "EMCAL_WSUC";
384 }
385
386 // check that we have a valid geometry name
387 if (!(contains(mGeoName, "EMCAL_PDC06") || contains(mGeoName, "EMCAL_WSUC") || contains(mGeoName, "EMCAL_COMPLETE") ||
388 contains(mGeoName, "EMCAL_COMPLETEV1") || contains(mGeoName, "EMCAL_COMPLETE12SMV1") ||
389 contains(mGeoName, "EMCAL_FIRSTYEAR") || contains(mGeoName, "EMCAL_FIRSTYEARV1"))) {
390 LOG(fatal) << "Init, " << mGeoName << " is an undefined geometry!\n";
391 }
392
393 // Option to know whether we have the "half" supermodule(s) or not
394 mKey110DEG = 0;
395 if (contains(mGeoName, "COMPLETE") || contains(mGeoName, "PDC06") || contains(mGeoName, "12SM")) {
396 mKey110DEG = 1; // for GetAbsCellId
397 }
398 if (contains(mGeoName, "COMPLETEV1")) {
399 mKey110DEG = 0;
400 }
401
402 mnSupModInDCAL = 0;
403 if (contains(mGeoName, "DCAL_DEV")) {
404 mnSupModInDCAL = 10;
405 } else if (contains(mGeoName, "DCAL_8SM")) {
406 mnSupModInDCAL = 8;
407 } else if (contains(mGeoName, "DCAL")) {
408 mnSupModInDCAL = 6;
409 }
410
411 // JLK 13-Apr-2008
412 // default parameters are those of EMCAL_COMPLETE geometry
413 // all others render variations from these at the end of
414 // geometry-name specific options
415
416 mNumberOfSuperModules = 12; // 12 = 6 * 2 (6 in phi, 2 in Z)
417 mNPhi = 12; // module granularity in phi within smod (azimuth)
418 mNZ = 24; // module granularity along Z within smod (eta)
419 mNPHIdiv = mNETAdiv = 2; // tower granularity within module
420 mArm1PhiMin = 80.0; // degrees, Starting EMCAL Phi position
421 mArm1PhiMax = 200.0; // degrees, Ending EMCAL Phi position
422 mArm1EtaMin = -0.7; // pseudorapidity, Starting EMCAL Eta position
423 mArm1EtaMax = +0.7; // pseudorapidity, Ending EMCAL Eta position
424 mIPDistance = 428.0; // cm, radial distance to front face from nominal vertex point
425 mPhiGapForSM = 2.; // cm, only for final TRD1 geometry
426 mFrontSteelStrip = 0.025; // 0.025cm = 0.25mm (13-may-05 from V.Petrov)
427 mPassiveScintThick = 0.8; // 0.8cm = 8mm (13-may-05 from V.Petrov)
428 mLateralSteelStrip = 0.01; // 0.01cm = 0.1mm (13-may-05 from V.Petrov) - was 0.025
429 mTrd1Angle = 1.5; // in degrees
430
431 mSampling = 1.; // should be calculated with call to DefineSamplingFraction()
432 mNECLayers = 77; // (13-may-05 from V.Petrov) - can be changed with additional options
433 mECScintThick = 0.176; // scintillator layer thickness
434 mECPbRadThickness = 0.144; // lead layer thickness
435
436 mPhiModuleSize = 12.26 - mPhiGapForSM / Float_t(mNPhi); // first assumption
438
439 mZLength = 700.; // Z coverage (cm)
440 mPhiSuperModule = 20.; // phi in degree
441 mDCALInnerEdge = mIPDistance * TMath::Tan(mTrd1Angle * 8. * TMath::DegToRad());
442
443 // modifications to the above for PDC06 geometry
444 if (contains(mGeoName, "PDC06")) { // 18-may-05 - about common structure
445 mECScintThick = mECPbRadThickness = 0.16; // (13-may-05 from V.Petrov)
446 }
447
448 // modifications to the above for WSUC geometry
449 if (contains(mGeoName, "WSUC")) { // 18-may-05 - about common structure
450 mNumberOfSuperModules = 2; // 27-may-05; Nov 24,2010 for TB
451 mNPhi = mNZ = 4;
452 mTrd1AlFrontThick = 1.0; // one cm
453 // Bond paper - two sheets around Sc tile
454 mTrd1BondPaperThick = 0.01; // 0.01cm = 0.1 mm
455
456 mPhiModuleSize = 12.0;
458 mLateralSteelStrip = 0.015; // 0.015cm = 0.15mm
459 }
460
461 // In 2009-2010 data taking runs only 4 SM, in the upper position.
462 if (contains(mGeoName, "FIRSTYEAR")) {
464 mArm1PhiMax = 120.0;
465 }
466
467 if (contains(mGeoName, "FIRSTYEARV1") || contains(mGeoName, "COMPLETEV1") || contains(mGeoName, "COMPLETE12SMV1")) {
468 // Oct 26,2010 : First module has tilt = 0.75 degree :
469 // look to AliEMCALShishKebabTrd1Module::DefineFirstModule(key)
470 // New sizes from production drawing, added Al front plate.
471 // The thickness of sampling is change due to existing two sheets of paper.
472
473 // Will replace fFrontSteelStrip
474 mTrd1AlFrontThick = 1.0; // one cm
475 // Bond paper - two sheets around Sc tile
476 mTrd1BondPaperThick = 0.01; // 0.01cm = 0.1 mm
477
478 mPhiModuleSize = 12.0;
480 mLateralSteelStrip = 0.015; // 0.015cm = 0.15mm
481
482 if (contains(mGeoName, "COMPLETEV1")) {
484 mArm1PhiMax = 180.0;
485 } else if (contains(mGeoName, "COMPLETE12SMV1")) {
487 mArm1PhiMax = 200.0;
488 }
489 if (contains(mGeoName, "DCAL")) {
491 mArm1PhiMax = 320.0;
492 if (contains(mGeoName, "DCAL_8SM")) {
493 mArm1PhiMax = 340.0; // degrees, End of DCAL Phi position
494 } else if (contains(mGeoName, "DCAL_DEV")) {
495 mArm1PhiMin = 40.0; // degrees, Starting EMCAL(shifted) Phi position
496 }
498 }
499 }
500
501 //
502 // Init EMCal/DCal SMs type array
503 mEMCSMSystem.clear();
505
506 for (Int_t i = 0; i < mNumberOfSuperModules; i++) {
508 }
509
510 Int_t iSM = 0;
511
512 //
513 // BASIC EMCAL SM
514 if (contains(mGeoName, "WSUC")) {
515 for (int i = 0; i < 2; i++) {
517 iSM++;
518 }
519 } else if (contains(mGeoName, "FIRSTYEAR")) {
520 for (int i = 0; i < 4; i++) {
522 iSM++;
523 }
524 } else if (contains(mGeoName, "PDC06") || contains(mGeoName, "COMPLETE")) {
525 for (int i = 0; i < 10; i++) {
527 iSM++;
528 }
529 }
530
531 //
532 // EMCAL 110SM
533 if (mKey110DEG && contains(mGeoName, "12SM")) {
534 for (int i = 0; i < 2; i++) {
536 if (contains(mGeoName, "12SMV1")) {
538 }
539 iSM++;
540 }
541 }
542
543 //
544 // DCAL SM
545 if (mnSupModInDCAL && contains(mGeoName, "DCAL")) {
546 if (contains(mGeoName, "8SM")) {
547 for (int i = 0; i < mnSupModInDCAL - 2; i++) {
549 iSM++;
550 }
551 for (int i = 0; i < 2; i++) {
552 mEMCSMSystem[iSM] = DCAL_EXT;
553 iSM++;
554 }
555 } else {
556 for (int i = 0; i < mnSupModInDCAL; i++) {
558 iSM++;
559 }
560 }
561 }
562
563 // constant for transition absid <--> indexes
566 mNCells = 0;
567 for (int i = 0; i < mNumberOfSuperModules; i++) {
568 if (GetSMType(i) == EMCAL_STANDARD) {
570 } else if (GetSMType(i) == EMCAL_HALF) {
572 } else if (GetSMType(i) == EMCAL_THIRD) {
574 } else if (GetSMType(i) == DCAL_STANDARD) {
575 mNCells += 2 * mNCellsInSupMod / 3;
576 } else if (GetSMType(i) == DCAL_EXT) {
578 } else {
579 LOG(error) << "Uknown SuperModule Type !!\n";
580 }
581 }
582
584 if (mNPhiSuperModule < 1) {
586 }
587
588 mPhiTileSize = mPhiModuleSize / double(mNPHIdiv) - mLateralSteelStrip; // 13-may-05
589 mEtaTileSize = mEtaModuleSize / double(mNETAdiv) - mLateralSteelStrip; // 13-may-05
590
592 if (contains(mGeoName, "V1")) {
593 Double_t ws = mECScintThick + mECPbRadThickness + 2. * mTrd1BondPaperThick; // sampling width
594 // Number of Pb tiles = Number of Sc tiles - 1
596 }
597 m2Trd1Dx2 = mEtaModuleSize + 2. * mLongModuleSize * TMath::Tan(mTrd1Angle * TMath::DegToRad() / 2.);
598
599 if (!contains(mGeoName, "WSUC")) {
601 }
602
603 // These parameters are used to create the mother volume to hold the supermodules
604 // 2cm padding added to allow for misalignments - JLK 30-May-2008
605 mEnvelop[0] = mIPDistance - 1.; // mother volume inner radius
606 mEnvelop[1] = mIPDistance + mShellThickness + 1.; // mother volume outer r.
607 mEnvelop[2] = mZLength + 2.; // mother volume length
608
609 // Local coordinates
610 mParSM[0] = GetShellThickness() / 2.;
611 mParSM[1] = GetPhiModuleSize() * GetNPhi() / 2.;
612 mParSM[2] = mZLength / 4.; // divide by 4 to get half-length of SM
613
614 // SM phi boundaries - (0,1),(2,3) ... - has the same boundaries;
618 Double_t kfSupermodulePhiWidth = mPhiSuperModule * TMath::DegToRad();
619 mPhiCentersOfSM[0] = (mArm1PhiMin + mPhiSuperModule / 2.) * TMath::DegToRad(); // Define from First SM
620 mPhiCentersOfSMSec[0] = mPhiCentersOfSM[0]; // the same in the First SM
621 mPhiBoundariesOfSM[0] = mPhiCentersOfSM[0] - TMath::ATan2(mParSM[1], mIPDistance); // 1th and 2th modules)
622 mPhiBoundariesOfSM[1] = mPhiCentersOfSM[0] + TMath::ATan2(mParSM[1], mIPDistance);
623
624 if (mNumberOfSuperModules > 2) { // 2 to Max
625 Int_t tmpSMType = GetSMType(2);
626 for (int i = 1; i < mNPhiSuperModule; i++) {
627 mPhiBoundariesOfSM[2 * i] += mPhiBoundariesOfSM[2 * i - 2] + kfSupermodulePhiWidth;
628 if (tmpSMType == GetSMType(2 * i)) {
629 mPhiBoundariesOfSM[2 * i + 1] += mPhiBoundariesOfSM[2 * i - 1] + kfSupermodulePhiWidth;
630 } else {
631 // changed SM Type, redefine the [2*i+1] Boundaries
632 tmpSMType = GetSMType(2 * i);
633 switch (GetSMType(2 * i)) {
634 case EMCAL_STANDARD:
635 mPhiBoundariesOfSM[2 * i + 1] = mPhiBoundariesOfSM[2 * i] + kfSupermodulePhiWidth;
636 break;
637 case EMCAL_HALF:
638 mPhiBoundariesOfSM[2 * i + 1] = mPhiBoundariesOfSM[2 * i] + 2. * TMath::ATan2((mParSM[1]) / 2, mIPDistance);
639 break;
640 case EMCAL_THIRD:
641 mPhiBoundariesOfSM[2 * i + 1] = mPhiBoundariesOfSM[2 * i] + 2. * TMath::ATan2((mParSM[1]) / 3, mIPDistance);
642 break;
643 case DCAL_STANDARD:
644 mPhiBoundariesOfSM[2 * i] = (mDCALPhiMin - mArm1PhiMin) * TMath::DegToRad() + mPhiBoundariesOfSM[0];
645 mPhiBoundariesOfSM[2 * i + 1] = (mDCALPhiMin - mArm1PhiMin) * TMath::DegToRad() + mPhiBoundariesOfSM[1];
646 break;
647 case DCAL_EXT:
648 mPhiBoundariesOfSM[2 * i + 1] = mPhiBoundariesOfSM[2 * i] + 2. * TMath::ATan2((mParSM[1]) / 3, mIPDistance);
649 break;
650 default:
651 break;
652 };
653 }
654 mPhiCentersOfSM[i] = (mPhiBoundariesOfSM[2 * i] + mPhiBoundariesOfSM[2 * i + 1]) / 2.;
655 mPhiCentersOfSMSec[i] = mPhiBoundariesOfSM[2 * i] + TMath::ATan2(mParSM[1], mIPDistance);
656 }
657 }
658
659 // inner extend in eta (same as outer part) for DCal (0.189917), //calculated from the smallest gap (1# cell to the
660 // 80-degree-edge),
661 const double INNNER_EXTENDED_PHI =
662 1.102840997; // calculated from the smallest gap (1# cell to the 80-degree-edge), too complicatd to explain...
663 mDCALInnerExtandedEta = -TMath::Log(
664 TMath::Tan((TMath::Pi() / 2. - 8 * mTrd1Angle * TMath::DegToRad() +
665 (TMath::Pi() / 2 - mNZ * mTrd1Angle * TMath::DegToRad() - TMath::ATan(TMath::Exp(mArm1EtaMin)) * 2)) /
666 2.));
667
669 mDCALPhiMax = mDCALPhiMin; // DCAl extention will not be included
670 for (Int_t i = 0; i < mNumberOfSuperModules; i += 2) {
671 switch (GetSMType(i)) {
672 case EMCAL_STANDARD:
673 mEMCALPhiMax += 20.;
674 break;
675 case EMCAL_HALF:
676 mEMCALPhiMax += mPhiSuperModule / 2. + INNNER_EXTENDED_PHI;
677 break;
678 case EMCAL_THIRD:
679 mEMCALPhiMax += mPhiSuperModule / 3. + 4.0 * INNNER_EXTENDED_PHI / 3.0;
680 break;
681 case DCAL_STANDARD:
682 mDCALPhiMax += 20.;
684 break;
685 case DCAL_EXT:
686 mDCALPhiMax += mPhiSuperModule / 3. + 4.0 * INNNER_EXTENDED_PHI / 3.0;
687 break;
688 default:
689 LOG(error) << "Unkown SM Type!!\n";
690 break;
691 };
692 }
693 // for compatible reason
694 // if(fNumberOfSuperModules == 4) {fEMCALPhiMax = fArm1PhiMax ;}
695 if (mNumberOfSuperModules == 12) {
697 }
698
699 // called after setting of scintillator and lead layer parameters
700 // called now in AliEMCALv0::CreateGeometry() - 15/03/16
701 // DefineSamplingFraction(mcname,mctitle);
702}
703
704void Geometry::GetGlobal(const Double_t* loc, Double_t* glob, int iSM) const
705{
706 const TGeoHMatrix* m = GetMatrixForSuperModule(iSM);
707 if (m) {
708 m->LocalToMaster(loc, glob);
709 } else {
710 LOG(fatal) << "Geo matrixes are not loaded \n";
711 }
712}
713
714void Geometry::GetGlobal(const TVector3& vloc, TVector3& vglob, int iSM) const
715{
716 Double_t tglob[3], tloc[3];
717 vloc.GetXYZ(tloc);
718 GetGlobal(tloc, tglob, iSM);
719 vglob.SetXYZ(tglob[0], tglob[1], tglob[2]);
720}
721
722void Geometry::GetGlobal(Int_t absId, Double_t glob[3]) const
723{
724 double loc[3];
725
726 memset(glob, 0, sizeof(Double_t) * 3);
727 try {
728 auto cellpos = RelPosCellInSModule(absId);
729 loc[0] = cellpos.X();
730 loc[1] = cellpos.Y();
731 loc[2] = cellpos.Z();
732 } catch (InvalidCellIDException& e) {
733 LOG(error) << e.what();
734 return;
735 }
736
737 Int_t nSupMod = std::get<0>(GetCellIndex(absId));
738 const TGeoHMatrix* m = GetMatrixForSuperModule(nSupMod);
739 if (m) {
740 m->LocalToMaster(loc, glob);
741 } else {
742 LOG(fatal) << "Geo matrixes are not loaded \n";
743 }
744}
745
746void Geometry::GetGlobal(Int_t absId, TVector3& vglob) const
747{
748 Double_t glob[3];
749
750 GetGlobal(absId, glob);
751 vglob.SetXYZ(glob[0], glob[1], glob[2]);
752}
753
754std::tuple<double, double> Geometry::EtaPhiFromIndex(Int_t absId) const
755{
756 TVector3 vglob;
757 GetGlobal(absId, vglob);
758 return std::make_tuple(vglob.Eta(), vglob.Phi());
759}
760
761int Geometry::GetAbsCellId(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const
762{
763 // 0 <= nSupMod < fNumberOfSuperModules
764 // 0 <= nModule < fNPHI * fNZ ( fNPHI * fNZ/2 for fKey110DEG=1)
765 // 0 <= nIphi < fNPHIdiv
766 // 0 <= nIeta < fNETAdiv
767 // 0 <= absid < fNCells
768 int cellid = 0; // have to change from 0 to fNCells-1
769 for (int i = 0; i < supermoduleID; i++) {
770 switch (GetSMType(i)) {
771 case EMCAL_STANDARD:
772 cellid += mNCellsInSupMod;
773 break;
774 case EMCAL_HALF:
775 cellid += mNCellsInSupMod / 2;
776 break;
777 case EMCAL_THIRD:
778 cellid += mNCellsInSupMod / 3;
779 break;
780 case DCAL_STANDARD:
781 cellid += 2 * mNCellsInSupMod / 3;
782 break;
783 case DCAL_EXT:
784 cellid += mNCellsInSupMod / 3;
785 break;
786 default:
788 };
789 }
790
791 cellid += mNCellsInModule * moduleID;
792 cellid += mNPHIdiv * phiInModule;
793 cellid += etaInModule;
794 if (!CheckAbsCellId(cellid)) {
795 throw InvalidCellIDException(cellid);
796 }
797
798 return cellid;
799}
800
801std::tuple<int, int, int> Geometry::GetModuleIndexesFromCellIndexesInSModule(int supermoduleID, int phiInSupermodule, int etaInSupermodule) const
802{
803 int nModulesInSMPhi = GetNumberOfModuleInPhiDirection(supermoduleID);
804
805 int moduleEta = etaInSupermodule / mNETAdiv,
806 modulePhi = phiInSupermodule / mNPHIdiv,
807 moduleID = moduleEta * nModulesInSMPhi + modulePhi;
808 int etaInModule = etaInSupermodule % mNETAdiv,
809 phiInModule = phiInSupermodule % mNPHIdiv;
810 // return std::make_tuple(modulePhi, moduleEta, moduleID);
811 return std::make_tuple(phiInModule, etaInModule, moduleID);
812}
813
814Int_t Geometry::GetAbsCellIdFromCellIndexes(Int_t nSupMod, Int_t iphi, Int_t ieta) const
815{
816 // Check if the indeces correspond to existing SM or tower indeces
817 if (iphi < 0 || iphi >= EMCAL_ROWS || ieta < 0 || ieta >= EMCAL_COLS || nSupMod < 0 ||
818 nSupMod >= GetNumberOfSuperModules()) {
819 LOG(debug) << "Wrong cell indexes : SM " << nSupMod << ", column (eta) " << ieta << ", row (phi) " << iphi;
820 return -1;
821 }
822
823 auto indexmod = GetModuleIndexesFromCellIndexesInSModule(nSupMod, iphi, ieta);
824
825 Int_t nIeta = ieta % mNETAdiv,
826 nIphi = iphi % mNPHIdiv;
827 nIeta = mNETAdiv - 1 - nIeta;
828 return GetAbsCellId(nSupMod, std::get<2>(indexmod), nIphi, nIeta);
829}
830
831std::tuple<int, int> Geometry::GlobalRowColFromIndex(int cellID) const
832{
833 if (!CheckAbsCellId(cellID)) {
834 throw InvalidCellIDException(cellID);
835 }
836 auto [supermodule, module, phiInModule, etaInModule] = GetCellIndex(cellID);
837 auto [row, col] = GetCellPhiEtaIndexInSModule(supermodule, module, phiInModule, etaInModule);
838 // add offsets (row / col per supermodule)
839 if (supermodule == 13 || supermodule == 15 || supermodule == 17) {
840 // DCal odd SMs need shift of the col. index in oder to get the global col. index
841 col += 16;
842 }
843 if (supermodule % 2) {
844 col += mNZ * 2;
845 }
846 int sector = supermodule / 2;
847 if (sector > 0) {
848 for (int isec = 0; isec < sector; isec++) {
849 auto smtype = GetSMType(isec * 2);
850 auto nphism = (smtype == EMCAL_THIRD || smtype == DCAL_EXT) ? GetNPhi() / 3 : GetNPhi();
851 row += 2 * nphism;
852 }
853 }
854 return std::make_tuple(row, col);
855}
856
857std::tuple<int, int, int> Geometry::GetPositionInSupermoduleFromGlobalRowCol(int row, int col) const
858{
859 if (col < 0 || col >= 4 * GetNEta()) {
860 throw RowColException(row, col);
861 }
862 int side = col < GetNEta() * 2 ? 0 : 1,
863 colSM = col % (GetNEta() * 2);
864 int sector = -1,
865 rowSM = row;
866 for (int isec = 0; isec < GetNPhiSuperModule(); isec++) {
867 auto smtype = GetSMType(isec * 2);
868 auto nphism = GetNPhi() * 2;
869 if (smtype == EMCAL_THIRD || smtype == DCAL_EXT) {
870 nphism /= 3;
871 }
872 if (rowSM < nphism) {
873 sector = isec;
874 break;
875 }
876 rowSM -= nphism;
877 }
878 if (sector < 0) {
879 throw RowColException(row, col);
880 }
881 int supermodule = sector * 2 + side;
882 if (supermodule == 13 || supermodule == 15 || supermodule == 17) {
883 // DCal odd SMs need shift of the col. index as global col index includes PHOS hole
884 colSM -= 16;
885 if (colSM < 0) {
886 throw RowColException(row, col); // Position inside PHOS hole specified
887 }
888 }
889 if (supermodule == 12 || supermodule == 14 || supermodule == 16) {
890 if (colSM > 32) {
891 throw RowColException(row, col); // Position inside PHOS hole specified
892 }
893 }
894 return std::make_tuple(supermodule, rowSM, colSM);
895}
896
898{
900 return GetAbsCellIdFromCellIndexes(supermodule, rowSM, colSM);
901}
902
903std::tuple<int, int, int, int> Geometry::GetCellIndexFromGlobalRowCol(int row, int col) const
904{
906 auto indexmod = GetModuleIndexesFromCellIndexesInSModule(supermodule, rowSM, colSM);
907
908 Int_t colInModule = colSM % mNETAdiv,
909 rowInMOdule = rowSM % mNPHIdiv;
910 colInModule = mNETAdiv - 1 - colInModule;
911 return std::make_tuple(supermodule, std::get<2>(indexmod), rowInMOdule, colInModule);
912}
913
914int Geometry::GlobalCol(int cellID) const
915{
916 return std::get<1>(GlobalRowColFromIndex(cellID));
917}
918
919int Geometry::GlobalRow(int cellID) const
920{
921 return std::get<0>(GlobalRowColFromIndex(cellID));
922}
923
924Int_t Geometry::SuperModuleNumberFromEtaPhi(Double_t eta, Double_t phi) const
925{
926 if (TMath::Abs(eta) > mEtaMaxOfTRD1) {
927 throw InvalidPositionException(eta, phi);
928 }
929
930 phi = TVector2::Phi_0_2pi(phi); // move phi to (0,2pi) boundaries
931 Int_t nphism = mNumberOfSuperModules / 2;
932 Int_t nSupMod = 0;
933 for (Int_t i = 0; i < nphism; i++) {
934 LOG(debug) << "Sec " << i << ": Min " << mPhiBoundariesOfSM[2 * i] << ", Max " << mPhiBoundariesOfSM[2 * i + 1];
935 if (phi >= mPhiBoundariesOfSM[2 * i] && phi <= mPhiBoundariesOfSM[2 * i + 1]) {
936 nSupMod = 2 * i;
937 if (eta < 0.0) {
938 nSupMod++;
939 }
940
941 if (GetSMType(nSupMod) == DCAL_STANDARD) { // Gap between DCAL
942 if (TMath::Abs(eta) < GetNEta() / 3 * mTrd1Angle * TMath::DegToRad()) {
943 throw InvalidPositionException(eta, phi);
944 }
945 }
946
947 LOG(debug) << "eta " << eta << " phi " << phi << " (" << std::setw(5) << std::setprecision(2)
948 << phi * TMath::RadToDeg() << ") : nSupMod " << nSupMod << ": #bound " << i;
949 return nSupMod;
950 }
951 }
952 throw InvalidPositionException(eta, phi);
953}
954
955Int_t Geometry::GetAbsCellIdFromEtaPhi(Double_t eta, Double_t phi) const
956{
957 Int_t nSupMod = SuperModuleNumberFromEtaPhi(eta, phi);
958
959 // phi index first
960 phi = TVector2::Phi_0_2pi(phi);
961 Double_t phiLoc = phi - mPhiCentersOfSMSec[nSupMod / 2];
962 Int_t nphi = mPhiCentersOfCells.size();
963 switch (GetSMType(nSupMod)) {
964 case EMCAL_HALF:
965 nphi /= 2;
966 case EMCAL_THIRD:
967 case DCAL_EXT:
968 nphi /= 3;
969 break;
970 default:
971 // All other supermodules have full number of cells in phi
972 break;
973 };
974
975 Double_t dmin = TMath::Abs(mPhiCentersOfCells[0] - phiLoc),
976 d = 0.;
977 Int_t iphi = 0;
978 for (Int_t i = 1; i < nphi; i++) {
979 d = TMath::Abs(mPhiCentersOfCells[i] - phiLoc);
980 if (d < dmin) {
981 dmin = d;
982 iphi = i;
983 }
984 // printf(" i %i : d %f : dmin %f : fPhiCentersOfCells[i] %f \n", i, d, dmin, fPhiCentersOfCells[i]);
985 }
986 // odd SM are turned with respect of even SM - reverse indexes
987 LOG(debug2) << " iphi " << iphi << " : dmin " << dmin << " (phi " << phi << ", phiLoc " << phiLoc << ")\n";
988
989 // eta index
990 Double_t absEta = TMath::Abs(eta);
991 Int_t neta = mCentersOfCellsEtaDir.size(),
992 etaShift = iphi * neta,
993 ieta = 0;
994 if (GetSMType(nSupMod) == DCAL_STANDARD) {
995 ieta += 16; // jump 16 cells for DCSM
996 }
997 dmin = TMath::Abs(mEtaCentersOfCells[etaShift + ieta] - absEta);
998 for (Int_t i = ieta + 1; i < neta; i++) {
999 d = TMath::Abs(mEtaCentersOfCells[i + etaShift] - absEta);
1000 if (d < dmin) {
1001 dmin = d;
1002 ieta = i;
1003 }
1004 }
1005
1006 if (GetSMType(nSupMod) == DCAL_STANDARD) {
1007 ieta -= 16; // jump 16 cells for DCSM
1008 }
1009
1010 LOG(debug2) << " ieta " << ieta << " : dmin " << dmin << " (eta=" << eta << ") : nSupMod " << nSupMod;
1011
1012 // patch for mapping following alice convention
1013 if (nSupMod % 2 ==
1014 0) { // 47 + 16 -ieta for DCSM, 47 - ieta for others, revert the ordering on A side in order to keep convention.
1015 ieta = (neta - 1) - ieta;
1016 if (GetSMType(nSupMod) == DCAL_STANDARD) {
1017 ieta -= 16; // recover cells for DCSM
1018 }
1019 }
1020
1021 return GetAbsCellIdFromCellIndexes(nSupMod, iphi, ieta);
1022}
1023
1024std::tuple<int, int, int, int> Geometry::CalculateCellIndex(Int_t absId) const
1025{
1026 if (!CheckAbsCellId(absId)) {
1027 throw InvalidCellIDException(absId);
1028 }
1029
1030 Int_t tmp = absId;
1031 Int_t test = absId;
1032
1033 Int_t nSupMod;
1034 for (nSupMod = -1; test >= 0;) {
1035 nSupMod++;
1036 tmp = test;
1037 switch (GetSMType(nSupMod)) {
1038 case EMCAL_STANDARD:
1040 break;
1041 case EMCAL_HALF:
1042 test -= mNCellsInSupMod / 2;
1043 break;
1044 case DCAL_STANDARD:
1045 test -= 2 * mNCellsInSupMod / 3;
1046 break;
1047 case EMCAL_THIRD:
1048 case DCAL_EXT:
1049 test -= mNCellsInSupMod / 3;
1050 break;
1051 default:
1053 };
1054 }
1055
1056 Int_t nModule = tmp / mNCellsInModule;
1057 tmp = tmp % mNCellsInModule;
1058 Int_t nIphi = tmp / mNPHIdiv, nIeta = tmp % mNETAdiv;
1059 return std::make_tuple(nSupMod, nModule, nIphi, nIeta);
1060}
1061
1062std::tuple<int, int, int, int> Geometry::GetCellIndex(Int_t absId) const
1063{
1064 if (!CheckAbsCellId(absId)) {
1065 throw InvalidCellIDException(absId);
1066 }
1067 return mCellIndexLookup[absId];
1068}
1069
1070Int_t Geometry::GetSuperModuleNumber(Int_t absId) const { return std::get<0>(GetCellIndex(absId)); }
1071
1072std::tuple<int, int> Geometry::GetModulePhiEtaIndexInSModule(int supermoduleID, int moduleID) const
1073{
1074 int nModulesInPhi = -1;
1075 switch (GetSMType(supermoduleID)) {
1076 case EMCAL_HALF:
1077 nModulesInPhi = mNPhi / 2; // halfSM
1078 break;
1079 case EMCAL_THIRD:
1080 case DCAL_EXT:
1081 nModulesInPhi = mNPhi / 3; // 1/3 SM
1082 break;
1083 default:
1084 nModulesInPhi = mNPhi; // full SM
1085 break;
1086 };
1087 return std::make_tuple(int(moduleID % nModulesInPhi), int(moduleID / nModulesInPhi));
1088}
1089
1090std::tuple<int, int> Geometry::GetCellPhiEtaIndexInSModule(int supermoduleID, int moduleID, int phiInModule,
1091 int etaInModule) const
1092{
1093 auto [phiOfModule, etaOfModule] = GetModulePhiEtaIndexInSModule(supermoduleID, moduleID);
1094
1095 // ieta = etaOfModule*fNETAdiv + (1-etaInModule); // x(module) = -z(SM)
1096 int etaInSupermodule = etaOfModule * mNETAdiv + (mNETAdiv - 1 - etaInModule); // x(module) = -z(SM)
1097 int phiInSupermodule = phiOfModule * mNPHIdiv + phiInModule; // y(module) = y(SM)
1098
1099 if (phiInSupermodule < 0 || etaInSupermodule < 0) {
1100 LOG(debug) << " Supermodule " << supermoduleID << ", Module " << moduleID << " (phi " << phiInModule << ", eta " << etaInModule << ")"
1101 << " => in Supermodule: eta " << etaInSupermodule << ", phi " << phiInSupermodule;
1102 }
1103 return std::make_tuple(phiInSupermodule, etaInSupermodule);
1104}
1105
1106std::tuple<short, short> Geometry::GetTopologicalRowColumn(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const
1107{
1108 auto [iphi, ieta] = GetCellPhiEtaIndexInSModule(supermoduleID, moduleID, phiInModule, etaInModule);
1109 int row = iphi;
1110 int column = ieta;
1111
1112 // Add shifts wrt. supermodule and type of calorimeter
1113 // NOTE:
1114 // * Rows (phi) are arranged that one space is left empty between supermodules in phi
1115 // This is due to the physical gap that forbids clustering
1116 // * For DCAL, there is an additional empty column between two supermodules in eta
1117 // Again, this is to account for the gap in DCAL
1118
1119 row += supermoduleID / 2 * (24 + 1);
1120 // In DCAL, leave a gap between two SMs with same phi
1121 if (!IsDCALSM(supermoduleID)) { // EMCAL
1122 column += supermoduleID % 2 * 48;
1123 } else {
1124 column += supermoduleID % 2 * (48 + 1);
1125 }
1126
1127 return std::make_tuple(static_cast<short>(row), static_cast<short>(column));
1128}
1129
1130std::tuple<int, int> Geometry::ShiftOnlineToOfflineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const
1131{
1132 if (supermoduleID == 13 || supermoduleID == 15 || supermoduleID == 17) {
1133 // DCal odd SMs
1134 ieta -= 16; // Same cabling mapping as for EMCal, not considered offline.
1135 } else if (supermoduleID == 18 || supermoduleID == 19) {
1136 // DCal 1/3 SMs
1137 iphi -= 16; // Needed due to cabling mistake.
1138 }
1139 return std::tuple<int, int>(iphi, ieta);
1140}
1141
1142std::tuple<int, int> Geometry::ShiftOfflineToOnlineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const
1143{
1144 if (supermoduleID == 13 || supermoduleID == 15 || supermoduleID == 17) {
1145 // DCal odd SMs
1146 ieta += 16; // Same cabling mapping as for EMCal, not considered offline.
1147 } else if (supermoduleID == 18 || supermoduleID == 19) {
1148 // DCal 1/3 SMs
1149 iphi += 16; // Needed due to cabling mistake.
1150 }
1151 return std::tuple<int, int>(iphi, ieta);
1152}
1153
1155{
1156 // Shift index taking into account the difference between standard SM
1157 // and SM of half (or one third) size in phi direction
1158
1159 Int_t phiindex = mCentersOfCellsPhiDir.size();
1160 Double_t zshift = 0.5 * GetDCALInnerEdge();
1161 Double_t xr, yr, zr;
1162
1163 if (!CheckAbsCellId(absId)) {
1164 throw InvalidCellIDException(absId);
1165 }
1166
1167 auto cellindex = GetCellIndex(absId);
1168 Int_t nSupMod = std::get<0>(cellindex), nModule = std::get<1>(cellindex), nIphi = std::get<2>(cellindex),
1169 nIeta = std::get<3>(cellindex);
1170 auto indexinsm = GetCellPhiEtaIndexInSModule(nSupMod, nModule, nIphi, nIeta);
1171 Int_t iphi = std::get<0>(indexinsm), ieta = std::get<1>(indexinsm);
1172
1173 // Get eta position. Careful with ALICE conventions (increase index decrease eta)
1174 Int_t ieta2 = ieta;
1175 if (nSupMod % 2 == 0) {
1176 ieta2 = (mCentersOfCellsEtaDir.size() - 1) -
1177 ieta; // 47-ieta, revert the ordering on A side in order to keep convention.
1178 }
1179
1180 if (GetSMType(nSupMod) == DCAL_STANDARD && nSupMod % 2) {
1181 ieta2 += 16; // DCAL revert the ordering on C side ...
1182 }
1183 zr = mCentersOfCellsEtaDir[ieta2];
1184 if (GetSMType(nSupMod) == DCAL_STANDARD) {
1185 zr -= zshift; // DCAL shift (SMALLER SM)
1186 }
1187 xr = mCentersOfCellsXDir[ieta2];
1188
1189 // Get phi position. Careful with ALICE conventions (increase index increase phi)
1190 Int_t iphi2 = iphi;
1191 if (GetSMType(nSupMod) == DCAL_EXT) {
1192 if (nSupMod % 2 != 0) {
1193 iphi2 = (phiindex / 3 - 1) - iphi; // 7-iphi [1/3SM], revert the ordering on C side in order to keep convention.
1194 }
1195 yr = mCentersOfCellsPhiDir[iphi2 + phiindex / 3];
1196 } else if (GetSMType(nSupMod) == EMCAL_HALF) {
1197 if (nSupMod % 2 != 0) {
1198 iphi2 = (phiindex / 2 - 1) - iphi; // 11-iphi [1/2SM], revert the ordering on C side in order to keep
1199 }
1200 // convention.
1201 yr = mCentersOfCellsPhiDir[iphi2 + phiindex / 4];
1202 } else if (GetSMType(nSupMod) == EMCAL_THIRD) {
1203 if (nSupMod % 2 != 0) {
1204 iphi2 = (phiindex / 3 - 1) - iphi; // 7-iphi [1/3SM], revert the ordering on C side in order to keep convention.
1205 }
1206 yr = mCentersOfCellsPhiDir[iphi2 + phiindex / 3];
1207 } else {
1208 if (nSupMod % 2 != 0) {
1209 iphi2 = (phiindex - 1) - iphi; // 23-iphi, revert the ordering on C side in order to keep conventi
1210 }
1211 yr = mCentersOfCellsPhiDir[iphi2];
1212 }
1213
1214 LOG(debug) << "absId " << absId << " nSupMod " << nSupMod << " iphi " << iphi << " ieta " << ieta << " xr " << xr
1215 << " yr " << yr << " zr " << zr;
1216 return o2::math_utils::Point3D<double>(xr, yr, zr);
1217}
1218
1220{
1221 // Shift index taking into account the difference between standard SM
1222 // and SM of half (or one third) size in phi direction
1223 Double_t xr, yr, zr;
1224 Int_t nphiIndex = mCentersOfCellsPhiDir.size();
1225 Double_t zshift = 0.5 * GetDCALInnerEdge();
1226 Int_t kDCalshift = 8; // wangml DCal cut first 8 modules(16 cells)
1227
1228 Int_t iphim = -1, ietam = -1;
1229 TVector2 v;
1230 if (!CheckAbsCellId(absId)) {
1231 throw InvalidCellIDException(absId);
1232 }
1233
1234 auto cellindex = GetCellIndex(absId);
1235 Int_t nSupMod = std::get<0>(cellindex), nModule = std::get<1>(cellindex), nIphi = std::get<2>(cellindex),
1236 nIeta = std::get<3>(cellindex);
1237 auto indmodep = GetModulePhiEtaIndexInSModule(nSupMod, nModule);
1238 iphim = std::get<0>(indmodep);
1239 ietam = std::get<1>(indmodep);
1240 auto indexinsm = GetCellPhiEtaIndexInSModule(nSupMod, nModule, nIphi, nIeta);
1241 Int_t iphi = std::get<0>(indexinsm), ieta = std::get<1>(indexinsm);
1242
1243 // Get eta position. Careful with ALICE conventions (increase index decrease eta)
1244 if (nSupMod % 2 == 0) {
1245 ietam = (mCentersOfCellsEtaDir.size() / 2 - 1) -
1246 ietam; // 24-ietam, revert the ordering on A side in order to keep convention.
1247 if (nIeta == 0) {
1248 nIeta = 1;
1249 } else {
1250 nIeta = 0;
1251 }
1252 }
1253
1254 if (GetSMType(nSupMod) == DCAL_STANDARD && nSupMod % 2) {
1255 ietam += kDCalshift; // DCAL revert the ordering on C side ....
1256 }
1257 const ShishKebabTrd1Module& mod = GetShishKebabModule(ietam);
1258 mod.GetPositionAtCenterCellLine(nIeta, distEff, v);
1259 xr = v.Y() - mParSM[0];
1260 zr = v.X() - mParSM[2];
1261 if (GetSMType(nSupMod) == DCAL_STANDARD) {
1262 zr -= zshift; // DCAL shift (SMALLER SM)
1263 }
1264
1265 // Get phi position. Careful with ALICE conventions (increase index increase phi)
1266 Int_t iphi2 = iphi;
1267 if (GetSMType(nSupMod) == DCAL_EXT) {
1268 if (nSupMod % 2 != 0) {
1269 iphi2 = (nphiIndex / 3 - 1) - iphi; // 7-iphi [1/3SM], revert the ordering on C side in order to keep convention.
1270 }
1271 yr = mCentersOfCellsPhiDir[iphi2 + nphiIndex / 3];
1272 } else if (GetSMType(nSupMod) == EMCAL_HALF) {
1273 if (nSupMod % 2 != 0) {
1274 iphi2 = (nphiIndex / 2 - 1) - iphi; // 11-iphi [1/2SM], revert the ordering on C side in order to keep
1275 }
1276 // convention.
1277 yr = mCentersOfCellsPhiDir[iphi2 + nphiIndex / 2];
1278 } else if (GetSMType(nSupMod) == EMCAL_THIRD) {
1279 if (nSupMod % 2 != 0) {
1280 iphi2 = (nphiIndex / 3 - 1) - iphi; // 7-iphi [1/3SM], revert the ordering on C side in order to keep convention.
1281 }
1282 yr = mCentersOfCellsPhiDir[iphi2 + nphiIndex / 3];
1283 } else {
1284 if (nSupMod % 2 != 0) {
1285 iphi2 = (nphiIndex - 1) - iphi; // 23-iphi, revert the ordering on C side in order to keep convention.
1286 }
1287 yr = mCentersOfCellsPhiDir[iphi2];
1288 }
1289
1290 LOG(debug) << "absId " << absId << " nSupMod " << nSupMod << " iphi " << iphi << " ieta " << ieta << " xr " << xr
1291 << " yr " << yr << " zr " << zr;
1292 return math_utils::Point3D<double>(xr, yr, zr);
1293}
1294
1296{
1297 LOG(debug2) << " o2::emcal::Geometry::CreateListOfTrd1Modules() started\n";
1298
1299 if (!mShishKebabTrd1Modules.size()) {
1300 for (int iz = 0; iz < mNZ; iz++) {
1301 if (iz == 0) {
1302 // mod = new AliEMCALShishKebabTrd1Module(TMath::Pi()/2.,this);
1303 mShishKebabTrd1Modules.emplace_back(ShishKebabTrd1Module(TMath::Pi() / 2., this));
1304 } else {
1306 }
1307 }
1308 } else {
1309 LOG(debug2) << " Already exits :\n";
1310 }
1311
1314 LOG(debug2) << " mShishKebabTrd1Modules has " << mShishKebabTrd1Modules.size() << " modules : max eta "
1315 << std::setw(5) << std::setprecision(4) << mEtaMaxOfTRD1;
1316
1317 // define grid for cells in eta(z) and x directions in local coordinates system of SM
1318 // Works just for 2x2 case only -- ?? start here
1319 //
1320 //
1321 // Define grid for cells in phi(y) direction in local coordinates system of SM
1322 // as for 2X2 as for 3X3 - Nov 8,2006
1323 //
1324 LOG(debug2) << " Cells grid in phi directions : size " << mCentersOfCellsPhiDir.size();
1325
1326 Int_t ind = 0; // this is phi index
1327 Int_t ieta = 0, nModule = 0;
1328 Double_t xr = 0., zr = 0., theta = 0., phi = 0., eta = 0., r = 0., x = 0., y = 0.;
1329 TVector3 vglob;
1330 Double_t ytCenterModule = 0.0, ytCenterCell = 0.0;
1331
1334
1335 Double_t r0 = mIPDistance + mLongModuleSize / 2.;
1336 for (Int_t it = 0; it < mNPhi; it++) { // cycle on modules
1337 ytCenterModule = -mParSM[1] + mPhiModuleSize * (2 * it + 1) / 2; // center of module
1338 for (Int_t ic = 0; ic < mNPHIdiv; ic++) { // cycle on cells in module
1339 if (mNPHIdiv == 2) {
1340 ytCenterCell = ytCenterModule + mPhiTileSize * (2 * ic - 1) / 2.;
1341 } else if (mNPHIdiv == 3) {
1342 ytCenterCell = ytCenterModule + mPhiTileSize * (ic - 1);
1343 } else if (mNPHIdiv == 1) {
1344 ytCenterCell = ytCenterModule;
1345 }
1346 mCentersOfCellsPhiDir[ind] = ytCenterCell;
1347 // Define grid on phi direction
1348 // Grid is not the same for different eta bin;
1349 // Effect is small but is still here
1350 phi = TMath::ATan2(ytCenterCell, r0);
1351 mPhiCentersOfCells[ind] = phi;
1352
1353 LOG(debug2) << " ind " << std::setw(2) << std::setprecision(2) << ind << " : y " << std::setw(8)
1354 << std::setprecision(3) << mCentersOfCellsPhiDir[ind];
1355 ind++;
1356 }
1357 }
1358
1362
1363 LOG(debug2) << " Cells grid in eta directions : size " << mCentersOfCellsEtaDir.size();
1364
1365 for (Int_t it = 0; it < mNZ; it++) {
1366 const ShishKebabTrd1Module& trd1 = GetShishKebabModule(it);
1367 nModule = mNPhi * it;
1368 for (Int_t ic = 0; ic < mNETAdiv; ic++) {
1369 if (mNPHIdiv == 2) {
1370 trd1.GetCenterOfCellInLocalCoordinateofSM(ic, xr, zr); // case of 2X2
1371 auto indexinsm = GetCellPhiEtaIndexInSModule(0, nModule, 0, ic);
1372 ieta = std::get<1>(indexinsm);
1373 }
1374 if (mNPHIdiv == 3) {
1375 trd1.GetCenterOfCellInLocalCoordinateofSM3X3(ic, xr, zr); // case of 3X3
1376 auto indexinsm = GetCellPhiEtaIndexInSModule(0, nModule, 0, ic);
1377 ieta = std::get<1>(indexinsm);
1378 }
1379 if (mNPHIdiv == 1) {
1380 trd1.GetCenterOfCellInLocalCoordinateofSM1X1(xr, zr); // case of 1X1
1381 auto indexinsm = GetCellPhiEtaIndexInSModule(0, nModule, 0, ic);
1382 ieta = std::get<1>(indexinsm);
1383 }
1384 mCentersOfCellsXDir[ieta] = float(xr) - mParSM[0];
1385 mCentersOfCellsEtaDir[ieta] = float(zr) - mParSM[2];
1386 // Define grid on eta direction for each bin in phi
1387 for (int iphi = 0; iphi < mCentersOfCellsPhiDir.size(); iphi++) {
1388 x = xr + trd1.GetRadius();
1389 y = mCentersOfCellsPhiDir[iphi];
1390 r = TMath::Sqrt(x * x + y * y + zr * zr);
1391 theta = TMath::ACos(zr / r);
1393 // ind = ieta*fCentersOfCellsPhiDir.GetSize() + iphi;
1394 ind = iphi * mCentersOfCellsEtaDir.size() + ieta;
1395 mEtaCentersOfCells[ind] = eta;
1396 }
1397 // printf(" ieta %i : xr + trd1->GetRadius() %f : zr %f : eta %f \n", ieta, xr + trd1->GetRadius(), zr, eta);
1398 }
1399 }
1400
1401 for (Int_t i = 0; i < mCentersOfCellsEtaDir.size(); i++) {
1402 LOG(debug2) << " ind " << std::setw(2) << std::setprecision(2) << i + 1 << " : z " << std::setw(8)
1403 << std::setprecision(3) << mCentersOfCellsEtaDir[i] << " : x " << std::setw(8)
1404 << std::setprecision(3) << mCentersOfCellsXDir[i];
1405 }
1406}
1407
1409{
1410 if (mShishKebabTrd1Modules.size() && neta >= 0 && neta < mShishKebabTrd1Modules.size()) {
1411 return mShishKebabTrd1Modules.at(neta);
1412 }
1414}
1415
1416Bool_t Geometry::Impact(const TParticle* particle) const
1417{
1418 Bool_t in = kFALSE;
1419 Int_t absID = 0;
1420 math_utils::Point3D<double> vimpact = {0, 0, 0};
1421
1422 ImpactOnEmcal({particle->Vx(), particle->Vy(), particle->Vz()}, particle->Theta(), particle->Phi(), absID, vimpact);
1423
1424 if (absID >= 0) {
1425 in = kTRUE;
1426 }
1427
1428 return in;
1429}
1430
1431void Geometry::ImpactOnEmcal(const math_utils::Point3D<double>& vtx, Double_t theta, Double_t phi, Int_t& absId, math_utils::Point3D<double>& vimpact) const
1432{
1433 math_utils::Vector3D<double> p(TMath::Sin(theta) * TMath::Cos(phi), TMath::Sin(theta) * TMath::Sin(phi), TMath::Cos(theta));
1434
1435 vimpact.SetXYZ(0, 0, 0);
1436 absId = -1;
1437 if (phi == 0 || theta == 0) {
1438 return;
1439 }
1440
1442 Double_t factor = (mIPDistance - vtx.Y()) / p.Y();
1443 direction = vtx + factor * p;
1444
1445 // from particle direction -> tower hitted
1446 absId = GetAbsCellIdFromEtaPhi(direction.Eta(), direction.Phi());
1447
1448 // tower absID hitted -> tower/module plane (evaluated at the center of the tower)
1449
1450 Double_t loc[3], loc2[3], loc3[3];
1451 Double_t glob[3] = {}, glob2[3] = {}, glob3[3] = {};
1452
1453 try {
1454 RelPosCellInSModule(absId).GetCoordinates(loc[0], loc[1], loc[2]);
1455 } catch (InvalidCellIDException& e) {
1456 LOG(error) << e.what();
1457 return;
1458 }
1459
1460 // loc is cell center of tower
1461 auto cellindex = GetCellIndex(absId);
1462 Int_t nSupMod = std::get<0>(cellindex), nModule = std::get<1>(cellindex), nIphi = std::get<2>(cellindex),
1463 nIeta = std::get<3>(cellindex);
1464 // look at 2 neighbours-s cell using nIphi={0,1} and nIeta={0,1}
1465 Int_t nIphi2 = -1, nIeta2 = -1, absId2 = -1, absId3 = -1;
1466 if (nIeta == 0) {
1467 nIeta2 = 1;
1468 } else {
1469 nIeta2 = 0;
1470 }
1471 absId2 = GetAbsCellId(nSupMod, nModule, nIphi, nIeta2);
1472 if (nIphi == 0) {
1473 nIphi2 = 1;
1474 } else {
1475 nIphi2 = 0;
1476 }
1477 absId3 = GetAbsCellId(nSupMod, nModule, nIphi2, nIeta);
1478
1479 // 2nd point on emcal cell plane
1480 try {
1481 RelPosCellInSModule(absId2).GetCoordinates(loc2[0], loc2[1], loc2[2]);
1482 } catch (InvalidCellIDException& e) {
1483 LOG(error) << e.what();
1484 return;
1485 }
1486
1487 // 3rd point on emcal cell plane
1488 try {
1489 RelPosCellInSModule(absId3).GetCoordinates(loc3[0], loc3[1], loc3[2]);
1490 } catch (InvalidCellIDException& e) {
1491 LOG(error) << e.what();
1492 return;
1493 }
1494
1495 // Get Matrix
1496 const TGeoHMatrix* m = GetMatrixForSuperModule(nSupMod);
1497 if (m) {
1498 m->LocalToMaster(loc, glob);
1499 m->LocalToMaster(loc2, glob2);
1500 m->LocalToMaster(loc3, glob3);
1501 } else {
1502 LOG(fatal) << "Geo matrixes are not loaded \n";
1503 }
1504
1505 // Equation of Plane from glob,glob2,glob3 (Ax+By+Cz+D=0)
1506 Double_t a = glob[1] * (glob2[2] - glob3[2]) + glob2[1] * (glob3[2] - glob[2]) + glob3[1] * (glob[2] - glob2[2]);
1507 Double_t b = glob[2] * (glob2[0] - glob3[0]) + glob2[2] * (glob3[0] - glob[0]) + glob3[2] * (glob[0] - glob2[0]);
1508 Double_t c = glob[0] * (glob2[1] - glob3[1]) + glob2[0] * (glob3[1] - glob[1]) + glob3[0] * (glob[1] - glob2[1]);
1509 Double_t d = glob[0] * (glob2[1] * glob3[2] - glob3[1] * glob2[2]) +
1510 glob2[0] * (glob3[1] * glob[2] - glob[1] * glob3[2]) +
1511 glob3[0] * (glob[1] * glob2[2] - glob2[1] * glob[2]);
1512 d = -d;
1513
1514 // shift equation of plane from tower/module center to surface along vector (A,B,C) normal to tower/module plane
1515 Double_t dist = mLongModuleSize / 2.;
1516 Double_t norm = TMath::Sqrt(a * a + b * b + c * c);
1517 Double_t glob4[3] = {};
1519 math_utils::Point3D<double> point = {glob[0], glob[1], glob[2]};
1520 if (point.Dot(dir) < 0) {
1521 dist *= -1;
1522 }
1523 glob4[0] = glob[0] - dist * a / norm;
1524 glob4[1] = glob[1] - dist * b / norm;
1525 glob4[2] = glob[2] - dist * c / norm;
1526 d = glob4[0] * a + glob4[1] * b + glob4[2] * c;
1527 d = -d;
1528
1529 // Line determination (2 points for equation of line : vtx and direction)
1530 // impact between line (particle) and plane (module/tower plane)
1531 Double_t den = a * (vtx.X() - direction.X()) + b * (vtx.Y() - direction.Y()) + c * (vtx.Z() - direction.Z());
1532 if (den == 0) {
1533 LOG(error) << "ImpactOnEmcal() No solution :\n";
1534 return;
1535 }
1536
1537 Double_t length = a * vtx.X() + b * vtx.Y() + c * vtx.Z() + d;
1538 length /= den;
1539
1540 vimpact.SetXYZ(vtx.X() + length * (direction.X() - vtx.X()), vtx.Y() + length * (direction.Y() - vtx.Y()),
1541 vtx.Z() + length * (direction.Z() - vtx.Z()));
1542
1543 // shift vimpact from tower/module surface to center along vector (A,B,C) normal to tower/module plane
1544 vimpact.SetXYZ(vimpact.Z() + dist * a / norm, vimpact.Y() + dist * b / norm, vimpact.Z() + dist * c / norm);
1545}
1546
1548{
1549 if (IsInEMCALOrDCAL(pnt) == EMCAL_ACCEPTANCE) {
1550 return kTRUE;
1551 } else {
1552 return kFALSE;
1553 }
1554}
1555
1557{
1558 if (IsInEMCALOrDCAL(pnt) == DCAL_ACCEPTANCE) {
1559 return kTRUE;
1560 } else {
1561 return kFALSE;
1562 }
1563}
1564
1566{
1567 Double_t r = sqrt(pnt.X() * pnt.X() + pnt.Y() * pnt.Y());
1568
1569 if (r <= mEnvelop[0]) {
1570 return NON_ACCEPTANCE;
1571 } else {
1572 Double_t theta = TMath::ATan2(r, pnt.Z());
1573 Double_t eta;
1574 if (theta == 0) {
1575 eta = 9999;
1576 } else {
1577 eta = -TMath::Log(TMath::Tan(theta / 2.));
1578 }
1579 if (eta < mArm1EtaMin || eta > mArm1EtaMax) {
1580 return NON_ACCEPTANCE;
1581 }
1582
1583 Double_t phi = TMath::ATan2(pnt.Y(), pnt.X()) * 180. / TMath::Pi();
1584 if (phi < 0) {
1585 phi += 360; // phi should go from 0 to 360 in this case
1586 }
1587
1588 if (phi >= mArm1PhiMin && phi <= mEMCALPhiMax) {
1589 return EMCAL_ACCEPTANCE;
1590 } else if (phi >= mDCALPhiMin && phi <= mDCALStandardPhiMax && TMath::Abs(eta) > mDCALInnerExtandedEta) {
1591 return DCAL_ACCEPTANCE;
1592 } else if (phi > mDCALStandardPhiMax && phi <= mDCALPhiMax) {
1593 return DCAL_ACCEPTANCE;
1594 }
1595 return NON_ACCEPTANCE;
1596 }
1597}
1598
1599const TGeoHMatrix* Geometry::GetMatrixForSuperModule(Int_t smod) const
1600{
1601 if (smod < 0 || smod > mNumberOfSuperModules) {
1602 LOG(fatal) << "Wrong supermodule index -> " << smod;
1603 }
1604
1605 if (!SMODULEMATRIX[smod]) {
1606 if (gGeoManager) {
1607 LOG(info) << "Loading EMCAL misalignment matrix for SM " << smod << " from GeoManager.";
1609 } else {
1610 LOG(fatal) << "Cannot find EMCAL misalignment matrices! Recover them either: \n"
1611 << "\t - importing TGeoManager from file geometry.root or \n"
1612 << "\t - from OADB in file OADB/EMCAL/EMCALlocal2master.root or \n"
1613 << "\t - from OCDB in directory OCDB/EMCAL/Align/Data/ or \n"
1614 << "\t - from AliESDs (not in AliAOD) via AliESDRun::GetEMCALMatrix(Int_t superModIndex). \n"
1615 << "Store them via AliEMCALGeometry::SetMisalMatrix(Int_t superModIndex)";
1616 }
1617 }
1618
1619 return SMODULEMATRIX[smod];
1620}
1621
1622const TGeoHMatrix* Geometry::GetMatrixForSuperModuleFromArray(Int_t smod) const
1623{
1624 if (smod < 0 || smod > mNumberOfSuperModules) {
1625 LOG(fatal) << "Wrong supermodule index -> " << smod;
1626 }
1627
1628 return SMODULEMATRIX[smod];
1629}
1630
1631const TGeoHMatrix* Geometry::GetMatrixForSuperModuleFromGeoManager(Int_t smod) const
1632{
1633 const Int_t buffersize = 255;
1634 char path[buffersize];
1635 Int_t tmpType = -1;
1636 Int_t smOrder = 0;
1637
1638 // Get the order for SM
1639 for (Int_t i = 0; i < smod + 1; i++) {
1640 if (GetSMType(i) == tmpType) {
1641 smOrder++;
1642 } else {
1643 tmpType = GetSMType(i);
1644 smOrder = 1;
1645 }
1646 }
1647
1648 Int_t smType = GetSMType(smod);
1649 TString smName = "";
1650
1651 if (smType == EMCAL_STANDARD) {
1652 smName = "SMOD";
1653 } else if (smType == EMCAL_HALF) {
1654 smName = "SM10";
1655 } else if (smType == EMCAL_THIRD) {
1656 smName = "SM3rd";
1657 } else if (smType == DCAL_STANDARD) {
1658 smName = "DCSM";
1659 } else if (smType == DCAL_EXT) {
1660 smName = "DCEXT";
1661 } else {
1662 LOG(error) << "Unkown SM Type!!\n";
1663 }
1664
1665 snprintf(path, buffersize, "/cave/barrel_1/%s_%d", smName.Data(), smOrder);
1666
1667 if (!gGeoManager->cd(path)) {
1668 LOG(fatal) << "Geo manager can not find path " << path << "!\n";
1669 }
1670
1671 return gGeoManager->GetCurrentMatrix();
1672}
1673
1674void Geometry::RecalculateTowerPosition(Float_t drow, Float_t dcol, const Int_t sm, const Float_t depth,
1675 const Float_t misaligTransShifts[15], const Float_t misaligRotShifts[15],
1676 Float_t global[3]) const
1677{
1678 // To use in a print later
1679 Float_t droworg = drow;
1680 Float_t dcolorg = dcol;
1681
1682 if (gGeoManager) {
1683 // Recover some stuff
1684
1685 const Int_t nSMod = mNumberOfSuperModules;
1686
1687 gGeoManager->cd("/cave/barrel_1/");
1688 TGeoNode* geoXEn1 = gGeoManager->GetCurrentNode();
1689 TGeoNodeMatrix* geoSM[nSMod];
1690 TGeoVolume* geoSMVol[nSMod];
1691 TGeoShape* geoSMShape[nSMod];
1692 TGeoBBox* geoBox[nSMod];
1693 TGeoMatrix* geoSMMatrix[nSMod];
1694
1695 for (int iSM = 0; iSM < nSMod; iSM++) {
1696 geoSM[iSM] = dynamic_cast<TGeoNodeMatrix*>(geoXEn1->GetDaughter(iSM));
1697 geoSMVol[iSM] = geoSM[iSM]->GetVolume();
1698 geoSMShape[iSM] = geoSMVol[iSM]->GetShape();
1699 geoBox[iSM] = dynamic_cast<TGeoBBox*>(geoSMShape[iSM]);
1700 geoSMMatrix[iSM] = geoSM[iSM]->GetMatrix();
1701 }
1702
1703 if (sm % 2 == 0) {
1704 dcol = 47. - dcol;
1705 drow = 23. - drow;
1706 }
1707
1708 Int_t istrip = 0;
1709 Float_t z0 = 0;
1710 Float_t zb = 0;
1711 Float_t zIs = 0;
1712
1713 Float_t x, y, z; // return variables in terry's RF
1714
1715 //***********************************************************
1716 // Do not like this: too many hardcoded values, is it not already stored somewhere else?
1717 // : need more comments in the code
1718 //***********************************************************
1719
1720 Float_t dz = 6.0; // base cell width in eta
1721 Float_t dx = 6.004; // base cell width in phi
1722
1723 // Float_t L = 26.04; // active tower length for hadron (lead+scint+paper)
1724 // we use the geant numbers 13.87*2=27.74
1725 Float_t teta1 = 0.;
1726
1727 // Do some basic checks
1728 if (dcol >= 47.5 || dcol < -0.5) {
1729 LOG(error) << "Bad tower coordinate dcol=" << dcol << ", where dcol >= 47.5 || dcol<-0.5; org: " << dcolorg;
1730 return;
1731 }
1732 if (drow >= 23.5 || drow < -0.5) {
1733 LOG(error) << "Bad tower coordinate drow=" << drow << ", where drow >= 23.5 || drow<-0.5; org: " << droworg;
1734 return;
1735 }
1736 if (sm >= nSMod || sm < 0) {
1737 LOG(error) << "Bad SM number sm=" << nSMod << ", where sm >= " << sm << " || sm < 0\n";
1738 return;
1739 }
1740
1741 istrip = int((dcol + 0.5) / 2);
1742
1743 // tapering angle
1744 teta1 = TMath::DegToRad() * istrip * 1.5;
1745
1746 // calculation of module corner along z
1747 // as a function of strip
1748
1749 for (int is = 0; is <= istrip; is++) {
1750 teta1 = TMath::DegToRad() * (is * 1.5 + 0.75);
1751 if (is == 0) {
1752 zIs = zIs + 2 * dz * TMath::Cos(teta1);
1753 } else {
1754 zIs =
1755 zIs + 2 * dz * TMath::Cos(teta1) + 2 * dz * TMath::Sin(teta1) * TMath::Tan(teta1 - 0.75 * TMath::DegToRad());
1756 }
1757 }
1758
1759 z0 = dz * (dcol - 2 * istrip + 0.5);
1760 zb = (2 * dz - z0 - depth * TMath::Tan(teta1));
1761
1762 z = zIs - zb * TMath::Cos(teta1);
1763 y = depth / TMath::Cos(teta1) + zb * TMath::Sin(teta1);
1764
1765 x = (drow + 0.5) * dx;
1766
1767 // moving the origin from terry's RF
1768 // to the GEANT one
1769
1770 double xx = y - geoBox[sm]->GetDX();
1771 double yy = -x + geoBox[sm]->GetDY();
1772 double zz = z - geoBox[sm]->GetDZ();
1773 const double localIn[3] = {xx, yy, zz};
1774 double dglobal[3];
1775 // geoSMMatrix[sm]->Print();
1776 // printf("TFF Local (row = %d, col = %d, x = %3.2f, y = %3.2f, z = %3.2f)\n", iroworg, icolorg, localIn[0],
1777 // localIn[1], localIn[2]);
1778 geoSMMatrix[sm]->LocalToMaster(localIn, dglobal);
1779 // printf("TFF Global (row = %2.0f, col = %2.0f, x = %3.2f, y = %3.2f, z = %3.2f)\n", drow, dcol, dglobal[0],
1780 // dglobal[1], dglobal[2]);
1781
1782 // apply global shifts
1783 if (sm == 2 || sm == 3) { // sector 1
1784 global[0] = dglobal[0] + misaligTransShifts[3] + misaligRotShifts[3] * TMath::Sin(TMath::DegToRad() * 20);
1785 global[1] = dglobal[1] + misaligTransShifts[4] + misaligRotShifts[4] * TMath::Cos(TMath::DegToRad() * 20);
1786 global[2] = dglobal[2] + misaligTransShifts[5];
1787 } else if (sm == 0 || sm == 1) { // sector 0
1788 global[0] = dglobal[0] + misaligTransShifts[0];
1789 global[1] = dglobal[1] + misaligTransShifts[1];
1790 global[2] = dglobal[2] + misaligTransShifts[2];
1791 } else {
1792 LOG(info) << "Careful, correction not implemented yet!\n";
1793 global[0] = dglobal[0];
1794 global[1] = dglobal[1];
1795 global[2] = dglobal[2];
1796 }
1797 } else {
1798 LOG(fatal) << "Geometry boxes information, check that geometry.root is loaded\n";
1799 }
1800}
1801
1802void Geometry::SetMisalMatrix(const TGeoHMatrix* m, Int_t smod) const
1803{
1804 if (smod >= 0 && smod < mNumberOfSuperModules) {
1805 if (!SMODULEMATRIX[smod]) {
1806 SMODULEMATRIX[smod] = new TGeoHMatrix(*m); // Set only if not set yet
1807 }
1808 } else {
1809 LOG(fatal) << "Wrong supermodule index -> " << smod << std::endl;
1810 }
1811}
1812
1813void Geometry::SetMisalMatrixFromCcdb(const char* path, int timestamp) const
1814{
1815 LOG(info) << "Using CCDB to obtain EMCal alignment.";
1817 std::map<std::string, std::string> metadata; // can be empty
1818 api.init("http://alice-ccdb.cern.ch");
1819 TObjArray* matrices = api.retrieveFromTFileAny<TObjArray>(path, metadata, timestamp);
1820
1821 for (int iSM = 0; iSM < mNumberOfSuperModules; ++iSM) {
1822 TGeoHMatrix* mat = reinterpret_cast<TGeoHMatrix*>(matrices->At(iSM));
1823 if (mat) {
1824
1825 SetMisalMatrix(mat, iSM);
1826 } else {
1827 LOG(info) << "Could not obtain Alignment Matrix for SM " << iSM;
1828 }
1829 }
1830}
1831
1832Bool_t Geometry::IsDCALSM(Int_t iSupMod) const
1833{
1834 if (mEMCSMSystem[iSupMod] == DCAL_STANDARD || mEMCSMSystem[iSupMod] == DCAL_EXT) {
1835 return kTRUE;
1836 }
1837
1838 return kFALSE;
1839}
1840
1841Bool_t Geometry::IsDCALExtSM(Int_t iSupMod) const
1842{
1843 if (mEMCSMSystem[iSupMod] == DCAL_EXT) {
1844 return kTRUE;
1845 }
1846
1847 return kFALSE;
1848}
1849
1850Double_t Geometry::GetPhiCenterOfSMSec(Int_t nsupmod) const
1851{
1852 int i = nsupmod / 2;
1853 return mPhiCentersOfSMSec[i];
1854}
1855
1856Double_t Geometry::GetPhiCenterOfSM(Int_t nsupmod) const
1857{
1858 int i = nsupmod / 2;
1859 return mPhiCentersOfSM[i];
1860}
1861
1862std::tuple<double, double> Geometry::GetPhiBoundariesOfSM(Int_t nSupMod) const
1863{
1864 int i;
1865 if (nSupMod < 0 || nSupMod > 12 + mnSupModInDCAL - 1) {
1866 throw InvalidModuleException(nSupMod, 12 + mnSupModInDCAL);
1867 }
1868 i = nSupMod / 2;
1869 return std::make_tuple((Double_t)mPhiBoundariesOfSM[2 * i], (Double_t)mPhiBoundariesOfSM[2 * i + 1]);
1870}
1871
1872std::tuple<double, double> Geometry::GetPhiBoundariesOfSMGap(Int_t nPhiSec) const
1873{
1874 if (nPhiSec < 0 || nPhiSec > 5 + mnSupModInDCAL / 2 - 1) {
1875 throw InvalidModuleException(nPhiSec, 5 + mnSupModInDCAL / 2);
1876 }
1877 return std::make_tuple(mPhiBoundariesOfSM[2 * nPhiSec + 1], mPhiBoundariesOfSM[2 * nPhiSec + 2]);
1878}
1879
1880std::tuple<int, int, int> Geometry::getOnlineID(int towerID)
1881{
1882 auto cellindex = GetCellIndex(towerID);
1883 auto supermoduleID = std::get<0>(cellindex);
1884 auto etaphi = GetCellPhiEtaIndexInSModule(supermoduleID, std::get<1>(cellindex), std::get<2>(cellindex), std::get<3>(cellindex));
1885 auto etaphishift = ShiftOfflineToOnlineCellIndexes(supermoduleID, std::get<0>(etaphi), std::get<1>(etaphi));
1886 int row = std::get<0>(etaphishift), col = std::get<1>(etaphishift);
1887
1888 int ddlInSupermoudel = -1;
1889 if (0 <= row && row < 8) {
1890 ddlInSupermoudel = 0; // first cable row
1891 } else if (8 <= row && row < 16 && 0 <= col && col < 24) {
1892 ddlInSupermoudel = 0; // first half;
1893 } else if (8 <= row && row < 16 && 24 <= col && col < 48) {
1894 ddlInSupermoudel = 1; // second half;
1895 } else if (16 <= row && row < 24) {
1896 ddlInSupermoudel = 1; // third cable row
1897 }
1898 if (supermoduleID % 2 == 1) {
1899 ddlInSupermoudel = 1 - ddlInSupermoudel; // swap for odd=C side, to allow us to cable both sides the same
1900 }
1901
1902 return std::make_tuple(supermoduleID * 2 + ddlInSupermoudel, row, col);
1903}
1904
1905std::tuple<bool, int, int> Geometry::areAbsIDsFromSameTCard(int absId1, int absId2) const
1906{
1907
1908 int rowDiff = -100;
1909 int colDiff = -100;
1910
1911 if (absId1 == absId2) {
1912 return {false, rowDiff, colDiff};
1913 }
1914
1915 // Check if in same SM, if not for sure not same TCard
1916 const int sm1 = GetSuperModuleNumber(absId1);
1917 const int sm2 = GetSuperModuleNumber(absId2);
1918 if (sm1 != sm2) {
1919 return {false, rowDiff, colDiff};
1920 }
1921
1922 // Get the column and row of each absId
1923 const auto [_, iTower1, iIphi1, iIeta1] = GetCellIndex(absId1);
1924 const auto [row1, col1] = GetCellPhiEtaIndexInSModule(sm1, iTower1, iIphi1, iIeta1);
1925
1926 const auto [__, iTower2, iIphi2, iIeta2] = GetCellIndex(absId2);
1927 const auto [row2, col2] = GetCellPhiEtaIndexInSModule(sm2, iTower2, iIphi2, iIeta2);
1928
1929 // Define corner of TCard for absId1
1930 const int tcardRow0 = row1 - row1 % 8;
1931 const int tcardCol0 = col1 - col1 % 2;
1932
1933 // Difference of absId2 from corner of absId1's TCard
1934 const int rowOffset = row2 - tcardRow0;
1935 const int colOffset = col2 - tcardCol0;
1936
1937 // Differences between the two cells directly
1938 rowDiff = row1 - row2;
1939 colDiff = col1 - col2;
1940
1941 const bool sameTCard = (rowOffset >= 0 && rowOffset < 8 &&
1942 colOffset >= 0 && colOffset < 2);
1943 return {sameTCard, rowDiff, colDiff};
1944}
std::ostringstream debug
int32_t i
uint32_t supermodule
Definition RawData.h:3
uint32_t side
Definition RawData.h:0
uint32_t col
Definition RawData.h:4
uint32_t c
Definition RawData.h:2
void init(std::string const &hosts)
Definition CcdbApi.cxx:160
std::enable_if<!std::is_base_of< o2::conf::ConfigurableParam, T >::value, T * >::type retrieveFromTFileAny(std::string const &path, std::map< std::string, std::string > const &metadata, long timestamp=-1, std::map< std::string, std::string > *headers=nullptr, std::string const &etag="", const std::string &createdNotAfter="", const std::string &createdNotBefore="") const
Definition CcdbApi.h:657
Error handling access to non-initialized geometry.
EMCAL geometry definition.
Definition Geometry.h:42
std::tuple< int, int, int > GetModuleIndexesFromCellIndexesInSModule(int supermoduleID, int phiInSupermodule, int etaInSupermodule) const
Transition from cell indexes (iphi, ieta) to module indexes (iphim, ietam, nModule)
Definition Geometry.cxx:801
Float_t mFrontSteelStrip
13-may-05
Definition Geometry.h:695
static Geometry * GetInstanceFromRunNumber(Int_t runNumber, const std::string_view="", const std::string_view mcname="TGeant3", const std::string_view mctitle="")
Instanciate geometry depending on the run number. Mostly used in analysis and MC anchors.
Definition Geometry.cxx:242
Float_t mEtaMaxOfTRD1
Max eta in case of TRD1 geometry (see AliEMCALShishKebabTrd1Module)
Definition Geometry.h:662
Int_t GetSuperModuleNumber(Int_t absId) const
Get cell SM, from absolute ID number.
void RecalculateTowerPosition(Float_t drow, Float_t dcol, const Int_t sm, const Float_t depth, const Float_t misaligTransShifts[15], const Float_t misaligRotShifts[15], Float_t global[3]) const
std::tuple< int, int, int, int > CalculateCellIndex(Int_t absId) const
Calculate cell SM, module numbers from absolute ID number.
Float_t mIPDistance
Radial Distance of the inner surface of the EMCAL.
Definition Geometry.h:676
Float_t mDCALPhiMin
Minimum angular position of DCAL in Phi (degrees)
Definition Geometry.h:663
Float_t mPassiveScintThick
13-may-05
Definition Geometry.h:697
Float_t mPhiTileSize
Size of phi tile.
Definition Geometry.h:673
int GlobalRow(int cellID) const
Get row number of cell in global numbering scheme.
Definition Geometry.cxx:919
Int_t mNZ
Number of Towers in the Z direction.
Definition Geometry.h:675
Float_t mECScintThick
cm, Thickness of the scintillators
Definition Geometry.h:686
const TGeoHMatrix * GetMatrixForSuperModuleFromGeoManager(Int_t smod) const
Provides shift-rotation matrix for EMCAL from the TGeoManager.
std::tuple< int, int > GetCellPhiEtaIndexInSModule(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const
Get eta-phi indexes of cell in SM.
std::tuple< int, int, int, int > GetCellIndex(Int_t absId) const
Get cell SM, module numbers from absolute ID number.
std::tuple< double, double > GetPhiBoundariesOfSM(Int_t nSupMod) const
Int_t mNECLayers
number of scintillator layers
Definition Geometry.h:687
std::vector< Double_t > mPhiCentersOfCells
[fNPhi*fNPHIdiv] from center of SM (-10. < phi < +10.)
Definition Geometry.h:649
const TGeoHMatrix * GetMatrixForSuperModule(Int_t smod) const
Provides shift-rotation matrix for EMCAL from externally set matrix or from TGeoManager.
Geometry & operator=(const Geometry &rvalue)
Assignment operator.
Definition Geometry.cxx:193
Float_t GetDCALInnerEdge() const
Definition Geometry.h:186
std::tuple< int, int, int > getOnlineID(int towerID)
Get link ID, row and column from cell ID, have a look here: https://alice.its.cern....
Int_t GetNEta() const
Get the number of modules in supermodule in #eta direction.
Definition Geometry.h:200
Float_t mEtaModuleSize
Eta -> Y.
Definition Geometry.h:672
void SetMisalMatrixFromCcdb(const char *path="Users/m/mhemmer/EMCAL/Config/GeometryAligned", int timestamp=10000) const
Double_t GetPhiCenterOfSMSec(Int_t nsupmod) const
Float_t GetPhiModuleSize() const
Definition Geometry.h:211
Bool_t CheckAbsCellId(Int_t absId) const
Check whether a cell number is valid.
Definition Geometry.h:726
Float_t mZLength
Total length in z direction.
Definition Geometry.h:681
Float_t mEnvelop[3]
The GEANT TUB for the detector.
Definition Geometry.h:657
Float_t mDCALInnerEdge
Inner edge for DCAL.
Definition Geometry.h:668
void CreateListOfTrd1Modules()
Float_t mParSM[3]
SM sizes as in GEANT (TRD1)
Definition Geometry.h:670
Float_t GetShellThickness() const
Definition Geometry.h:184
Int_t mNCellsInModule
Number cell in module.
Definition Geometry.h:643
std::vector< EMCALSMType > mEMCSMSystem
geometry structure
Definition Geometry.h:693
Bool_t IsDCALExtSM(Int_t nSupMod) const
Check if iSupMod is a valid DCal 1/3rd SM.
std::tuple< int, int, int, int > GetCellIndexFromGlobalRowCol(int row, int col) const
Get the cell indices from global position in the EMCAL.
Definition Geometry.cxx:903
Double_t GetPhiCenterOfSM(Int_t nsupmod) const
const TGeoHMatrix * SMODULEMATRIX[EMCAL_MODULES]
Orientations of EMCAL super modules.
Definition Geometry.h:719
std::tuple< int, int > GetModulePhiEtaIndexInSModule(int supermoduleID, int moduleID) const
Get eta-phi indexes of module in SM.
std::tuple< int, int, int > GetPositionInSupermoduleFromGlobalRowCol(int row, int col) const
Get the posision (row, col) of a global row-col position.
Definition Geometry.cxx:857
const TGeoHMatrix * GetMatrixForSuperModuleFromArray(Int_t smod) const
Provides shift-rotation matrix for EMCAL from fkSModuleMatrix[smod].
Float_t mTrd1AlFrontThick
Thickness of the Al front plate.
Definition Geometry.h:708
Float_t mSampling
Sampling factor.
Definition Geometry.h:682
Float_t mDCALPhiMax
Maximum angular position of DCAL in Phi (degrees)
Definition Geometry.h:664
Float_t mTrd1BondPaperThick
Thickness of the Bond Paper sheet.
Definition Geometry.h:709
Int_t mNCells
Number of cells in calo.
Definition Geometry.h:654
std::tuple< bool, int, int > areAbsIDsFromSameTCard(int absId1, int absId2) const
Check if 2 cells belong to the same T-Card.
Float_t mArm1PhiMin
Minimum angular position of EMCAL in Phi (degrees)
Definition Geometry.h:660
void SetMisalMatrix(const TGeoHMatrix *m, Int_t smod) const
Int_t mNETAdiv
Number eta division of module.
Definition Geometry.h:641
std::tuple< short, short > GetTopologicalRowColumn(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const
Get topological row and column of cell in SM (same as for clusteriser with artifical gaps)
Float_t mLateralSteelStrip
13-may-05
Definition Geometry.h:696
std::tuple< int, int > GlobalRowColFromIndex(int cellID) const
get (Column,Row) pair of cell in global numbering scheme
Definition Geometry.cxx:831
Float_t mTrd1Angle
angle in x-z plane (in degree)
Definition Geometry.h:703
Float_t mECPbRadThickness
cm, Thickness of the Pb radiators
Definition Geometry.h:685
void GetGlobal(const Double_t *loc, Double_t *glob, int ind) const
Figure out the global coordinates from local coordinates on a supermodule.
Definition Geometry.cxx:704
std::vector< Double_t > mEtaCentersOfCells
[fNEta*fNETAdiv*fNPhi*fNPHIdiv], positive direction (eta>0); eta depend from phi position;
Definition Geometry.h:653
int GetAbsCellId(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const
Get cell absolute ID number from location module (2 times 2 cells) of a super module.
Definition Geometry.cxx:761
Float_t mArm1PhiMax
Maximum angular position of EMCAL in Phi (degrees)
Definition Geometry.h:661
static Geometry * GetInstance()
Get geometry instance. It should have been set before.
Definition Geometry.cxx:213
std::tuple< double, double > EtaPhiFromIndex(Int_t absId) const
Figure out the eta/phi coordinates of a cell.
Definition Geometry.cxx:754
Int_t mKey110DEG
For calculation abs cell id; 19-oct-05.
Definition Geometry.h:638
std::vector< Double_t > mPhiCentersOfSMSec
Phi of centers of section where SM lies; size is fNumberOfSuperModules/2.
Definition Geometry.h:646
std::vector< Double_t > mCentersOfCellsEtaDir
Size fNEta*fNETAdiv (for TRD1 only) (eta or z in SM, in cm)
Definition Geometry.h:650
std::string mGeoName
Geometry name string.
Definition Geometry.h:637
Int_t GetNPhiSuperModule() const
Definition Geometry.h:219
int GlobalCol(int cellID) const
Get column number of cell in global numbering scheme.
Definition Geometry.cxx:914
Float_t mShellThickness
Total thickness in (x,y) direction.
Definition Geometry.h:680
EMCALSMType GetSMType(Int_t nSupMod) const
Definition Geometry.h:244
std::vector< Double_t > mCentersOfCellsPhiDir
Size fNPhi*fNPHIdiv (for TRD1 only) (phi or y in SM, in cm)
Definition Geometry.h:651
Int_t mNPHIdiv
Number phi division of module.
Definition Geometry.h:642
Int_t mnSupModInDCAL
For calculation abs cell id; 06-nov-12.
Definition Geometry.h:639
Float_t mEtaTileSize
Size of eta tile.
Definition Geometry.h:674
int SuperModuleNumberFromEtaPhi(Double_t eta, Double_t phi) const
Given a global eta/phi point check if it belongs to a supermodule covered region.
Definition Geometry.cxx:924
void DefineEMC(std::string_view mcname, std::string_view mctitle)
Init function of previous class EMCGeometry.
Definition Geometry.cxx:366
Bool_t IsDCALSM(Int_t nSupMod) const
Check if iSupMod is a valid DCal standard SM.
Bool_t Impact(const TParticle *particle) const
Check if particle falls in the EMCal/DCal geometry.
std::vector< ShishKebabTrd1Module > mShishKebabTrd1Modules
List of modules.
Definition Geometry.h:669
Bool_t IsInEMCAL(const math_utils::Point3D< double > &pnt) const
Checks whether point is inside the EMCal volume.
math_utils::Point3D< double > RelPosCellInSModule(Int_t absId, Double_t distEf) const
Look to see what the relative position inside a given cell is for a recpoint.
const std::string & GetName() const
Definition Geometry.h:111
Float_t mEMCALPhiMax
Maximum angular position of EMCAL in Phi (degrees)
Definition Geometry.h:665
Bool_t IsInDCAL(const math_utils::Point3D< double > &pnt) const
Checks whether point is inside the DCal volume.
Int_t mNPhiSuperModule
9 - number supermodule in phi direction
Definition Geometry.h:700
Float_t mPhiModuleSize
Phi -> X.
Definition Geometry.h:671
Float_t mPhiGapForSM
Gap betweeen supermodules in phi direction.
Definition Geometry.h:705
Geometry()=default
Default constructor. It must be kept public for root persistency purposes, but should never be called...
Float_t m2Trd1Dx2
2*dx2 for TRD1
Definition Geometry.h:704
AcceptanceType_t IsInEMCALOrDCAL(const math_utils::Point3D< double > &pnt) const
Checks whether point is inside the EMCal volume (included DCal)
Int_t GetNumberOfSuperModules() const
Definition Geometry.h:209
~Geometry()
Destructor.
Definition Geometry.cxx:199
Int_t GetNumberOfModuleInPhiDirection(Int_t nSupMod) const
Definition Geometry.h:472
std::tuple< int, int > ShiftOnlineToOfflineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const
Adapt cell indices in supermodule to online indexing.
Float_t mPhiSuperModule
Phi of normal supermodule (20, in degree)
Definition Geometry.h:699
Float_t mArm1EtaMin
Minimum pseudorapidity position of EMCAL in Eta.
Definition Geometry.h:658
void ImpactOnEmcal(const math_utils::Point3D< double > &vtx, Double_t theta, Double_t phi, Int_t &absId, math_utils::Point3D< double > &vimpact) const
Get the impact coordinates on EMCAL.
Int_t GetNPhi() const
Get the number of modules in supermodule in #phi direction.
Definition Geometry.h:204
Float_t mArm1EtaMax
Maximum pseudorapidity position of EMCAL in Eta.
Definition Geometry.h:659
Int_t GetAbsCellIdFromCellIndexes(Int_t nSupMod, Int_t iphi, Int_t ieta) const
Transition from super module number (nSupMod) and cell indexes (ieta,iphi) to cell absolute ID number...
Definition Geometry.cxx:814
std::vector< Double_t > mCentersOfCellsXDir
Size fNEta*fNETAdiv (for TRD1 only) ( x in SM, in cm)
Definition Geometry.h:656
Int_t mNPhi
Number of Towers in the PHI direction.
Definition Geometry.h:655
const ShishKebabTrd1Module & GetShishKebabModule(Int_t neta) const
Get the Module parameters for a eta.
Int_t mNumberOfSuperModules
default is 12 = 6 * 2
Definition Geometry.h:690
std::tuple< double, double > GetPhiBoundariesOfSMGap(Int_t nPhiSec) const
std::vector< Double_t > mPhiBoundariesOfSM
Phi boundaries of SM in rad; size is fNumberOfSuperModules;.
Definition Geometry.h:644
std::vector< std::tuple< int, int, int, int > > mCellIndexLookup
Lookup table for cell indices.
Definition Geometry.h:720
void DefineSamplingFraction(const std::string_view mcname="", const std::string_view mctitle="")
Set the value of the Sampling used to calibrate the MC hits energy (check)
Definition Geometry.cxx:308
Float_t mDCALStandardPhiMax
Special edge for the case that DCAL contian extension.
Definition Geometry.h:666
Float_t mLongModuleSize
Size of long module.
Definition Geometry.h:677
int GetAbsCellIdFromEtaPhi(Double_t eta, Double_t phi) const
Get cell absolute ID number from eta and phi location.
Definition Geometry.cxx:955
Int_t mNCellsInSupMod
Number cell in super module.
Definition Geometry.h:640
int GetCellAbsIDFromGlobalRowCol(int row, int col) const
Get the absolute cell ID from global position in the EMCAL.
Definition Geometry.cxx:897
std::tuple< int, int > ShiftOfflineToOnlineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const
Adapt cell indices in supermodule to offline indexing.
std::vector< Double_t > mPhiCentersOfSM
Phi of centers of SM; size is fNumberOfSuperModules/2.
Definition Geometry.h:645
Float_t mDCALInnerExtandedEta
DCAL inner edge in Eta (with some extension)
Definition Geometry.h:667
Exception handling non-existing cell IDs.
const char * what() const noexcept final
Access to error message of the exception.
Error Handling when an invalid module ID (outside the limits) is called.
Exception handling errors due to positions not in the EMCAL area.
Exception handling improper or uninitialized supermodule types.
Handling error for invalid positions in row-column space.
Main class for TRD1 geometry of Shish-Kebab case.
void GetCenterOfCellInLocalCoordinateofSM3X3(Int_t ieta, Double_t &xr, Double_t &zr) const
void GetCenterOfCellInLocalCoordinateofSM1X1(Double_t &xr, Double_t &zr) const
static Double_t ThetaToEta(Double_t theta)
void GetPositionAtCenterCellLine(Int_t ieta, Double_t dist, TVector2 &v) const
const TVector2 & GetCenterOfCellInLocalCoordinateofSM(Int_t ieta) const
GLint GLenum GLint x
Definition glcorearb.h:403
const GLfloat * m
Definition glcorearb.h:4066
const GLdouble * v
Definition glcorearb.h:832
GLuint const GLchar * name
Definition glcorearb.h:781
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLint y
Definition glcorearb.h:270
GLuint GLsizei GLsizei * length
Definition glcorearb.h:790
GLint GLint GLsizei GLsizei GLsizei depth
Definition glcorearb.h:470
GLsizei const GLchar *const * path
Definition glcorearb.h:3591
GLboolean r
Definition glcorearb.h:1233
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
const std::string DEFAULT_GEOMETRY
@ EMCAL_MODULES
Number of modules, 12 for EMCal + 8 for DCAL.
Definition Constants.h:24
@ EMCAL_ROWS
Number of rows per module for EMCAL.
Definition Constants.h:25
@ EMCAL_COLS
Number of columns per module for EMCAL.
Definition Constants.h:26
FIXME: do not use data model tables.
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< int > row