Project
Loading...
Searching...
No Matches
Detector.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
14
16
17#include "DetectorsBase/Stack.h"
20
24
25// FairRoot includes
26#include "FairDetector.h" // for FairDetector
27#include "FairRootManager.h" // for FairRootManager
28#include "FairRootManager.h"
29#include "FairRun.h" // for FairRun
30#include "FairRuntimeDb.h" // for FairRuntimeDb
31#include "FairVolume.h" // for FairVolume
32
33#include "TGeoManager.h" // for TGeoManager, gGeoManager
34#include "TGeoPcon.h" // for TGeoPcon
35#include "TGeoTube.h" // for TGeoTube
36#include "TGeoVolume.h" // for TGeoVolume, TGeoVolumeAssembly
37#include "TString.h" // for TString, operator+
38#include "TVirtualMC.h" // for gMC, TVirtualMC
39#include "TVirtualMCStack.h" // for TVirtualMCStack
40
41#include <fairlogger/Logger.h> // for LOG, LOG_IF
42
43#include <cstdio> // for NULL, snprintf
44
45#define MAX_SENSORS 2000
46
47class FairModule;
48
49class TGeoMedium;
50
51class TParticle;
52
53using namespace o2::ft3;
54using o2::itsmft::Hit;
55
56//_________________________________________________________________________________________________
58 : o2::base::DetImpl<Detector>("FT3", kTRUE),
59 mTrackData(),
60 mHits(o2::utils::createSimVector<o2::itsmft::Hit>())
61{
62}
63
64//_________________________________________________________________________________________________
66{
67 // Build a basic parametrized FT3 detector with nLayers equally spaced between z_first and z_first+z_length
68 // Covering pseudo rapidity [etaIn,etaOut]. Silicon thinkness computed to match layer x/X0
69
70 LOG(info) << "Building FT3 Detector: Conical Telescope";
71
72 const int numberOfLayers = param.nLayers;
73 const auto z_first = param.z0;
74 const auto z_length = param.zLength;
75 const auto etaIn = param.etaIn;
76 const auto etaOut = param.etaOut;
77 const auto Layerx2X0 = param.Layerx2X0;
78 mLayerName[IdxBackwardDisks].resize(numberOfLayers);
79 mLayerName[IdxForwardDisks].resize(numberOfLayers);
80
81 for (int direction : {IdxBackwardDisks, IdxForwardDisks}) {
82 for (int layerNumber = 0; layerNumber < numberOfLayers; layerNumber++) {
83 std::string layerName = GeometryTGeo::getFT3LayerPattern() + std::to_string(layerNumber + numberOfLayers * direction);
84 mLayerName[direction][layerNumber] = layerName;
85
86 // Adds evenly spaced layers
87 const float layerZ = z_first + (layerNumber * z_length / numberOfLayers) * std::copysign(1, z_first);
88 const float rIn = std::abs(layerZ * std::tan(2.f * std::atan(std::exp(-etaIn))));
89 const float rOut = std::abs(layerZ * std::tan(2.f * std::atan(std::exp(-etaOut))));
90 const bool isMiddleLayer = layerNumber < 3;
91 auto& thisLayer = mLayers[direction].emplace_back(direction, layerNumber, layerName, layerZ, rIn, rOut, Layerx2X0, isMiddleLayer);
92 }
93 }
94}
95
96//_________________________________________________________________________________________________
98{
99 // Build FT3 detector according to
100 // https://indico.cern.ch/event/992488/contributions/4174473/attachments/2168881/3661331/tracker_parameters_werner_jan_11_2021.pdf
101
102 LOG(info) << "Building FT3 Detector: V1";
103
104 const int numberOfLayers = 10;
105 const float sensorThickness = 30.e-4;
106 const float layersx2X0 = 1.e-2;
107 const std::vector<std::array<float, 4>> layersConfig{
108 {26., .5, 3., 0.1f * layersx2X0}, // {z_layer, r_in, r_out, Layerx2X0}
109 {30., .5, 3., 0.1f * layersx2X0},
110 {34., .5, 3., 0.1f * layersx2X0},
111 {77., 3.5, 35., layersx2X0},
112 {100., 3.5, 35., layersx2X0},
113 {122., 3.5, 35., layersx2X0},
114 {150., 3.5, 80.f, layersx2X0},
115 {180., 3.5, 80.f, layersx2X0},
116 {220., 3.5, 80.f, layersx2X0},
117 {279., 3.5, 80.f, layersx2X0}};
118
119 mLayerName[IdxBackwardDisks].resize(numberOfLayers);
120 mLayerName[IdxForwardDisks].resize(numberOfLayers);
121
122 for (auto direction : {IdxBackwardDisks, IdxForwardDisks}) {
123 for (int layerNumber = 0; layerNumber < numberOfLayers; layerNumber++) {
124 std::string directionName = std::to_string(direction);
125 std::string layerName = GeometryTGeo::getFT3LayerPattern() + directionName + std::string("_") + std::to_string(layerNumber);
126 mLayerName[direction][layerNumber] = layerName;
127 auto& z = layersConfig[layerNumber][0];
128
129 auto& rIn = layersConfig[layerNumber][1];
130 auto& rOut = layersConfig[layerNumber][2];
131 auto& x0 = layersConfig[layerNumber][3];
132
133 LOG(info) << "Adding Layer " << layerName << " at z = " << z;
134 // Add layers
135 const bool isMiddleLayer = layerNumber < 3;
136 auto& thisLayer = mLayers[direction].emplace_back(direction, layerNumber, layerName, z, rIn, rOut, x0, isMiddleLayer);
137 }
138 }
139}
140
141//_________________________________________________________________________________________________
143{
144 // Build FT3 detector according to
145 // https://www.overleaf.com/project/6051acc870e39aaeb4653621
146
147 LOG(info) << "Building FT3 Detector: V3b";
148
149 const int numberOfLayers = 12;
150 float sensorThickness = 30.e-4;
151 float layersx2X0 = 1.e-2;
152 std::vector<std::array<float, 4>> layersConfig{
153 {26., .5, 3., 0.1f * layersx2X0}, // {z_layer, r_in, r_out, Layerx2X0}
154 {30., .5, 3., 0.1f * layersx2X0},
155 {34., .5, 3., 0.1f * layersx2X0},
156 {77., 5.0, 35., layersx2X0},
157 {100., 5.0, 35., layersx2X0},
158 {122., 5.0, 35., layersx2X0},
159 {150., 5.5, 80.f, layersx2X0},
160 {180., 6.6, 80.f, layersx2X0},
161 {220., 8.1, 80.f, layersx2X0},
162 {279., 10.2, 80.f, layersx2X0},
163 {340., 12.5, 80.f, layersx2X0},
164 {400., 14.7, 80.f, layersx2X0}};
165
166 mLayerName[IdxBackwardDisks].resize(numberOfLayers);
167 mLayerName[IdxForwardDisks].resize(numberOfLayers);
168
169 for (auto direction : {IdxBackwardDisks, IdxForwardDisks}) {
170 for (int layerNumber = 0; layerNumber < numberOfLayers; layerNumber++) {
171 std::string directionName = std::to_string(direction);
172 std::string layerName = GeometryTGeo::getFT3LayerPattern() + directionName + std::string("_") + std::to_string(layerNumber);
173 mLayerName[direction][layerNumber] = layerName;
174 auto& z = layersConfig[layerNumber][0];
175
176 auto& rIn = layersConfig[layerNumber][1];
177 auto& rOut = layersConfig[layerNumber][2];
178 auto& x0 = layersConfig[layerNumber][3];
179
180 LOG(info) << "Adding Layer " << layerName << " at z = " << z;
181 // Add layers
182 const bool isMiddleLayer = layerNumber < 3;
183 auto& thisLayer = mLayers[direction].emplace_back(direction, layerNumber, layerName, z, rIn, rOut, x0, isMiddleLayer);
184 }
185 }
186}
187
189{
190 // Build the FT3 detector according to changes proposed during
191 // https://indico.cern.ch/event/1407704/
192 // to adhere to the changes that were presented at the ALICE 3 Upgrade days in March 2024
193 // Inner radius at C-side to 7 cm
194 // Inner radius at A-side stays at 5 cm
195 // 06.02.2025 update: IRIS layers are now in TRK
196
197 LOG(info) << "Building FT3 Detector: After Upgrade Days March 2024 version";
198
199 const int numberOfLayers = 9;
200 const float sensorThickness = 30.e-4;
201 const float layersx2X0 = 1.e-2;
202 const std::vector<std::array<float, 4>> layersConfigCSide{
203 {77., 7.0, 35., layersx2X0}, // {z_layer, r_in, r_out, Layerx2X0}
204 {100., 7.0, 35., layersx2X0},
205 {122., 7.0, 35., layersx2X0},
206 {150., 7.0, 68.f, layersx2X0},
207 {180., 7.0, 68.f, layersx2X0},
208 {220., 7.0, 68.f, layersx2X0},
209 {260., 7.0, 68.f, layersx2X0},
210 {300., 7.0, 68.f, layersx2X0},
211 {350., 7.0, 68.f, layersx2X0}};
212
213 const std::vector<std::array<float, 4>> layersConfigASide{
214 {77., 5.0, 35., layersx2X0}, // {z_layer, r_in, r_out, Layerx2X0}
215 {100., 5.0, 35., layersx2X0},
216 {122., 5.0, 35., layersx2X0},
217 {150., 5.0, 68.f, layersx2X0},
218 {180., 5.0, 68.f, layersx2X0},
219 {220., 5.0, 68.f, layersx2X0},
220 {260., 5.0, 68.f, layersx2X0},
221 {300., 5.0, 68.f, layersx2X0},
222 {350., 5.0, 68.f, layersx2X0}};
223
224 mLayerName[IdxBackwardDisks].resize(numberOfLayers);
225 mLayerName[IdxForwardDisks].resize(numberOfLayers);
226
227 for (auto direction : {IdxBackwardDisks, IdxForwardDisks}) {
228 for (int layerNumber = 0; layerNumber < numberOfLayers; layerNumber++) {
229 std::string directionName = std::to_string(direction);
230 std::string layerName = GeometryTGeo::getFT3LayerPattern() + directionName + std::string("_") + std::to_string(layerNumber);
231 mLayerName[direction][layerNumber] = layerName;
232 float z, rIn, rOut, x0;
233 if (direction == 0) { // C-Side
234 z = layersConfigCSide[layerNumber][0];
235 rIn = layersConfigCSide[layerNumber][1];
236 rOut = layersConfigCSide[layerNumber][2];
237 x0 = layersConfigCSide[layerNumber][3];
238 } else if (direction == 1) { // A-Side
239 z = layersConfigASide[layerNumber][0];
240 rIn = layersConfigASide[layerNumber][1];
241 rOut = layersConfigASide[layerNumber][2];
242 x0 = layersConfigASide[layerNumber][3];
243 }
244
245 LOG(info) << "Adding Layer " << layerName << " at z = " << z;
246 // Add layers
247 const bool isMiddleLayer = layerNumber < 3;
248 auto& thisLayer = mLayers[direction].emplace_back(direction, layerNumber, layerName, z, rIn, rOut, x0, isMiddleLayer);
249 }
250 }
251}
252
254{
255 // Build the FT3 detector according to v3 layout
256 // https://indico.cern.ch/event/1596309/contributions/6728167/attachments/3190117/5677220/2025-12-10-AW-ALICE3planning.pdf
257 // Middle disks inner radius 10 cm
258 // Outer disks inner radius 20 cm
259
260 LOG(info) << "Building FT3 Detector: v3 scoping version";
261
262 const int numberOfLayers = 6;
263 const float sensorThickness = 30.e-4;
264 const float layersx2X0 = 1.e-2;
265 using LayerConfig = std::array<float, 4>; // {z_layer, r_in, r_out, Layerx2X0}
266 const std::array<LayerConfig, numberOfLayers> layersConfigCSide{LayerConfig{77., 10.0, 35., layersx2X0},
267 LayerConfig{100., 10.0, 35., layersx2X0},
268 LayerConfig{122., 10.0, 35., layersx2X0},
269 LayerConfig{150., 20.0, 68.f, layersx2X0},
270 LayerConfig{180., 20.0, 68.f, layersx2X0},
271 LayerConfig{220., 20.0, 68.f, layersx2X0}};
272
273 const std::array<LayerConfig, numberOfLayers> layersConfigASide{LayerConfig{77., 10.0, 35., layersx2X0},
274 LayerConfig{100., 10.0, 35., layersx2X0},
275 LayerConfig{122., 10.0, 35., layersx2X0},
276 LayerConfig{150., 20.0, 68.f, layersx2X0},
277 LayerConfig{180., 20.0, 68.f, layersx2X0},
278 LayerConfig{220., 20.0, 68.f, layersx2X0}};
279 const std::array<bool, numberOfLayers> enabled{true, true, true, true, true, true}; // To enable or disable layers for debug purpose
280
281 for (int direction : {IdxBackwardDisks, IdxForwardDisks}) {
282 mLayerName[direction].clear();
283 const std::array<LayerConfig, numberOfLayers>& layerConfig = (direction == IdxBackwardDisks) ? layersConfigCSide : layersConfigASide;
284 for (int layerNumber = 0; layerNumber < numberOfLayers; layerNumber++) {
285 if (!enabled[layerNumber]) {
286 continue;
287 }
288 const std::string directionName = std::to_string(direction);
289 const std::string layerName = GeometryTGeo::getFT3LayerPattern() + directionName + std::string("_") + std::to_string(layerNumber);
290 mLayerName[direction].push_back(layerName.c_str());
291 const float z = layerConfig[layerNumber][0];
292 const float rIn = layerConfig[layerNumber][1];
293 const float rOut = layerConfig[layerNumber][2];
294 const float x0 = layerConfig[layerNumber][3];
295 LOG(info) << "buildFT3ScopingV3 -> Adding Layer " << layerNumber << "/" << numberOfLayers << " " << layerName << " at z = " << z;
296 // Add layers
297 const bool isMiddleLayer = layerNumber < 3;
298 auto& thisLayer = mLayers[direction].emplace_back(direction, layerNumber, layerName, z, rIn, rOut, x0, isMiddleLayer);
299 }
300 }
301}
302
303//_________________________________________________________________________________________________
305{
306 // Build FT3 detector according to the scoping document
307
308 LOG(info) << "Building FT3 Detector: Scoping document version";
309
310 const int numberOfLayers = 12;
311 const float sensorThickness = 30.e-4;
312 const float layersx2X0 = 1.e-2;
313 const std::vector<std::array<float, 4>> layersConfig{
314 {26., .5, 2.5, 0.1f * layersx2X0}, // {z_layer, r_in, r_out, Layerx2X0}
315 {30., .5, 2.5, 0.1f * layersx2X0},
316 {34., .5, 2.5, 0.1f * layersx2X0},
317 {77., 5.0, 35., layersx2X0},
318 {100., 5.0, 35., layersx2X0},
319 {122., 5.0, 35., layersx2X0},
320 {150., 5.0, 68.f, layersx2X0},
321 {180., 5.0, 68.f, layersx2X0},
322 {220., 5.0, 68.f, layersx2X0},
323 {260., 5.0, 68.f, layersx2X0},
324 {300., 5.0, 68.f, layersx2X0},
325 {350., 5.0, 68.f, layersx2X0}};
326
327 mLayerName[IdxBackwardDisks].resize(numberOfLayers);
328 mLayerName[IdxForwardDisks].resize(numberOfLayers);
329
330 for (auto direction : {IdxBackwardDisks, IdxForwardDisks}) {
331 for (int layerNumber = 0; layerNumber < numberOfLayers; layerNumber++) {
332 std::string directionName = std::to_string(direction);
333 std::string layerName = GeometryTGeo::getFT3LayerPattern() + directionName + std::string("_") + std::to_string(layerNumber);
334 mLayerName[direction][layerNumber] = layerName;
335 auto& z = layersConfig[layerNumber][0];
336 auto& rIn = layersConfig[layerNumber][1];
337 auto& rOut = layersConfig[layerNumber][2];
338 auto& x0 = layersConfig[layerNumber][3];
339
340 LOG(info) << "Adding Layer " << layerName << " at z = " << z;
341 // Add layers
342 const bool isMiddleLayer = layerNumber < 3;
343 auto& thisLayer = mLayers[direction].emplace_back(direction, layerNumber, layerName, z, rIn, rOut, x0, isMiddleLayer);
344 }
345 }
346}
347
348//_________________________________________________________________________________________________
350 : o2::base::DetImpl<Detector>("FT3", active),
351 mTrackData(),
352 mHits(o2::utils::createSimVector<o2::itsmft::Hit>())
353{
354 buildFT3ScopingV3(); // v3 Dec 25
355}
356
357//_________________________________________________________________________________________________
359 : o2::base::DetImpl<Detector>(rhs),
360 mTrackData(),
361
363 mHits(o2::utils::createSimVector<o2::itsmft::Hit>())
364{
365 mLayerName = rhs.mLayerName;
366 mActiveSensorMap = rhs.mActiveSensorMap;
367}
368
369//_________________________________________________________________________________________________
371{
372
373 if (mHits) {
374 // delete mHits;
376 }
377}
378
379//_________________________________________________________________________________________________
381{
382 // The standard = operator
383 // Inputs:
384 // Detector &h the sourse of this copy
385 // Outputs:
386 // none.
387 // Return:
388 // A copy of the sourse hit h
389
390 if (this == &rhs) {
391 return *this;
392 }
393
394 // base class assignment
396
397 mLayerName = rhs.mLayerName;
398 mActiveSensorMap = rhs.mActiveSensorMap;
399 mLayers = rhs.mLayers;
400 mTrackData = rhs.mTrackData;
401
403 mHits = nullptr;
404
405 return *this;
406}
407
408//_________________________________________________________________________________________________
410{
411 // Define the list of sensitive volumes
412 LOG(info) << "Initialize FT3 O2Detector";
413
414 defineSensitiveVolumes();
415}
416
417//_________________________________________________________________________________________________
418bool Detector::ProcessHits(FairVolume* vol)
419{
420 // This method is called from the MC stepping
421 if (!(fMC->TrackCharge())) {
422 return kFALSE;
423 }
424
425 int volID = vol->getMCid();
426
427 auto it = mActiveSensorMap.find(volID);
428 if (it == mActiveSensorMap.end()) {
429 return kFALSE; // Not a sensitive volume
430 }
431
432 int lay = it->second;
433
434 auto stack = (o2::data::Stack*)fMC->GetStack();
435
436 bool startHit = false, stopHit = false;
437 unsigned char status = 0;
438 if (fMC->IsTrackEntering()) {
439 status |= Hit::kTrackEntering;
440 }
441 if (fMC->IsTrackInside()) {
442 status |= Hit::kTrackInside;
443 }
444 if (fMC->IsTrackExiting()) {
445 status |= Hit::kTrackExiting;
446 }
447 if (fMC->IsTrackOut()) {
448 status |= Hit::kTrackOut;
449 }
450 if (fMC->IsTrackStop()) {
451 status |= Hit::kTrackStopped;
452 }
453 if (fMC->IsTrackAlive()) {
454 status |= Hit::kTrackAlive;
455 }
456
457 // track is entering or created in the volume
458 if ((status & Hit::kTrackEntering) || (status & Hit::kTrackInside && !mTrackData.mHitStarted)) {
459 startHit = true;
460 } else if ((status & (Hit::kTrackExiting | Hit::kTrackOut | Hit::kTrackStopped))) {
461 stopHit = true;
462 }
463
464 // increment energy loss at all steps except entrance
465 if (!startHit) {
466 mTrackData.mEnergyLoss += fMC->Edep();
467 }
468 if (!(startHit | stopHit)) {
469 return kFALSE; // do noting
470 }
471 if (startHit) {
472 mTrackData.mEnergyLoss = 0.;
473 fMC->TrackMomentum(mTrackData.mMomentumStart);
474 fMC->TrackPosition(mTrackData.mPositionStart);
475 mTrackData.mTrkStatusStart = status;
476 mTrackData.mHitStarted = true;
477 }
478 if (stopHit) {
479 TLorentzVector positionStop;
480 fMC->TrackPosition(positionStop);
481 // Retrieve the indices with the volume path
482 int chipindex = lay;
483
484 Hit* p = addHit(stack->GetCurrentTrackNumber(), chipindex, mTrackData.mPositionStart.Vect(), positionStop.Vect(),
485 mTrackData.mMomentumStart.Vect(), mTrackData.mMomentumStart.E(), positionStop.T(),
486 mTrackData.mEnergyLoss, mTrackData.mTrkStatusStart, status);
487 // p->SetTotalEnergy(vmc->Etot());
488
489 // RS: not sure this is needed
490 // Increment number of Detector det points in TParticle
491 stack->addHit(GetDetId());
492 }
493
494 return kTRUE;
495}
496
497//_________________________________________________________________________________________________
498void Detector::createMaterials()
499{
500 int ifield = 2;
501 float fieldm = 10.0;
503
504 float tmaxfdSi = 0.1; // .10000E+01; // Degree
505 float stemaxSi = 0.0075; // .10000E+01; // cm
506 float deemaxSi = 0.1; // 0.30000E-02; // Fraction of particle's energy 0<deemax<=1
507 float epsilSi = 1.0E-4; // .10000E+01;
508 float stminSi = 0.0; // cm "Default value used"
509
510 float tmaxfdAir = 0.1; // .10000E+01; // Degree
511 float stemaxAir = .10000E+01; // cm
512 float deemaxAir = 0.1; // 0.30000E-02; // Fraction of particle's energy 0<deemax<=1
513 float epsilAir = 1.0E-4; // .10000E+01;
514 float stminAir = 0.0; // cm "Default value used"
515
516 // AIR
517 float aAir[4] = {12.0107, 14.0067, 15.9994, 39.948};
518 float zAir[4] = {6., 7., 8., 18.};
519 float wAir[4] = {0.000124, 0.755267, 0.231781, 0.012827};
520 float dAir = 1.20479E-3;
521
522 o2::base::Detector::Mixture(1, "AIR$", aAir, zAir, dAir, 4, wAir);
523 o2::base::Detector::Medium(1, "AIR$", 1, 0, ifield, fieldm, tmaxfdAir, stemaxAir, deemaxAir, epsilAir, stminAir);
524
525 o2::base::Detector::Material(3, "SILICON$", 0.28086E+02, 0.14000E+02, 0.23300E+01, 0.93600E+01, 0.99900E+03);
526 o2::base::Detector::Medium(3, "SILICON$", 3, 0, ifield, fieldm, tmaxfdSi, stemaxSi, deemaxSi, epsilSi, stminSi);
527}
528
529//_________________________________________________________________________________________________
530void Detector::EndOfEvent() { Reset(); }
531
532//_________________________________________________________________________________________________
534{
535 // This will create a branch in the output tree called Hit, setting the last
536 // parameter to kFALSE means that this collection will not be written to the file,
537 // it will exist only during the simulation
538
539 if (FairRootManager::Instance()) {
540 FairRootManager::Instance()->RegisterAny(addNameTo("Hit").data(), mHits, kTRUE);
541 }
542}
543
544//_________________________________________________________________________________________________
546{
547 if (!o2::utils::ShmManager::Instance().isOperational()) {
548 mHits->clear();
549 }
550}
551
552//_________________________________________________________________________________________________
554{
555 // Create detector materials
556 createMaterials();
557
558 // Construct the detector geometry
559 createGeometry();
560}
561
562//_________________________________________________________________________________________________
563void Detector::createGeometry()
564{
565
566 TGeoVolume* volFT3 = new TGeoVolumeAssembly(GeometryTGeo::getFT3VolPattern());
567 TGeoVolume* volIFT3 = new TGeoVolumeAssembly(GeometryTGeo::getFT3InnerVolPattern());
568
569 LOG(info) << "FT3: createGeometry volume name = " << GeometryTGeo::getFT3VolPattern();
570
571 TGeoVolume* vALIC = gGeoManager->GetVolume("barrel");
572 if (!vALIC) {
573 LOG(fatal) << "Could not find the top volume";
574 }
575
576 TGeoVolume* A3IPvac = gGeoManager->GetVolume("OUT_PIPEVACUUM");
577 if (!A3IPvac) {
578 LOG(info) << "Running simulation with no beam pipe.";
579 }
580
581 // This will need to adapt to the new scheme
582 if (!A3IPvac) {
583 for (int direction : {IdxBackwardDisks, IdxForwardDisks}) { // Backward layers at mLayers[0]; Forward layers at mLayers[1]
584 const std::string directionString = direction ? "Forward" : "Backward";
585 LOG(info) << " Creating FT3 without beampipe " << directionString << " layers:";
586 for (int iLayer = 0; iLayer < mLayers[direction].size(); iLayer++) {
587 mLayers[direction][iLayer].createLayer(volFT3);
588 }
589 }
590 vALIC->AddNode(volFT3, 2, new TGeoTranslation(0., 30., 0.));
591 } else { // If beampipe is enabled append inner disks to beampipe filling volume, this should be temporary.
592 for (int direction : {IdxBackwardDisks, IdxForwardDisks}) {
593 const std::string directionString = direction ? "Forward" : "Backward";
594 LOG(info) << " Creating FT3 " << directionString << " layers:";
595 for (int iLayer = 0; iLayer < mLayers[direction].size(); iLayer++) {
596 LOG(info) << " Creating " << directionString << " layer " << iLayer;
597 if (mLayers[direction][iLayer].getIsInMiddleLayer()) { // ML disks
598 mLayers[direction][iLayer].createLayer(volIFT3);
599 } else {
600 mLayers[direction][iLayer].createLayer(volFT3);
601 }
602 }
603 }
604 A3IPvac->AddNode(volIFT3, 2, new TGeoTranslation(0., 0., 0.));
605 vALIC->AddNode(volFT3, 2, new TGeoTranslation(0., 30., 0.));
606 }
607}
608
609//_________________________________________________________________________________________________
610void Detector::defineSensitiveVolumes()
611{
612 TGeoManager* geoManager = gGeoManager;
613
614 // Get the flat list of ALL volumes present in the geometry
615 TObjArray* allVolumes = geoManager->GetListOfVolumes();
616 int nVolumes = allVolumes->GetEntriesFast();
617
618 LOG(info) << "Adding FT3 Sensitive Volumes by iterating over all geometry volumes...";
619
620 for (int direction : {IdxBackwardDisks, IdxForwardDisks}) {
621 for (int iLayer = 0; iLayer < getNumberOfLayers(); iLayer++) {
622 int iSens = 0;
623
624 // Build the "signatures" (prefixes) of the names for the various layouts for this specific layer and direction:
625
626 // 1. Trapezoidal/Cylindrical (format: FT3Sensor_<dir>_<layer>)
627 std::string sig1 = Form("%s_%d_%d", GeometryTGeo::getFT3SensorPattern(), direction, iLayer);
628
629 // 2. Segmented front/back (format: FT3Sensor_front_<layer>_<dir>_...)
630 std::string sig2 = "FT3Sensor_front_" + std::to_string(iLayer) + "_" + std::to_string(direction);
631 std::string sig3 = "FT3Sensor_back_" + std::to_string(iLayer) + "_" + std::to_string(direction);
632
633 // 3. SegmentedStave (format: FT3Sensor_<layer>_<dir>_...)
634 // Add the trailing underscore to avoid confusing it with sig1
635 std::string sig4 = "FT3Sensor_" + std::to_string(iLayer) + "_" + std::to_string(direction) + "_";
636
637 // Iterate over all existing volumes to find matches
638 for (int i = 0; i < nVolumes; ++i) {
639 TGeoVolume* v = (TGeoVolume*)allVolumes->At(i);
640 std::string vName = v->GetName();
641
642 // Explicitly exclude the inactive silicon regions created in FT3Module
643 if (vName.find("Inactive") != std::string::npos || vName.find("inactive") != std::string::npos) {
644 continue;
645 }
646
647 // Check if the volume name matches one of our active sensors
648 bool isMatch = false;
649 if (vName == sig1) {
650 isMatch = true; // Exact match for Trapezoidal/Cylindrical layouts
651 } else if (vName.find(sig2) == 0 || vName.find(sig3) == 0 || vName.find(sig4) == 0) {
652 isMatch = true; // Prefix match for Segmented and SegmentedStave layouts
653 }
654
655 if (isMatch) {
656 AddSensitiveVolume(v);
657 int volID = gMC ? TVirtualMC::GetMC()->VolId(vName.c_str()) : 0;
658 if (volID > 0) {
659 mActiveSensorMap[volID] = iLayer;
660 }
661 iSens++;
662 }
663 }
664
665 if (iSens == 0) {
666 LOG(error) << "NO sensitive volume found for direction " << direction << ", layer " << iLayer;
667 } else {
668 LOG(info) << iSens << " sensitive volume(s) added for direction " << direction << " layer " << iLayer;
669 }
670 }
671 }
672}
673
674//_________________________________________________________________________________________________
675Hit* Detector::addHit(int trackID, int detID, const TVector3& startPos, const TVector3& endPos,
676 const TVector3& startMom, double startE, double endTime, double eLoss, unsigned char startStatus,
677 unsigned char endStatus)
678{
679 mHits->emplace_back(trackID, detID, startPos, endPos, startMom, startE, endTime, eLoss, startStatus, endStatus);
680 return &(mHits->back());
681}
682
Definition of the Stack class.
Definition of the ITSMFT Hit class.
Definition of the FT3Layer class.
int32_t i
ClassImp(IdPath)
uint32_t stack
Definition RawData.h:1
Definition of the GeometryTGeo class.
Definition of the Detector class.
Detector & operator=(const Detector &)
Definition Detector.cxx:46
void Mixture(Int_t imat, const char *name, Float_t *a, Float_t *z, Float_t dens, Int_t nlmat, Float_t *wmat)
Definition Detector.cxx:66
void Medium(Int_t numed, const char *name, Int_t nmat, Int_t isvol, Int_t ifield, Float_t fieldm, Float_t tmaxfd, Float_t stemax, Float_t deemax, Float_t epsil, Float_t stmin, Float_t *ubuf=nullptr, Int_t nbuf=0)
Definition Detector.cxx:72
static void initFieldTrackingParams(int &mode, float &maxfield)
Definition Detector.cxx:143
virtual void InitializeO2Detector()=0
Definition Detector.cxx:98
void Material(Int_t imat, const char *name, Float_t a, Float_t z, Float_t dens, Float_t radl, Float_t absl, Float_t *buf=nullptr, Int_t nwbuf=0)
Definition Detector.cxx:59
std::string addNameTo(const char *ext) const
Definition Detector.h:150
std::array< std::vector< TString >, 2 > mLayerName
Definition Detector.h:118
void buildFT3ScopingV3()
Definition Detector.cxx:253
static constexpr int IdxBackwardDisks
Definition Detector.h:100
void EndOfEvent() override
Definition Detector.cxx:530
static constexpr int IdxForwardDisks
Definition Detector.h:99
o2::itsmft::Hit * addHit(int trackID, int detID, const TVector3 &startPos, const TVector3 &endPos, const TVector3 &startMom, double startE, double endTime, double eLoss, unsigned char startStatus, unsigned char endStatus)
This method is an example of how to add your own point of type Hit to the clones array.
Definition Detector.cxx:675
void Register() override
Registers the produced collections in FAIRRootManager.
Definition Detector.cxx:533
void ConstructGeometry() override
Base class to create the detector geometry.
Definition Detector.cxx:553
Detector()
Default constructor.
Definition Detector.cxx:57
Bool_t ProcessHits(FairVolume *v=nullptr) override
This method is called for each step during simulation (see FairMCApplication::Stepping())
Definition Detector.cxx:418
void buildBasicFT3(const FT3BaseParam &param)
Definition Detector.cxx:65
void buildFT3NewVacuumVessel()
Definition Detector.cxx:188
void Reset() override
Has to be called after each event to reset the containers.
Definition Detector.cxx:545
static const char * getFT3SensorPattern()
static const char * getFT3LayerPattern()
static const char * getFT3VolPattern()
static const char * getFT3InnerVolPattern()
static ShmManager & Instance()
Definition ShmManager.h:61
const GLdouble * v
Definition glcorearb.h:832
GLenum GLenum GLsizei const GLuint GLboolean enabled
Definition glcorearb.h:2513
GLboolean * data
Definition glcorearb.h:298
GLuint GLfloat x0
Definition glcorearb.h:5034
GLenum GLfloat param
Definition glcorearb.h:271
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
void freeSimVector(std::vector< T > *ptr)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::string to_string(gsl::span< T, Size > span)
Definition common.h:52
Common utility functions.
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"