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#include "TPCSimulation/Point.h"
17
18#include "DetectorsBase/Stack.h"
20
21#include "FairVolume.h" // for FairVolume
22
23#include "TVirtualMC.h" // for TVirtualMC, gMC
24
25#include <cstddef> // for NULL
26
27#include "FairGeoVolume.h"
28#include "FairGeoNode.h"
29#include "FairGeoLoader.h"
30#include "FairGeoInterface.h"
31#include "FairRun.h"
32#include "FairRuntimeDb.h"
33#include <fairlogger/Logger.h>
34#include "FairRootManager.h"
35
36#include "TSystem.h"
37#include "TVirtualMC.h"
38
39#include "TFile.h"
40
41// geo stuff
42#include "TGeoManager.h"
43#include "TGeoVolume.h"
44#include "TGeoPcon.h"
45#include "TGeoTube.h"
46#include "TGeoCone.h"
47#include "TGeoPgon.h"
48#include "TGeoTrd1.h"
49#include "TGeoCompositeShape.h"
50#include "TGeoPara.h"
51#include "TGeoPhysicalNode.h"
52#include "TGeoHalfSpace.h"
53#include "TGeoArb8.h"
54#include "TGeoMatrix.h"
55
56#include <iostream>
57#include <cmath>
58
59using std::cout;
60using std::endl;
61using std::ifstream;
62using std::ios_base;
63using namespace o2::tpc;
64
65Detector::Detector(Bool_t active) : o2::base::DetImpl<Detector>("TPC", active), mGeoFileName()
66{
67 for (int i = 0; i < Sector::MAXSECTOR; ++i) {
68 mHitsPerSectorCollection[i] = o2::utils::createSimVector<o2::tpc::HitGroup>(); // new std::vector<o2::tpc::HitGroup>;
69 }
70}
71
72// forward default constructor
73Detector::Detector() : o2::tpc::Detector(kTRUE) {}
74
76 : o2::base::DetImpl<Detector>(rhs),
77 mHitCounter(0),
78 mElectronCounter(0),
79 mStepCounter(0),
80 mGeoFileName(rhs.mGeoFileName)
81{
82 for (int i = 0; i < Sector::MAXSECTOR; ++i) {
83 mHitsPerSectorCollection[i] = o2::utils::createSimVector<o2::tpc::HitGroup>(); // new std::vector<o2::tpc::HitGroup>;new std::vector<o2::tpc::HitGroup>;
84 }
85}
86
88{
89 for (int i = 0; i < Sector::MAXSECTOR; ++i) {
90 // mHitsPerSectorCollection[i]->clear();
91 o2::utils::freeSimVector(mHitsPerSectorCollection[i]);
92 }
93 std::cout << "Produced hits " << mHitCounter << "\n";
94 std::cout << "Produced electrons " << mElectronCounter << "\n";
95 std::cout << "Stepping called " << mStepCounter << "\n";
96}
97
99{
100 // Define the list of sensitive volumes
101 defineSensitiveVolumes();
102}
103
104Bool_t Detector::ProcessHits(FairVolume* vol)
105{
106 mStepCounter++;
107 auto& gasParam = ParameterGas::Instance();
108 auto& detParam = ParameterDetector::Instance();
109 const Int_t kMaxDistRef = 15; // maximal difference between 2 stored references - the parameter should be 15 cm as default
110 static Double_t lastReferenceR = 0; // keeps last reference point in radius (cm)
111
112 /* This method is called from the MC stepping for the sensitive volume only */
113 // LOG(info) << "tpc::ProcessHits";
114 const double trackCharge = fMC->TrackCharge();
115 if (static_cast<int>(trackCharge) == 0) {
116
117 // set a very large step size for neutral particles
118 fMC->SetMaxStep(1.e10);
119 return kFALSE; // take only charged particles
120 }
121
122 // ===| SET THE LENGTH OF THE NEXT ENERGY LOSS STEP |=========================
123 // (We do this first so we can skip out of the method in the following
124 // Eloss->Ionization part.)
125 //
126 // In all cases we have multiple collisions and we use 2mm (+ some
127 // random shift to avoid binning effects), which was tuned for GEANT4, see
128 // https://indico.cern.ch/event/316891/contributions/732168/
129
130 const double rnd = fMC->GetRandom()->Rndm();
131 fMC->SetMaxStep(0.2 + (2. * rnd - 1.) * 0.05); // 2 mm +- rndm*0.5mm step
132
133 // ===| check active sector |=================================================
134 //
135 // Get the sector ID and check if the sector is active
136 static thread_local TLorentzVector position;
137 fMC->TrackPosition(position);
138 // for processing reasons in the digitizer, the sectors are shifted by -10deg, so sector 0 will be
139 // from -10 - +10 deg instead of 0-20 and so on
140 const int sectorID = static_cast<int>(Sector::ToShiftedSector(position.X(), position.Y(), position.Z()));
141 // const int sectorID = static_cast<int>(Sector::ToSector(position.X(), position.Y(), position.Z()));
142 // TODO: Temporary hack to process only one sector
143 // if (sectorID != 0) return kFALSE;
144
145 // ---| momentum and beta gamma |---
146 static TLorentzVector momentum; // static to make avoid creation/deletion of this expensive object
147 fMC->TrackMomentum(momentum);
148
149 const float time = fMC->TrackTime() * 1.0e9;
150 const int trackID = fMC->GetStack()->GetCurrentTrackNumber();
151 const int detID = vol->getMCid();
152 o2::data::Stack* stack = (o2::data::Stack*)fMC->GetStack();
153 if (fMC->IsTrackEntering() || fMC->IsTrackExiting()) {
154 stack->addTrackReference(o2::TrackReference(position.X(), position.Y(), position.Z(), momentum.X(), momentum.Y(),
155 momentum.Z(), fMC->TrackLength(), time, trackID, GetDetId()));
156 lastReferenceR = fMC->TrackLength();
157 }
158 if (TMath::Abs(lastReferenceR - fMC->TrackLength()) > kMaxDistRef) {
159 stack->addTrackReference(o2::TrackReference(position.X(), position.Y(), position.Z(), momentum.X(), momentum.Y(),
160 momentum.Z(), fMC->TrackLength(), time, trackID, GetDetId()));
161 lastReferenceR = fMC->TrackLength();
162 }
163
164 // ---| remove clusters between the IFC and the FC strips |---
165 // those should not enter the active readout area
166 // do coarse selection before, to limit number of transformations
167 if (detParam.ExcludeFCGap) {
168 const auto rCluster = std::sqrt(position.X() * position.X() + position.Y() * position.Y());
169 const float rodRin = 81.5 + 2.2; // radial position of the inner field cage rods + radial size of the field cage rods
170 const float rodRout = 254.25 + 2.2; // radial position of the outer field cage rods + radial size of the field cage rods
171 const float fcLxIn = 82.428409; // position of the inner FC strips in local x = cos(10 deg) * rodRin;
172 const float fcLxOut = 252.55395; // position of the outer FC strips in local x = cos(10 deg) * rodRin;
173
174 if (rCluster < rodRin || rCluster > fcLxOut) {
175 const int sectorIDnonShift = static_cast<int>(Sector::ToSector(position.X(), position.Y(), position.Z()));
176 const double alpha = TMath::DegToRad() * (10. + sectorIDnonShift * 20.);
177 const double cs = std::cos(-alpha), sn = std::sin(-alpha);
178 const auto localX = position.X() * cs - position.Y() * sn;
179 // fine cut
180 if (localX < fcLxIn || localX > fcLxOut) {
181 return kFALSE;
182 }
183 }
184 }
185
186 // ===| CONVERT THE ENERGY LOSS TO IONIZATION ELECTRONS |=====================
187 //
188 // The energy loss is implemented directly below and taken GEANT3,
189 // ILOSS model 5 (in gfluct.F), which gives
190 // the energy loss in a single collision (NA49 model).
191 // TODO: Add discussion about drawback
192
193 Int_t numberOfElectrons = 0;
194 // I.H. - the type expected in addHit is short
195
196 // ---| Stepsize in cm |---
197 const double stepSize = fMC->TrackStep();
198
199 double betaGamma = momentum.P() / fMC->TrackMass();
200 betaGamma = TMath::Max(betaGamma, 7.e-3); // protection against too small bg
201
202 // ---| number of primary ionisations per cm |---
203 const double primaryElectronsPerCM =
204 gasParam.Nprim * BetheBlochAleph(static_cast<float>(betaGamma), gasParam.BetheBlochParam[0],
205 gasParam.BetheBlochParam[1], gasParam.BetheBlochParam[2],
206 gasParam.BetheBlochParam[3], gasParam.BetheBlochParam[4]);
207
208 // ---| mean number of collisions and random for this event |---
209 const double meanNcoll = stepSize * trackCharge * trackCharge * primaryElectronsPerCM;
210 const int nColl = static_cast<int>(fMC->GetRandom()->Poisson(meanNcoll));
211
212 // Variables needed to generate random powerlaw distributed energy loss
213 const double alpha_p1 = 1. - gasParam.Exp; // NA49/G3 value
214 const double oneOverAlpha_p1 = 1. / alpha_p1;
215 const double eMin = gasParam.Ipot;
216 const double eMax = gasParam.Eend;
217 const double kMin = TMath::Power(eMin, alpha_p1);
218 const double kMax = TMath::Power(eMax, alpha_p1);
219 const double wIon = gasParam.Wion;
220
221 for (Int_t n = 0; n < nColl; n++) {
222 // Use GEANT3 / NA49 expression:
223 // P(eDep) ~ k * edep^-gasParam.getExp()
224 // eMin(~I) < eDep < eMax(300 electrons)
225 // k fixed so that Int_Emin^EMax P(Edep) = 1.
226 const double rndm = fMC->GetRandom()->Rndm();
227 const double eDep = TMath::Power((kMax - kMin) * rndm + kMin, oneOverAlpha_p1);
228 int nel_step = static_cast<int>(((eDep - eMin) / wIon) + 1);
229 nel_step = TMath::Min(nel_step, 300); // 300 electrons corresponds to 10 keV
230 numberOfElectrons += nel_step;
231 }
232
233 // LOG(info) << "tpc::AddHit" << FairLogger::endl << "Eloss: "
234 //<< fMC->Edep() << ", Nelectrons: "
235 //<< numberOfElectrons;
236
237 if (numberOfElectrons <= 0) { // Could maybe be smaller than 0 due to the Gamma function
238 return kFALSE;
239 }
240
241 // ADD HIT
242 static thread_local int oldTrackId = trackID;
243 static thread_local int oldDetId = detID;
244 static thread_local int groupCounter = 0;
245 static thread_local int oldSectorId = sectorID;
246
247 // a new group is starting -> put it into the container
248 static thread_local HitGroup* currentgroup = nullptr;
249 if (groupCounter == 0) {
250 mHitsPerSectorCollection[sectorID]->emplace_back(trackID);
251 currentgroup = &(mHitsPerSectorCollection[sectorID]->back());
252 }
253 if (trackID == oldTrackId && oldSectorId == sectorID) {
254 groupCounter++;
255 mHitCounter++;
256 mElectronCounter += numberOfElectrons;
257 currentgroup->addHit(position.X(), position.Y(), position.Z(), time, numberOfElectrons);
258
259 // add last buffered hit, which was not yet added to the currentgroup
260 if (mHitLast.GetEnergyLoss() >= 0) {
261 currentgroup->addHit(mHitLast.GetX(), mHitLast.GetY(), mHitLast.GetZ(), mHitLast.GetTime(), mHitLast.GetEnergyLoss());
262 mHitLast.mELoss = -1;
263 groupCounter++;
264 mHitCounter++;
265 mElectronCounter += mHitLast.GetEnergyLoss();
266 }
267 }
268 // finish group
269 else {
270 oldTrackId = trackID;
271 oldSectorId = sectorID;
272 groupCounter = 0;
273
274 // buffer this hit, otherwise it wouldnt be stored in the HitGroup
275 mHitLast = ElementalHit(position.X(), position.Y(), position.Z(), time, numberOfElectrons);
276 }
277
278 // LOG(info) << "tpc::AddHit" << FairLogger::endl
279 //<< " -- " << trackNumberID <<"," << volumeID << " " << vol->GetName()
280 //<< ", Pos: (" << position.X() << ", " << position.Y() <<", "<< position.Z()<< ", " << r << ") "
281 //<< ", Mom: (" << momentum.Px() << ", " << momentum.Py() << ", " << momentum.Pz() << ") "
282 //<< " Time: "<< time <<", Len: " << length << ", Nelectrons: " <<
283 // numberOfElectrons;
284 // I.H. - the code above does not compile if uncommented
285
286 // Increment number of Detector det points in TParticle
287 stack->addHit(GetDetId());
288
289 return kTRUE;
290}
291
293{
294 if (!o2::utils::ShmManager::Instance().isOperational()) {
295 for (int i = 0; i < Sector::MAXSECTOR; ++i) {
296 mHitsPerSectorCollection[i]->clear();
297 }
298 }
299}
300
302{
309 auto* mgr = FairRootManager::Instance();
310 for (int i = 0; i < Sector::MAXSECTOR; ++i) {
311 mgr->RegisterAny(getHitBranchNames(i).c_str(), mHitsPerSectorCollection[i], kTRUE);
312 }
313}
314
316{
317 for (int i = 0; i < Sector::MAXSECTOR; ++i) {
318 // cout << "CLEARED (reset) hitsCollection " << i << " " << mHitsPerSectorCollection[i] << endl;
319 mHitsPerSectorCollection[i]->clear();
320 }
321}
322
324{
325 // Create the detector materials
326 CreateMaterials();
327
328 // Load geometry
329 // LoadGeometryFromFile();
330 ConstructTPCGeometry();
331}
332
333void Detector::CreateMaterials()
334{
335 //-----------------------------------------------
336 // Create Materials for for TPC simulations
337 //-----------------------------------------------
338
339 //-----------------------------------------------------------------
340 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
341 //-----------------------------------------------------------------
342
343 const auto& gasParam = ParameterGas::Instance();
344
345 Int_t iSXFLD = 2;
346 Float_t sXMGMX = 10.0;
347 // init the field tracking params
349
350 Float_t amat[7]; // atomic numbers
351 Float_t zmat[7]; // z
352 Float_t wmat[7]; // proportions
353
354 Float_t density;
355
356 // TODO: load pressure and temperature values from CCDB
357 const Double_t pressure = gasParam.Pressure; // in mbar
358 const Double_t temperature = gasParam.Temperature + 273.15; // in K
359
360 // densities were taken for these values
361 const Double_t t1 = 293.15; // 20°C in K
362 const Double_t p1 = 1013.25; // 1 atm in mbars
363
364 // sanity check - temperature between 10 and 30 deg, pressure between 800 and 1200 mbar
365 Double_t ptCorr = 1.;
366 if (TMath::Abs(temperature - 293.15) > 10. || TMath::Abs(pressure - 1000.) > 200.) {
367 ptCorr = 1.;
368 } else {
369 ptCorr = (pressure * t1) / (p1 * temperature);
370 }
371 LOG(info) << "Setting gas density correction to: " << ptCorr;
372
373 //***************** Gases *************************
374
375 //--------------------------------------------------------------
376 // gases - air and CO2
377 //--------------------------------------------------------------
378
379 // CO2
380
381 amat[0] = 12.011;
382 amat[1] = 15.9994;
383
384 zmat[0] = 6.;
385 zmat[1] = 8.;
386
387 wmat[0] = 0.2729;
388 wmat[1] = 0.7271;
389
390 density = 1.842e-3;
391
392 o2::base::Detector::Mixture(10, "CO2", amat, zmat, density * ptCorr, 2, wmat);
393 //
394 // Air
395 //
396 amat[0] = 15.9994;
397 amat[1] = 14.007;
398 //
399 zmat[0] = 8.;
400 zmat[1] = 7.;
401 //
402 wmat[0] = 0.233;
403 wmat[1] = 0.767;
404 //
405 density = 0.001205;
406
407 o2::base::Detector::Mixture(11, "Air", amat, zmat, density * ptCorr, 2, wmat);
408
409 //----------------------------------------------------------------
410 // drift gases 5 mixtures, 5 materials
411 //----------------------------------------------------------------
412 //
413 // Drift gases 1 - nonsensitive, 2 - sensitive, 3 - for Kr
414 // Composition by % of volume, values at 20deg and 1 atm.
415 //
416 // get the geometry title - defined in Config.C
417 //
418 //--------------------------------------------------------------
419 // predefined gases, composition taken from param file
420 //--------------------------------------------------------------
421 TString names[6] = {"Ne", "Ar", "CO2", "N", "CF4", "CH4"};
422 TString gname;
423
425 // container in the future
426 Float_t comp[6] = {90. / 105., 0., 10. / 105., 5. / 105., 0., 0.};
427 // indices:
428 // 0-Ne, 1-Ar, 2-CO2, 3-N, 4-CF4, 5-CH4
429 //
430 // elements' masses
431 //
432 amat[0] = 20.18; // Ne
433 amat[1] = 39.95; // Ar
434 amat[2] = 12.011; // C
435 amat[3] = 15.9994; // O
436 amat[4] = 14.007; // N
437 amat[5] = 18.998; // F
438 amat[6] = 1.; // H
439 //
440 // elements' atomic numbers
441 //
442 //
443 zmat[0] = 10.; // Ne
444 zmat[1] = 18.; // Ar
445 zmat[2] = 6.; // C
446 zmat[3] = 8.; // O
447 zmat[4] = 7.; // N
448 zmat[5] = 9.; // F
449 zmat[6] = 1.; // H
450 //
451 // Mol masses
452 //
453 Float_t wmol[6];
454 wmol[0] = 20.18; // Ne
455 wmol[1] = 39.948; // Ar
456 wmol[2] = 44.0098; // CO2
457 wmol[3] = 2. * 14.0067; // N2
458 wmol[4] = 88.0046; // CF4
459 wmol[5] = 16.011; // CH4
460 //
461 Float_t wtot = 0.; // total mass of the mixture
462 for (Int_t i = 0; i < 6; i++) {
463 wtot += *(comp + i) * wmol[i];
464 }
465 wmat[0] = comp[0] * amat[0] / wtot; // Ne
466 wmat[1] = comp[1] * amat[1] / wtot; // Ar
467 wmat[2] = (comp[2] * amat[2] + comp[4] * amat[2] + comp[5] * amat[2]) / wtot; // C
468 wmat[3] = comp[2] * amat[3] * 2. / wtot; // O
469 wmat[4] = comp[3] * amat[4] * 2. / wtot; // N
470 wmat[5] = comp[4] * amat[5] * 4. / wtot; // F
471 wmat[6] = comp[5] * amat[6] * 4. / wtot; // H
472 //
473 // densities (NTP)
474 //
475 Float_t dens[6] = {0.839e-3, 1.661e-3, 1.842e-3, 1.165e-3, 3.466e-3, 0.668e-3};
476 //
477 density = 0.;
478 for (Int_t i = 0; i < 6; i++) {
479 density += comp[i] * dens[i];
480 }
481 //
482 // names
483 //
484 Int_t cnt = 0;
485 for (Int_t i = 0; i < 6; i++) {
486 if (comp[i]) {
487 if (cnt) {
488 gname += "-";
489 }
490 gname += names[i];
491 cnt++;
492 }
493 }
494 TString gname1, gname2, gname3;
495 gname1 = gname + "-1";
496 gname2 = gname + "-2";
497 gname3 = gname + "-3";
498 //
499 // take only elements with nonzero weights
500 //
501 Float_t amat1[6], zmat1[6], wmat1[6];
502 cnt = 0;
503 for (Int_t i = 0; i < 7; i++) {
504 if (wmat[i]) {
505 zmat1[cnt] = zmat[i];
506 amat1[cnt] = amat[i];
507 wmat1[cnt] = wmat[i];
508 cnt++;
509 }
510 }
511
512 //
513 o2::base::Detector::Mixture(12, gname1.Data(), amat1, zmat1, density * ptCorr, cnt, wmat1); // nonsensitive
514 o2::base::Detector::Mixture(13, gname2.Data(), amat1, zmat1, density * ptCorr, cnt, wmat1); // sensitive
515 o2::base::Detector::Mixture(40, gname3.Data(), amat1, zmat1, density * ptCorr, cnt, wmat1); // sensitive Kr
516
517 //----------------------------------------------------------------------
518 // solid materials
519 //----------------------------------------------------------------------
520
521 // Kevlar C14H22O2N2
522
523 amat[0] = 12.011;
524 amat[1] = 1.;
525 amat[2] = 15.999;
526 amat[3] = 14.006;
527
528 zmat[0] = 6.;
529 zmat[1] = 1.;
530 zmat[2] = 8.;
531 zmat[3] = 7.;
532
533 wmat[0] = 14.;
534 wmat[1] = 22.;
535 wmat[2] = 2.;
536 wmat[3] = 2.;
537
538 density = 1.45;
539
540 o2::base::Detector::Mixture(14, "Kevlar", amat, zmat, density, -4, wmat);
541
542 // NOMEX
543
544 amat[0] = 12.011;
545 amat[1] = 1.;
546 amat[2] = 15.999;
547 amat[3] = 14.006;
548
549 zmat[0] = 6.;
550 zmat[1] = 1.;
551 zmat[2] = 8.;
552 zmat[3] = 7.;
553
554 wmat[0] = 14.;
555 wmat[1] = 22.;
556 wmat[2] = 2.;
557 wmat[3] = 2.;
558
559 density = 0.029;
560
561 o2::base::Detector::Mixture(15, "NOMEX", amat, zmat, density, -4, wmat);
562
563 // Makrolon C16H18O3
564
565 amat[0] = 12.011;
566 amat[1] = 1.;
567 amat[2] = 15.999;
568
569 zmat[0] = 6.;
570 zmat[1] = 1.;
571 zmat[2] = 8.;
572
573 wmat[0] = 16.;
574 wmat[1] = 18.;
575 wmat[2] = 3.;
576
577 density = 1.2;
578
579 o2::base::Detector::Mixture(16, "Makrolon", amat, zmat, density, -3, wmat);
580
581 // Tedlar C2H3F
582
583 amat[0] = 12.011;
584 amat[1] = 1.;
585 amat[2] = 18.998;
586
587 zmat[0] = 6.;
588 zmat[1] = 1.;
589 zmat[2] = 9.;
590
591 wmat[0] = 2.;
592 wmat[1] = 3.;
593 wmat[2] = 1.;
594
595 density = 1.71;
596
597 o2::base::Detector::Mixture(17, "Tedlar", amat, zmat, density, -3, wmat);
598
599 // Mylar C5H4O2
600
601 amat[0] = 12.011;
602 amat[1] = 1.;
603 amat[2] = 15.9994;
604
605 zmat[0] = 6.;
606 zmat[1] = 1.;
607 zmat[2] = 8.;
608
609 wmat[0] = 5.;
610 wmat[1] = 4.;
611 wmat[2] = 2.;
612
613 density = 1.39;
614
615 o2::base::Detector::Mixture(18, "Mylar", amat, zmat, density, -3, wmat);
616 // material for "prepregs"
617 // Epoxy - C14 H20 O3
618 // Quartz SiO2
619 // Carbon C
620 // prepreg1 60% C-fiber, 40% epoxy (vol)
621 amat[0] = 12.011;
622 amat[1] = 1.;
623 amat[2] = 15.994;
624
625 zmat[0] = 6.;
626 zmat[1] = 1.;
627 zmat[2] = 8.;
628
629 wmat[0] = 0.923;
630 wmat[1] = 0.023;
631 wmat[2] = 0.054;
632
633 density = 1.859;
634
635 o2::base::Detector::Mixture(19, "Prepreg1", amat, zmat, density, 3, wmat);
636
637 // prepreg2 60% glass-fiber, 40% epoxy
638
639 amat[0] = 12.01;
640 amat[1] = 1.;
641 amat[2] = 15.994;
642 amat[3] = 28.086;
643
644 zmat[0] = 6.;
645 zmat[1] = 1.;
646 zmat[2] = 8.;
647 zmat[3] = 14.;
648
649 wmat[0] = 0.194;
650 wmat[1] = 0.023;
651 wmat[2] = 0.443;
652 wmat[3] = 0.34;
653
654 density = 1.82;
655
656 o2::base::Detector::Mixture(20, "Prepreg2", amat, zmat, density, 4, wmat);
657
658 // prepreg3 50% glass-fiber, 50% epoxy
659
660 amat[0] = 12.01;
661 amat[1] = 1.;
662 amat[2] = 15.994;
663 amat[3] = 28.086;
664
665 zmat[0] = 6.;
666 zmat[1] = 1.;
667 zmat[2] = 8.;
668 zmat[3] = 14.;
669
670 wmat[0] = 0.257;
671 wmat[1] = 0.03;
672 wmat[2] = 0.412;
673 wmat[3] = 0.3;
674
675 density = 1.725;
676
677 o2::base::Detector::Mixture(21, "Prepreg3", amat, zmat, density, 4, wmat);
678
679 // G10 60% SiO2 40% epoxy
680
681 amat[0] = 12.01;
682 amat[1] = 1.;
683 amat[2] = 15.994;
684 amat[3] = 28.086;
685
686 zmat[0] = 6.;
687 zmat[1] = 1.;
688 zmat[2] = 8.;
689 zmat[3] = 14.;
690
691 wmat[0] = 0.194;
692 wmat[1] = 0.023;
693 wmat[2] = 0.443;
694 wmat[3] = 0.340;
695
696 density = 1.7;
697
698 o2::base::Detector::Mixture(22, "G10", amat, zmat, density, 4, wmat);
699
700 // Al
701
702 amat[0] = 26.98;
703 zmat[0] = 13.;
704
705 density = 2.7;
706
707 o2::base::Detector::Material(23, "Al", amat[0], zmat[0], density, 999., 999.);
708
709 // Si (for electronics
710
711 amat[0] = 28.086;
712 zmat[0] = 14.;
713
714 density = 2.33;
715
716 o2::base::Detector::Material(24, "Si", amat[0], zmat[0], density, 999., 999.);
717
718 // Cu
719
720 amat[0] = 63.546;
721 zmat[0] = 29.;
722
723 density = 8.96;
724
725 o2::base::Detector::Material(25, "Cu", amat[0], zmat[0], density, 999., 999.);
726
727 // brass
728
729 amat[0] = 63.546;
730 zmat[0] = 29.;
731 //
732 amat[1] = 65.409;
733 zmat[1] = 30.;
734 //
735 wmat[0] = 0.6;
736 wmat[1] = 0.4;
737
738 //
739 density = 8.23;
740
741 //
742 o2::base::Detector::Mixture(33, "Brass", amat, zmat, density, 2, wmat);
743
744 // Epoxy - C14 H20 O3
745
746 amat[0] = 12.011;
747 amat[1] = 1.;
748 amat[2] = 15.9994;
749
750 zmat[0] = 6.;
751 zmat[1] = 1.;
752 zmat[2] = 8.;
753
754 wmat[0] = 14.;
755 wmat[1] = 20.;
756 wmat[2] = 3.;
757
758 density = 1.25;
759
760 o2::base::Detector::Mixture(26, "Epoxy", amat, zmat, density, -3, wmat);
761
762 // Epoxy - C14 H20 O3 for glue
763
764 amat[0] = 12.011;
765 amat[1] = 1.;
766 amat[2] = 15.9994;
767
768 zmat[0] = 6.;
769 zmat[1] = 1.;
770 zmat[2] = 8.;
771
772 wmat[0] = 14.;
773 wmat[1] = 20.;
774 wmat[2] = 3.;
775
776 density = 1.25;
777
778 density *= 1.25;
779
780 o2::base::Detector::Mixture(35, "Epoxy1", amat, zmat, density, -3, wmat);
781 //
782 // epoxy film - 90% epoxy, 10% glass fiber
783 //
784 amat[0] = 12.01;
785 amat[1] = 1.;
786 amat[2] = 15.994;
787 amat[3] = 28.086;
788
789 zmat[0] = 6.;
790 zmat[1] = 1.;
791 zmat[2] = 8.;
792 zmat[3] = 14.;
793
794 wmat[0] = 0.596;
795 wmat[1] = 0.071;
796 wmat[2] = 0.257;
797 wmat[3] = 0.076;
798
799 density = 1.345;
800
801 o2::base::Detector::Mixture(34, "Epoxy-film", amat, zmat, density, 4, wmat);
802
803 // Plexiglas C5H8O2
804
805 amat[0] = 12.011;
806 amat[1] = 1.;
807 amat[2] = 15.9994;
808
809 zmat[0] = 6.;
810 zmat[1] = 1.;
811 zmat[2] = 8.;
812
813 wmat[0] = 5.;
814 wmat[1] = 8.;
815 wmat[2] = 2.;
816
817 density = 1.18;
818
819 o2::base::Detector::Mixture(27, "Plexiglas", amat, zmat, density, -3, wmat);
820
821 // Carbon
822
823 amat[0] = 12.011;
824 zmat[0] = 6.;
825 density = 2.265;
826
827 o2::base::Detector::Material(28, "C", amat[0], zmat[0], density, 999., 999.);
828
829 // Fe (steel for the inner heat screen)
830
831 amat[0] = 55.845;
832
833 zmat[0] = 26.;
834
835 density = 7.87;
836
837 o2::base::Detector::Material(29, "Fe", amat[0], zmat[0], density, 999., 999.);
838 //
839 // Peek - (C6H4-O-OC6H4-O-C6H4-CO)n
840 amat[0] = 12.011;
841 amat[1] = 1.;
842 amat[2] = 15.9994;
843
844 zmat[0] = 6.;
845 zmat[1] = 1.;
846 zmat[2] = 8.;
847
848 wmat[0] = 19.;
849 wmat[1] = 12.;
850 wmat[2] = 3.;
851 //
852 density = 1.3;
853 //
854 o2::base::Detector::Mixture(30, "Peek", amat, zmat, density, -3, wmat);
855 //
856 // Ceramics - Al2O3
857 //
858 amat[0] = 26.98;
859 amat[1] = 15.9994;
860
861 zmat[0] = 13.;
862 zmat[1] = 8.;
863
864 wmat[0] = 2.;
865 wmat[1] = 3.;
866
867 density = 3.97;
868
869 o2::base::Detector::Mixture(31, "Alumina", amat, zmat, density, -2, wmat);
870 //
871 // Ceramics for resistors
872 //
873 amat[0] = 26.98;
874 amat[1] = 15.9994;
875
876 zmat[0] = 13.;
877 zmat[1] = 8.;
878
879 wmat[0] = 2.;
880 wmat[1] = 3.;
881
882 density = 3.97;
883 //
884 density *= 1.25;
885
886 o2::base::Detector::Mixture(36, "Alumina1", amat, zmat, density, -2, wmat);
887 //
888 // liquids
889 //
890
891 // water
892
893 amat[0] = 1.;
894 amat[1] = 15.9994;
895
896 zmat[0] = 1.;
897 zmat[1] = 8.;
898
899 wmat[0] = 2.;
900 wmat[1] = 1.;
901
902 density = 1.;
903
904 o2::base::Detector::Mixture(32, "Water", amat, zmat, density, -2, wmat);
905
906 //----------------------------------------------------------
907 // tracking media for gases
908 //----------------------------------------------------------
909
910 o2::base::Detector::Medium(kAir, "Air", 11, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
911 o2::base::Detector::Medium(kDriftGas1, "DriftGas1", 12, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
912 o2::base::Detector::Medium(kDriftGas2, "DriftGas2", 13, 1, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
913 o2::base::Detector::Medium(kCO2, "CO2", 10, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
914 o2::base::Detector::Medium(kDriftGas3, "DriftGas3", 40, 1, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
915 //-----------------------------------------------------------
916 // tracking media for solids
917 //-----------------------------------------------------------
918
919 o2::base::Detector::Medium(kAl, "Al", 23, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
920 o2::base::Detector::Medium(kKevlar, "Kevlar", 14, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
921 o2::base::Detector::Medium(kNomex, "Nomex", 15, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
922 o2::base::Detector::Medium(kMakrolon, "Makrolon", 16, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
923 o2::base::Detector::Medium(kMylar, "Mylar", 18, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
924 o2::base::Detector::Medium(kTedlar, "Tedlar", 17, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
925 //
926 o2::base::Detector::Medium(kPrepreg1, "Prepreg1", 19, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
927 o2::base::Detector::Medium(kPrepreg2, "Prepreg2", 20, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
928 o2::base::Detector::Medium(kPrepreg3, "Prepreg3", 21, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
929 o2::base::Detector::Medium(kEpoxy, "Epoxy", 26, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
930
931 o2::base::Detector::Medium(kCu, "Cu", 25, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
932 o2::base::Detector::Medium(kSi, "Si", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
933 o2::base::Detector::Medium(kG10, "G10", 22, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
934 o2::base::Detector::Medium(kPlexiglas, "Plexiglas", 27, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
935 o2::base::Detector::Medium(kSteel, "Steel", 29, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
936 o2::base::Detector::Medium(kPeek, "Peek", 30, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
937 o2::base::Detector::Medium(kAlumina, "Alumina", 31, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
938 o2::base::Detector::Medium(kWater, "Water", 32, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
939 o2::base::Detector::Medium(kBrass, "Brass", 33, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
940 o2::base::Detector::Medium(kEpoxyfm, "Epoxyfm", 34, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
941 o2::base::Detector::Medium(kEpoxy1, "Epoxy1", 35, 0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
942 o2::base::Detector::Medium(kAlumina1, "Alumina1", 36, 0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
943}
944
945void Detector::ConstructTPCGeometry()
946{
947 //
948 // Create the geometry of Time Projection Chamber version 2
949 //
950 // Begin_Html
951 /*
952 * <img src="picts/AliTPC.gif">
953 */
954 // End_Html
955 // Begin_Html
956 /*
957 * <img src="picts/AliTPCv2Tree.gif">
958 */
959 // End_Html
960
961 //----------------------------------------------------------
962 // This geometry is written using TGeo class
963 // Firstly the shapes are defined, and only then the volumes
964 // What is recognized by the MC are volumes
965 //----------------------------------------------------------
966 //
967 // tpc - this will be the mother volume
968 //
969
970 // if (!mParam) {
971 // LOG(error) << "TPC Parameters not available, cannot create Geometry";
972 // return;
973 // }
974
975 //
976 // here I define a volume TPC
977 // retrive the medium name with "TPC_" as a leading string
978 //
979 auto* tpc = new TGeoPcon(0., 360., 30); // 30 sections
980 //
981 tpc->DefineSection(0, -289.6, 77., 278.);
982 tpc->DefineSection(1, -262.1, 77., 278.);
983 //
984 tpc->DefineSection(2, -262.1, 83.1, 278.);
985 tpc->DefineSection(3, -260., 83.1, 278.);
986 //
987 tpc->DefineSection(4, -260., 70., 278.);
988 tpc->DefineSection(5, -259.6, 70., 278.);
989 //
990 tpc->DefineSection(6, -259.6, 68.1, 278.);
991 tpc->DefineSection(7, -253.6, 68.1, 278.);
992 //
993 tpc->DefineSection(8, -253.6, 67.88, 278.); // hs
994 tpc->DefineSection(9, -74.0, 60.68, 278.); // hs
995 //
996 tpc->DefineSection(10, -74.0, 60.1, 278.);
997 tpc->DefineSection(11, -73.3, 60.1, 278.);
998 //
999 tpc->DefineSection(12, -73.3, 56.9, 278.);
1000 tpc->DefineSection(13, -68.5, 56.9, 278.);
1001 //
1002 tpc->DefineSection(14, -68.5, 60., 278.);
1003 tpc->DefineSection(15, -64.7, 60., 278.);
1004 //
1005 tpc->DefineSection(16, -64.7, 56.9, 278.);
1006 tpc->DefineSection(17, 73.3, 56.9, 278.);
1007 //
1008 tpc->DefineSection(18, 73.3, 60.1, 278.);
1009 tpc->DefineSection(19, 74.0, 60.1, 278.);
1010 //
1011 tpc->DefineSection(20, 74.0, 60.68, 278.); // hs
1012 tpc->DefineSection(21, 253.6, 65.38, 278.); // hs
1013 //
1014 tpc->DefineSection(22, 253.6, 65.6, 278.);
1015 tpc->DefineSection(23, 259.6, 65.6, 278.);
1016 //
1017 tpc->DefineSection(24, 259.6, 70.0, 278.);
1018 tpc->DefineSection(25, 260., 70.0, 278.);
1019 //
1020 tpc->DefineSection(26, 260., 83.1, 278.);
1021 tpc->DefineSection(27, 262.1, 83.1, 278.);
1022 //
1023 tpc->DefineSection(28, 262.1, 77., 278);
1024 tpc->DefineSection(29, 289.6, 77., 278.);
1025 //
1026 TGeoMedium* m1 = gGeoManager->GetMedium("TPC_Air");
1027 auto* v1 = new TGeoVolume("TPC_M", tpc, m1);
1028 //
1029 // drift volume - sensitive volume, extended beyond the
1030 // endcaps, because of the alignment
1031 //
1032 auto* dvol = new TGeoPcon(0., 360., 6);
1033 dvol->DefineSection(0, -260., 74.5, 264.4);
1034 dvol->DefineSection(1, -253.6, 74.5, 264.4);
1035 //
1036 dvol->DefineSection(2, -253.6, 76.6774, 258.);
1037 dvol->DefineSection(3, 253.6, 76.6774, 258.);
1038 //
1039 dvol->DefineSection(4, 253.6, 74.5, 264.4);
1040 dvol->DefineSection(5, 260., 74.5, 264.4);
1041 //
1042 TGeoMedium* m5 = gGeoManager->GetMedium("TPC_DriftGas2");
1043 auto* v9 = new TGeoVolume("TPC_Drift", dvol, m5);
1044 //
1045 v1->AddNode(v9, 1);
1046 //
1047 // outer insulator
1048 //
1049 auto* tpco = new TGeoPcon(0., 360., 6); // insulator
1050 //
1051 tpco->DefineSection(0, -256.6, 264.8, 278.);
1052 tpco->DefineSection(1, -253.6, 264.8, 278.);
1053 //
1054 tpco->DefineSection(2, -253.6, 258., 278.);
1055 tpco->DefineSection(3, 250.6, 258., 278.);
1056 //
1057 tpco->DefineSection(4, 250.6, 258., 275.5);
1058 tpco->DefineSection(5, 253.6, 258., 275.5);
1059 //
1060 TGeoMedium* m2 = gGeoManager->GetMedium("TPC_CO2");
1061 auto* v2 = new TGeoVolume("TPC_OI", tpco, m2);
1062 //
1063 TGeoRotation* segrot; // segment rotations
1064 //
1065 // outer containment vessel
1066 //
1067 auto* tocv = new TGeoPcon(0., 360., 6); // containment vessel
1068 //
1069 tocv->DefineSection(0, -256.6, 264.8, 278.);
1070 tocv->DefineSection(1, -253.6, 264.8, 278.);
1071 //
1072 tocv->DefineSection(2, -253.6, 274.8124, 278.);
1073 tocv->DefineSection(3, 247.6, 274.8124, 278.);
1074 //
1075 tocv->DefineSection(4, 247.6, 270.4, 278.);
1076 tocv->DefineSection(5, 250.6, 270.4, 278.);
1077 //
1078 TGeoMedium* m3 = gGeoManager->GetMedium("TPC_Al");
1079 auto* v3 = new TGeoVolume("TPC_OCV", tocv, m3);
1080 //
1081 auto* to1 = new TGeoTubeSeg(274.8174, 277.995, 252.1, 0., 59.9); // epoxy
1082 auto* to2 = new TGeoTubeSeg(274.8274, 277.985, 252.1, 0., 59.9); // tedlar
1083 auto* to3 = new TGeoTubeSeg(274.8312, 277.9812, 252.1, 0., 59.9); // prepreg2
1084 auto* to4 = new TGeoTubeSeg(274.9062, 277.9062, 252.1, 0., 59.9); // nomex
1085 auto* tog5 = new TGeoTubeSeg(274.8174, 277.995, 252.1, 59.9, 60.); // epoxy
1086 //
1087 TGeoMedium* sm1 = gGeoManager->GetMedium("TPC_Epoxy");
1088 TGeoMedium* sm2 = gGeoManager->GetMedium("TPC_Tedlar");
1089 TGeoMedium* sm3 = gGeoManager->GetMedium("TPC_Prepreg2");
1090 TGeoMedium* sm4 = gGeoManager->GetMedium("TPC_Nomex");
1091 //
1092 TGeoMedium* smep = gGeoManager->GetMedium("TPC_Epoxy1");
1093 //
1094 auto* tov1 = new TGeoVolume("TPC_OCV1", to1, sm1);
1095 auto* tov2 = new TGeoVolume("TPC_OCV2", to2, sm2);
1096 auto* tov3 = new TGeoVolume("TPC_OCV3", to3, sm3);
1097 auto* tov4 = new TGeoVolume("TPC_OCV4", to4, sm4);
1098 auto* togv5 = new TGeoVolume("TPC_OCVG5", tog5, sm1);
1099 //
1100 TGeoMedium* mhs = gGeoManager->GetMedium("TPC_Steel");
1101 TGeoMedium* m12 = gGeoManager->GetMedium("TPC_Water");
1102 //-------------------------------------------------------
1103 // Tpc Outer Field Cage
1104 // daughters - composite (sandwich)
1105 //-------------------------------------------------------
1106
1107 auto* tofc = new TGeoPcon(0., 360., 6);
1108 //
1109 tofc->DefineSection(0, -253.6, 258., 269.6);
1110 tofc->DefineSection(1, -250.6, 258., 269.6);
1111 //
1112 tofc->DefineSection(2, -250.6, 258., 260.0676);
1113 tofc->DefineSection(3, 250.6, 258., 260.0676);
1114 //
1115 tofc->DefineSection(4, 250.6, 258., 275.5);
1116 tofc->DefineSection(5, 253.6, 258., 275.5);
1117 //
1118 auto* v4 = new TGeoVolume("TPC_TOFC", tofc, m3);
1119 // sandwich
1120 auto* tf1 = new TGeoTubeSeg(258.0, 260.0676, 252.1, 0., 59.9); // tedlar
1121 auto* tf2 = new TGeoTubeSeg(258.0038, 260.0638, 252.1, 0., 59.9); // prepreg3
1122 auto* tf3 = new TGeoTubeSeg(258.0338, 260.0338, 252.1, 0., 59.9); // nomex
1123 auto* tfg4 = new TGeoTubeSeg(258.0, 260.0676, 252.1, 59.9, 60.); // epoxy glue
1124 //
1125 TGeoMedium* sm5 = gGeoManager->GetMedium("TPC_Prepreg3");
1126 //
1127 auto* tf1v = new TGeoVolume("TPC_OFC1", tf1, sm2);
1128 auto* tf2v = new TGeoVolume("TPC_OFC2", tf2, sm5);
1129 auto* tf3v = new TGeoVolume("TPC_OFC3", tf3, sm4);
1130 auto* tfg4v = new TGeoVolume("TPC_OFCG4", tfg4, smep);
1131 //
1132 // outer part - positioning
1133 //
1134 tov1->AddNode(tov2, 1);
1135 tov2->AddNode(tov3, 1);
1136 tov3->AddNode(tov4, 1); // ocv
1137 //
1138 tf1v->AddNode(tf2v, 1);
1139 tf2v->AddNode(tf3v, 1); // ofc
1140 //
1141 auto* t200 = new TGeoVolumeAssembly("TPC_OCVSEG");
1142 auto* t300 = new TGeoVolumeAssembly("TPC_OFCSEG");
1143 //
1144 // assembly OCV and OFC
1145 //
1146 // 1st - no rotation
1147 t200->AddNode(tov1, 1);
1148 t200->AddNode(togv5, 1);
1149 t300->AddNode(tf1v, 1);
1150 t300->AddNode(tfg4v, 1);
1151 // 2nd - rotation 60 deg
1152 segrot = new TGeoRotation();
1153 segrot->RotateZ(60.);
1154 t200->AddNode(tov1, 2, segrot);
1155 t200->AddNode(togv5, 2, segrot);
1156 t300->AddNode(tf1v, 2, segrot);
1157 t300->AddNode(tfg4v, 2, segrot);
1158 // 3rd rotation 120 deg
1159 segrot = new TGeoRotation();
1160 segrot->RotateZ(120.);
1161 t200->AddNode(tov1, 3, segrot);
1162 t200->AddNode(togv5, 3, segrot);
1163 t300->AddNode(tf1v, 3, segrot);
1164 t300->AddNode(tfg4v, 3, segrot);
1165 // 4th rotation 180 deg
1166 segrot = new TGeoRotation();
1167 segrot->RotateZ(180.);
1168 t200->AddNode(tov1, 4, segrot);
1169 t200->AddNode(togv5, 4, segrot);
1170 t300->AddNode(tf1v, 4, segrot);
1171 t300->AddNode(tfg4v, 4, segrot);
1172 // 5th rotation 240 deg
1173 segrot = new TGeoRotation();
1174 segrot->RotateZ(240.);
1175 t200->AddNode(tov1, 5, segrot);
1176 t200->AddNode(togv5, 5, segrot);
1177 t300->AddNode(tf1v, 5, segrot);
1178 t300->AddNode(tfg4v, 5, segrot);
1179 // 6th rotation 300 deg
1180 segrot = new TGeoRotation();
1181 segrot->RotateZ(300.);
1182 t200->AddNode(tov1, 6, segrot);
1183 t200->AddNode(togv5, 6, segrot);
1184 t300->AddNode(tf1v, 6, segrot);
1185 t300->AddNode(tfg4v, 6, segrot);
1186 //
1187 v3->AddNode(t200, 1, new TGeoTranslation(0., 0., -1.5));
1188 v4->AddNode(t300, 1);
1189 //
1190 v2->AddNode(v3, 1);
1191 v2->AddNode(v4, 1);
1192 //
1193 v1->AddNode(v2, 1);
1194 //
1195 // Outer field cage guard rings. Inner placed in the drift gas, outer placed in the outer insulator (CO2)
1196 //
1197 auto* ogri = new TGeoTube(257.985, 258., 0.6); // placed in the drift volume
1198 auto* ogro = new TGeoTube(260.0676, 260.0826, 0.6); // placed in the outer insulator
1199 //
1200 auto* ogriv = new TGeoVolume("TPC_OGRI", ogri, m3);
1201 auto* ogrov = new TGeoVolume("TPC_OGRO", ogro, m3);
1202 //
1203 for (Int_t i = 0; i < 24; i++) {
1204 v9->AddNode(ogriv, (i + 1), new TGeoTranslation(0., 0., (i + 1) * 10));
1205 v9->AddNode(ogriv, (i + 25), new TGeoTranslation(0., 0., -(i + 1) * 10));
1206 v2->AddNode(ogrov, (i + 1), new TGeoTranslation(0., 0., (i + 1) * 10));
1207 v2->AddNode(ogrov, (i + 25), new TGeoTranslation(0., 0., -(i + 1) * 10));
1208 }
1209 //--------------------------------------------------------------------
1210 // Tpc Inner INsulator (CO2)
1211 // the cones, the central drum and the inner f.c. sandwich with a piece
1212 // of the flane will be placed in the TPC
1213 //--------------------------------------------------------------------
1214 auto* tpci = new TGeoPcon(0., 360., 4);
1215 //
1216 tpci->DefineSection(0, -253.6, 68.4, 76.6774);
1217 tpci->DefineSection(1, -74.0, 61.2, 76.6774);
1218 //
1219 tpci->DefineSection(2, 74.0, 61.2, 76.6774);
1220 //
1221 tpci->DefineSection(3, 253.6, 65.9, 76.6774);
1222 //
1223 auto* v5 = new TGeoVolume("TPC_INI", tpci, m2);
1224 //
1225 // now the inner field cage - only part of flanges (2 copies)
1226 //
1227 auto* tif1 = new TGeoTube(69.9, 76.6774, 1.5);
1228 auto* v6 = new TGeoVolume("TPC_IFC1", tif1, m3);
1229 //
1230 //---------------------------------------------------------
1231 // Tpc Inner Containment vessel - Muon side
1232 //---------------------------------------------------------
1233 auto* tcms = new TGeoPcon(0., 360., 10);
1234 //
1235 tcms->DefineSection(0, -259.1, 68.1, 74.2);
1236 tcms->DefineSection(1, -253.6, 68.1, 74.2);
1237 //
1238 tcms->DefineSection(2, -253.6, 68.1, 68.4);
1239 tcms->DefineSection(3, -74.0, 60.9, 61.2);
1240 //
1241 tcms->DefineSection(4, -74.0, 60.1, 61.2);
1242 tcms->DefineSection(5, -73.3, 60.1, 61.2);
1243 //
1244 tcms->DefineSection(6, -73.3, 56.9, 61.2);
1245 tcms->DefineSection(7, -73.0, 56.9, 61.2);
1246 //
1247 tcms->DefineSection(8, -73.0, 56.9, 58.8);
1248 tcms->DefineSection(9, -71.3, 56.9, 58.8);
1249 //
1250 auto* v7 = new TGeoVolume("TPC_ICVM", tcms, m3);
1251 //------------------------------------------------
1252 // Heat screen muon side
1253 //------------------------------------------------
1254
1255 auto* thsm = new TGeoCone(89.8, 67.88, 68.1, 60.68, 60.9);
1256 auto* thsmw = new TGeoCone(89.8, 67.94, 68.04, 60.74, 60.84);
1257 auto* hvsm = new TGeoVolume("TPC_HSM", thsm, mhs); // steel
1258 auto* hvsmw = new TGeoVolume("TPC_HSMW", thsmw, m12); // water
1259 // assembly heat screen muon
1260 hvsm->AddNode(hvsmw, 1);
1261 //-----------------------------------------------
1262 // inner containment vessel - shaft side
1263 //-----------------------------------------------
1264 auto* tcss = new TGeoPcon(0., 360., 10);
1265 //
1266 tcss->DefineSection(0, 71.3, 56.9, 58.8);
1267 tcss->DefineSection(1, 73.0, 56.9, 58.8);
1268 //
1269 tcss->DefineSection(2, 73.0, 56.9, 61.2);
1270 tcss->DefineSection(3, 73.3, 56.9, 61.2);
1271 //
1272 tcss->DefineSection(4, 73.3, 60.1, 61.2);
1273 tcss->DefineSection(5, 74.0, 60.1, 61.2);
1274 //
1275 tcss->DefineSection(6, 74.0, 60.9, 61.2);
1276 tcss->DefineSection(7, 253.6, 65.6, 65.9);
1277 //
1278 tcss->DefineSection(8, 253.6, 65.6, 74.2);
1279 tcss->DefineSection(9, 258.1, 65.6, 74.2);
1280 //
1281 auto* v8 = new TGeoVolume("TPC_ICVS", tcss, m3);
1282 //-------------------------------------------------
1283 // Heat screen shaft side
1284 //--------------------------------------------------
1285 auto* thss = new TGeoCone(89.8, 60.68, 60.9, 65.38, 65.6);
1286 auto* thssw = new TGeoCone(89.8, 60.74, 60.84, 65.44, 65.54);
1287 auto* hvss = new TGeoVolume("TPC_HSS", thss, mhs); // steel
1288 auto* hvssw = new TGeoVolume("TPC_HSSW", thssw, m12); // water
1289 // assembly heat screen shaft
1290 hvss->AddNode(hvssw, 1);
1291 //-----------------------------------------------
1292 // Inner field cage
1293 // define 4 parts and make an assembly
1294 //-----------------------------------------------
1295 // part1 - Al - 2 copies
1296 auto* t1 = new TGeoTube(76.6774, 78.845, 0.75);
1297 auto* tv1 = new TGeoVolume("TPC_IFC2", t1, m3);
1298 // sandwich - outermost parts - 2 copies
1299 //
1300 // segment outermost
1301 //
1302 auto* t2 = new TGeoTubeSeg(76.6774, 78.845, 74.175, 350., 109.4); // tedlar 38 microns
1303 auto* t3 = new TGeoTubeSeg(76.6812, 78.8412, 74.175, 350., 109.4); // prepreg2 500 microns
1304 auto* t4 = new TGeoTubeSeg(76.7312, 78.7912, 74.175, 350., 109.4); // prepreg3 300 microns
1305 auto* t5 = new TGeoTubeSeg(76.7612, 78.7612, 74.175, 350., 109.4); // nomex 2 cm
1306 auto* tepox1 = new TGeoTubeSeg(76.6774, 78.845, 74.175, 109.4, 110.); // epoxy
1307 auto* tpr1 = new TGeoTubeSeg(78.845, 78.885, 74.175, 109., 111.);
1308
1309 // volumes for the outer part
1310 auto* tv2 = new TGeoVolume("TPC_IFC3", t2, sm2);
1311 auto* tv3 = new TGeoVolume("TPC_IFC4", t3, sm3);
1312 auto* tv4 = new TGeoVolume("TPC_IFC5", t4, sm5);
1313 auto* tv5 = new TGeoVolume("TPC_IFC6", t5, sm4);
1314 auto* tvep1 = new TGeoVolume("TPC_IFEPOX1", tepox1, smep);
1315 auto* tvpr1 = new TGeoVolume("TPC_PRSTR1", tpr1, sm2);
1316 //
1317 // middle parts - 2 copies
1318 //
1319 // segment middle
1320 //
1321 auto* t6 = new TGeoTubeSeg(76.6774, 78.795, 5., 350., 109.4); // tedlar 38 microns
1322 auto* t7 = new TGeoTubeSeg(76.6812, 78.7912, 5., 350., 109.4); // prepreg2 250 microns
1323 auto* t8 = new TGeoTubeSeg(76.7062, 78.7662, 5., 350., 109.4); // prepreg3 300 microns
1324 auto* t9 = new TGeoTubeSeg(76.7362, 78.7362, 5., 350., 109.4); // nomex 2 cm
1325 auto* tepox2 = new TGeoTubeSeg(76.6774, 78.795, 5., 109.4, 110.); // epoxy
1326 auto* tpr2 = new TGeoTubeSeg(78.795, 78.835, 5., 109., 111.);
1327 // volumes for the middle part
1328 auto* tv6 = new TGeoVolume("TPC_IFC7", t6, sm2);
1329 auto* tv7 = new TGeoVolume("TPC_IFC8", t7, sm3);
1330 auto* tv8 = new TGeoVolume("TPC_IFC9", t8, sm5);
1331 auto* tv9 = new TGeoVolume("TPC_IFC10", t9, sm4);
1332 auto* tvep2 = new TGeoVolume("TPC_IFEPOX2", tepox2, smep);
1333 auto* tvpr2 = new TGeoVolume("TPC_PRSTR2", tpr2, sm2);
1334 // central part - 1 copy
1335 //
1336 // segment central part
1337 //
1338 auto* t10 = new TGeoTubeSeg(76.6774, 78.785, 93.75, 350., 109.4); // tedlar 38 microns
1339 auto* t11 = new TGeoTubeSeg(76.6812, 78.7812, 93.75, 350., 109.4); // prepreg3 500 microns
1340 auto* t12 = new TGeoTubeSeg(76.7312, 78.7312, 93.75, 350., 109.4); // nomex 2 cm
1341 auto* tepox3 = new TGeoTubeSeg(76.6774, 78.785, 93.75, 109.4, 110.); // epoxy
1342 auto* tpr3 = new TGeoTubeSeg(78.785, 78.825, 93.75, 109., 111.);
1343 // volumes for the central part
1344 auto* tv10 = new TGeoVolume("TPC_IFC11", t10, sm2);
1345 auto* tv11 = new TGeoVolume("TPC_IFC12", t11, sm5);
1346 auto* tv12 = new TGeoVolume("TPC_IFC13", t12, sm4);
1347 auto* tvep3 = new TGeoVolume("TPC_IFEPOX3", tepox3, smep);
1348 auto* tvpr3 = new TGeoVolume("TPC_PRSTR3", tpr3, sm2);
1349 //
1350 // creating a sandwich for the outer par,t tv2 is the mother
1351 //
1352 tv2->AddNode(tv3, 1);
1353 tv3->AddNode(tv4, 1);
1354 tv4->AddNode(tv5, 1);
1355 //
1356 // creating a sandwich for the middle part, tv6 is the mother
1357 //
1358 tv6->AddNode(tv7, 1);
1359 tv7->AddNode(tv8, 1);
1360 tv8->AddNode(tv9, 1);
1361 //
1362 // creating a sandwich for the central part, tv10 is the mother
1363 //
1364 tv10->AddNode(tv11, 1);
1365 tv11->AddNode(tv12, 1);
1366 //
1367 auto* tv100 = new TGeoVolumeAssembly("TPC_IFC"); // ifc itself - 3 segments
1368
1369 //
1370 // first segment - no rotation
1371 //
1372 // central
1373 tv100->AddNode(tv10, 1); // sandwich
1374 tv100->AddNode(tvep3, 1); // epoxy
1375 tv100->AddNode(tvpr3, 1); // prepreg strip
1376 // middle
1377 tv100->AddNode(tv6, 1, new TGeoTranslation(0., 0., -98.75)); // sandwich1
1378 tv100->AddNode(tv6, 2, new TGeoTranslation(0., 0., 98.75)); // sandwich2
1379 tv100->AddNode(tvep2, 1, new TGeoTranslation(0., 0., -98.75)); // epoxy
1380 tv100->AddNode(tvep2, 2, new TGeoTranslation(0., 0., 98.75)); // epoxy
1381 tv100->AddNode(tvpr2, 1, new TGeoTranslation(0., 0., -98.75)); // prepreg strip
1382 tv100->AddNode(tvpr2, 2, new TGeoTranslation(0., 0., 98.75));
1383 // outer
1384 tv100->AddNode(tv2, 1, new TGeoTranslation(0., 0., -177.925)); // sandwich
1385 tv100->AddNode(tv2, 2, new TGeoTranslation(0., 0., 177.925));
1386 tv100->AddNode(tvep1, 1, new TGeoTranslation(0., 0., -177.925)); // epoxy
1387 tv100->AddNode(tvep1, 2, new TGeoTranslation(0., 0., 177.925));
1388 tv100->AddNode(tvpr1, 1, new TGeoTranslation(0., 0., -177.925)); // prepreg strip
1389 tv100->AddNode(tvpr1, 2, new TGeoTranslation(0., 0., -177.925));
1390 //
1391 // second segment - rotation 120 deg.
1392 //
1393 segrot = new TGeoRotation();
1394 segrot->RotateZ(120.);
1395 //
1396 // central
1397 tv100->AddNode(tv10, 2, segrot); // sandwich
1398 tv100->AddNode(tvep3, 2, segrot); // epoxy
1399 tv100->AddNode(tvpr3, 2, segrot); // prepreg strip
1400 // middle
1401 tv100->AddNode(tv6, 3, new TGeoCombiTrans(0., 0., -98.75, segrot)); // sandwich1
1402 tv100->AddNode(tv6, 4, new TGeoCombiTrans(0., 0., 98.75, segrot)); // sandwich2
1403 tv100->AddNode(tvep2, 3, new TGeoCombiTrans(0., 0., -98.75, segrot)); // epoxy
1404 tv100->AddNode(tvep2, 4, new TGeoCombiTrans(0., 0., 98.75, segrot)); // epoxy
1405 tv100->AddNode(tvpr2, 3, new TGeoCombiTrans(0., 0., -98.75, segrot)); // prepreg strip
1406 tv100->AddNode(tvpr2, 4, new TGeoCombiTrans(0., 0., 98.75, segrot));
1407 // outer
1408 tv100->AddNode(tv2, 3, new TGeoCombiTrans(0., 0., -177.925, segrot)); // sandwich
1409 tv100->AddNode(tv2, 4, new TGeoCombiTrans(0., 0., 177.925, segrot));
1410 tv100->AddNode(tvep1, 3, new TGeoCombiTrans(0., 0., -177.925, segrot)); // epoxy
1411 tv100->AddNode(tvep1, 4, new TGeoCombiTrans(0., 0., 177.925, segrot));
1412 tv100->AddNode(tvpr1, 3, new TGeoCombiTrans(0., 0., -177.925, segrot)); // prepreg strip
1413 tv100->AddNode(tvpr1, 4, new TGeoCombiTrans(0., 0., 177.925, segrot));
1414 //
1415 // third segment - rotation 240 deg.
1416 //
1417 segrot = new TGeoRotation();
1418 segrot->RotateZ(240.);
1419 //
1420 // central
1421 tv100->AddNode(tv10, 3, segrot); // sandwich
1422 tv100->AddNode(tvep3, 3, segrot); // epoxy
1423 tv100->AddNode(tvpr3, 3, segrot); // prepreg strip
1424 // middle
1425 tv100->AddNode(tv6, 5, new TGeoCombiTrans(0., 0., -98.75, segrot)); // sandwich1
1426 tv100->AddNode(tv6, 6, new TGeoCombiTrans(0., 0., 98.75, segrot)); // sandwich2
1427 tv100->AddNode(tvep2, 5, new TGeoCombiTrans(0., 0., -98.75, segrot)); // epoxy
1428 tv100->AddNode(tvep2, 6, new TGeoCombiTrans(0., 0., 98.75, segrot)); // epoxy
1429 tv100->AddNode(tvpr2, 5, new TGeoCombiTrans(0., 0., -98.75, segrot)); // prepreg strip
1430 tv100->AddNode(tvpr2, 6, new TGeoCombiTrans(0., 0., 98.75, segrot));
1431 // outer
1432 tv100->AddNode(tv2, 5, new TGeoCombiTrans(0., 0., -177.925, segrot)); // sandwich
1433 tv100->AddNode(tv2, 6, new TGeoCombiTrans(0., 0., 177.925, segrot));
1434 tv100->AddNode(tvep1, 5, new TGeoCombiTrans(0., 0., -177.925, segrot)); // epoxy
1435 tv100->AddNode(tvep1, 6, new TGeoCombiTrans(0., 0., 177.925, segrot));
1436 tv100->AddNode(tvpr1, 5, new TGeoCombiTrans(0., 0., -177.925, segrot)); // prepreg strip
1437 tv100->AddNode(tvpr1, 6, new TGeoCombiTrans(0., 0., 177.925, segrot));
1438 // Al parts - rings
1439 tv100->AddNode(tv1, 1, new TGeoTranslation(0., 0., -252.85));
1440 tv100->AddNode(tv1, 2, new TGeoTranslation(0., 0., 252.85));
1441 //
1442 v5->AddNode(v6, 1, new TGeoTranslation(0., 0., -252.1));
1443 v5->AddNode(v6, 2, new TGeoTranslation(0., 0., 252.1));
1444 v1->AddNode(v5, 1);
1445 v1->AddNode(v7, 1);
1446 v1->AddNode(v8, 1);
1447 v1->AddNode(hvsm, 1, new TGeoTranslation(0., 0., -163.8));
1448 v1->AddNode(hvss, 1, new TGeoTranslation(0., 0., 163.8));
1449 v9->AddNode(tv100, 1);
1450 //
1451 // guard rings for IFC - outer placed in inner insulator, inner placed in the drift gas (3 different radii)
1452 // AL, 1.2 cm wide, 0.015 cm thick, volumes TPC_IGR1 - outer, TPC_IGR2-4 - inner
1453 //
1454 auto* igro = new TGeoTube(76.6624, 76.6774, 0.6); // inner part, ends at inner radius of the IFC
1455 auto* igrio = new TGeoTube(78.845, 78.86, 0.6); // outer part
1456 auto* igrim = new TGeoTube(78.795, 78.81, 0.6);
1457 auto* igric = new TGeoTube(78.785, 78.8, 0.6);
1458 //
1459 // volumes
1460 //
1461 auto* igrov = new TGeoVolume("TPC_IGR1", igro, m3);
1462 auto* igriov = new TGeoVolume("TPC_IGR2", igrio, m3);
1463 auto* igrimv = new TGeoVolume("TPC_IGR3", igrim, m3);
1464 auto* igricv = new TGeoVolume("TPC_IGR4", igric, m3);
1465 //
1466 // outer guard rings for IFC placement - every 10 cm
1467 //
1468 for (Int_t i = 0; i < 24; i++) {
1469 v5->AddNode(igrov, (i + 1), new TGeoTranslation(0., 0., (i + 1) * 10));
1470 v5->AddNode(igrov, (i + 25), new TGeoTranslation(0., 0., -(i + 1) * 10));
1471 }
1472 //
1473 // inner guard rings for IFC placement
1474 //
1475 for (Int_t i = 0; i < 9; i++) {
1476 v9->AddNode(igricv, (i + 1), new TGeoTranslation(0., 0., (i + 1) * 10));
1477 v9->AddNode(igricv, (i + 10), new TGeoTranslation(0., 0., -(i + 1) * 10));
1478 }
1479 v9->AddNode(igrimv, 1, new TGeoTranslation(0., 0., 100.));
1480 v9->AddNode(igrimv, 2, new TGeoTranslation(0., 0., -100.));
1481 //
1482 for (Int_t i = 0; i < 13; i++) {
1483 v9->AddNode(igriov, i + 1, new TGeoTranslation(0., 0., 100 + (i + 1) * 10));
1484 v9->AddNode(igriov, i + 14, new TGeoTranslation(0., 0., -(100 + (i + 1) * 10)));
1485 }
1486 //
1487 // central drum
1488 //
1489 // flange + sandwich
1490 //
1491 auto* cfl = new TGeoPcon(0., 360., 6);
1492 cfl->DefineSection(0, -71.1, 59.7, 61.2);
1493 cfl->DefineSection(1, -68.6, 59.7, 61.2);
1494 //
1495 cfl->DefineSection(2, -68.6, 60.6124, 61.2);
1496 cfl->DefineSection(3, 68.6, 60.6124, 61.2);
1497 //
1498 cfl->DefineSection(4, 68.6, 59.7, 61.2);
1499 cfl->DefineSection(5, 71.1, 59.7, 61.2);
1500 //
1501 auto* cflv = new TGeoVolume("TPC_CDR", cfl, m3);
1502 // sandwich
1503 auto* cd1 = new TGeoTubeSeg(60.6224, 61.19, 69.8, 0.05, 119.95);
1504 auto* cd2 = new TGeoTubeSeg(60.6262, 61.1862, 69.8, 0.05, 119.95);
1505 auto* cd3 = new TGeoTubeSeg(60.6462, 61.1662, 69.8, 0.05, 119.95);
1506 auto* cd4 = new TGeoTubeSeg(60.6562, 61.1562, 69.8, 0.05, 119.95);
1507 auto* tepox4 = new TGeoTubeSeg(60.6224, 61.19, 69.8, 359.95, 0.05); // epoxy glue 0.01 deg
1508 //
1509 TGeoMedium* sm6 = gGeoManager->GetMedium("TPC_Prepreg1");
1510 TGeoMedium* sm8 = gGeoManager->GetMedium("TPC_Epoxyfm");
1511 auto* cd1v = new TGeoVolume("TPC_CDR1", cd1, sm2); // tedlar
1512 auto* cd2v = new TGeoVolume("TPC_CDR2", cd2, sm6); // prepreg1
1513 auto* cd3v = new TGeoVolume("TPC_CDR3", cd3, sm8); // epoxy film
1514 auto* cd4v = new TGeoVolume("TPC_CDR4", cd4, sm4); // nomex
1515 auto* tvep4 = new TGeoVolume("TPC_IFEPOX4", tepox4, smep); // epoxy glue
1516 //
1517 // joints between sections 1 deg prepreg1 placed in nomex at lower and upper radius + 0.1 deg of glue (epoxy)
1518 //
1519 auto* cdjl = new TGeoTubeSeg(60.6562, 60.6762, 69.8, 0., 1.0); // lower, to be rotated when positioned
1520 auto* cdju = new TGeoTubeSeg(61.1362, 61.1562, 69.8, 0., 1.0); // upper, to be rotated when positioned
1521 //
1522 auto* cdjlv = new TGeoVolume("TPC_CDJL", cdjl, sm6);
1523 auto* cdjuv = new TGeoVolume("TPC_CDJU", cdju, sm6);
1524 //
1525 // seals for central drum 2 copies
1526 //
1527 auto* cs = new TGeoTube(56.9, 61.2, 0.1);
1528 TGeoMedium* sm7 = gGeoManager->GetMedium("TPC_Mylar");
1529 auto* csv = new TGeoVolume("TPC_CDRS", cs, sm7);
1530 v1->AddNode(csv, 1, new TGeoTranslation(0., 0., -71.2));
1531 v1->AddNode(csv, 2, new TGeoTranslation(0., 0., 71.2));
1532 //
1533 // seal collars
1534 auto* se = new TGeoPcon(0., 360., 6);
1535 se->DefineSection(0, -72.8, 59.7, 61.2);
1536 se->DefineSection(1, -72.3, 59.7, 61.2);
1537 //
1538 se->DefineSection(2, -72.3, 58.85, 61.2);
1539 se->DefineSection(3, -71.6, 58.85, 61.2);
1540 //
1541 se->DefineSection(4, -71.6, 59.7, 61.2);
1542 se->DefineSection(5, -71.3, 59.7, 61.2);
1543 //
1544 auto* sev = new TGeoVolume("TPC_CDCE", se, m3);
1545 //
1546 auto* si = new TGeoTube(56.9, 58.8, 1.);
1547 auto* siv = new TGeoVolume("TPC_CDCI", si, m3);
1548 //
1549 // define reflection matrix
1550 //
1551 auto* ref = new TGeoRotation("ref", 90., 0., 90., 90., 180., 0.);
1552 //
1553 cd1v->AddNode(cd2v, 1);
1554 cd2v->AddNode(cd3v, 1);
1555 cd3v->AddNode(cd4v, 1); // sandwich
1556 //
1557 // joints, lower part and upper parts, placed in nomex
1558 //
1559 segrot = new TGeoRotation();
1560 segrot->RotateZ(0.05);
1561 cd4v->AddNode(cdjlv, 1, segrot);
1562 cd4v->AddNode(cdjuv, 1, segrot);
1563 segrot = new TGeoRotation();
1564 segrot->RotateZ(118.95);
1565 cd4v->AddNode(cdjlv, 2, segrot);
1566 cd4v->AddNode(cdjuv, 2, segrot);
1567 //
1568 // according to the conversion data, segments have an additionaly rotated by 4.6 deg
1569 //
1570 // first segment
1571 segrot = new TGeoRotation();
1572 segrot->RotateZ(4.6);
1573 cflv->AddNode(cd1v, 1, segrot);
1574 cflv->AddNode(tvep4, 1, segrot);
1575 // second segment
1576 segrot = new TGeoRotation();
1577 segrot->RotateZ(124.6);
1578 cflv->AddNode(cd1v, 2, segrot);
1579 cflv->AddNode(tvep4, 2, segrot);
1580 // third segment
1581 segrot = new TGeoRotation();
1582 segrot->RotateZ(244.6);
1583 cflv->AddNode(cd1v, 3, segrot);
1584 cflv->AddNode(tvep4, 3, segrot);
1585 //
1586 // heating strips
1587 //
1588 auto* hstr = new TGeoTubeSeg(60.6124, 60.6224, 68.5, 0., 1.25);
1589 auto* hstrv = new TGeoVolume("TPC_HSTR", hstr, m1); // air, caved out from cflv, first strip starts at 0 deg
1590 for (Int_t i = 0; i < 144; i++) {
1591 Double_t alpha = 1.25 + i * 2.5;
1592 segrot = new TGeoRotation();
1593 segrot->RotateZ(alpha);
1594 cflv->AddNode(hstrv, i + 1, segrot);
1595 }
1596 //
1597 v1->AddNode(siv, 1, new TGeoTranslation(0., 0., -69.9));
1598 v1->AddNode(siv, 2, new TGeoTranslation(0., 0., 69.9));
1599 v1->AddNode(sev, 1);
1600 v1->AddNode(sev, 2, ref);
1601 v1->AddNode(cflv, 1);
1602 //
1603 // central membrane - 2 rings and a mylar membrane - assembly
1604 //
1605 auto* ih = new TGeoTube(81.05, 84.05, 0.3);
1606 auto* oh = new TGeoTube(250., 256., 0.5);
1607 auto* mem = new TGeoTube(84.05, 250., 0.00115);
1608
1609 //
1610 TGeoMedium* m4 = gGeoManager->GetMedium("TPC_G10");
1611 //
1612 auto* ihv = new TGeoVolume("TPC_IHVH", ih, m3);
1613 auto* ohv = new TGeoVolume("TPC_OHVH", oh, m3);
1614
1615 auto* memv = new TGeoVolume("TPC_HV", mem, sm7);
1616 //
1617 auto* cm = new TGeoVolumeAssembly("TPC_HVMEM");
1618 cm->AddNode(ihv, 1);
1619 cm->AddNode(ohv, 1);
1620 cm->AddNode(memv, 1);
1621
1622 v9->AddNode(cm, 1);
1623 //
1624 // end caps - they are make as an assembly of single segments
1625 // containing both readout chambers
1626 //
1627 Double_t openingAngle = 10. * TMath::DegToRad();
1628 Double_t thick = 1.5; // rib
1629 Double_t shift = thick / TMath::Sin(openingAngle);
1630 //
1631 Double_t lowEdge = 86.3; // hole in the wheel
1632 Double_t upEdge = 240.4; // hole in the wheel
1633 //
1634 new TGeoTubeSeg("tpc_ssec", 74.5, 264.4, 3., 0., 20.);
1635 //
1636 auto* tpc_hole = new TGeoPgon("tpc_hole", 0., 20., 1, 4);
1637 //
1638 tpc_hole->DefineSection(0, -3.5, lowEdge - shift, upEdge - shift);
1639 tpc_hole->DefineSection(1, -1.5, lowEdge - shift, upEdge - shift);
1640 //
1641 tpc_hole->DefineSection(2, -1.5, lowEdge - shift, upEdge + 3. - shift);
1642 tpc_hole->DefineSection(3, 3.5, lowEdge - shift, upEdge + 3. - shift);
1643 //
1644 Double_t ys = shift * TMath::Sin(openingAngle);
1645 Double_t xs = shift * TMath::Cos(openingAngle);
1646 auto* tr = new TGeoTranslation("tr", xs, ys, 0.);
1647 tr->RegisterYourself();
1648 auto* chamber = new TGeoCompositeShape("tpc_ssec-tpc_hole:tr");
1649 auto* sv = new TGeoVolume("TPC_WSEG", chamber, m3);
1650 auto* bar = new TGeoPgon("bar", 0., 20., 1, 2);
1651 bar->DefineSection(0, -3., 131.5 - shift, 136.5 - shift);
1652 bar->DefineSection(1, 1.5, 131.5 - shift, 136.5 - shift);
1653 auto* barv = new TGeoVolume("TPC_WBAR", bar, m3);
1654 auto* ch = new TGeoVolumeAssembly("TPC_WCH"); // empty segment
1655 //
1656 ch->AddNode(sv, 1);
1657 ch->AddNode(barv, 1, tr);
1658 //
1659 // readout chambers
1660 //
1661 // IROC first
1662 //
1663 TGeoMedium* m6 = gGeoManager->GetMedium("TPC_Makrolon");
1664 //
1665 auto* ibody = new TGeoTrd1(13.8742, 21.3328, 4.29, 21.15);
1666 auto* ibdv = new TGeoVolume("TPC_IROCB", ibody, m3);
1667 // empty space
1668 auto* emp = new TGeoTrd1(12.3742, 19.8328, 4.05, 19.65);
1669 auto* empv = new TGeoVolume("TPC_IROCE", emp, m1);
1670 ibdv->AddNode(empv, 1, new TGeoTranslation(0., -0.24, 0.));
1671 // bars
1672 Double_t tga = (19.8328 - 12.3742) / 39.3;
1673 Double_t xmin, xmax;
1674 // 1
1675 xmin = 9.65 * tga + 12.3742;
1676 xmax = 10.05 * tga + 12.3742;
1677 //
1678 auto* ib1 = new TGeoTrd1(xmin, xmax, 2.06, 0.2);
1679 auto* ib1v = new TGeoVolume("TPC_IRB1", ib1, m3);
1680 empv->AddNode(ib1v, 1, new TGeoTranslation("tt1", 0., 1.99, -9.8));
1681 // 2
1682 xmin = 19.5 * tga + 12.3742;
1683 xmax = 19.9 * tga + 12.3742;
1684 //
1685 auto* ib2 = new TGeoTrd1(xmin, xmax, 2.06, 0.2);
1686 auto* ib2v = new TGeoVolume("TPC_TRB2", ib2, m3);
1687 empv->AddNode(ib2v, 1, new TGeoTranslation(0., 1.99, 0.05));
1688 // 3
1689 xmin = 29.35 * tga + 12.3742;
1690 xmax = 29.75 * tga + 12.3742;
1691 //
1692 auto* ib3 = new TGeoTrd1(xmin, xmax, 2.06, 0.2);
1693 auto* ib3v = new TGeoVolume("TPC_IRB3", ib3, m3);
1694 empv->AddNode(ib3v, 1, new TGeoTranslation(0., 1.99, 9.9));
1695 //
1696 // connectors alu body & strong back
1697 //
1698 auto* conn = new TGeoBBox(0.4, 0.24, 1.937); // connectors alu body
1699 auto* connv = new TGeoVolume("TPC_RCCON", conn, m6); // makrolon
1700 //
1701 auto* connb = new TGeoBBox(0.5, 0.25, 2.375); // connectors strong back
1702 auto* connbv = new TGeoVolume("TPC_RCCONB", connb, m6);
1703 //
1704 // strong back
1705 //
1706 auto* icsb = new TGeoTrd1(14.5974, 23.3521, 0.25, 24.825);
1707 auto* icsbv = new TGeoVolume("TPC_ISB", icsb, m4);
1708 //
1709 // file with positions of connectors
1710 //
1711 TString fileName(gSystem->Getenv("VMCWORKDIR"));
1712 fileName += "/Detectors/Geometry/TPC/conn_iroc.dat";
1713 ifstream in;
1714 in.open(fileName.Data(), ios_base::in); // asci file
1715 if (!in.is_open()) {
1716 LOG(fatal) << "Cannot open input file : " << fileName.Data();
1717 }
1718 for (Int_t i = 0; i < 132; i++) {
1719 Double_t x, z, ang;
1720 in >> ang >> x >> z;
1721 //
1722 ang = -ang;
1723 z -= 1102.; // changing the reference frame from the beam to the volume
1724 z *= 0.1;
1725 x *= 0.1;
1726 //
1727 auto* rrr = new TGeoRotation();
1728 rrr->RotateY(ang);
1729 //
1730 ibdv->AddNode(connv, i + 1, new TGeoCombiTrans(x, 4.05, z, rrr)); // connectors alu body
1731 icsbv->AddNode(connbv, i + 1, new TGeoCombiTrans(x, 0., z + 1.725, rrr)); // connectors strong back
1732 }
1733 in.close();
1734 //
1735 // "cap"
1736 new TGeoTrd1("icap", 14.5974, 23.3521, 1.19, 24.825);
1737 // "hole"
1738 new TGeoTrd1("ihole", 13.8742, 21.3328, 1.2, 21.15);
1739 auto* tr1 = new TGeoTranslation("tr1", 0., 0., 1.725);
1740 tr1->RegisterYourself();
1741 auto* ic = new TGeoCompositeShape("icap-ihole:tr1");
1742 auto* icv = new TGeoVolume("TPC_IRCAP", ic, m3);
1743 //
1744 // pad plane G10 3.2 mm thick
1745 //
1746 auto* icpp = new TGeoTrd1(14.5974, 23.3521, 0.16, 24.825);
1747 auto* icppv = new TGeoVolume("TPC_IPP", icpp, m4); // pad plane
1748 //
1749 // gem
1750 //
1751 new TGeoTrd1("igem", 14.5974, 23.3521, 0.1, 24.825);
1752 new TGeoTrd1("igemh", 14.5974 - .5, 23.3521 - .5, 0.11, 24.825 - .5);
1753 auto* icgem = new TGeoCompositeShape("igem-igemh");
1754 auto* icgemv = new TGeoVolume("TPC_ICGEM", icgem, m4);
1755 //
1756 // assembly of the iroc
1757 //
1758 auto* iroc = new TGeoVolumeAssembly("TPC_IROC");
1759 //
1760 iroc->AddNode(ibdv, 1); // main body
1761 iroc->AddNode(icv, 1, new TGeoTranslation(0., 3.1, -1.725)); // cap
1762 iroc->AddNode(icsbv, 1, new TGeoTranslation(0., 4.54, -1.725)); // strong back
1763 iroc->AddNode(icppv, 1, new TGeoTranslation(0., 4.95, -1.725)); // pad plane
1764 iroc->AddNode(icgemv, 1, new TGeoTranslation(0., 5.21, -1.725)); // gem
1765 //
1766 // OROC
1767 //
1768 auto* obody = new TGeoTrd1(22.2938, 40.5084, 4.29, 51.65);
1769 auto* obdv = new TGeoVolume("TPC_OROCB", obody, m3);
1770 auto* oemp = new TGeoTrd1(20.7938, 39.0084, 3.89, 50.15);
1771 auto* oempv = new TGeoVolume("TPC_OROCE", oemp, m1);
1772 obdv->AddNode(oempv, 1, new TGeoTranslation(0., -0.4, 0.));
1773 //
1774 // horizontal bars
1775 //
1776 tga = (39.0084 - 20.7938) / 100.3;
1777 xmin = tga * 14.2 + 20.7938;
1778 xmax = tga * 14.6 + 20.7938;
1779 auto* ob1 = new TGeoTrd1(xmin, xmax, 2.94, 0.2);
1780 auto* ob1v = new TGeoVolume("TPC_ORB1", ob1, m3);
1781 //
1782 xmin = 30.4 * tga + 20.7938;
1783 xmax = 32.1 * tga + 20.7938;
1784 auto* ob2 = new TGeoTrd1(xmin, xmax, 2.94, 1.05);
1785 auto* ob2v = new TGeoVolume("TPC_ORB2", ob2, m3);
1786 //
1787 xmin = 51.5 * tga + 20.7938;
1788 xmax = 51.9 * tga + 20.7938;
1789 auto* ob3 = new TGeoTrd1(xmin, xmax, 2.94, 0.2);
1790 auto* ob3v = new TGeoVolume("TPC_ORB3", ob3, m3);
1791 //
1792 xmin = 68.5 * tga + 20.7938;
1793 xmax = 70.6 * tga + 20.7938;
1794 auto* ob4 = new TGeoTrd1(xmin, xmax, 2.94, 1.05);
1795 auto* ob4v = new TGeoVolume("TPC_ORB4", ob4, m3);
1796 //
1797 xmin = 89.9 * tga + 20.7938;
1798 xmax = 90.3 * tga + 20.7938;
1799 auto* ob5 = new TGeoTrd1(xmin, xmax, 2.94, 0.2);
1800 auto* ob5v = new TGeoVolume("TPC_ORB5", ob5, m3);
1801 //
1802 oempv->AddNode(ob1v, 1, new TGeoTranslation(0., 0.59, -35.75));
1803 oempv->AddNode(ob2v, 1, new TGeoTranslation(0., 0.59, -18.7));
1804 oempv->AddNode(ob3v, 1, new TGeoTranslation(0., 0.59, 1.55));
1805 oempv->AddNode(ob4v, 1, new TGeoTranslation(0., 0.59, 19.4));
1806 oempv->AddNode(ob5v, 1, new TGeoTranslation(0., 0.59, 39.95));
1807 //
1808 // connectors, identical as for iroc, but I prefer to have separate volumes for better control
1809 //
1810 auto* conno = new TGeoBBox(0.4, 0.4, 1.937); // alu body
1811 auto* connov = new TGeoVolume("TPC_RCCONO", conno, m6); // makrolon
1812 //
1813 auto* connob = new TGeoBBox(0.5, 0.25, 2.375);
1814 auto* connobv = new TGeoVolume("TPC_RCCONOB", connob, m6); // strong back
1815 //
1816 // cap
1817 //
1818 new TGeoTrd1("ocap", 23.3875, 43.524, 1.19, 57.1);
1819 new TGeoTrd1("ohole", 22.2938, 40.5084, 1.19, 51.65);
1820 auto* tr5 = new TGeoTranslation("tr5", 0., 0., -2.15);
1821 tr5->RegisterYourself();
1822 auto* oc = new TGeoCompositeShape("ocap-ohole:tr5");
1823 auto* ocv = new TGeoVolume("TPC_ORCAP", oc, m3);
1824 //
1825 // stron back 5 mm
1826 //
1827 auto* osb = new TGeoTrd1(23.3874, 43.524, 0.25, 57.1);
1828 auto* osbv = new TGeoVolume("TPC_OSB", osb, m4);
1829 //
1830 // pad plane 3.2 mm
1831 //
1832 auto* opp = new TGeoTrd1(23.3874, 43.524, 0.16, 57.1);
1833 auto* oppv = new TGeoVolume("TPC_OPP", opp, m4);
1834 //
1835 // gem
1836 //
1837 new TGeoTrd1("ogem", 23.3874, 43.524, 0.1, 57.1);
1838 //
1839 // ogemh1 - first "hole" ends at z = 36.25, ogemh2 - second "hole" starts at z = 74.15
1840 //
1841 new TGeoTrd1("ogemh1", 22.548, 28.579, 0.1, 17.625); // ogemh1 - lower hole
1842 new TGeoTrd1("ogemh2", 28.949, 35.297, 0.1, 18.45); // ogemh2 - middle hole
1843 new TGeoTrd1("ogemh3", 35.667, 42.332, 0.1, 19.425); // ogemh3 - upper hole
1844 //
1845 auto* tr2 = new TGeoTranslation("tr2", 0., 0., 18.125 - 57.1);
1846 auto* tr3 = new TGeoTranslation("tr3", 0., 0., 55.2 - 57.1);
1847 auto* tr4 = new TGeoTranslation("tr4", 0., 0., 94.175 - 57.1);
1848 tr2->RegisterYourself();
1849 tr3->RegisterYourself();
1850 tr4->RegisterYourself();
1851 auto* ocgem = new TGeoCompositeShape("ogem-ogemh1:tr2-ogemh2:tr3-ogemh3:tr4");
1852 auto* ocgemv = new TGeoVolume("TPC_OCGEM", ocgem, m4);
1853 //
1854 // holes for connectors
1855 //
1856 fileName = gSystem->Getenv("VMCWORKDIR");
1857 fileName += "/Detectors/Geometry/TPC/conn_oroc.dat";
1858 in.open(fileName.Data(), ios_base::in); // asci file
1859 if (!in.is_open()) {
1860 LOG(fatal) << "Cannot open input file : " << fileName.Data();
1861 }
1862 for (Int_t i = 0; i < 232; i++) {
1863 Double_t x, z, ang;
1864 in >> ang >> x >> z;
1865 ang = -ang;
1866 z -= 1884.5;
1867 z *= 0.1;
1868 x *= 0.1;
1869 //
1870 auto* rrr = new TGeoRotation();
1871 rrr->RotateY(ang);
1872 obdv->AddNode(connov, i + 1, new TGeoCombiTrans(x, 3.89, z, rrr));
1873 osbv->AddNode(connobv, i + 1, new TGeoCombiTrans(x, 0., z - 2.15, rrr));
1874 }
1875 in.close();
1876 //
1877 auto* oroc = new TGeoVolumeAssembly("TPC_OROC");
1878 //
1879 oroc->AddNode(obdv, 1);
1880 oroc->AddNode(ocv, 1, new TGeoTranslation(0., 3.1, 2.15));
1881 oroc->AddNode(osbv, 1, new TGeoTranslation(0., 4.54, 2.15));
1882 oroc->AddNode(oppv, 1, new TGeoTranslation(0., 4.95, 2.15));
1883 oroc->AddNode(ocgemv, 1, new TGeoTranslation(0., 5.21, 2.15));
1884 //
1885 // now iroc and oroc are placed into a sector...
1886 //
1887 auto* secta = new TGeoVolumeAssembly("TPC_SECT"); // a-side
1888 auto* sectc = new TGeoVolumeAssembly("TPC_SECT"); // c-side
1889 TGeoRotation rot1("rot1", 90., 90., 0.);
1890 TGeoRotation rot2("rot2");
1891 rot2.RotateY(10.);
1892 auto* rot = new TGeoRotation("rot");
1893 *rot = rot1 * rot2;
1894 //
1895 Double_t x0, y0;
1896 //
1897 // a and c sides separately because of possible different z-coordinate of first gem
1898 //
1899 x0 = 110.2 * TMath::Cos(openingAngle);
1900 y0 = 110.2 * TMath::Sin(openingAngle);
1901 auto* combi1a = new TGeoCombiTrans("combi1", x0, y0, 1.09, rot); // a-side
1902 auto* combi1c = new TGeoCombiTrans("combi1", x0, y0, 1.09, rot); // c-side
1903 x0 = 188.45 * TMath::Cos(openingAngle);
1904 y0 = 188.45 * TMath::Sin(openingAngle);
1905 auto* combi2a = new TGeoCombiTrans("combi2", x0, y0, 1.09, rot); // a-side
1906 auto* combi2c = new TGeoCombiTrans("combi2", x0, y0, 1.09, rot); // c-side
1907 //
1908 //
1909 // A-side
1910 //
1911 secta->AddNode(ch, 1);
1912 secta->AddNode(iroc, 1, combi1a);
1913 secta->AddNode(oroc, 1, combi2a);
1914 //
1915 // C-side
1916 //
1917 sectc->AddNode(ch, 1);
1918 sectc->AddNode(iroc, 1, combi1c);
1919 sectc->AddNode(oroc, 1, combi2c);
1920 //
1921 // now I try to make wheels...
1922 //
1923 auto* wheela = new TGeoVolumeAssembly("TPC_ENDCAP");
1924 auto* wheelc = new TGeoVolumeAssembly("TPC_ENDCAP");
1925 //
1926 for (Int_t i = 0; i < 18; i++) {
1927 Double_t phi = (20. * i);
1928 auto* rwh = new TGeoRotation();
1929 rwh->RotateZ(phi);
1930 wheela->AddNode(secta, i + 1, rwh);
1931 wheelc->AddNode(sectc, i + 1, rwh);
1932 }
1933 // wheels in the drift volume!
1934
1935 auto* combi3 = new TGeoCombiTrans("combi3", 0., 0., 256.6, ref);
1936 v9->AddNode(wheela, 1, combi3);
1937 v9->AddNode(wheelc, 2, new TGeoTranslation(0., 0., -256.6));
1938 //_____________________________________________________________
1939 // service support wheel
1940 //_____________________________________________________________
1941 auto* sw = new TGeoPgon(0., 20., 1, 2);
1942 sw->DefineSection(0, -4., 80.5, 251.75);
1943 sw->DefineSection(1, 4., 80.5, 251.75);
1944 auto* swv = new TGeoVolume("TPC_SWSEG", sw, m3); // Al
1945 //
1946 thick = 1.;
1947 shift = thick / TMath::Sin(openingAngle);
1948 auto* sh = new TGeoPgon(0., 20., 1, 2);
1949 sh->DefineSection(0, -4., 81.5 - shift, 250.75 - shift);
1950 sh->DefineSection(1, 4., 81.5 - shift, 250.75 - shift);
1951 auto* shv = new TGeoVolume("TPC_SWS1", sh, m1); // Air
1952 //
1953 TGeoMedium* m9 = gGeoManager->GetMedium("TPC_Si");
1954 auto* el = new TGeoPgon(0., 20., 1, 2);
1955 el->DefineSection(0, -1.872, 81.5 - shift, 250.75 - shift);
1956 el->DefineSection(1, 1.872, 81.5 - shift, 250.75 - shift);
1957 auto* elv = new TGeoVolume("TPC_ELEC", el, m9); // Si
1958 //
1959 shv->AddNode(elv, 1);
1960 //
1961 //
1962 ys = shift * TMath::Sin(openingAngle);
1963 xs = shift * TMath::Cos(openingAngle);
1964 swv->AddNode(shv, 1, new TGeoTranslation(xs, ys, 0.));
1965 // cover
1966 auto* co = new TGeoPgon(0., 20., 1, 2);
1967 co->DefineSection(0, -0.5, 77., 255.25);
1968 co->DefineSection(1, 0.5, 77., 255.25);
1969 auto* cov = new TGeoVolume("TPC_SWC1", co, m3); // Al
1970 // hole in a cover
1971 auto* coh = new TGeoPgon(0., 20., 1, 2);
1972 shift = 4. / TMath::Sin(openingAngle);
1973 coh->DefineSection(0, -0.5, 85. - shift, 247.25 - shift);
1974 coh->DefineSection(1, 0.5, 85. - shift, 247.25 - shift);
1975 //
1976 auto* cohv = new TGeoVolume("TPC_SWC2", coh, m1);
1977 //
1978 ys = shift * TMath::Sin(openingAngle);
1979 xs = shift * TMath::Cos(openingAngle);
1980 cov->AddNode(cohv, 1, new TGeoTranslation(xs, ys, 0.));
1981 //
1982 // Sector as an Assembly
1983 //
1984 auto* swhs = new TGeoVolumeAssembly("TPC_SSWSEC");
1985 swhs->AddNode(swv, 1);
1986 swhs->AddNode(cov, 1, new TGeoTranslation(0., 0., -4.5));
1987 swhs->AddNode(cov, 2, new TGeoTranslation(0., 0., 4.5));
1988 //
1989 // SSW as an Assembly of sectors
1990 //
1991 TGeoRotation* rsw[18];
1992 auto* swheel = new TGeoVolumeAssembly("TPC_SSWHEEL");
1993 for (Int_t i = 0; i < 18; i++) {
1994 Double_t phi = (20. * i);
1995 rsw[i] = new TGeoRotation();
1996 rsw[i]->RotateZ(phi);
1997 swheel->AddNode(swhs, i + 1, rsw[i]);
1998 }
1999 v1->AddNode(swheel, 1, new TGeoTranslation(0., 0., -284.6));
2000 v1->AddNode(swheel, 2, new TGeoTranslation(0., 0., 284.6));
2001
2002 // sensitive strips - strip "0" is always set
2003 // conditional
2005 Int_t totrows = 159;
2006 // totrows = mParam->GetNRowLow() + mParam->GetNRowUp();
2007 Double_t* upar;
2008 upar = nullptr;
2009 gGeoManager->Volume("TPC_Strip", "PGON", m5->GetId(), upar);
2010 upar = new Double_t[10];
2011 upar[0] = 0.;
2012 upar[1] = 360.;
2013 upar[2] = 18.;
2014 upar[3] = 2.;
2015 //
2016 upar[4] = -124.8;
2017 upar[7] = 124.8;
2018
2020 // Double_t rlow=mParam->GetPadRowRadiiLow(0);
2021 Double_t rlow = 85.225; // cm
2022
2023 upar[5] = rlow;
2024 upar[6] = rlow + .01;
2025 upar[8] = upar[5];
2026 upar[9] = upar[6];
2027 //
2028 gGeoManager->Node("TPC_Strip", 1, "TPC_Drift", 0., 0., 124.82, 0, kTRUE, upar, 10);
2029 gGeoManager->Node("TPC_Strip", totrows + 1, "TPC_Drift", 0., 0., -124.82, 0, kTRUE, upar, 10);
2030 //
2031 // now, strips optionally
2032 //
2033 // if(mSens){
2034 // //lower sectors
2035 // for(Int_t i=2;i<mParam->GetNRowLow()+1;i++){
2036 // rlow=mParam->GetPadRowRadiiLow(i-1);
2037 // upar[5]=rlow;
2038 // upar[6]=rlow+.01;
2039 // upar[8]=upar[5];
2040 // upar[9]=upar[6];
2041 // gGeoManager->Node("TPC_Strip",i,
2042 // "TPC_Drift",0.,0.,124.82,0,kTRUE,upar,10);
2043 // gGeoManager->Node("TPC_Strip",totrows+i,
2044 // "TPC_Drift",0.,0.,-124.82,0,kTRUE,upar,10);
2045 // }
2046 // //upper sectors
2047 // for(Int_t i=1;i<mParam->GetNRowUp()+1;i++){
2048 // rlow=mParam->GetPadRowRadiiUp(i-1);
2049 // upar[5]=rlow;
2050 // upar[6]=rlow+.01;
2051 // upar[8]=upar[5];
2052 // upar[9]=upar[6];
2053 // gGeoManager->Node("TPC_Strip",i+mParam->GetNRowLow(),
2054 // "TPC_Drift",0.,0.,124.82,0,kTRUE,upar,10);
2055 // gGeoManager->Node("TPC_Strip",totrows+i+mParam->GetNRowLow(),
2056 // "TPC_Drift",0.,0.,-124.82,0,kTRUE,upar,10);
2057 // }
2058 // }//strips
2059 //----------------------------------------------------------
2060 // TPC Support Rods - MAKROLON
2061 //----------------------------------------------------------
2062 //
2063 TGeoMedium* m7 = gGeoManager->GetMedium("TPC_Cu");
2064 TGeoMedium* m10 = gGeoManager->GetMedium("TPC_Alumina");
2065 TGeoMedium* m11 = gGeoManager->GetMedium("TPC_Peek");
2066 TGeoMedium* m13 = gGeoManager->GetMedium("TPC_Brass");
2067 TGeoMedium* m14 = gGeoManager->GetMedium("TPC_Alumina1");
2068 //
2069 // tpc rod is an assembly of 10 long parts and 2 short parts
2070 // connected with alu rings and plagged on both sides.
2071 //
2072 //
2073 // tpc rod long
2074 //
2075 auto* rod = new TGeoPcon("rod", 0., 360., 6);
2076 rod->DefineSection(0, -10.43, 1.92, 2.08);
2077 rod->DefineSection(1, -9.75, 1.92, 2.08);
2078
2079 rod->DefineSection(2, -9.75, 1.8, 2.2);
2080 rod->DefineSection(3, 9.75, 1.8, 2.2);
2081
2082 rod->DefineSection(4, 9.75, 1.92, 2.08);
2083 rod->DefineSection(5, 10.43, 1.92, 2.08);
2084 //
2085 auto* mrodl = new TGeoVolume("TPC_mrodl", rod, m6);
2086 //
2087 // tpc rod short
2088 //
2089 auto* rod1 = new TGeoPcon("rod1", 0., 360., 6);
2090 rod1->DefineSection(0, -8.93, 1.92, 2.08);
2091 rod1->DefineSection(1, -8.25, 1.92, 2.08);
2092
2093 rod1->DefineSection(2, -8.25, 1.8, 2.2);
2094 rod1->DefineSection(3, 8.25, 1.8, 2.2);
2095
2096 rod1->DefineSection(4, 8.25, 1.92, 2.08);
2097 rod1->DefineSection(5, 8.93, 1.92, 2.08);
2098 //
2099 auto* mrods = new TGeoVolume("TPC_mrods", rod1, m6);
2100 //
2101 // below is for the resistor rod
2102 //
2103 // hole for the brass connectors
2104 //
2105
2106 new TGeoTube("hhole", 0., 0.3, 0.3);
2107 //
2108 // transformations for holes - initialy they
2109 // are placed at x=0 and negative y
2110 //
2111 auto* rhole = new TGeoRotation();
2112 rhole->RotateX(90.);
2113 TGeoCombiTrans* transf[13];
2114 Char_t name[30];
2115 for (Int_t i = 0; i < 13; i++) {
2116 snprintf(name, 30, "transf%d", i);
2117 transf[i] = new TGeoCombiTrans(name, 0., -2., -9. + i * 1.5, rhole);
2118 transf[i]->RegisterYourself();
2119 }
2120 // union expression for holes
2121 TString operl("hhole:transf0");
2122 for (Int_t i = 1; i < 13; i++) {
2123 snprintf(name, 30, "+hhole:transf%d", i);
2124 operl.Append(name);
2125 }
2126 //
2127 TString opers("hhole:transf1");
2128 for (Int_t i = 2; i < 12; i++) {
2129 snprintf(name, 30, "+hhole:transf%d", i);
2130 opers.Append(name);
2131 }
2132 // union of holes
2133 new TGeoCompositeShape("hlv", operl.Data());
2134 new TGeoCompositeShape("hsv", opers.Data());
2135 //
2136 auto* rodl = new TGeoCompositeShape("rodl", "rod-hlv");
2137 auto* rods = new TGeoCompositeShape("rods", "rod1-hsv");
2138 // rods - volumes - makrolon rods with holes
2139 auto* rodlv = new TGeoVolume("TPC_rodl", rodl, m6);
2140 auto* rodsv = new TGeoVolume("TPC_rods", rods, m6);
2141 // brass connectors
2142 // connectors
2143 auto* bcon = new TGeoTube(0., 0.3, 0.3); // connectors
2144 auto* bconv = new TGeoVolume("TPC_bcon", bcon, m13);
2145 //
2146 // hooks holding strips
2147 //
2148 new TGeoBBox("hk1", 0.625, 0.015, 0.75);
2149 new TGeoBBox("hk2", 0.625, 0.015, 0.15);
2150 auto* tr21 = new TGeoTranslation("tr21", 0., -0.03, -0.6);
2151 auto* tr12 = new TGeoTranslation("tr12", 0., -0.03, 0.6);
2152 tr21->RegisterYourself();
2153 tr12->RegisterYourself();
2154
2155 auto* hook = new TGeoCompositeShape("hook", "hk1+hk2:tr21+hk2:tr12");
2156 auto* hookv = new TGeoVolume("TPC_hook", hook, m13);
2157 //
2158 // assembly of the short rod with connectors and hooks
2159 //
2160 //
2161 // short rod
2162 //
2163 auto* spart = new TGeoVolumeAssembly("TPC_spart");
2164 //
2165 spart->AddNode(rodsv, 1);
2166 for (Int_t i = 1; i < 12; i++) {
2167 spart->AddNode(bconv, i, transf[i]);
2168 }
2169 for (Int_t i = 0; i < 11; i++) {
2170 spart->AddNode(hookv, i + 1, new TGeoTranslation(0., -2.315, -7.5 + i * 1.5));
2171 }
2172 //
2173 // long rod
2174 //
2175 auto* lpart = new TGeoVolumeAssembly("TPC_lpart");
2176 //
2177 lpart->AddNode(rodlv, 1);
2178 for (Int_t i = 0; i < 13; i++) {
2179 lpart->AddNode(bconv, i + 12, transf[i]);
2180 }
2181 for (Int_t i = 0; i < 13; i++) {
2182 lpart->AddNode(hookv, i + 12, new TGeoTranslation(0., -2.315, -9. + i * 1.5));
2183 }
2184 //
2185 // alu ring
2186 //
2187 new TGeoTube("ring1", 2.1075, 2.235, 0.53);
2188 new TGeoTube("ring2", 1.7925, 1.89, 0.43);
2189 new TGeoTube("ring3", 1.89, 2.1075, 0.05);
2190 auto* ring = new TGeoCompositeShape("ring", "ring1+ring2+ring3");
2191 auto* ringv = new TGeoVolume("TPC_ring", ring, m3);
2192 //
2193 // rod assembly
2194 //
2195 auto* tpcrrod = new TGeoVolumeAssembly("TPC_rrod"); // rrod
2196 auto* tpcmrod = new TGeoVolumeAssembly("TPC_mrod"); // makrolon rod
2197 // long pieces
2198 for (Int_t i = 0; i < 11; i++) {
2199 tpcrrod->AddNode(ringv, i + 1, new TGeoTranslation(0., 0., -105. + i * 21));
2200 tpcmrod->AddNode(ringv, i + 12, new TGeoTranslation(0., 0., -105. + i * 21));
2201 }
2202 for (Int_t i = 0; i < 10; i++) {
2203 tpcrrod->AddNode(lpart, i + 1, new TGeoTranslation(0., 0., -94.5 + i * 21)); // resistor rod
2204 tpcmrod->AddNode(mrodl, i + 1, new TGeoTranslation(0., 0., -94.5 + i * 21)); // makrolon rod
2205 }
2206 //
2207 // right plug - identical for all rods
2208 //
2209 auto* tpcrp = new TGeoPcon(0., 360., 6);
2210 //
2211 tpcrp->DefineSection(0, 123.05, 1.89, 2.1075);
2212 tpcrp->DefineSection(1, 123.59, 1.89, 2.1075);
2213 //
2214 tpcrp->DefineSection(2, 123.59, 1.8, 2.2);
2215 tpcrp->DefineSection(3, 127., 1.8, 2.2);
2216 //
2217 tpcrp->DefineSection(4, 127., 0., 2.2);
2218 tpcrp->DefineSection(5, 127.5, 0., 2.2);
2219 //
2220 auto* tpcrpv = new TGeoVolume("TPC_RP", tpcrp, m6);
2221 //
2222 // adding short pieces and right plug
2223 //
2224 tpcrrod->AddNode(spart, 1, new TGeoTranslation(0., 0., -114.));
2225 tpcrrod->AddNode(spart, 2, new TGeoTranslation(0., 0., 114.));
2226 tpcrrod->AddNode(ringv, 23, new TGeoTranslation(0., 0., -123.));
2227 tpcrrod->AddNode(ringv, 24, new TGeoTranslation(0., 0., 123.));
2228 tpcrrod->AddNode(tpcrpv, 1);
2229 //
2230 tpcmrod->AddNode(mrods, 1, new TGeoTranslation(0., 0., -114.));
2231 tpcmrod->AddNode(mrods, 2, new TGeoTranslation(0., 0., 114.));
2232 tpcmrod->AddNode(ringv, 25, new TGeoTranslation(0., 0., -123.));
2233 tpcmrod->AddNode(ringv, 26, new TGeoTranslation(0., 0., 123.));
2234 tpcmrod->AddNode(tpcrpv, 2);
2235 //
2236 // from the ringv position to the CM is 3.0 cm!
2237 //----------------------------------------
2238 //
2239 //
2240 // HV rods - makrolon + 0.58cm (diameter) Cu ->check the length
2241 auto* hvr = new TGeoTube(0., 1.465, 123.);
2242 auto* hvc = new TGeoTube(0., 0.29, 123.);
2243 //
2244 auto* hvrv = new TGeoVolume("TPC_HV_Rod", hvr, m6);
2245 auto* hvcv = new TGeoVolume("TPC_HV_Cable", hvc, m7);
2246 hvrv->AddNode(hvcv, 1);
2247 //
2248 // resistor rod
2249 //
2250 auto* cr = new TGeoTube(0., 0.45, 123.);
2251 auto* cw = new TGeoTube(0., 0.15, 123.);
2252 auto* crv = new TGeoVolume("TPC_CR", cr, m10);
2253 auto* cwv = new TGeoVolume("TPC_W", cw, m12);
2254 //
2255 // ceramic rod with water
2256 //
2257 crv->AddNode(cwv, 1);
2258 //
2259 // peek rod
2260 //
2261 auto* pr = new TGeoTube(0.2, 0.35, 123.);
2262 auto* prv = new TGeoVolume("TPC_PR", pr, m11);
2263 //
2264 // copper plates with connectors
2265 //
2266 new TGeoTube("tub", 0., 1.7, 0.025);
2267 //
2268 // half space - points on the plane and a normal vector
2269 //
2270 Double_t n[3], p[3];
2271 Double_t slope = TMath::Tan(22. * TMath::DegToRad());
2272 Double_t intp = 1.245;
2273 //
2274 Double_t b = slope * slope + 1.;
2275 p[0] = intp * slope / b;
2276 p[1] = -intp / b;
2277 p[2] = 0.;
2278 //
2279 n[0] = -p[0];
2280 n[1] = -p[1];
2281 n[2] = 0.;
2282 Double_t norm;
2283 norm = TMath::Sqrt(n[0] * n[0] + n[1] * n[1]);
2284 n[0] /= norm;
2285 n[1] /= norm;
2286 //
2287 new TGeoHalfSpace("sp1", p, n);
2288 //
2289 slope = -slope;
2290 //
2291 p[0] = intp * slope / b;
2292 p[1] = -intp / b;
2293 //
2294 n[0] = -p[0];
2295 n[1] = -p[1];
2296 norm = TMath::Sqrt(n[0] * n[0] + n[1] * n[1]);
2297 n[0] /= norm;
2298 n[1] /= norm;
2299 //
2300 new TGeoHalfSpace("sp2", p, n);
2301 // holes for rods
2302 // holes
2303 new TGeoTube("h1", 0., 0.5, 0.025);
2304 new TGeoTube("h2", 0., 0.35, 0.025);
2305 // translations:
2306 auto* ttr11 = new TGeoTranslation("ttr11", -0.866, 0.5, 0.);
2307 auto* ttr22 = new TGeoTranslation("ttr22", 0.866, 0.5, 0.);
2308 ttr11->RegisterYourself();
2309 ttr22->RegisterYourself();
2310 // elastic connector
2311 new TGeoBBox("elcon", 0.72, 0.005, 0.3);
2312 auto* crr1 = new TGeoRotation();
2313 crr1->RotateZ(-22.);
2314 auto* ctr1 = new TGeoCombiTrans("ctr1", -0.36011, -1.09951, -0.325, crr1);
2315 ctr1->RegisterYourself();
2316 auto* cs1 = new TGeoCompositeShape("cs1", "(((((tub-h1:ttr11)-h1:ttr22)-sp1)-sp2)-h2)+elcon:ctr1");
2317 //
2318 auto* csvv = new TGeoVolume("TPC_RR_CU", cs1, m7);
2319 //
2320 // resistor rod assembly 2 ceramic rods, peak rod, Cu plates
2321 // and resistors
2322 //
2323 auto* rrod = new TGeoVolumeAssembly("TPC_RRIN");
2324 // rods
2325 rrod->AddNode(crv, 1, ttr11);
2326 rrod->AddNode(crv, 2, ttr22);
2327 rrod->AddNode(prv, 1);
2328 // Cu plates
2329 for (Int_t i = 0; i < 165; i++) {
2330 rrod->AddNode(csvv, i + 1, new TGeoTranslation(0., 0., -122.675 + i * 1.5));
2331 }
2332 // resistors
2333 auto* res = new TGeoTube(0., 0.15, 0.5);
2334 auto* resv = new TGeoVolume("TPC_RES", res, m14);
2335 auto* ress = new TGeoVolumeAssembly("TPC_RES_CH");
2336 ress->AddNode(resv, 1, new TGeoTranslation(0.2, 0., 0.));
2337 ress->AddNode(resv, 2, new TGeoTranslation(-0.2, 0., 0.));
2338 //
2339 auto* crr2 = new TGeoRotation();
2340 crr2->RotateY(30.);
2341 auto* crr3 = new TGeoRotation();
2342 crr3->RotateY(-30.);
2343 //
2344 for (Int_t i = 0; i < 164; i += 2) {
2345 rrod->AddNode(ress, i + 1, new TGeoCombiTrans(0., 1.2, -121.925 + i * 1.5, crr2));
2346 rrod->AddNode(ress, i + 2, new TGeoCombiTrans(0., 1.2, -121.925 + (i + 1) * 1.5, crr3));
2347 }
2348
2349 tpcrrod->AddNode(rrod, 1, new TGeoCombiTrans(0., 0., 0.5, crr1));
2350 //
2351 // rod left head with holders - inner
2352 //
2353 // first element - support for inner holder TPC_IHS
2354 Double_t shift1[3] = {0.0, -0.175, 0.0};
2355
2356 new TGeoBBox("tpcihs1", 4.7, 0.66, 2.35);
2357 new TGeoBBox("tpcihs2", 4.7, 0.485, 1.0, shift1);
2358 new TGeoBBox("tpcihs3", 1.5, 0.485, 2.35, shift1);
2359 new TGeoTube("tpcihs4", 0.0, 2.38, 0.1);
2360 //
2361 Double_t pointstrap[16];
2362 pointstrap[0] = 0.0;
2363 pointstrap[1] = 0.0;
2364 pointstrap[2] = 0.0;
2365 pointstrap[3] = 1.08;
2366 pointstrap[4] = 2.3;
2367 pointstrap[5] = 1.08;
2368 pointstrap[6] = 3.38;
2369 pointstrap[7] = 0.0;
2370 pointstrap[8] = 0.0;
2371 pointstrap[9] = 0.0;
2372 pointstrap[10] = 0.0;
2373 pointstrap[11] = 1.08;
2374 pointstrap[12] = 2.3;
2375 pointstrap[13] = 1.08;
2376 pointstrap[14] = 3.38;
2377 pointstrap[15] = 0.0;
2378 //
2379 auto* tpcihs5 = new TGeoArb8("tpcihs5", 0.6, pointstrap);
2380 //
2381 // half space - cutting "legs"
2382 //
2383 p[0] = 0.0;
2384 p[1] = 0.105;
2385 p[2] = 0.0;
2386 //
2387 n[0] = 0.0;
2388 n[1] = 1.0;
2389 n[2] = 0.0;
2390
2391 new TGeoHalfSpace("cutil1", p, n);
2392
2393 //
2394 // transformations
2395 //
2396 auto* trans2 = new TGeoTranslation("trans2", 0.0, 2.84, 2.25);
2397 trans2->RegisterYourself();
2398 auto* trans3 = new TGeoTranslation("trans3", 0.0, 2.84, -2.25);
2399 trans3->RegisterYourself();
2400 // support - composite volume
2401 //
2402 auto* tpcihs6 =
2403 new TGeoCompositeShape("tpcihs6", "tpcihs1-(tpcihs2+tpcihs3)-(tpcihs4:trans2)-(tpcihs4:trans3)-cutil1");
2404 //
2405 // volumes - all makrolon
2406 //
2407 auto* tpcihss = new TGeoVolume("TPC_IHSS", tpcihs6, m6); // support
2408 auto* tpcihst = new TGeoVolume("TPC_IHSTR", tpcihs5, m6); // trapesoid
2409 // now assembly
2410 auto* rot111 = new TGeoRotation();
2411 rot111->RotateY(180.0);
2412 //
2413 auto* tpcihs = new TGeoVolumeAssembly("TPC_IHS"); // assembly of the support
2414 tpcihs->AddNode(tpcihss, 1);
2415 tpcihs->AddNode(tpcihst, 1, new TGeoTranslation(-4.7, 0.66, 0.0));
2416 tpcihs->AddNode(tpcihst, 2, new TGeoCombiTrans(4.7, 0.66, 0.0, rot111));
2417 //
2418 // two rod holders (TPC_IRH) assembled with the support
2419 //
2420 new TGeoBBox("tpcirh1", 4.7, 1.33, 0.5);
2421 shift1[0] = -3.65;
2422 shift1[1] = 0.53;
2423 shift1[2] = 0.;
2424 new TGeoBBox("tpcirh2", 1.05, 0.8, 0.5, shift1);
2425 shift1[0] = 3.65;
2426 shift1[1] = 0.53;
2427 shift1[2] = 0.;
2428 new TGeoBBox("tpcirh3", 1.05, 0.8, 0.5, shift1);
2429 shift1[0] = 0.0;
2430 shift1[1] = 1.08;
2431 shift1[2] = 0.;
2432 new TGeoBBox("tpcirh4", 1.9, 0.25, 0.5, shift1);
2433 new TGeoTube("tpcirh5", 0, 1.9, 5);
2434 //
2435 auto* trans4 = new TGeoTranslation("trans4", 0, 0.83, 0.0);
2436 trans4->RegisterYourself();
2437 //
2438 auto* tpcirh6 = new TGeoCompositeShape("tpcirh6", "tpcirh1-tpcirh2-tpcirh3-(tpcirh5:trans4)-tpcirh4");
2439 //
2440 // now volume
2441 //
2442 auto* tpcirh = new TGeoVolume("TPC_IRH", tpcirh6, m6);
2443 //
2444 // and all together...
2445 //
2446 TGeoVolume* tpciclamp = new TGeoVolumeAssembly("TPC_ICLP");
2447 tpciclamp->AddNode(tpcihs, 1);
2448 tpciclamp->AddNode(tpcirh, 1, new TGeoTranslation(0, 1.99, 1.1));
2449 tpciclamp->AddNode(tpcirh, 2, new TGeoTranslation(0, 1.99, -1.1));
2450 //
2451 // and now left inner "head"
2452 //
2453 auto* inplug = new TGeoPcon("inplug", 0.0, 360.0, 13);
2454
2455 inplug->DefineSection(0, 0.3, 0.0, 2.2);
2456 inplug->DefineSection(1, 0.6, 0.0, 2.2);
2457
2458 inplug->DefineSection(2, 0.6, 0.0, 1.75);
2459 inplug->DefineSection(3, 0.7, 0.0, 1.75);
2460
2461 inplug->DefineSection(4, 0.7, 1.55, 1.75);
2462 inplug->DefineSection(5, 1.6, 1.55, 1.75);
2463
2464 inplug->DefineSection(6, 1.6, 1.55, 2.2);
2465 inplug->DefineSection(7, 1.875, 1.55, 2.2);
2466
2467 inplug->DefineSection(8, 2.47, 1.75, 2.2);
2468
2469 inplug->DefineSection(9, 2.47, 1.75, 2.08);
2470 inplug->DefineSection(10, 2.57, 1.8, 2.08);
2471
2472 inplug->DefineSection(11, 2.57, 1.92, 2.08);
2473 inplug->DefineSection(12, 2.95, 1.92, 2.08);
2474 //
2475 shift1[0] = 0.0;
2476 shift1[1] = -2.09;
2477 shift1[2] = 1.075;
2478 //
2479 new TGeoBBox("pcuti", 1.5, 0.11, 1.075, shift1);
2480 //
2481 auto* inplleft = new TGeoCompositeShape("inplleft", "inplug-pcuti");
2482 auto* tpcinlplug = new TGeoVolume("TPC_INPLL", inplleft, m6);
2483 //
2484 // holder + plugs
2485 //
2486 TGeoVolume* tpcihpl = new TGeoVolumeAssembly("TPC_IHPL"); // holder+2 plugs (reflected)
2487 tpcihpl->AddNode(tpcinlplug, 1);
2488 tpcihpl->AddNode(tpcinlplug, 2, ref);
2489 tpcihpl->AddNode(tpciclamp, 1, new TGeoTranslation(0.0, -2.765, 0.0));
2490 //
2491 // outer holders and clamps
2492 //
2493
2494 // outer membrane holder (between rods)
2495 pointstrap[0] = 0.0;
2496 pointstrap[1] = 0.0;
2497 pointstrap[2] = 0.0;
2498 pointstrap[3] = 2.8;
2499 pointstrap[4] = 3.1;
2500 pointstrap[5] = 2.8 - 3.1 * TMath::Tan(15. * TMath::DegToRad());
2501 pointstrap[6] = 3.1;
2502 pointstrap[7] = 0.0;
2503 pointstrap[8] = 0.0;
2504 pointstrap[9] = 0.0;
2505 pointstrap[10] = 0.0;
2506 pointstrap[11] = 2.8;
2507 pointstrap[12] = 3.1;
2508 pointstrap[13] = 2.8 - 3.1 * TMath::Tan(15. * TMath::DegToRad());
2509 pointstrap[14] = 3.1;
2510 pointstrap[15] = 0.0;
2511 //
2512 auto* tpcomh1 = new TGeoArb8("tpcomh1", 1.05, pointstrap);
2513 auto* tpcomh2 = new TGeoBBox("tpcomh2", 0.8, 1.4, 6);
2514 //
2515 auto* tpcomh1v = new TGeoVolume("TPC_OMH1", tpcomh1, m7);
2516 auto* tpcomh2v = new TGeoVolume("TPC_OMH2", tpcomh2, m7);
2517 //
2518 TGeoVolume* tpcomh3v = new TGeoVolumeAssembly("TPC_OMH3"); // assembly1
2519 tpcomh3v->AddNode(tpcomh1v, 1, new TGeoTranslation(0.8, -1.4, 4.95));
2520 tpcomh3v->AddNode(tpcomh1v, 2, new TGeoTranslation(0.8, -1.4, -4.95));
2521 tpcomh3v->AddNode(tpcomh2v, 1);
2522 //
2523 shift1[0] = 0.9;
2524 shift1[1] = -1.85;
2525 shift1[2] = 0.0;
2526 //
2527 new TGeoBBox("tpcomh3", 1.65, 1.15, 3.4);
2528 auto* tpcomh4 = new TGeoBBox("tpcomh4", 0.75, 0.7, 3.4, shift1);
2529 //
2530 // halfspace 1
2531 //
2532 p[0] = 0.0;
2533 p[1] = -1.05;
2534 p[2] = -3.4;
2535 //
2536 n[0] = 0.0;
2537 n[1] = -1.0 * TMath::Tan(30. * TMath::DegToRad());
2538 n[2] = 1.0;
2539 //
2540 new TGeoHalfSpace("cutomh1", p, n);
2541 //
2542 // halfspace 2
2543 //
2544 p[0] = 0.0;
2545 p[1] = -1.05;
2546 p[2] = 3.4;
2547 //
2548 n[0] = 0.0;
2549 n[1] = -1.0 * TMath::Tan(30. * TMath::DegToRad());
2550 n[2] = -1.0;
2551 //
2552 new TGeoHalfSpace("cutomh2", p, n);
2553 //
2554 // halfspace 3
2555 //
2556 p[0] = -1.65;
2557 p[1] = 0.0;
2558 p[2] = -0.9;
2559 //
2560 n[0] = 1.0 * TMath::Tan(75. * TMath::DegToRad());
2561 n[1] = 0.0;
2562 n[2] = 1.0;
2563 //
2564 new TGeoHalfSpace("cutomh3", p, n);
2565 //
2566 // halfspace 4
2567 //
2568 p[0] = -1.65;
2569 p[1] = 0.0;
2570 p[2] = 0.9;
2571 //
2572 n[0] = 1.0 * TMath::Tan(75 * TMath::DegToRad());
2573 n[1] = 0.0;
2574 n[2] = -1.0;
2575 //
2576 new TGeoHalfSpace("cutomh4", p, n);
2577 //
2578 // halsfspace 5
2579 //
2580 p[0] = 1.65;
2581 p[1] = -1.05;
2582 p[2] = 0.0;
2583 //
2584 n[0] = -1.0;
2585 n[1] = -1.0 * TMath::Tan(20. * TMath::DegToRad());
2586 n[2] = 0.0;
2587 //
2588 new TGeoHalfSpace("cutomh5", p, n);
2589 //
2590 auto* tpcomh5 = new TGeoCompositeShape("tpcomh5", "tpcomh3-cutomh1-cutomh2-cutomh3-cutomh4-cutomh5");
2591 //
2592 auto* tpcomh5v = new TGeoVolume("TPC_OMH5", tpcomh5, m6);
2593 auto* tpcomh4v = new TGeoVolume("TPC_OMH6", tpcomh4, m6);
2594 //
2595 auto* tpcomh7v = new TGeoVolumeAssembly("TPC_OMH7");
2596 tpcomh7v->AddNode(tpcomh5v, 1);
2597 tpcomh7v->AddNode(tpcomh4v, 1);
2598 //
2599 // full membrane holder - tpcomh3v + tpcomh7v
2600 //
2601 auto* tpcomh = new TGeoVolumeAssembly("TPC_OMH");
2602 tpcomh->AddNode(tpcomh3v, 1, new TGeoTranslation(1.5, 0., 0.));
2603 tpcomh->AddNode(tpcomh3v, 2, new TGeoCombiTrans(-1.5, 0., 0., rot111));
2604 tpcomh->AddNode(tpcomh7v, 1, new TGeoTranslation(0.65 + 1.5, 2.55, 0.0));
2605 tpcomh->AddNode(tpcomh7v, 2, new TGeoCombiTrans(-0.65 - 1.5, 2.55, 0.0, rot111));
2606 //
2607 // outer rod holder support
2608 //
2609 new TGeoBBox("tpcohs1", 3.8, 0.675, 2.35);
2610 //
2611 shift1[0] = 0.0;
2612 shift1[1] = 0.175;
2613 shift1[2] = 0.0;
2614 //
2615 new TGeoBBox("tpcohs2", 1.5, 0.5, 2.35, shift1);
2616 new TGeoBBox("tpcohs3", 3.8, 0.5, 0.85, shift1);
2617 //
2618 shift1[0] = 0.0;
2619 shift1[1] = -1.175;
2620 shift1[2] = 0.0;
2621 //
2622 auto* tpcohs4 = new TGeoBBox("tpsohs4", 3.1, 0.5, 0.7, shift1);
2623 //
2624 auto* tpcohs4v = new TGeoVolume("TPC_OHS4", tpcohs4, m6);
2625 //
2626 p[0] = 0.0;
2627 p[1] = -0.186;
2628 p[2] = 0.0;
2629 //
2630 n[0] = 0.0;
2631 n[1] = -1.0;
2632 n[2] = 0.0;
2633 //
2634 new TGeoHalfSpace("cutohs1", p, n);
2635 //
2636 auto* tpcohs5 = new TGeoCompositeShape("tpcohs5", "tpcohs1-tpcohs2-tpcohs3-cutohs1");
2637 auto* tpcohs5v = new TGeoVolume("TPC_OHS5", tpcohs5, m6);
2638 //
2639 auto* tpcohs = new TGeoVolumeAssembly("TPC_OHS");
2640 tpcohs->AddNode(tpcohs5v, 1);
2641 tpcohs->AddNode(tpcohs4v, 1);
2642 //
2643 // outer rod holder itself
2644 //
2645 shift1[0] = 0.0;
2646 shift1[1] = 1.325;
2647 shift1[2] = 0.0;
2648 new TGeoBBox("tpcorh1", 3.1, 1.825, 0.55); // from this box we cut pieces...
2649 //
2650 shift1[0] = -3.1;
2651 shift1[1] = -0.5;
2652 shift1[2] = 0.0;
2653 //
2654 new TGeoBBox("tpcorh2", 0.5, 2.75, 1.1, shift1);
2655 //
2656 shift1[0] = 3.1;
2657 shift1[1] = -0.5;
2658 shift1[2] = 0.0;
2659 //
2660 new TGeoBBox("tpcorh3", 0.5, 2.75, 1.1, shift1);
2661 //
2662 shift1[0] = 0.0;
2663 shift1[1] = -0.5;
2664 shift1[2] = -0.95;
2665 //
2666 new TGeoBBox("tpcorh4", 3.9, 2.75, 0.5, shift1);
2667 //
2668 shift1[0] = 0.0;
2669 shift1[1] = -0.5;
2670 shift1[2] = 0.0;
2671 //
2672 new TGeoBBox("tpcorh5", 1.95, 0.5, 1.1, shift1);
2673 //
2674 shift1[0] = 0.0;
2675 shift1[1] = -0.5;
2676 shift1[2] = 0.55;
2677 //
2678 new TGeoBBox("tpcorh6", 2.4, 0.5, 0.6, shift1);
2679 //
2680 new TGeoTube("tpcorh7", 0, 1.95, 0.85);
2681 new TGeoTube("tpcorh8", 0, 2.4, 0.6);
2682 //
2683 auto* trans33 = new TGeoTranslation("trans33", 0.0, 0.0, 0.55);
2684 trans33->RegisterYourself();
2685 //
2686 auto* tpcorh9 =
2687 new TGeoCompositeShape("tpcorh9", "tpcorh1-tpcorh2-tpcorh3-tpcorh4-tpcorh5-tpcorh6-(tpcorh8:trans33)-tpcorh7");
2688 //
2689 auto* tpcorh9v = new TGeoVolume("TPC_ORH", tpcorh9, m6); // outer rod holder
2690 //
2691 // now 2 holders together
2692 //
2693 auto* tpcorh = new TGeoVolumeAssembly("TPC_ORH2");
2694 //
2695 tpcorh->AddNode(tpcorh9v, 1, new TGeoTranslation(0.0, 0.0, 1.25));
2696 tpcorh->AddNode(tpcorh9v, 2, new TGeoCombiTrans(0.0, 0.0, -1.25, rot111));
2697 //
2698 // outer rod plug left
2699 //
2700 auto* outplug = new TGeoPcon("outplug", 0.0, 360.0, 13);
2701
2702 outplug->DefineSection(0, 0.5, 0.0, 2.2);
2703 outplug->DefineSection(1, 0.7, 0.0, 2.2);
2704 outplug->DefineSection(2, 0.7, 1.55, 2.2);
2705 outplug->DefineSection(3, 0.8, 1.55, 2.2);
2706 outplug->DefineSection(4, 0.8, 1.55, 1.75);
2707 outplug->DefineSection(5, 1.2, 1.55, 1.75);
2708 outplug->DefineSection(6, 1.2, 1.55, 2.2);
2709 outplug->DefineSection(7, 1.875, 1.55, 2.2);
2710 outplug->DefineSection(8, 2.47, 1.75, 2.2);
2711 outplug->DefineSection(9, 2.47, 1.75, 2.08);
2712 outplug->DefineSection(10, 2.57, 1.8, 2.08);
2713 outplug->DefineSection(11, 2.57, 1.92, 2.08);
2714 outplug->DefineSection(12, 2.95, 1.92, 2.08);
2715 //
2716 shift1[0] = 0.0;
2717 shift1[1] = 2.09;
2718 shift1[2] = 1.01;
2719
2720 new TGeoBBox("cutout", 2.5, 0.11, 1.01, shift1);
2721 //
2722
2723 auto* outplleft = new TGeoCompositeShape("outplleft", "outplug-cutout");
2724 auto* outplleftv = new TGeoVolume("TPC_OPLL", outplleft, m6);
2725 //
2726 // support + holder + plug
2727 //
2728
2729 auto* tpcohpl = new TGeoVolumeAssembly("TPC_OHPL");
2730 //
2731 tpcohpl->AddNode(outplleftv, 1); // plug
2732 tpcohpl->AddNode(outplleftv, 2, ref); // plug reflected
2733 tpcohpl->AddNode(tpcorh, 1); // rod holder
2734 tpcohpl->AddNode(tpcohs, 1, new TGeoTranslation(0.0, 3.925, 0)); // support
2735 //
2736
2737 //
2738 // main membrane holder
2739 //
2740 pointstrap[0] = 0.0;
2741 pointstrap[1] = 0.0;
2742 pointstrap[2] = 0.0;
2743 pointstrap[3] = 2.8;
2744 pointstrap[4] = 3.1;
2745 pointstrap[5] = 1.96;
2746 pointstrap[6] = 3.1;
2747 pointstrap[7] = 0.0;
2748 pointstrap[8] = 0.0;
2749 pointstrap[9] = 0.0;
2750 pointstrap[10] = 0.0;
2751 pointstrap[11] = 2.8;
2752 pointstrap[12] = 3.1;
2753 pointstrap[13] = 1.96;
2754 pointstrap[14] = 3.1;
2755 pointstrap[15] = 0.0;
2756 //
2757 auto* tpcmmh1 = new TGeoArb8("tpcmmh1", 1.75, pointstrap);
2758 auto* tpcmmh2 = new TGeoBBox("tpcmmh2", 0.8, 1.4, 12.5);
2759 //
2760 auto* tpcmmh1v = new TGeoVolume("TPC_MMH1", tpcmmh1, m6);
2761 auto* tpcmmh2v = new TGeoVolume("TPC_MMH2", tpcmmh2, m6);
2762 //
2763 auto* tpcmmhs = new TGeoVolumeAssembly("TPC_MMHS");
2764 tpcmmhs->AddNode(tpcmmh1v, 1, new TGeoTranslation(0.8, -1.4, 10.75));
2765 tpcmmhs->AddNode(tpcmmh1v, 2, new TGeoTranslation(0.8, -1.4, -10.75));
2766 tpcmmhs->AddNode(tpcmmh2v, 1);
2767 //
2768 // main membrahe holder clamp
2769 //
2770 shift1[0] = -0.75;
2771 shift1[1] = -1.15;
2772 shift1[2] = 0.0;
2773 //
2774 new TGeoBBox("tpcmmhc1", 1.65, 1.85, 8.9);
2775 new TGeoBBox("tpcmmhc2", 0.9, 0.7, 8.9, shift1);
2776 //
2777 // half spaces - cuts
2778 //
2779 p[0] = -1.65;
2780 p[1] = 0.0;
2781 p[2] = -0.9;
2782 //
2783 n[0] = 8.0;
2784 n[1] = 0.0;
2785 n[2] = 8.0 * TMath::Tan(13. * TMath::DegToRad());
2786 //
2787 new TGeoHalfSpace("cutmmh1", p, n);
2788 //
2789 p[0] = -1.65;
2790 p[1] = 0.0;
2791 p[2] = 0.9;
2792 //
2793 n[0] = 8.0;
2794 n[1] = 0.0;
2795 n[2] = -8.0 * TMath::Tan(13. * TMath::DegToRad());
2796 //
2797 new TGeoHalfSpace("cutmmh2", p, n);
2798 //
2799 p[0] = 0.0;
2800 p[1] = 1.85;
2801 p[2] = -2.8;
2802 //
2803 n[0] = 0.0;
2804 n[1] = -6.1;
2805 n[2] = 6.1 * TMath::Tan(20. * TMath::DegToRad());
2806 //
2807 new TGeoHalfSpace("cutmmh3", p, n);
2808 //
2809 p[0] = 0.0;
2810 p[1] = 1.85;
2811 p[2] = 2.8;
2812 //
2813 n[0] = 0.0;
2814 n[1] = -6.1;
2815 n[2] = -6.1 * TMath::Tan(20 * TMath::DegToRad());
2816 //
2817 new TGeoHalfSpace("cutmmh4", p, n);
2818 //
2819 p[0] = 0.75;
2820 p[1] = 0.0;
2821 p[2] = -8.9;
2822 //
2823 n[0] = 2.4 * TMath::Tan(30 * TMath::DegToRad());
2824 n[1] = 0.0;
2825 n[2] = 2.4;
2826 //
2827 new TGeoHalfSpace("cutmmh5", p, n);
2828 //
2829 p[0] = 0.75;
2830 p[1] = 0.0;
2831 p[2] = 8.9;
2832 //
2833 n[0] = 2.4 * TMath::Tan(30 * TMath::DegToRad());
2834 n[1] = 0.0;
2835 n[2] = -2.4;
2836 //
2837 new TGeoHalfSpace("cutmmh6", p, n);
2838
2839 auto* tpcmmhc =
2840 new TGeoCompositeShape("TPC_MMHC", "tpcmmhc1-tpcmmhc2-cutmmh1-cutmmh2-cutmmh3-cutmmh4-cutmmh5-cutmmh6");
2841
2842 auto* tpcmmhcv = new TGeoVolume("TPC_MMHC", tpcmmhc, m6);
2843 //
2844 TGeoVolume* tpcmmh = new TGeoVolumeAssembly("TPC_MMH");
2845 //
2846 tpcmmh->AddNode(tpcmmhcv, 1, new TGeoTranslation(0.65 + 1.5, 1.85, 0.0));
2847 tpcmmh->AddNode(tpcmmhcv, 2, new TGeoCombiTrans(-0.65 - 1.5, 1.85, 0.0, rot111));
2848 tpcmmh->AddNode(tpcmmhs, 1, new TGeoTranslation(1.5, 0.0, 0.0));
2849 tpcmmh->AddNode(tpcmmhs, 2, new TGeoCombiTrans(-1.5, 0.0, 0.0, rot111));
2850 //
2851
2852 //
2853
2854 //--------------------------------------------
2855 //
2856 // guard ring resistor chain
2857 //
2858
2859 auto* gres1 = new TGeoTube(0., 0.375, 125.); // inside ifc
2860 //
2861 auto* vgres1 = new TGeoVolume("TPC_GRES1", gres1, m14);
2862
2863 //
2864 Double_t xrc, yrc;
2865 //
2866 xrc = 79.3 * TMath::Cos(350. * TMath::DegToRad());
2867 yrc = 79.3 * TMath::Sin(350. * TMath::DegToRad());
2868 //
2869 v9->AddNode(vgres1, 1, new TGeoTranslation(xrc, yrc, 126.9));
2870 v9->AddNode(vgres1, 2, new TGeoTranslation(xrc, yrc, -126.9));
2871 //
2872 xrc = 79.3 * TMath::Cos(190. * TMath::DegToRad());
2873 yrc = 79.3 * TMath::Sin(190. * TMath::DegToRad());
2874 //
2875 v9->AddNode(vgres1, 3, new TGeoTranslation(xrc, yrc, 126.9));
2876 v9->AddNode(vgres1, 4, new TGeoTranslation(xrc, yrc, -126.9));
2877 //------------------------------------------------------------------
2878 TGeoRotation refl("refl", 90., 0., 90., 90., 180., 0.);
2879 TGeoRotation rotrod("rotrod");
2880 //
2881 TGeoRotation* rotpos[2];
2882 //
2883 TGeoRotation* rotrod1[2];
2884 //
2885 // clamps holding rods
2886 //
2887 auto* clampi1 = new TGeoBBox("clampi1", 0.2, 3.1, 0.8);
2888 auto* clampi1v = new TGeoVolume("TPC_clampi1v", clampi1, m6);
2889 //
2890 pointstrap[0] = 0.49;
2891 pointstrap[1] = 0.375;
2892 //
2893 pointstrap[2] = 0.49;
2894 pointstrap[3] = -0.375;
2895 //
2896 pointstrap[4] = -0.49;
2897 pointstrap[5] = -0.375;
2898 //
2899 pointstrap[6] = -0.49;
2900 pointstrap[7] = 1.225;
2901 //
2902 pointstrap[8] = 0.49;
2903 pointstrap[9] = 0.375;
2904 //
2905 pointstrap[10] = 0.49;
2906 pointstrap[11] = -0.375;
2907 //
2908 pointstrap[12] = -0.49;
2909 pointstrap[13] = -0.375;
2910 //
2911 pointstrap[14] = -0.49;
2912 pointstrap[15] = 1.225;
2913 //
2914 auto* clitrap = new TGeoArb8("clitrap", 0.25, pointstrap);
2915 auto* clitrapv = new TGeoVolume("TPC_clitrapv", clitrap, m6);
2916 //
2917 auto* clamprot = new TGeoRotation();
2918 clamprot->RotateX(180.);
2919 //
2920 new TGeoBBox("clibox", 1.125, 3.1, .1);
2921 new TGeoTube("clitub", 0., 2.2, 0.1);
2922 //
2923 // copmisite shape for the clamp holder
2924 //
2925 auto* clitr1 = new TGeoTranslation("clitr1", 1.125, 0., 0.);
2926 clitr1->RegisterYourself();
2927 auto* clihold = new TGeoCompositeShape("clihold", "clibox-clitub:clitr1");
2928 auto* cliholdv = new TGeoVolume("TPC_cliholdv", clihold, m6);
2929 //
2930 // now assembly the whole inner clamp
2931 //
2932 TGeoVolume* iclamp = new TGeoVolumeAssembly("TPC_iclamp");
2933 //
2934 iclamp->AddNode(clampi1v, 1); // main box
2935 iclamp->AddNode(clitrapv, 1, new TGeoTranslation(0.69, -2.725, 0.35)); // trapezoids
2936 iclamp->AddNode(clitrapv, 2, new TGeoTranslation(0.69, -2.725, -0.35));
2937 iclamp->AddNode(clitrapv, 3, new TGeoCombiTrans(0.69, 2.725, 0.35, clamprot));
2938 iclamp->AddNode(clitrapv, 4, new TGeoCombiTrans(0.69, 2.725, -0.35, clamprot));
2939 iclamp->AddNode(cliholdv, 1, new TGeoTranslation(1.325, 0., 0.)); // holder
2940 //
2941 // outer clamps
2942 //
2943 auto* clampo1 = new TGeoBBox("clampo1", 0.25, 3.1, 1.);
2944 auto* clampo2 = new TGeoBBox("clampo2", 0.4, 0.85, 1.);
2945 //
2946 auto* clampo1v = new TGeoVolume("TPC_clampo1v", clampo1, m6);
2947 auto* clampo2v = new TGeoVolume("TPC_clampo2v", clampo2, m6);
2948 //
2949 auto* oclamp = new TGeoVolumeAssembly("TPC_oclamp");
2950 //
2951 oclamp->AddNode(clampo1v, 1);
2952 //
2953 oclamp->AddNode(clampo2v, 1, new TGeoTranslation(0.65, -2.25, 0));
2954 oclamp->AddNode(clampo2v, 2, new TGeoTranslation(0.65, 2.25, 0));
2955
2956 //
2957 pointstrap[0] = 0.375;
2958 pointstrap[1] = 0.75;
2959 pointstrap[2] = 0.375;
2960 pointstrap[3] = -0.35;
2961 pointstrap[4] = -0.375;
2962 pointstrap[5] = -0.35;
2963 pointstrap[6] = -0.375;
2964 pointstrap[7] = 0.35;
2965 //
2966 pointstrap[8] = 0.375;
2967 pointstrap[9] = 0.75;
2968 pointstrap[10] = 0.375;
2969 pointstrap[11] = -0.35;
2970 pointstrap[12] = -0.375;
2971 pointstrap[13] = -0.35;
2972 pointstrap[14] = -0.375;
2973 pointstrap[15] = 0.35;
2974 //
2975 auto* clotrap = new TGeoArb8("clotrap", 0.25, pointstrap);
2976 auto* clotrapv = new TGeoVolume("TPC_clotrapv", clotrap, m6);
2977 //
2978 oclamp->AddNode(clotrapv, 1, new TGeoTranslation(-0.625, -2.75, 0.35));
2979 oclamp->AddNode(clotrapv, 2, new TGeoTranslation(-0.625, -2.75, -0.35));
2980 oclamp->AddNode(clotrapv, 3, new TGeoCombiTrans(-0.625, 2.75, 0.35, clamprot));
2981 oclamp->AddNode(clotrapv, 4, new TGeoCombiTrans(-0.625, 2.75, -0.35, clamprot));
2982 //
2983 auto* clampo3 = new TGeoBBox("clampo3", 1.6, 0.45, .1);
2984 auto* clampo3v = new TGeoVolume("TPC_clampo3v", clampo3, m6);
2985 //
2986 oclamp->AddNode(clampo3v, 1, new TGeoTranslation(-1.85, 2.625, 0.));
2987 oclamp->AddNode(clampo3v, 2, new TGeoTranslation(-1.85, -2.625, 0));
2988 //
2989 auto* clampo4 = new TGeoTubeSeg("clampo4", 2.2, 3.1, 0.1, 90., 270.);
2990 auto* clampo4v = new TGeoVolume("TPC_clampo4v", clampo4, m6);
2991 //
2992 oclamp->AddNode(clampo4v, 1, new TGeoTranslation(-3.45, 0., 0.));
2993
2994 // v9 - drift gas
2995
2996 TGeoRotation rot102("rot102");
2997 rot102.RotateY(-90.);
2998
2999 for (Int_t i = 0; i < 18; i++) {
3000 Double_t angle, x, y;
3001 Double_t z, r;
3002 angle = TMath::DegToRad() * 20. * (Double_t)i;
3003 // inner rods
3004 r = 81.5;
3005 x = r * TMath::Cos(angle);
3006 y = r * TMath::Sin(angle);
3007 z = 126.;
3008 auto* rot12 = new TGeoRotation();
3009 rot12->RotateZ(-90.0 + i * 20.);
3010 v9->AddNode(tpcihpl, i + 1, new TGeoCombiTrans(x, y, 0., rot12));
3011 //
3012 if (i == 11) { // resistor rod inner
3013 rotrod.RotateZ(-90. + i * 20.);
3014 rotrod1[0] = new TGeoRotation();
3015 rotpos[0] = new TGeoRotation();
3016 //
3017 rotrod1[0]->RotateZ(90. + i * 20.);
3018 *rotpos[0] = refl * rotrod; // rotation+reflection
3019 v9->AddNode(tpcrrod, 1, new TGeoCombiTrans(x, y, z, rotrod1[0])); // A
3020 v9->AddNode(tpcrrod, 2, new TGeoCombiTrans(x, y, -z, rotpos[0])); // C
3021 } else {
3022 v9->AddNode(tpcmrod, i + 1, new TGeoTranslation(x, y, z)); // shaft
3023 v9->AddNode(tpcmrod, i + 19, new TGeoCombiTrans(x, y, -z, ref)); // muon
3024 }
3025 //
3026 // inner clamps positioning
3027 //
3028 r = 79.05;
3029 x = r * TMath::Cos(angle);
3030 y = r * TMath::Sin(angle);
3031 rot12 = new TGeoRotation();
3032 rot12->RotateZ(i * 20.);
3033 //
3034 // A-side
3035 v9->AddNode(iclamp, 7 * i + 1, new TGeoCombiTrans(x, y, 5.25, rot12));
3036 v9->AddNode(iclamp, 7 * i + 2, new TGeoCombiTrans(x, y, 38.25, rot12));
3037 v9->AddNode(iclamp, 7 * i + 3, new TGeoCombiTrans(x, y, 80.25, rot12));
3038 v9->AddNode(iclamp, 7 * i + 4, new TGeoCombiTrans(x, y, 122.25, rot12));
3039 v9->AddNode(iclamp, 7 * i + 5, new TGeoCombiTrans(x, y, 164.25, rot12));
3040 v9->AddNode(iclamp, 7 * i + 6, new TGeoCombiTrans(x, y, 206.25, rot12));
3041 v9->AddNode(iclamp, 7 * i + 7, new TGeoCombiTrans(x, y, 246.75, rot12));
3042 // C-side
3043 v9->AddNode(iclamp, 7 * i + 127, new TGeoCombiTrans(x, y, -5.25, rot12));
3044 v9->AddNode(iclamp, 7 * i + 128, new TGeoCombiTrans(x, y, -38.25, rot12));
3045 v9->AddNode(iclamp, 7 * i + 129, new TGeoCombiTrans(x, y, -80.25, rot12));
3046 v9->AddNode(iclamp, 7 * i + 130, new TGeoCombiTrans(x, y, -122.25, rot12));
3047 v9->AddNode(iclamp, 7 * i + 131, new TGeoCombiTrans(x, y, -164.25, rot12));
3048 v9->AddNode(iclamp, 7 * i + 132, new TGeoCombiTrans(x, y, -206.25, rot12));
3049 v9->AddNode(iclamp, 7 * i + 133, new TGeoCombiTrans(x, y, -246.75, rot12));
3050 //
3051 //--------------------------
3052 // outer rods
3053 r = 254.25;
3054 x = r * TMath::Cos(angle);
3055 y = r * TMath::Sin(angle);
3056 z = 126.;
3057 //
3058 // outer rod holder + outer left plug
3059 //
3060
3061 auto* rot33 = new TGeoRotation();
3062 rot33->RotateZ(-90 + i * 20.);
3063 //
3064 v9->AddNode(tpcohpl, i + 1, new TGeoCombiTrans(x, y, 0., rot33));
3065 //
3066 Double_t xxx = 256.297 * TMath::Cos((i * 20. + 10.) * TMath::DegToRad());
3067 Double_t yyy = 256.297 * TMath::Sin((i * 20. + 10.) * TMath::DegToRad());
3068 //
3069 TGeoRotation rot101("rot101");
3070 rot101.RotateZ(90. + i * 20. + 10.);
3071 auto* rot103 = new TGeoRotation("rot103");
3072 *rot103 = rot101 * rot102;
3073 //
3074 auto* trh100 = new TGeoCombiTrans(xxx, yyy, 0., rot103);
3075 //
3076 if (i == 2) {
3077 // main membrane holder
3078 v9->AddNode(tpcmmh, 1, trh100);
3079 } else {
3080 // "normal" membrane holder
3081 v9->AddNode(tpcomh, i + 1, trh100);
3082 }
3083
3084 //
3085 if (i == 3) { // resistor rod outer
3086 rotrod.RotateZ(90. + i * 20.);
3087 rotrod1[1] = new TGeoRotation();
3088 rotpos[1] = new TGeoRotation();
3089 rotrod1[1]->RotateZ(90. + i * 20.);
3090 *rotpos[1] = refl * rotrod; // rotation+reflection
3091 v9->AddNode(tpcrrod, 3, new TGeoCombiTrans(x, y, z, rotrod1[1])); // A
3092 v9->AddNode(tpcrrod, 4, new TGeoCombiTrans(x, y, -z, rotpos[1])); // C
3093 } else {
3094 v9->AddNode(tpcmrod, i + 37, new TGeoTranslation(x, y, z)); // shaft
3095 v9->AddNode(tpcmrod, i + 55, new TGeoCombiTrans(x, y, -z, ref)); // muon
3096 }
3097 if (i == 15) {
3098 v9->AddNode(hvrv, 1, new TGeoTranslation(x, y, z + 0.7)); // hv->A-side only
3099 }
3100 //
3101 // outer clamps
3102 //
3103 r = 256.9;
3104 x = r * TMath::Cos(angle);
3105 y = r * TMath::Sin(angle);
3106 rot12 = new TGeoRotation();
3107 rot12->RotateZ(i * 20.);
3108 //
3109 // A-side
3110 v9->AddNode(oclamp, 7 * i + 1, new TGeoCombiTrans(x, y, 5.25, rot12));
3111 v9->AddNode(oclamp, 7 * i + 2, new TGeoCombiTrans(x, y, 38.25, rot12));
3112 v9->AddNode(oclamp, 7 * i + 3, new TGeoCombiTrans(x, y, 80.25, rot12));
3113 v9->AddNode(oclamp, 7 * i + 4, new TGeoCombiTrans(x, y, 122.25, rot12));
3114 v9->AddNode(oclamp, 7 * i + 5, new TGeoCombiTrans(x, y, 164.25, rot12));
3115 v9->AddNode(oclamp, 7 * i + 6, new TGeoCombiTrans(x, y, 206.25, rot12));
3116 v9->AddNode(oclamp, 7 * i + 7, new TGeoCombiTrans(x, y, 246.75, rot12));
3117 // C-side
3118 v9->AddNode(oclamp, 7 * i + 127, new TGeoCombiTrans(x, y, -5.25, rot12));
3119 v9->AddNode(oclamp, 7 * i + 128, new TGeoCombiTrans(x, y, -38.25, rot12));
3120 v9->AddNode(oclamp, 7 * i + 129, new TGeoCombiTrans(x, y, -80.25, rot12));
3121 v9->AddNode(oclamp, 7 * i + 130, new TGeoCombiTrans(x, y, -122.25, rot12));
3122 v9->AddNode(oclamp, 7 * i + 131, new TGeoCombiTrans(x, y, -164.25, rot12));
3123 v9->AddNode(oclamp, 7 * i + 132, new TGeoCombiTrans(x, y, -206.25, rot12));
3124 v9->AddNode(oclamp, 7 * i + 133, new TGeoCombiTrans(x, y, -246.75, rot12));
3125
3126 } // end of rods positioning
3127
3128 TGeoVolume* alice = gGeoManager->GetVolume("barrel");
3129 alice->AddNode(v1, 1, new TGeoTranslation(0., 30., 0.));
3130
3131} // end of function
3132
3133void Detector::LoadGeometryFromFile()
3134{
3135 // ===| Read the TPC geometry from file |=====================================
3136 if (mGeoFileName.IsNull()) {
3137 LOG(fatal) << "TPC geometry file name not set";
3138 return;
3139 }
3140
3141 TFile* fGeoFile = TFile::Open(mGeoFileName);
3142 if (!fGeoFile || !fGeoFile->IsOpen() || fGeoFile->IsZombie()) {
3143 LOG(fatal) << "Could not open TPC geometry file '" << mGeoFileName << "'";
3144 return;
3145 }
3146
3147 TGeoVolume* tpcVolume = dynamic_cast<TGeoVolume*>(fGeoFile->Get("TPC_M"));
3148 if (!tpcVolume) {
3149 LOG(fatal) << "Could not retrieve TPC geometry from file '" << mGeoFileName << "'";
3150 return;
3151 }
3152
3153 LOG(info) << "Loaded TPC geometry from file '" << mGeoFileName << "'";
3154 TGeoVolume* alice = gGeoManager->GetVolume("barrel");
3155 alice->AddNode(tpcVolume, 1, new TGeoTranslation(0., 30., 0.));
3156}
3157
3158void Detector::defineSensitiveVolumes()
3159{
3160 TGeoManager* geoManager = gGeoManager;
3161 TGeoVolume* v = nullptr;
3162
3163 // const Int_t nSensitive=2;
3164 // const char* volumeNames[nSensitive]={"TPC_Drift","TPC_Strip"};
3165 const Int_t nSensitive = 1;
3166 const char* volumeNames[nSensitive] = {"TPC_Drift"};
3167
3168 // The names of the ITS sensitive volumes have the format: ITSUSensor(0...mNumberLayers-1)
3169 for (Int_t ivol = 0; ivol < nSensitive; ++ivol) {
3170 TString volumeName = volumeNames[ivol];
3171 v = geoManager->GetVolume(volumeName.Data());
3172 if (!v) {
3173 LOG(error) << "Could not find volume '" << volumeName << "'";
3174 continue;
3175 }
3176
3177 // set volume sentive
3178 AddSensitiveVolume(v);
3179 }
3180
3181 // Special sensitive volume parameters in case FLUKA is used as transport engine
3182 auto vmc = TVirtualMC::GetMC();
3183 if (strcmp(vmc->GetName(), "TFluka") == 0) {
3184 LOG(info) << "Setting special FLUKA parameters for TPC Driftgas";
3186 Int_t index = mgr.getMediumID("TPC", kDriftGas2);
3187 vmc->Gstpar(index, "PRIMIO_E", 20.77);
3188 vmc->Gstpar(index, "PRIMIO_N", 14.35);
3189 vmc->Gstpar(index, "LOSS", 14);
3190 vmc->Gstpar(index, "STRA", 4);
3191 }
3192}
3193
3194Double_t Detector::Gamma(Double_t k)
3195{
3196 static thread_local Double_t n = 0;
3197 static thread_local Double_t c1 = 0;
3198 static thread_local Double_t c2 = 0;
3199 static thread_local Double_t b1 = 0;
3200 static thread_local Double_t b2 = 0;
3201 if (k > 0) {
3202 if (k < 0.4) {
3203 n = 1. / k;
3204 } else if (k >= 0.4 && k < 4) {
3205 n = 1. / k + (k - 0.4) / k / 3.6;
3206 } else if (k >= 4.) {
3207 n = 1. / TMath::Sqrt(k);
3208 }
3209 b1 = k - 1. / n;
3210 b2 = k + 1. / n;
3211 c1 = (k < 0.4) ? 0 : b1 * (TMath::Log(b1) - 1.) / 2.;
3212 c2 = b2 * (TMath::Log(b2) - 1.) / 2.;
3213 }
3214 Double_t x;
3215 Double_t y = -1.;
3216 while (1) {
3217 Double_t nu1 = gRandom->Rndm();
3218 Double_t nu2 = gRandom->Rndm();
3219 Double_t w1 = c1 + TMath::Log(nu1);
3220 Double_t w2 = c2 + TMath::Log(nu2);
3221 y = n * (b1 * w2 - b2 * w1);
3222 if (y < 0) {
3223 continue;
3224 }
3225 x = n * (w2 - w1);
3226 if (TMath::Log(y) >= x) {
3227 break;
3228 }
3229 }
3230 return TMath::Exp(x);
3231}
3232
3233std::string Detector::getHitBranchNames(int probe) const
3234{
3235 if (probe >= 0 && probe < Sector::MAXSECTOR) {
3236 TString name;
3237 name.Form("%sHitsShiftedSector%d", GetName(), probe);
3238 return std::string(name.Data());
3239 }
3240 return std::string();
3241}
3242
3244
3245// Define Factory method for calling from the outside
3246extern "C" {
3251}
Definition of the Stack class.
int16_t time
Definition RawEventData.h:4
int32_t i
const GPUTPCGMMerger::trackCluster & b1
bool const GPUTPCGMMerger::trackCluster * c1
bool const GPUTPCGMMerger::trackCluster const clcomparestruct * c2
constexpr int p1()
constexpr to accelerate the coordinates changing
ClassImp(IdPath)
Definition of the parameter class for the detector.
Definition of the parameter class for the detector gas.
Class for TPC Point.
uint16_t slope
Definition RawData.h:1
uint32_t res
Definition RawData.h:0
uint32_t stack
Definition RawData.h:1
o2::base::Detector * create_detector_tpc(bool active)
virtual std::string getHitBranchNames(int probe) const =0
~Detector() override
Default Destructor.
Definition Detector.cxx:87
Detector()
Default Constructor.
Definition Detector.cxx:36
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
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
virtual void InitializeO2Detector()=0
Definition Detector.cxx:98
static MaterialManager & Instance()
void addTrackReference(const o2::TrackReference &p)
Definition Stack.h:332
void Reset() override
void ConstructGeometry() override
Double_t Gamma(Double_t k)
Bool_t ProcessHits(FairVolume *v=nullptr) override
void Register() override
void EndOfEvent() override
static o2::base::Detector * create(bool active)
Definition Detector.h:72
void addHit(float x, float y, float z, float time, short e)
Definition Point.h:120
static int ToShiftedSector(T x, T y, T z)
Definition Sector.h:126
static int ToSector(T x, T y, T z)
Definition Sector.h:107
static constexpr int MAXSECTOR
Definition Sector.h:44
static ShmManager & Instance()
Definition ShmManager.h:61
GLdouble n
Definition glcorearb.h:1982
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLint GLenum GLint x
Definition glcorearb.h:403
const GLdouble * v
Definition glcorearb.h:832
GLuint index
Definition glcorearb.h:781
GLuint const GLchar * name
Definition glcorearb.h:781
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLint y
Definition glcorearb.h:270
GLuint GLfloat x0
Definition glcorearb.h:5034
GLfloat angle
Definition glcorearb.h:4071
GLfloat GLfloat v1
Definition glcorearb.h:812
GLboolean r
Definition glcorearb.h:1233
GLfloat GLfloat GLfloat GLfloat v3
Definition glcorearb.h:814
GLfloat GLfloat GLfloat v2
Definition glcorearb.h:813
GLuint GLfloat GLfloat y0
Definition glcorearb.h:5034
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat t1
Definition glcorearb.h:5034
constexpr float cm
Definition SpecsV2.h:26
std::optional< Chamber > chamber(int chamberId)
Definition Chamber.cxx:17
std::array< uint64_t, 23 > m2
std::array< uint64_t, 17 > m5
std::array< uint64_t, 17 > m4
std::array< uint64_t, 23 > m3
std::array< uint64_t, 23 > m1
Global TPC definitions and constants.
Definition SimTraits.h:167
void freeSimVector(std::vector< T > *ptr)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
TStopwatch sw