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