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