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