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