]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDRing.cxx
Adding missing include
[u/mrichter/AliRoot.git] / FMD / AliFMDRing.cxx
CommitLineData
4347b38f 1/**************************************************************************
c4bcfd14 2 * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
4347b38f 3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/* $Id$ */
17
c4bcfd14 18//__________________________________________________________________
4347b38f 19//
20// Utility class to help implement collection of FMD modules into
21// rings. This is used by AliFMDSubDetector and AliFMD.
22//
23// The AliFMD object owns the AliFMDRing objects, and the
24// AliFMDSubDetector objects reference these. That is, the AliFMDRing
25// objects are share amoung the AliFMDSubDetector objects.
26//
27// Latest changes by Christian Holm Christensen
28//
5edc364e 29
30#include <math.h>
31
e802be3e 32#include "AliFMDRing.h" // ALIFMDRING_H
33#include "AliLog.h" // ALILOG_H
34#include "TMath.h" // ROOT_TMath
35#include "TH2.h" // ROOT_TH2
36#include "TVirtualMC.h" // ROOT_TVirtualMC
37#include "TVector2.h" // ROOT_TVector2
38#include "TBrowser.h" // ROOT_TBrowser
39#include "TString.h" // ROOT_TString
40#include "TArc.h" // ROOT_TArc
41#include "TObjArray.h" // ROOT_TObjArray
42#include "TXTRU.h" // ROOT_TXTRU
43#include "TNode.h" // ROOT_TNode
44#include "TRotMatrix.h" // ROOT_TRotMatrix
45#include "TList.h" // ROOT_TList
4347b38f 46
c4bcfd14 47const Char_t* AliFMDRing::fgkRingFormat = "FRG%c";
48const Char_t* AliFMDRing::fgkVirtualFormat = "FV%c%c";
49const Char_t* AliFMDRing::fgkActiveFormat = "FAC%c";
50const Char_t* AliFMDRing::fgkSectorFormat = "FSE%c";
51const Char_t* AliFMDRing::fgkStripFormat = "FST%c";
52const Char_t* AliFMDRing::fgkPrintboardFormat = "FP%c%c";
42403906 53
54
55//____________________________________________________________________
4347b38f 56ClassImp(AliFMDRing);
57
58//____________________________________________________________________
4347b38f 59AliFMDRing::AliFMDRing(Char_t id, Bool_t detailed)
60 : fId(id),
61 fDetailed(detailed),
c4bcfd14 62 fActiveId(0),
63 fPrintboardBottomId(0),
64 fPrintboardTopId(0),
65 fRingId(0),
66 fSectionId(0),
67 fStripId(0),
68 fVirtualBackId(0),
69 fVirtualFrontId(0),
70 fBondingWidth(0),
4347b38f 71 fWaferRadius(0),
72 fSiThickness(0),
73 fLowR(0),
74 fHighR(0),
75 fTheta(0),
76 fNStrips(0),
c4bcfd14 77 fRingDepth(0),
78 fLegRadius(0),
79 fLegLength(0),
80 fLegOffset(0),
81 fModuleSpacing(0),
82 fPrintboardThickness(0),
4347b38f 83 fShape(0),
84 fRotMatricies(0)
42403906 85{
86 // Construct a alifmdring.
87 //
88 // id Id of the ring (either 'i' or 'o').
c4bcfd14 89 // detailed Whether the strips are made or not.
90 //
91}
92
93//____________________________________________________________________
94AliFMDRing::AliFMDRing(const AliFMDRing& other)
95 : TObject(other),
96 fId(other.fId),
97 fDetailed(other.fDetailed),
98 fActiveId(other.fActiveId),
99 fPrintboardBottomId(other.fPrintboardBottomId),
100 fPrintboardTopId(other.fPrintboardTopId),
101 fRingId(other.fRingId),
102 fSectionId(other.fSectionId),
103 fStripId(other.fStripId),
104 fVirtualBackId(other.fVirtualBackId),
105 fVirtualFrontId(other.fVirtualFrontId),
106 fBondingWidth(other.fBondingWidth),
107 fWaferRadius(other.fWaferRadius),
108 fSiThickness(other.fSiThickness),
109 fLowR(other.fLowR),
110 fHighR(other.fHighR),
111 fTheta(other.fTheta),
112 fNStrips(other.fNStrips),
113 fRingDepth(other.fRingDepth),
114 fLegRadius(other.fLegRadius),
115 fLegLength(other.fLegLength),
116 fLegOffset(other.fLegOffset),
117 fModuleSpacing(other.fModuleSpacing),
118 fPrintboardThickness(other.fPrintboardThickness),
119 fRotations(other.fRotations),
120 fShape(other.fShape),
121 fRotMatricies(other.fRotMatricies)
122{
123 // Copy constructor of a AliFMDRing.
124}
125
126//____________________________________________________________________
127AliFMDRing&
128AliFMDRing::operator=(const AliFMDRing& other)
129{
130 // Assignment operator
131 //
132 fId = other.fId;
133 fDetailed = other.fDetailed;
134 fActiveId = other.fActiveId;
135 fPrintboardBottomId = other.fPrintboardBottomId;
136 fPrintboardTopId = other.fPrintboardTopId;
137 fRingId = other.fRingId;
138 fSectionId = other.fSectionId;
139 fStripId = other.fStripId;
140 fVirtualBackId = other.fVirtualBackId;
141 fVirtualFrontId = other.fVirtualFrontId;
142 fBondingWidth = other.fBondingWidth;
143 fWaferRadius = other.fWaferRadius;
144 fSiThickness = other.fSiThickness;
145 fLowR = other.fLowR;
146 fHighR = other.fHighR;
147 fTheta = other.fTheta;
148 fNStrips = other.fNStrips;
149 fRingDepth = other.fRingDepth;
150 fLegRadius = other.fLegRadius;
151 fLegLength = other.fLegLength;
152 fLegOffset = other.fLegOffset;
153 fModuleSpacing = other.fModuleSpacing;
154 fPrintboardThickness = other.fPrintboardThickness;
155 fRotations = other.fRotations;
156 if (other.fShape) {
157 if (other.fShape->IsA() == TXTRU::Class())
158 ((TXTRU*)other.fShape)->Copy(*fShape);
159 else
160 fShape = 0;
161 }
162 if (other.fRotMatricies) {
163 Int_t n = other.fRotMatricies->GetEntries();
164 if (!fRotMatricies) fRotMatricies = new TObjArray(n);
165 else fRotMatricies->Expand(n);
166 TIter next(other.fRotMatricies);
167 TObject* o = 0;
168 while ((o = next())) fRotMatricies->Add(o);
169 }
170 return *this;
42403906 171}
172
c4bcfd14 173
4347b38f 174
175//____________________________________________________________________
176void
177AliFMDRing::Init()
178{
179 // Initialize the ring object.
180 // DebugGuard guard("AliFMDRing::Init");
181 AliDebug(10, "AliFMDRing::Init");
182 fPolygon.Clear();
183 SetupCoordinates();
184}
185
186//____________________________________________________________________
187AliFMDRing::~AliFMDRing()
188{
42403906 189 // Destructor - deletes shape and rotation matricies
4347b38f 190 if (fShape) delete fShape;
191 if (fRotMatricies) delete fRotMatricies;
192}
193
194
195//____________________________________________________________________
196void
197AliFMDRing::Browse(TBrowser* /* b */)
198{
199 // DebugGuard guard("AliFMDRing::Browse");
200 AliDebug(10, "AliFMDRing::Browse");
201}
202
203
204//____________________________________________________________________
205void
206AliFMDRing::SetupCoordinates()
207{
208 // Calculates the parameters of the polygon shape.
209 //
210 // DebugGuard guard("AliFMDRing::SetupCoordinates");
211 AliDebug(10, "AliFMDRing::SetupCoordinates");
212 // Get out immediately if we have already done all this
213 if (fPolygon.GetNVerticies() > 1) return;
214
c4bcfd14 215 double tanTheta = TMath::Tan(fTheta * TMath::Pi() / 180.);
216 double tanTheta2 = TMath::Power(tanTheta,2);
4347b38f 217 double r2 = TMath::Power(fWaferRadius,2);
c4bcfd14 218 double yA = tanTheta * fLowR;
4347b38f 219 double lr2 = TMath::Power(fLowR, 2);
220 double hr2 = TMath::Power(fHighR,2);
c4bcfd14 221 double xD = fLowR + TMath::Sqrt(r2 - tanTheta2 * lr2);
222 double xD2 = TMath::Power(xD,2);
223 //double xD_2 = fLowR - TMath::Sqrt(r2 - tanTheta2 * lr2);
224 double yB = sqrt(r2 - hr2 + 2 * fHighR * xD - xD2);
225 double xC = ((xD + TMath::Sqrt(-tanTheta2 * xD2 + r2
226 + r2 * tanTheta2))
227 / (1 + tanTheta2));
228 double yC = tanTheta * xC;
229
230 fPolygon.AddVertex(fLowR, -yA);
231 fPolygon.AddVertex(xC, -yC);
232 fPolygon.AddVertex(fHighR, -yB);
233 fPolygon.AddVertex(fHighR, yB);
234 fPolygon.AddVertex(xC, yC);
235 fPolygon.AddVertex(fLowR, yA);
4347b38f 236}
237
238//____________________________________________________________________
239bool
240AliFMDRing::IsWithin(size_t moduleNo, double x, double y) const
241{
242 // Checks if a point (x,y) is inside the module with number moduleNo
243 //
244 // DebugGuard guard("AliFMDRing::IsWithin");
245 AliDebug(10, "AliFMDRing::IsWithin");
246 bool ret = false;
247 double r2 = x * x + y * y;
248 if (r2 < fHighR * fHighR && r2 > fLowR * fLowR) {
249 // double point_angle = TMath::ATan2(y, x);
250 // int n_modules = 360 / Int_t(fTheta * 2);
251 double m_angle = (.5 + moduleNo) * 2 * fTheta;
252 double m_radians = TMath::Pi() * m_angle / 180.;
253
254 // Rotate the point.
255 double xr = x * TMath::Cos(-m_radians) - y * TMath::Sin(-m_radians);
256 double yr = x * TMath::Sin(-m_radians) + y * TMath::Cos(-m_radians);
257
258 ret = fPolygon.Contains(xr,yr);
259 }
260 return ret;
261}
262
263
264
265
266//____________________________________________________________________
267void
268AliFMDRing::Draw(Option_t* option) const
269{
270 // Draw a the shape of the ring into a 2D histogram. Useful for
271 // superimposing the actual shape of the ring onto a scatter plot of
272 // hits in the detector.
273 //
274 // DebugGuard guard("AliFMDRing::Draw");
275 AliDebug(10, "AliFMDRing::Draw");
276 // The unrotated coordinates of the polygon verticies
277 if (fPolygon.GetNVerticies() < 1) return;
278
279 TVector2 v[6];
280 for (size_t i = 0; i < fPolygon.GetNVerticies(); i++)
281 v[i] = fPolygon.GetVertex(i);
282
283 Int_t nModules = 360 / Int_t(fTheta * 2);
284 Double_t dTheta = fTheta * 2;
285
286 TString opt(option);
287 if (opt.Contains("B", TString::kIgnoreCase)) {
288 opt.Remove(opt.Index("B", 1, TString::kIgnoreCase),1);
289 TH1* null = new TH2F("null", "Null",
290 100, -fHighR * 1.1, fHighR * 1.1,
291 100, -fHighR * 1.1, fHighR * 1.1);
292 null->SetStats(0);
293 null->Draw(opt.Data());
294 }
295
296 for (int i = 0; i < nModules; i++) {
297 Double_t theta = (i + .5) * dTheta;
298 AliFMDPolygon p;
299 for (int j = 0; j < 6; j++) {
300 TVector2 vr(v[j].Rotate(TMath::Pi() * theta / 180.));
301 if (!p.AddVertex(vr.X(),vr.Y())) {
302 // std::cerr << "Draw of polygon " << i << " failed" << std::endl;
303 break;
304 }
305 }
306 p.Draw(opt.Data(), Form("MOD%c_%d", fId, i));
307 }
308 if (opt.Contains("0", TString::kIgnoreCase)) {
309 TArc* arcH = new TArc(0,0, fHighR);
310 arcH->SetLineStyle(2);
311 arcH->SetLineColor(4);
312 arcH->Draw();
313
314 TArc* arcL = new TArc(0,0, fLowR);
315 arcL->SetLineStyle(2);
316 arcL->SetLineColor(4);
317 arcL->Draw();
318 }
319}
320
321//____________________________________________________________________
322void
323AliFMDRing::SetupGeometry(Int_t vacuumId, Int_t siId, Int_t pcbId,
324 Int_t pbRotId, Int_t idRotId)
325{
326 // Setup the geometry of the ring. It defines the volumes
327 // RNGI or RNGO which can later be positioned in a sub-detector
328 // volume.
329 //
330 // The hieracy of the RNGx volume is
331 //
42403906 332 // FRGx // Ring volume
333 // FVFx // Container of hybrid + legs
334 // FACx // Active volume (si sensor approx)
335 // FSEx // Section division
336 // FSTx // Strip division
337 // FPTx // Print board (bottom)
338 // FPBx // Print board (top)
339 // FLL // Support leg (long version)
340 // FVBx // Container of hybrid + legs
341 // FACx // Active volume (si sensor approx)
342 // FSEx // Section division
343 // FSTx // Strip division
344 // FPTx // Print board (bottom)
345 // FPBx // Print board (top)
346 // FSL // Support leg (long version)
4347b38f 347 //
348 // Parameters:
349 //
350 // vacuumId Medium of inactive virtual volumes
351 // siId Medium of Silicon sensor (active)
352 // pcbId Medium of print boards
353 // pbRotId Print board rotation matrix
354 // idRotId Identity rotation matrix
355 //
356 // DebugGuard guard("AliFMDRing::SetupGeometry");
357 AliDebug(10, "AliFMDRing::SetupGeometry");
358
359 const TVector2& bCorner = fPolygon.GetVertex(3); // Third corner
360 const TVector2& aCorner = fPolygon.GetVertex(5); // First corner
361 const TVector2& cCorner = fPolygon.GetVertex(4); // Second corner
362 TString name;
363 TString name2;
364 Double_t dStrip = (bCorner.Mod() - aCorner.Mod()) / fNStrips;
365 Double_t stripOff = aCorner.Mod();
366 Double_t rmin = fLowR;
367 Double_t rmax = bCorner.Mod();
368 Double_t pars[10];
369 fRingDepth = (fSiThickness
370 + fPrintboardThickness
371 + fLegLength
372 + fModuleSpacing);
373
374 // Ring virtual volume
375 pars[0] = rmin;
c4bcfd14 376 pars[1] = rmax;
4347b38f 377 pars[2] = fRingDepth / 2;
c4bcfd14 378 name = Form(fgkRingFormat, fId);
4347b38f 379 fRingId = gMC->Gsvolu(name.Data(), "TUBE", vacuumId, pars, 3);
380
381 // Virtual volume for modules with long legs
382 pars[1] = rmax;
383 pars[3] = -fTheta;
384 pars[4] = fTheta;
c4bcfd14 385 name = Form(fgkVirtualFormat, 'F', fId);
4347b38f 386 fVirtualFrontId = gMC->Gsvolu(name.Data(), "TUBS", vacuumId, pars, 5);
387
388 // Virtual volume for modules with long legs
389 pars[2] = (fRingDepth - fModuleSpacing) / 2;
c4bcfd14 390 name = Form(fgkVirtualFormat, 'B', fId);
4347b38f 391 fVirtualBackId = gMC->Gsvolu(name.Data(), "TUBS", vacuumId, pars, 5);
392
393 // Virtual mother volume for silicon
394 pars[2] = fSiThickness/2;
395 name2 = name;
c4bcfd14 396 name = Form(fgkActiveFormat, fId);
4347b38f 397 fActiveId = gMC->Gsvolu(name.Data(), "TUBS", vacuumId , pars, 5);
398
399 if (fDetailed) {
400 // Virtual sector volumes
401 name2 = name;
c4bcfd14 402 name = Form(fgkSectorFormat, fId);
4347b38f 403 gMC->Gsdvn2(name.Data(), name2.Data(), 2, 2, -fTheta, vacuumId);
404 fSectionId = gMC->VolId(name.Data());
405
406 // Active strip volumes
407 name2 = name;
c4bcfd14 408 name = Form(fgkStripFormat, fId);
42403906 409 gMC->Gsdvt2(name.Data(), name2.Data(), dStrip, 1,stripOff, siId, fNStrips);
4347b38f 410 fStripId = gMC->VolId(name.Data());
411 }
412
413 // Print-board on back of module
414 pars[4] = TMath::Tan(TMath::Pi() * fTheta / 180) * fBondingWidth;
415 // Top of the print board
416 pars[0] = cCorner.Y() - pars[4];
417 pars[1] = bCorner.Y() - pars[4];
418 pars[2] = fPrintboardThickness / 2; // PCB half thickness
419 pars[3] = (bCorner.X() - cCorner.X()) / 2;
c4bcfd14 420 name = Form(fgkPrintboardFormat, 'T', fId);
4347b38f 421 fPrintboardTopId = gMC->Gsvolu(name.Data(), "TRD1", pcbId, pars, 4);
422
423 // Bottom of the print board
424 pars[0] = aCorner.Y() - pars[4];
425 pars[1] = cCorner.Y() - pars[4];
426 pars[3] = (cCorner.X() - aCorner.X()) / 2;
c4bcfd14 427 name = Form(fgkPrintboardFormat, 'B', fId);
4347b38f 428 fPrintboardBottomId = gMC->Gsvolu(name.Data(), "TRD1", pcbId, pars, 4);
429
430 // Define rotation matricies
431 Int_t nModules = 360 / Int_t(fTheta * 2);
432 Double_t dTheta = fTheta * 2;
433 fRotations.Set(nModules);
434 for (int i = 0; i < nModules; i++) {
435 Double_t theta = (i + .5) * dTheta;
436 Int_t idrot = 0;
437 // Rotation matrix for virtual module volumes
438 gMC->Matrix(idrot, 90, theta, 90, fmod(90 + theta, 360), 0, 0);
439 fRotations[i] = idrot;
440 }
441
442
443 // Int_t nModules = 360 / Int_t(fTheta * 2);
444 // Double_t dTheta = fTheta * 2;
445 Double_t pbTopL = (bCorner.X() - cCorner.X());
446 Double_t pbBotL = (cCorner.X() - aCorner.X());
447 Double_t yoffset = ((TMath::Tan(TMath::Pi() * fTheta / 180)
448 * fBondingWidth));
449
450 for (int i = 0; i < nModules; i++) {
c4bcfd14 451 TString name2 = Form(fgkRingFormat, fId);
4347b38f 452
453 Int_t id = i;
454 // Double_t theta = (i + .5) * dTheta;
455 Bool_t isFront = (i % 2 == 1);
456 Double_t dz = 0;
457 Double_t w = fRingDepth - (isFront ? 0 : fModuleSpacing);
458
459 // Place virtual module volume
c4bcfd14 460 name = Form(fgkVirtualFormat, (isFront ? 'F' : 'B'), fId);
4347b38f 461 dz = (w - fRingDepth) / 2;
462 gMC->Gspos(name.Data(), id, name2.Data(), 0., 0., dz,fRotations[i]);
463
464 // We only need to place the children once, they are copied when
465 // we place the other virtual volumes.
466 if (i > 1) continue;
467 name2 = name;
468
469 // Place active silicon wafer - this is put so that the front of
470 // the silicon is on the edge of the virtual volume.
c4bcfd14 471 name = Form(fgkActiveFormat, fId);
4347b38f 472 dz = (w - fSiThickness) / 2;
473 gMC->Gspos(name.Data(), id, name2.Data(),0.,0.,dz,idRotId);
474
475 // Place print board. This is put immediately behind the silicon
c4bcfd14 476 name = Form(fgkPrintboardFormat, 'T', fId);
4347b38f 477 dz = w / 2 - fSiThickness - fPrintboardThickness / 2;
478 gMC->Gspos(name.Data(), id, name2.Data(),
479 fLowR + pbBotL + pbTopL / 2, 0, dz, pbRotId, "ONLY");
c4bcfd14 480 name = Form(fgkPrintboardFormat, 'B', fId);
4347b38f 481 gMC->Gspos(name.Data(), id, name2.Data(),
482 fLowR + pbBotL / 2, 0, dz, pbRotId, "ONLY");
483
484 // Support legs
485 // This is put immediately behind the pringboard.
486 dz = (w / 2 - fSiThickness - fPrintboardThickness
487 - (fLegLength + (isFront ? fModuleSpacing : 0)) /2);
42403906 488 name = (isFront ? "FLL" : "FSL");
4347b38f 489 gMC->Gspos(name.Data(), id*10 + 1, name2.Data(),
490 aCorner.X() + fLegOffset + fLegRadius, 0., dz, idRotId, "");
491 Double_t y = cCorner.Y() - yoffset - fLegOffset - fLegRadius;
492 gMC->Gspos(name.Data(),id*10+2,name2.Data(),cCorner.X(), y,dz,idRotId,"");
c4bcfd14 493 gMC->Gspos(name.Data(),id*10+3,name2.Data(),cCorner.X(), -y,dz,idRotId,"");
4347b38f 494 }
495}
496//____________________________________________________________________
497void
498AliFMDRing::Geometry(const char* mother, Int_t baseId, Double_t z,
499 Int_t /* pbRotId */, Int_t idRotId)
500{
501 // Positions a RNGx volume inside a mother.
502 //
503 // Parameters
504 //
505 // mother Mother volume to position the RNGx volume in
506 // baseId Base copy number
507 // z Z coordinate where the front of the active silicon
508 // should be in the mother volume, so we need to
509 // subtract half the ring width.
510 // idRotId Identity rotation matrix
511 //
512 // DebugGuard guard("AliFMDRing::Geometry");
513 AliDebug(10, "AliFMDRing::Geometry");
514 TString name;
515 Double_t offsetZ = (fSiThickness
516 + fPrintboardThickness
517 + fLegLength + fModuleSpacing) / 2;
c4bcfd14 518 name = Form(fgkRingFormat, fId);
4347b38f 519 gMC->Gspos(name.Data(), baseId, mother, 0., 0., z - offsetZ, idRotId, "");
520}
521
522//____________________________________________________________________
523void
524AliFMDRing::SimpleGeometry(TList* nodes,
525 TNode* mother,
526 Int_t colour,
527 Double_t z,
528 Int_t n)
529{
530 // Make a simple geometry of the ring for event display.
531 //
532 // The simple geometry is made from ROOT TNode and TShape objects.
533 // Note, that we cache the TShape and TRotMatrix objects used for
534 // this.
535 //
536 // Parameters
537 //
538 // nodes List of nodes to register all create nodes in
539 // mother Mother node to put the ring in.
540 // colour Colour of the nodes
541 // z Z position of the node in the mother volume
542 // n Detector number
543 //
544 // DebugGuard guard("AliFMDRing::SimpleGeometry");
545 AliDebug(10, "AliFMDRing::SimpleGeometry");
546 SetupCoordinates();
547
548 // If the shape hasn't been defined yet, we define it here.
549 if (!fShape) {
c4bcfd14 550 TString name(Form(fgkActiveFormat, fId));
4347b38f 551 TString title(Form("Shape of modules in %c Rings", fId));
552 Int_t n = fPolygon.GetNVerticies();
553 TXTRU* shape = new TXTRU(name.Data(), title.Data(), "void", n, 2);
554 for (Int_t i = 0; i < n; i++) {
555 const TVector2& v = fPolygon.GetVertex(i);
556 shape->DefineVertex(i, v.X(), v.Y());
557 }
558 shape->DefineSection(0, - fSiThickness / 2, 1, 0, 0);
559 shape->DefineSection(1, + fSiThickness / 2, 1, 0, 0);
560 fShape = shape;
561 fShape->SetLineColor(colour);
562 }
563
564 Int_t nModules = 360 / Int_t(fTheta * 2);
565 Double_t dTheta = fTheta * 2;
566
567 // If the roation matricies hasn't been defined yet, we do so here
568 if (!fRotMatricies) {
569 fRotMatricies = new TObjArray(nModules);
570 for (int i = 0; i < nModules; i++) {
571 Double_t theta = (i + .5) * dTheta;
572 TString name(Form("FMD_ring_%c_rot", fId));
573 TString title(Form("FMD Ring %c Rotation", fId));
574 TRotMatrix* rot =
575 new TRotMatrix(name.Data(), title.Data(),
576 90, theta, 90, fmod(90 + theta, 360), 0, 0);
577 fRotMatricies->AddAt(rot, i);
578 }
579 }
580
581 Double_t offsetZ = (fSiThickness
582 + fPrintboardThickness
583 + fLegLength + fModuleSpacing) / 2;
584
585 // Make all the nodes
586 for (int i = 0; i < nModules; i++) {
587 Bool_t isFront = (i % 2 == 1);
588 mother->cd();
589 TRotMatrix* rot = static_cast<TRotMatrix*>(fRotMatricies->At(i));
42403906 590 TString name(Form("FAC%c_%d_%d", fId, n, i));
4347b38f 591 TString title(Form("Active FMD%d volume in %c Ring", n, fId));
592 TNode* node = new TNode(name.Data(), title.Data(), fShape,
593 0, 0,
594 z - offsetZ + (isFront ? fModuleSpacing : 0),
595 rot);
596 node->SetLineColor(colour);
597 nodes->Add(node);
598 }
599}
600
601
602
603//____________________________________________________________________
604void
605AliFMDRing::Gsatt()
606{
607 // Set drawing attributes for the RING
608 //
609 // DebugGuard guard("AliFMDRing::Gsatt");
610 AliDebug(10, "AliFMDRing::Gsatt");
611 TString name;
c4bcfd14 612 name = Form(fgkRingFormat,fId);
4347b38f 613 gMC->Gsatt(name.Data(), "SEEN", 0);
614
c4bcfd14 615 name = Form(fgkVirtualFormat, 'T', fId);
4347b38f 616 gMC->Gsatt(name.Data(), "SEEN", 0);
617
c4bcfd14 618 name = Form(fgkVirtualFormat, 'B', fId);
4347b38f 619 gMC->Gsatt(name.Data(), "SEEN", 0);
620
c4bcfd14 621 name = Form(fgkActiveFormat,fId);
4347b38f 622 gMC->Gsatt(name.Data(), "SEEN", 1);
623
c4bcfd14 624 name = Form(fgkPrintboardFormat, 'T', fId);
4347b38f 625 gMC->Gsatt(name.Data(), "SEEN", 1);
626
c4bcfd14 627 name = Form(fgkPrintboardFormat, 'B',fId);
4347b38f 628 gMC->Gsatt(name.Data(), "SEEN", 1);
629}
630
631//
632// EOF
633//