]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDGeometryBuilder.cxx
Adding the new AliRawReaderChain to the Create() mactory method. The URI syntax is...
[u/mrichter/AliRoot.git] / FMD / AliFMDGeometryBuilder.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 /* $Id$ */
16 /** @file    AliFMDGeometryBuilder.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:41:17 2006
19     @brief   Class to build the FMD geometry 
20 */
21 //____________________________________________________________________
22 //                                                                          
23 // Builder of FMD geometry. 
24 //
25 // This class takes care of actually building the geometry using the 
26 // TGeo classes.  Various parameters are fecthed from the
27 // AliFMDGeometry manager.  
28 // Forward Multiplicity Detector based on Silicon wafers. This class
29 // contains the base procedures for the Forward Multiplicity detector
30 // Detector consists of 3 sub-detectors FMD1, FMD2, and FMD3, each of
31 // which has 1 or 2 rings of silicon sensors. 
32 //                                                       
33 // 
34
35 #include <TArrayD.h>            // ROOT_TArrayD
36 #include <TGeoManager.h>        // ROOT_TGeoManager
37 #include <TGeoMatrix.h>         // ROOT_TGeoMatrix
38 #include <TGeoTube.h>           // ROOT_TGeoTube
39 #include <TGeoTrd1.h>           // ROOT_TGeoTrd1
40 #include <TGeoCone.h>           // ROOT_TGeoTrd1
41 #include <TGeoVolume.h>         // ROOT_TGeoVolume
42 #include <TGeoXtru.h>           // ROOT_TGeoXtru
43 #include <TGeoPcon.h>           // ROOT_TGeoPcon
44 #include <TGeoCompositeShape.h>
45 #include <TMath.h>
46 #include <TVector2.h>           // ROOT_TVector2
47 #include <TVector3.h>           // ROOT_TVector3
48 //#include <TGeoMaterial.h>     // ROOT_TGeoMaterial
49 //#include <TGeoMedium.h>               // ROOT_TGeoMedium
50 //#include <TGeoPcon.h>         // ROOT_TGeoPcon
51 //#include <TGeoPolygon.h>      // ROOT_TGeoPolygon
52
53 #include "AliFMDGeometryBuilder.h"      // ALIFMDGEOSIMULATOR_H
54 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
55 #include "AliFMDDetector.h"     // ALIFMDDETECTOR_H
56 #include "AliFMDRing.h"         // ALIFMDRING_H
57 #include "AliFMD1.h"            // ALIFMD1_H
58 #include "AliFMD2.h"            // ALIFMD2_H
59 #include "AliFMD3.h"            // ALIFMD3_H
60 // #include "AliFMD.h"          // ALIFMD_H
61 #include "AliFMDDebug.h"                // ALILOG_H
62 #include <iostream>
63
64 //====================================================================
65 ClassImp(AliFMDGeometryBuilder)
66 #if 0
67   ; // This is here to keep Emacs for indenting the next line
68 #endif
69
70 //____________________________________________________________________
71 const Char_t* AliFMDGeometryBuilder::fgkActiveName      = "F%cAC";
72 const Char_t* AliFMDGeometryBuilder::fgkSectorName      = "F%cSC";
73 const Char_t* AliFMDGeometryBuilder::fgkStripName       = "F%cST";
74 const Char_t* AliFMDGeometryBuilder::fgkSensorName      = "F%cSE";
75 const Char_t* AliFMDGeometryBuilder::fgkPCBName         = "F%cPB";
76 const Char_t* AliFMDGeometryBuilder::fgkCuName          = "F%cCU";
77 const Char_t* AliFMDGeometryBuilder::fgkChipName        = "F%cCH";
78 const Char_t* AliFMDGeometryBuilder::fgkLongLegName     = "F%cLL";
79 const Char_t* AliFMDGeometryBuilder::fgkShortLegName    = "F%cSL";
80 const Char_t* AliFMDGeometryBuilder::fgkFrontVName      = "F%cFH";
81 const Char_t* AliFMDGeometryBuilder::fgkBackVName       = "F%cBH";
82 const Char_t* AliFMDGeometryBuilder::fgkRingTopName     = "F%cTV";
83 const Char_t* AliFMDGeometryBuilder::fgkRingBotName     = "F%cBV";
84 const Char_t* AliFMDGeometryBuilder::fgkHCName          = "F%dH%c";
85 const Char_t* AliFMDGeometryBuilder::fgkIHCName         = "F%dI%c";
86 const Char_t* AliFMDGeometryBuilder::fgkNoseName        = "F3SN";
87 const Char_t* AliFMDGeometryBuilder::fgkBackName        = "F%dSB";
88 const Char_t* AliFMDGeometryBuilder::fgkTopName         = "F%dSU";
89 const Char_t* AliFMDGeometryBuilder::fgkBeamName        = "F%dSL";
90 const Char_t* AliFMDGeometryBuilder::fgkFlangeName      = "F%dSF";
91 const Char_t* AliFMDGeometryBuilder::fgkFMDDCuName      = "F%cDC";
92 const Char_t* AliFMDGeometryBuilder::fgkFMDDPCBName     = "F%cDP";
93 const Char_t* AliFMDGeometryBuilder::fgkFMDDChipName    = "F%cDI";
94 const Char_t* AliFMDGeometryBuilder::fgkFMDDName        = "F%cDD";
95 const Char_t* AliFMDGeometryBuilder::fgkFMDName         = "F%dM%c";
96
97 //____________________________________________________________________
98 AliFMDGeometryBuilder::AliFMDGeometryBuilder() 
99   : TTask("FMD", "Geomtry builder"),
100     fActiveId(0),
101     fDetailed(kTRUE),
102     fUseAssembly(kTRUE),
103     fSectorOff(0),
104     fModuleOff(0),
105     fRingOff(0),
106     fDetectorOff(0),
107     fSi(0),
108     fC(0),
109     fAl(0),
110     fPCB(0),
111     fChip(0),
112     fAir(0),
113     fPlastic(0),
114     fCopper(0),
115     fSteel(0)
116 {
117   // Default constructor
118   fActiveId.Set(2);
119 }
120
121 //____________________________________________________________________
122 AliFMDGeometryBuilder::AliFMDGeometryBuilder(Bool_t detailed) 
123   : TTask("FMD", "Geometry builder"),
124     fActiveId(0),
125     fDetailed(detailed),
126     fUseAssembly(kTRUE),
127     fSectorOff(0),
128     fModuleOff(0),
129     fRingOff(0),
130     fDetectorOff(0),
131     fSi(0),
132     fC(0),
133     fAl(0),
134     fPCB(0),
135     fChip(0),
136     fAir(0),
137     fPlastic(0),
138     fCopper(0),
139     fSteel(0)
140 {
141   // Normal constructor
142   // 
143   // Parameters: 
144   // 
145   //      fmd           Pointer to AliFMD object 
146   //      detailed      Whether to make a detailed simulation or not 
147   // 
148   fActiveId.Set(2);
149 }
150
151
152 //____________________________________________________________________
153 TGeoVolume*
154 AliFMDGeometryBuilder::RingGeometry(AliFMDRing* r) 
155 {
156   // Setup the geometry of a ring.    The defined TGeoVolume is
157   // returned, and should be used when setting up the rest of the
158   // volumes. 
159   // 
160   // 
161   // Parameters:
162   //
163   //     r              Pointer to ring geometry object 
164   // 
165   // Returns:
166   //    pointer to ring volume 
167   //
168   if (!r) { 
169     AliError("Didn't get a ring object");
170     return 0;
171   }
172   Char_t        id       = r->GetId();
173   const Char_t* lName    = (id == 'i' || id == 'I' ? "inner" : "outer");
174   Double_t      siThick  = r->GetSiThickness();
175   const Int_t   knv      = r->GetNVerticies();
176   TVector2*     a        = r->GetVertex(5);
177   TVector2*     b        = r->GetVertex(3);
178   TVector2*     c        = r->GetVertex(4);
179   Double_t      theta    = r->GetTheta();
180   Double_t      off      = (TMath::Tan(TMath::Pi() * theta / 180) 
181                             * r->GetBondingWidth());
182   Double_t      rmax     = b->Mod();
183   Double_t      rmin     = r->GetLowR();
184   Double_t      pcbThick = r->GetPrintboardThickness();
185   Double_t      cuThick  = r->GetCopperThickness();
186   Double_t      chipThick= r->GetChipThickness();
187   Double_t      modSpace = r->GetModuleSpacing();
188   Double_t      legr     = r->GetLegRadius();
189   Double_t      legl     = r->GetLegLength();
190   Double_t      legoff   = r->GetLegOffset();
191   Int_t         ns       = r->GetNStrips();
192   Double_t      stripoff = a->Mod();
193   Double_t      dstrip   = (rmax - stripoff) / ns;
194   Double_t      space    = r->GetSpacing();
195   TArrayD       xs(knv);
196   TArrayD       ys(knv);
197   for (Int_t i = 0; i < knv; i++) {
198     // Reverse the order 
199     TVector2* vv = r->GetVertex(knv - 1 - i);
200     if (!vv) {
201       AliError(Form("Failed to get vertex # %d", knv - 1 - i));
202       continue;
203     }
204     xs[i] = vv->X();
205     ys[i] = vv->Y();
206   }
207   
208   // Shape of actual sensor 
209   TGeoXtru* sensorShape = new TGeoXtru(2);
210   sensorShape->DefinePolygon(knv, xs.fArray, ys.fArray);
211   sensorShape->DefineSection(0, - siThick/2);
212   sensorShape->DefineSection(1, siThick/2);
213   sensorShape->SetName(Form(fgkSensorName, id));
214   sensorShape->SetTitle(Form("FMD %s Sensor", lName));
215   TGeoVolume* sensorVolume = new TGeoVolume(Form(fgkSensorName, id), 
216                                             sensorShape, fSi);
217   sensorVolume->SetTitle(Form("FMD %s Sensor", lName));
218   sensorVolume->VisibleDaughters(kFALSE);
219   Int_t sid = sensorVolume->GetNumber();
220   fSectorOff   = -1;
221   fModuleOff   = 1;
222   fRingOff     = 2;
223   fDetectorOff = 3;
224   if (fDetailed) {
225     fSectorOff   = 1;
226     fModuleOff   = 4;
227     fRingOff     = 5;
228     fDetectorOff = 6;
229     // Virtual volume shape to divide - This volume is only defined if
230     // the geometry is set to be detailed. 
231     TGeoTubeSeg* activeShape = new TGeoTubeSeg(rmin, rmax, siThick/2, 
232                                                - theta, theta);
233     activeShape->SetName(Form(fgkActiveName, id));
234     activeShape->SetTitle(Form("FMD %s active area", lName));
235     TGeoVolume* activeVolume = new TGeoVolume(Form(fgkActiveName, id),
236                                               activeShape,fSi);
237     activeVolume->SetTitle(Form("FMD %s active area", lName));
238     TGeoVolume* sectorVolume = activeVolume->Divide(Form(fgkSectorName,id), 
239                                                       2, 2, -theta,0,0,"N");
240     sectorVolume->SetTitle(Form("FMD %s sector", lName));
241     TGeoVolume* stripVolume  = sectorVolume->Divide(Form(fgkStripName, id), 
242                                                     1, ns, stripoff, dstrip, 
243                                                     0, "SX");
244     stripVolume->SetTitle(Form("FMD %s strip", lName));
245     sid = stripVolume->GetNumber();
246     sensorVolume->AddNodeOverlap(activeVolume, 0);
247   }
248   
249   switch (id) {
250   case 'i': case 'I': fActiveId[0] = sid; break;
251   case 'o': case 'O': fActiveId[1] = sid; break;
252   }
253
254   // Shape of Printed circuit Board 
255   for (Int_t i = 0;       i < knv / 2; i++) ys[i] -= off;
256   for (Int_t i = knv / 2; i < knv;     i++) ys[i] += off;
257   TGeoXtru* pcbShape         = new TGeoXtru(2);
258   pcbShape->DefinePolygon(knv, xs.fArray, ys.fArray);
259   pcbShape->DefineSection(0, - pcbThick/2);
260   pcbShape->DefineSection(1, pcbThick/2);
261   pcbShape->SetName(Form(fgkPCBName, id));
262   pcbShape->SetTitle(Form("FMD %s hybrid PCB", lName));
263   TGeoVolume* pcbVolume      = new TGeoVolume(Form(fgkPCBName, id), 
264                                               pcbShape, fPCB);
265   pcbVolume->SetTitle(Form("FMD %s hybrid PCB", lName));
266
267   // Copper layer
268   TGeoXtru* cuShape       = new TGeoXtru(2);
269   cuShape->DefinePolygon(6, xs.fArray, ys.fArray);
270   cuShape->DefineSection(0, - cuThick/2);
271   cuShape->DefineSection(1, cuThick/2);
272   cuShape->SetTitle(Form("FMD %s hybrid copper", lName));
273   TGeoVolume* cuVolume    = new TGeoVolume(Form(fgkCuName,id),cuShape,fCopper);
274   cuVolume->SetTitle(Form("FMD %s hybrid copper", lName));
275
276   // Chip layer
277   TGeoXtru*   chipShape   = new TGeoXtru(2);
278   chipShape->DefinePolygon(6, xs.fArray, ys.fArray);
279   chipShape->DefineSection(0, - chipThick/2);
280   chipShape->DefineSection(1, chipThick/2);
281   chipShape->SetTitle(Form("FMD %s hybrid chip", lName));
282   TGeoVolume* chipVolume = new TGeoVolume(Form(fgkChipName,id),
283                                           chipShape,fChip);
284   chipVolume->SetTitle(Form("FMD %s hybrid chip", lName));
285
286   // Short leg shape 
287   TGeoTube*   shortLegShape  = new TGeoTube(0, legr, legl / 2);
288   shortLegShape->SetName(Form(fgkShortLegName, id));
289   shortLegShape->SetTitle(Form("FMD %s short support foot", lName));
290   TGeoVolume* shortLegVolume = new TGeoVolume(Form(fgkShortLegName, id), 
291                                               shortLegShape, fCopper);
292   shortLegVolume->SetTitle(Form("FMD %s short support foot", lName));
293   // Long leg shape
294   TGeoTube*   longLegShape   = new TGeoTube(0, legr, (legl + modSpace) / 2);
295   longLegShape->SetName(Form(fgkLongLegName, id));
296   longLegShape->SetTitle(Form("FMD %s long support foot", lName));
297   TGeoVolume* longLegVolume  = new TGeoVolume(Form(fgkLongLegName, id), 
298                                               longLegShape, fCopper);
299   longLegVolume->SetTitle(Form("FMD %s long support foot", lName));
300   
301   
302   // Back container volume 
303   TGeoVolume* backVolume     = new TGeoVolumeAssembly(Form(fgkBackVName, id));
304   backVolume->SetTitle(Form("FMD %s back module", lName));
305   Double_t x = 0;
306   Double_t y = 0;
307   Double_t z = siThick / 2;
308   backVolume->AddNode(sensorVolume, 0, new TGeoTranslation(x, y, z));
309   z          += siThick / 2 + space + pcbThick / 2;
310   backVolume->AddNode(pcbVolume, 0, new TGeoTranslation(x,y,z));
311   z          += (pcbThick + cuThick) / 2;
312   backVolume->AddNode(cuVolume, 0, new TGeoTranslation(0, 0, z));
313   z          += (cuThick + chipThick) / 2;
314   backVolume->AddNode(chipVolume, 0, new TGeoTranslation(0, 0, z));
315   x          =  a->X() + legoff + legr;
316   y          =  0;
317   z          += pcbThick / 2 + legl / 2;
318   backVolume->AddNode(shortLegVolume, 0, new TGeoTranslation(x,y,z));
319   x          =  c->X();
320   y          =  c->Y() - legoff - legr - off;
321   backVolume->AddNode(shortLegVolume, 1, new TGeoTranslation(x,y,z));
322   y          =  -y;
323   backVolume->AddNode(shortLegVolume, 2, new TGeoTranslation(x,y,z));
324
325   // Front container volume 
326   TGeoVolume* frontVolume    = new TGeoVolumeAssembly(Form(fgkFrontVName, id));
327   frontVolume->SetTitle(Form("FMD %s front module", lName));
328   x         =  0;
329   y         =  0;
330   z         = siThick / 2;
331   frontVolume->AddNode(sensorVolume, 0, new TGeoTranslation(x, y, z));
332   z          += siThick / 2 + space + pcbThick / 2;
333   frontVolume->AddNode(pcbVolume, 0, new TGeoTranslation(x,y,z));
334   z          += (pcbThick + cuThick) / 2;
335   frontVolume->AddNode(cuVolume, 0, new TGeoTranslation(0, 0, z));
336   z          += (cuThick + chipThick) / 2;
337   frontVolume->AddNode(chipVolume, 0, new TGeoTranslation(0, 0, z));
338   x         =  a->X() + legoff + legr;
339   y         =  0;
340   z         += pcbThick / 2 + (legl + modSpace)/ 2;
341   frontVolume->AddNode(longLegVolume, 0, new TGeoTranslation(x,y,z));
342   x         =  c->X();
343   y         =  c->Y() - legoff - legr - off;
344   frontVolume->AddNode(longLegVolume, 1, new TGeoTranslation(x,y,z));
345   y         =  -y;
346   frontVolume->AddNode(longLegVolume, 2, new TGeoTranslation(x,y,z));
347
348
349   // FMDD 
350   Double_t ddlr = r->GetFMDDLowR();
351   Double_t ddhr = r->GetFMDDHighR();
352   Double_t ddpt = r->GetFMDDPrintboardThickness();
353   Double_t ddct = r->GetFMDDCopperThickness();
354   Double_t ddit = r->GetFMDDChipThickness();
355   Double_t ddt  = ddpt + ddct + ddit;
356   
357   TGeoShape* fmddPcbShape  = new TGeoTubeSeg(ddlr, ddhr, ddpt/2,0,180);
358   TGeoShape* fmddCuShape   = new TGeoTubeSeg(ddlr, ddhr, ddct/2,0,180);
359   TGeoShape* fmddChipShape = new TGeoTubeSeg(ddlr, ddhr, ddit/2,0,180);
360   fmddPcbShape->SetName(Form(fgkFMDDPCBName, id));
361   fmddCuShape->SetName(Form(fgkFMDDCuName, id));
362   fmddChipShape->SetName(Form(fgkFMDDChipName, id));
363   if (id == 'O' || id == 'o') { 
364     TString pcbName(fmddPcbShape->GetName());
365     TString cuName(fmddCuShape->GetName());
366     TString chipName(fmddChipShape->GetName());
367     
368     fmddPcbShape->SetName(Form("%s_inner",  pcbName.Data()));
369     fmddCuShape->SetName(Form("%s_inner",   cuName.Data()));
370     fmddChipShape->SetName(Form("%s_inner", chipName.Data()));
371     new TGeoBBox(Form("%s_clip",  pcbName.Data()), ddlr+3, ddhr/2, ddpt);
372     new TGeoBBox(Form("%s_clip",  cuName.Data()),  ddlr+3, ddhr/2, ddpt);
373     new TGeoBBox(Form("%s_clip",  chipName.Data()),ddlr+3, ddhr/2, ddpt);
374     TGeoTranslation* trans = new TGeoTranslation(Form("%s_trans",
375                                                       pcbName.Data()), 
376                                                  0, ddhr/2, 0);
377     trans->RegisterYourself();
378     fmddPcbShape = new TGeoCompositeShape(pcbName.Data(), 
379                                           Form("%s_inner*%s_clip:%s_trans",
380                                                pcbName.Data(), 
381                                                pcbName.Data(), 
382                                                pcbName.Data())); 
383     fmddCuShape = new TGeoCompositeShape(cuName.Data(), 
384                                          Form("%s_inner*%s_clip:%s_trans",
385                                               cuName.Data(), 
386                                               cuName.Data(), 
387                                               pcbName.Data()));
388     fmddChipShape = new TGeoCompositeShape(chipName.Data(), 
389                                            Form("%s_inner*%s_clip:%s_trans",
390                                                 chipName.Data(), 
391                                                 chipName.Data(), 
392                                                 pcbName.Data()));
393   }
394   fmddPcbShape->SetTitle(Form("FMD %s digitiser PCB", lName));
395   fmddCuShape->SetTitle(Form("FMD %s digitiser copper", lName));
396   fmddChipShape->SetTitle(Form("FMD %s digitiser chip", lName));
397
398   TGeoVolume*  fmddPcbVolume = new TGeoVolume(Form(fgkFMDDPCBName, id),
399                                               fmddPcbShape, fPCB);
400   TGeoVolume*  fmddCuVolume  = new TGeoVolume(Form(fgkFMDDCuName, id),
401                                               fmddCuShape, fCopper);
402   TGeoVolume*  fmddChipVolume= new TGeoVolume(Form(fgkFMDDChipName, id),
403                                               fmddChipShape, fChip);
404   fmddPcbVolume->SetTitle(Form("FMD %s digitiser PCB", lName));
405   fmddCuVolume->SetTitle(Form("FMD %s digitiser copper", lName));
406   fmddChipVolume->SetTitle(Form("FMD %s digitiser chip", lName));
407
408   // Half ring mother volumes. 
409   TGeoVolume* ringTopVolume = new TGeoVolumeAssembly(Form(fgkRingTopName,id));
410   TGeoVolume* ringBotVolume = new TGeoVolumeAssembly(Form(fgkRingBotName,id));
411   TGeoVolume* halfRing      = ringTopVolume;
412   ringTopVolume->SetTitle(Form("FMD %s top half-ring", lName));
413   ringBotVolume->SetTitle(Form("FMD %s bottom half-ring", lName));
414   
415   // Adding modules to half-rings
416   Int_t    nmod =  r->GetNModules();
417   AliFMDDebug(10, ("making %d modules in ring %c", nmod, id));
418   for (Int_t i = 0; i < nmod; i++) {
419     if (i == nmod / 2) halfRing = ringBotVolume;
420     Bool_t      front =  (i % 2 == 0);
421     TGeoVolume* vol   =  (front ? frontVolume : backVolume);
422     // vol->AddNode(sensorVolume, i, new TGeoTranslation(0,0,siThick/2));
423     Double_t    z1    =  (i % 2) * modSpace;
424     Double_t    th    =  (2 * i + 1) * theta;
425     TGeoMatrix* mat1  =  new TGeoCombiTrans(0,0,z1,0); 
426     mat1->RotateZ(th);
427     mat1->SetName(Form("FMD%c_module_%02d", id, i));
428     mat1->SetTitle(Form("FMD %s module %2d matrix", lName, i));
429     halfRing->AddNode(vol, i, mat1);
430 #if 0
431     Double_t    z2    =  z1 + siThick / 2 + space;
432     Double_t    th    =  (2 * i + 1) * theta;
433     AliFMDDebug(20, ("Placing copy %d of %s and %s in %s at z=%f and %f, "
434                       "and theta=%f", i, sensorVolume->GetName(), 
435                       vol->GetName(), halfRing->GetName(), z1, z2, th));
436     TGeoMatrix* mat1  =  new TGeoCombiTrans(0,0,z1,0); 
437     mat1->RotateZ(th);
438     halfRing->AddNode(sensorVolume, i, mat1);
439     TGeoMatrix* mat2  =  new TGeoCombiTrans(0,0,z2,0); 
440     mat2->RotateZ(th);
441     halfRing->AddNode(vol, i, mat2);
442 #endif
443   }
444
445   // Add the FMDD 
446   Double_t zi = r->GetFullDepth() - ddt;
447   Int_t    n  = 2;
448   for (Int_t i = 0; i  < n; i++) {
449     halfRing             = (i == 0 ? ringTopVolume : ringBotVolume);
450     Double_t      phi    = 360. / n * i;
451     TGeoRotation* rot    = new TGeoRotation(Form("FMDD%c rotation %d", id, i));
452     rot->RotateZ(phi);
453     rot->SetTitle(Form("FMD %s digitiser rotation %2d", lName, i));
454     z         =  zi + ddpt / 2;
455     halfRing->AddNode(fmddPcbVolume, i, new TGeoCombiTrans(0,0,z,rot));
456     z          += (ddpt + ddct) / 2;
457     halfRing->AddNode(fmddCuVolume, i, new TGeoCombiTrans(0,0,z,rot));
458     z          += (ddct + ddit) / 2;
459     halfRing->AddNode(fmddChipVolume, i, new TGeoCombiTrans(0,0,z,rot));
460   }
461   
462
463   return 0;
464 }
465
466 //____________________________________________________________________
467 TGeoShape*
468 AliFMDGeometryBuilder::HoneycombShape(Int_t id, Char_t ring,
469                                       double r1, double r2, 
470                                       double w, double t, double c)
471 {
472   // Make a honey comb shape from passed parameters.
473   // Parameters: 
474   //   id       Detector identifier (1,2, or 3)
475   //   ring     Ring identifier ('I' or 'O')
476   //   r1       Inner radius
477   //   r2       Outer radius
478   //   w        width 
479   //   t        Thickness of material 
480   //   c        Clearing from horizontal. 
481   // Return 
482   //   Pointer to newly allocated composite shape. 
483   TString      form  = Form("FMD%d%c_%%c_%%c", id, ring);
484   double       a1    = TMath::ATan2(c, r1) * 180  / TMath::Pi();
485
486   TString      fn    = Form(form.Data(),'F','1');
487   TString      bn    = Form(form.Data(),'B','1');
488   TString      cn    = Form(form.Data(),'C','O');
489   TString      in    = Form(form.Data(),'R','I');
490   TString      on    = Form(form.Data(),'R','O');
491   TString      en    = Form(form.Data(),'E','X');
492   double       y     = c;
493   double       x     = r1 * TMath::Cos(TMath::Pi()*a1/180);
494   new TGeoTubeSeg(fn.Data(),r1,r2,t/2,0,180);
495   new TGeoTubeSeg(bn.Data(),r1,r2,t/2,0,180);
496   new TGeoBBox(cn.Data(),(r2-r1)/2,t/2,w/2);
497   new TGeoTubeSeg(in.Data(),r1,r1+t,w/2,0,180);
498   new TGeoTubeSeg(on.Data(),r2-t,r2,w/2,0,180);
499   new TGeoBBox(en.Data(),r2+.005,c/2+.005,w/2+.005);
500     
501   TString          ftn = Form(form.Data(),'F','T');
502   TString          btn = Form(form.Data(),'F','B');
503   TString          ltn = Form(form.Data(),'C','L');
504   TString          rtn = Form(form.Data(),'C','R');
505   TString          etn = Form(form.Data(),'E','X');
506   (new TGeoTranslation(ftn.Data(),0,0,+w/2-t/2))->RegisterYourself();
507   (new TGeoTranslation(btn.Data(),0,0,-w/2+t/2))->RegisterYourself();
508   (new TGeoTranslation(ltn.Data(),-(x+(r2-r1)/2), y+t/2,0))->RegisterYourself();
509   (new TGeoTranslation(rtn.Data(),(x+(r2-r1)/2), y+t/2,0))->RegisterYourself();
510   (new TGeoTranslation(etn.Data(),0, c/2,0))->RegisterYourself();
511   
512   TString comp(Form("(%s:%s+%s:%s+%s+%s+%s:%s+%s:%s)-%s:%s", 
513                     fn.Data(),ftn.Data(),
514                     bn.Data(),btn.Data(),
515                     in.Data(),on.Data(),
516                     cn.Data(),ltn.Data(),
517                     cn.Data(),rtn.Data(),
518                     en.Data(),etn.Data()));
519   TGeoCompositeShape* shape = new TGeoCompositeShape(comp.Data());
520   shape->SetName(Form(fgkHCName,id,ring));
521   shape->SetTitle(Form("FMD%d%c Honeycomb shape", id, ring));
522   return shape;
523 }
524
525
526 //____________________________________________________________________
527 TGeoVolume*
528 AliFMDGeometryBuilder::DetectorGeometry(AliFMDDetector* d, 
529                                         TGeoVolume* topMother, 
530                                         TGeoVolume* botMother, 
531                                         Double_t    zMother, 
532                                         TGeoVolume* innerTop, 
533                                         TGeoVolume* innerBot, 
534                                         TGeoVolume* outerTop, 
535                                         TGeoVolume* outerBot) 
536 {
537   // Common stuff for setting up the FMD1, FMD2, and FMD3 geometries.
538   // This includes putting the Honeycomb support plates and the rings
539   // into the mother volumes.   
540   // 
541   // Parameeters:
542   //    d         The detector geometry to use 
543   //    mother    The mother volume of the detector 
544   //    zmother   The midpoint in global coordinates of detector vol.
545   //    inner     Pointer to inner ring volume 
546   //    outer     Pointer to outer ring volume
547   //
548   // Returns:
549   //    Pointer to mother (detector volume) 
550   // 
551   if (!d) return 0;
552   // Loop over the defined rings 
553   for (int i = 0; i < 2; i++) {
554     AliFMDRing* r     = 0;
555     Double_t    lowr  = 0;
556     Double_t    highr = 0;
557     Double_t    rz    = 0;
558     TGeoVolume* tvol  = 0;
559     TGeoVolume* bvol  = 0;
560     switch (i) {
561     case 0: 
562       r      = d->GetInner();
563       lowr   = d->GetInnerHoneyLowR();
564       highr  = d->GetInnerHoneyHighR();
565       rz     = d->GetInnerZ();
566       tvol   = innerTop;
567       bvol   = innerBot;
568       break;
569     case 1: 
570       r      = d->GetOuter();
571       lowr   = d->GetOuterHoneyLowR();
572       highr  = d->GetOuterHoneyHighR();
573       rz     = d->GetOuterZ();
574       tvol   = outerTop;
575       bvol   = outerBot;
576       break;
577     }
578     if (!r) continue;
579     Char_t   c       = r->GetId();
580     Int_t    id      = d->GetId();
581     Double_t hcThick = r->GetHoneycombThickness();
582     Double_t alThick = r->GetAlThickness();
583     Double_t z       = TMath::Abs(rz - zMother);
584
585     // Place ring in mother volume
586     // TGeoMatrix*matrix=new TGeoTranslation(Form("FMD%d%c trans",id,c),0,0,0);
587     AliFMDDebug(1, ("Placing volumes %s and %s in %s and %s at z=%f", 
588                      tvol->GetName(), bvol->GetName(), 
589                      topMother->GetName(), botMother->GetName(), z));
590     topMother->AddNode(tvol, Int_t(c), new TGeoTranslation(0,0,z));
591     botMother->AddNode(bvol, Int_t(c), new TGeoTranslation(0,0,z));
592
593     // Honeycomp 
594     TGeoShape*   hcSha = HoneycombShape(id, c, lowr, highr, hcThick, alThick);
595     TGeoVolume*  hcVol = new TGeoVolume(Form(fgkHCName,id,c),hcSha,fAl);
596     hcVol->SetTitle(Form("FMD%d%c honeycomb shell", id, c));
597     
598     z += (r->GetSiThickness() + 
599           r->GetSpacing() + 
600           r->GetPrintboardThickness() + 
601           r->GetCopperThickness() + 
602           r->GetChipThickness() + 
603           r->GetModuleSpacing() +
604           r->GetLegLength() + 
605           r->GetHoneycombThickness() + 
606           r->GetFMDDPrintboardThickness() - 
607           hcThick / 2); 
608
609     AliFMDDebug(15, ("Placing a copy of %s in %s and %s at z=%f", 
610                       hcVol->GetName(), topMother->GetName(), 
611                       botMother->GetName(), z));
612     // Add to top 
613     topMother->AddNode(hcVol, 0, new TGeoTranslation(0, 0, z));
614
615     // Add to bottom
616     TGeoMatrix*   bhcMatrix = new TGeoCombiTrans(0,0,z,0);
617     bhcMatrix->SetName(Form("FMD%d%c_honeycomp", id, c));
618     bhcMatrix->SetTitle(Form("FMD%d%c honeycomp", id, c));
619     bhcMatrix->RotateZ(180);
620     botMother->AddNode(hcVol, 1, bhcMatrix);
621   }
622   return 0;
623 }
624
625 //____________________________________________________________________
626 TGeoVolume*
627 AliFMDGeometryBuilder::FMD1Geometry(AliFMD1* fmd1, 
628                                     TGeoVolume* innerTop, 
629                                     TGeoVolume* innerBot) 
630 {
631   // Setup the FMD1 geometry.  The FMD1 only has one ring, and no
632   // special support as it is at the momement. 
633   // 
634   // See also AliFMDGeometryBuilder::DetectorGeometry 
635   // 
636   if (!fmd1 || !innerTop || !innerBot) return 0;
637   AliFMDRing* r             = fmd1->GetInner();
638   Double_t    z             = fmd1->GetInnerZ();  
639   Double_t    disce         = 2;
640   Double_t    backlr        = fmd1->GetInnerHoneyHighR();
641   Double_t    backhr        = fmd1->GetInnerHoneyHighR()+5;
642   Double_t    backth        = 0.2;
643   Double_t    toplr         = r->GetLowR();
644   Double_t    tophr         = fmd1->GetInnerHoneyHighR()+disce;
645   Double_t    wallbh        = (r->GetFullDepth() + disce);
646   Double_t    wallth        = wallbh+0.1;
647   
648   TGeoVolume* fmd1TopVolume = new TGeoVolumeAssembly(Form(fgkFMDName, 
649                                                           fmd1->GetId(), 'T'));
650   fmd1TopVolume->SetTitle("FMD1 top half");
651   TGeoVolume* fmd1BotVolume = new TGeoVolumeAssembly(Form(fgkFMDName, 
652                                                           fmd1->GetId(), 'B'));
653   fmd1BotVolume->SetTitle("FMD1 bottom half");
654   
655   // Basic detector geometry 
656   DetectorGeometry(fmd1, fmd1TopVolume, fmd1BotVolume, z, 
657                    innerTop, innerBot, 0, 0);
658
659
660   // Back
661   TGeoTubeSeg* backShape  = new TGeoTubeSeg(backlr, backhr, backth / 2, 0, 180);
662   TGeoTubeSeg* wallbShape = new TGeoTubeSeg(backlr, backlr + backth, 
663                                             wallbh/2, 0, 180);
664   TGeoTubeSeg* topShape   = new TGeoTubeSeg(toplr, tophr, backth / 2, 0, 180);
665   TGeoTubeSeg* walltShape = new TGeoTubeSeg(tophr, tophr + backth, 
666                                             wallth/2, 0, 180);
667   TGeoVolume*  backVolume = new TGeoVolume(Form(fgkBackName, fmd1->GetId()), 
668                                            backShape, fC);
669   TGeoVolume*  wallbVolume= new TGeoVolume(Form(fgkFlangeName, fmd1->GetId()), 
670                                            wallbShape, fC);
671   TGeoVolume*  topVolume  = new TGeoVolume(Form(fgkTopName, fmd1->GetId()), 
672                                            topShape, fC);
673   TGeoVolume*  walltVolume= new TGeoVolume(Form(fgkBeamName, fmd1->GetId()), 
674                                            walltShape, fC);
675   backShape->SetName(Form(fgkBackName, fmd1->GetId()));
676   wallbShape->SetName(Form(fgkFlangeName, fmd1->GetId()));
677   topShape->SetName(Form(fgkTopName, fmd1->GetId()));
678   walltShape->SetName(Form(fgkBeamName, fmd1->GetId()));
679   backShape->SetTitle("FMD1 back saucer rim");
680   wallbShape->SetTitle("FMD1 back saucer wall");
681   topShape->SetTitle("FMD1 top lid");
682   walltShape->SetTitle("FMD1 top lid wall");
683   backVolume->SetFillColor(kGray);
684   topVolume->SetFillColor(kGray);
685   wallbVolume->SetFillColor(kGray);
686   walltVolume->SetFillColor(kGray);
687   backVolume->SetTitle("FMD1 back saucer rim");
688   wallbVolume->SetTitle("FMD1 back saucer wall");
689   topVolume->SetTitle("FMD1 top lid");
690   walltVolume->SetTitle("FMD1 top lid wall");
691   
692   // Place volumes
693   Double_t zb = TMath::Abs(fmd1->GetInnerZ() - z);
694   Double_t zi = zb;
695   Int_t    n  = 2;
696   
697   // Place top cover
698   zi -= disce / 2 + backth / 2;
699   zb =  zi;
700   for (Int_t i = 0; i  < 2; i++) {
701     TGeoVolume*   mother = (i == 0 ? fmd1TopVolume : fmd1BotVolume);
702     Double_t      phi    = 360. / n * i;
703     TGeoRotation* rot    = new TGeoRotation(Form("FMD1 top rotation %d",i));
704     rot->RotateZ(phi);
705     TGeoMatrix* matrix   = new TGeoCombiTrans(Form("FMD1 top wall trans %d", i),
706                                               0, 0, zi, rot);
707     mother->AddNode(topVolume, i, matrix);    
708   }
709   // Place outer wall
710   zi += wallth / 2 + backth / 2;
711   for (Int_t i = 0; i  < 2; i++) {
712     TGeoVolume*   mother = (i == 0 ? fmd1TopVolume : fmd1BotVolume);
713     Double_t      phi    = 360. / n * i;
714     TGeoRotation* rot    = new TGeoRotation(Form("FMD1 outer wall rotation %d",
715                                                  i));
716     rot->RotateZ(phi);
717     TGeoMatrix* matrix   = new TGeoCombiTrans(Form("FMD1 outer wall trans %d",
718                                                    i), 0, 0, zi, rot);
719     mother->AddNode(walltVolume, i, matrix);    
720   }
721   // Place back
722   zi += wallth / 2 + backth / 2; // + disce / 2;
723   for (Int_t i = 0; i  < 2; i++) {
724     TGeoVolume*   mother = (i == 0 ? fmd1TopVolume : fmd1BotVolume);
725     Double_t      phi    = 360. / n * i;
726     TGeoRotation* rot    = new TGeoRotation(Form("FMD1 back rotation %d", i));
727     rot->RotateZ(phi);
728     TGeoMatrix* matrix   = new TGeoCombiTrans(Form("FMD1 back trans %d", i),
729                                              0, 0, zi, rot);
730     mother->AddNode(backVolume, i, matrix);    
731   }
732   // Place inner wall
733   zi -= wallbh / 2 + backth / 2; // + disce / 2;
734   for (Int_t i = 0; i  < 2; i++) {
735     TGeoVolume*   mother = (i == 0 ? fmd1TopVolume : fmd1BotVolume);
736     Double_t      phi    = 360. / n * i;
737     TGeoRotation* rot    = new TGeoRotation(Form("FMD1 inner wall rotation %d",
738                                                  i)); 
739     rot->RotateZ(phi);
740     TGeoMatrix*   matrix = new TGeoCombiTrans(Form("FMD1 inner wall trans %d", 
741                                                    i), 0, 0, zi, rot);
742     mother->AddNode(wallbVolume, i, matrix);    
743   }
744
745
746   // Must add this after filling the assembly.
747   TGeoVolume* top    = gGeoManager->GetVolume("ALIC");
748   // TGeoMatrix* matrix = new TGeoTranslation("FMD1 trans", 0, 0, z);
749   TGeoRotation* rot = new TGeoRotation("FMD1 rotatation");
750   rot->RotateZ(-90);
751   TGeoMatrix* matrix = new TGeoCombiTrans("FMD1 trans", 0, 0, z, rot);
752   AliFMDDebug(5, ("Placing volumes %s and %s in ALIC at z=%f", 
753                    fmd1TopVolume->GetName(), fmd1BotVolume->GetName(), z));
754   top->AddNode(fmd1TopVolume, fmd1->GetId(), matrix);
755   top->AddNode(fmd1BotVolume, fmd1->GetId(), matrix);
756
757   return 0;
758 }
759
760 //____________________________________________________________________
761 TGeoVolume*
762 AliFMDGeometryBuilder::FMD2Geometry(AliFMD2* fmd2, 
763                                     TGeoVolume* innerTop, 
764                                     TGeoVolume* innerBot, 
765                                     TGeoVolume* outerTop,
766                                     TGeoVolume* outerBot) 
767 {
768   // Setup the FMD2 geometry.  The FMD2 has no
769   // special support as it is at the momement. 
770   // 
771   // See also AliFMDGeometryBuilder::DetectorGeometry 
772   // 
773   if (!fmd2 || !innerTop || !innerBot || !outerTop || !outerBot) return 0;
774   AliFMDRing* r             = fmd2->GetOuter();
775   Double_t    z             = fmd2->GetOuterZ();  
776   Double_t    framelr       = fmd2->GetOuterHoneyHighR()+0.5;
777   Double_t    framehr       = fmd2->GetOuterHoneyHighR()+1.8;
778   Double_t    framelz       = -.5;
779   Double_t    framehz       = (fmd2->GetInnerZ()-z) + r->GetFullDepth() + .5;
780   Double_t    framel        = framehz - framelz;
781   Double_t    coverlr       = fmd2->GetInner()->GetLowR()+1;
782   Double_t    backth        = 0.05;
783
784   TGeoVolume* fmd2TopVolume = new TGeoVolumeAssembly(Form(fgkFMDName, 
785                                                           fmd2->GetId(), 'T'));
786   TGeoVolume* fmd2BotVolume = new TGeoVolumeAssembly(Form(fgkFMDName, 
787                                                           fmd2->GetId(), 'B'));
788   fmd2TopVolume->SetTitle("FMD2 top half");
789   fmd2BotVolume->SetTitle("FMD2 bottom half");
790   
791   DetectorGeometry(fmd2, fmd2TopVolume, fmd2BotVolume, z, 
792                    innerTop, innerBot, outerTop, outerBot);
793
794   TGeoShape*  cylinderShape   = new TGeoTubeSeg(framelr,framehr,framel/2,0,180);
795   TGeoVolume* cylinderVolume  = new TGeoVolume(Form(fgkBackName, fmd2->GetId()),
796                                                cylinderShape, fC);
797   TGeoShape*  coverShape      = new TGeoTubeSeg(coverlr,framehr,backth/2,0,180);
798   TGeoVolume* coverVolume     = new TGeoVolume(Form(fgkTopName, fmd2->GetId()), 
799                                                coverShape, fC);
800   cylinderShape->SetName(Form(fgkBackName, fmd2->GetId()));
801   cylinderShape->SetTitle("FMD2 cylinder");
802   cylinderVolume->SetTitle("FMD2 cylinder");
803   cylinderVolume->SetTransparency(63);
804   coverShape->SetName(Form(fgkTopName, fmd2->GetId()));
805   coverShape->SetTitle("FMD2 cover");
806   coverVolume->SetTitle("FMD2 cover");
807   coverVolume->SetTransparency(63);
808   
809   for (Int_t i = 0; i  < 2; i++) {
810     TGeoVolume*   mother = (i == 0 ? fmd2TopVolume : fmd2BotVolume);
811     
812     Double_t      phi    = 360. / 2 * i;
813     TGeoRotation* rot    = new TGeoRotation(Form("FMD2 support rot %d",i)); 
814     rot->RotateZ(phi);
815     TGeoMatrix*   matrix = new TGeoCombiTrans(Form("FMD2 cyl trans %d", i),
816                                               0, 0, framelz+framel/2, rot);
817     mother->AddNode(cylinderVolume, i, matrix);    
818     matrix               = new TGeoCombiTrans(Form("FMD2 fcov trans %d", i),
819                                               0, 0, framelz-backth/2, rot);
820     mother->AddNode(coverVolume, 2*i+0, matrix);    
821     matrix               = new TGeoCombiTrans(Form("FMD2 bcov trans %d", i),
822                                               0, 0, framelz+framel+backth/2, 
823                                               rot);
824     mother->AddNode(coverVolume, 2*i+1, matrix);    
825   }
826
827
828   Double_t    f1l           = 10;
829   Double_t    f1w           = 6;
830   Double_t    f1d           = 1;
831   
832   TGeoBBox*   flange1Shape  = new TGeoBBox(f1l/2, f1w/2, f1d/2);
833   TGeoVolume* flange1Volume = new TGeoVolume(Form(fgkFlangeName, fmd2->GetId()),
834                                              flange1Shape, fAl);
835   TGeoBBox*   flange2Shape  = new TGeoBBox(f1w/2, f1d/2, (framel+backth)/2);
836   TGeoVolume* flange2Volume = new TGeoVolume(Form("F%dSG", fmd2->GetId()),
837                                              flange2Shape, fAl);
838   flange1Shape->SetName(Form(fgkFlangeName, fmd2->GetId()));
839   flange1Shape->SetTitle("FMD2 vertical flange");
840   flange1Volume->SetTitle("FMD2 vertical flange");
841   flange2Shape->SetName(Form("F%dSG", fmd2->GetId()));
842   flange2Shape->SetTitle("FMD2 horizontal flange");
843   flange2Volume->SetTitle("FMD2 horizontal flange ");
844   
845   flange1Volume->SetTransparency(42);
846   for (Int_t i = 0; i  < 4; i++) {
847     TGeoVolume*   mother = (i < 2 ? fmd2TopVolume : fmd2BotVolume);
848     
849     Double_t      phi    = 360. / 4 * i - 45;
850     Double_t      rphi   = TMath::Pi()*phi/180;
851     Double_t      x      = (framelr + f1l/2) * TMath::Sin(rphi);
852     Double_t      y      = (framelr + f1l/2) * TMath::Cos(rphi);
853     TGeoRotation* rot    = new TGeoRotation(Form("FMD2 support rot %d",i)); 
854     rot->RotateZ(phi);
855     TGeoMatrix*   matrix = new TGeoCombiTrans(Form("FMD2 flange 1 trans %d", i),
856                                               x,y, framelz-backth-f1d/2, rot);
857     mother->AddNode(flange1Volume, 2*i+0, matrix);    
858     matrix               = new TGeoCombiTrans(Form("FMD2 flange 2 trans %d", i),
859                                               x,y,framelz+framel+backth+f1d/2, 
860                                               rot);
861     mother->AddNode(flange1Volume, 2*i+1, matrix);    
862     Double_t x1 = x - (f1w-f1d) / 2 * TMath::Cos(rphi); 
863     Double_t y1 = y + (f1w-f1d) / 2 * TMath::Sin(rphi);
864     matrix               = new TGeoCombiTrans(Form("FMD2 flange 3 trans %d", i),
865                                               x1,y1,framelz+framel/2, rot);
866     mother->AddNode(flange2Volume, 2*i+0, matrix);    
867     Double_t x2 = x + (f1w-f1d) / 2 * TMath::Cos(rphi); 
868     Double_t y2 = y - (f1w-f1d) / 2 * TMath::Sin(rphi);
869     matrix               = new TGeoCombiTrans(Form("FMD2 flange 4 trans %d", i),
870                                               x2,y2,framelz+framel/2, rot);
871     mother->AddNode(flange2Volume, 2*i+1, matrix);    
872   }
873   
874   // Must be done after filling the assemblies 
875   TGeoVolume* top = gGeoManager->GetVolume("ALIC");
876   TGeoMatrix* matrix = new TGeoTranslation("FMD2 trans", 0, 0, z);
877   AliFMDDebug(5, ("Placing volumes %s and %s in ALIC at z=%f", 
878                    fmd2TopVolume->GetName(), fmd2BotVolume->GetName(), z));
879   top->AddNode(fmd2TopVolume, fmd2->GetId(), matrix);
880   top->AddNode(fmd2BotVolume, fmd2->GetId(), matrix);
881
882
883   return 0;
884 }
885   
886 #if 1
887 //____________________________________________________________________
888 TGeoVolume*
889 AliFMDGeometryBuilder::FMD3Geometry(AliFMD3* fmd3, 
890                                     TGeoVolume* innerTop, 
891                                     TGeoVolume* innerBot, 
892                                     TGeoVolume* outerTop,
893                                     TGeoVolume* outerBot) 
894 {
895   // Setup the FMD3 geometry.  The FMD2 has a rather elaborate support
896   // structure, as the support will also support the vacuum
897   // beam-pipe. 
898   // 
899   // See also AliFMDGeometryBuilder::DetectorGeometry 
900   // 
901   if (!fmd3 || !innerTop || !innerBot || !outerTop || !outerBot) return 0;
902
903   //__________________________________________________________________
904   // Basic detector set-up.
905   TGeoVolume* fmd3TopVolume = new TGeoVolumeAssembly(Form(fgkFMDName, 
906                                                           fmd3->GetId(), 'T'));
907   TGeoVolume* fmd3BotVolume = new TGeoVolumeAssembly(Form(fgkFMDName, 
908                                                           fmd3->GetId(), 'B'));
909   fmd3TopVolume->SetTitle("FMD3 top half");
910   fmd3BotVolume->SetTitle("FMD3 bottom half");
911   DetectorGeometry(fmd3, fmd3TopVolume, fmd3BotVolume, fmd3->GetInnerZ(), 
912                    innerTop, innerBot, outerTop, outerBot);
913
914   //__________________________________________________________________
915   // Mother for all support material
916   TGeoVolumeAssembly* support = new TGeoVolumeAssembly("F3SU");
917   support->SetTitle("FMD3 support");
918
919   //__________________________________________________________________
920   // Base of cone
921   const TObjArray& radii    = fmd3->ConeRadii();
922   Int_t            nRadii   = radii.GetEntriesFast();
923   TGeoPcon*        coneBase = new TGeoPcon("FMD3_cone_base", 0., 180., nRadii);
924   TVector3*        r5       = 0;
925   TVector3*        r4       = 0;
926   for (Int_t i = 0; i < nRadii; i++) { 
927     TVector3* v = static_cast<TVector3*>(radii.At(i));
928     coneBase->DefineSection(i, v->X(), v->Y(), v->Z());
929     if      (i == 5) r5 = v;
930     else if (i == 4) r4 = v;
931   }
932   TString          coneComb("(FMD3_cone_base");
933
934   //__________________________________________________________________
935   // Flanges 
936   double    flangeDepth    = fmd3->GetFlangeDepth() / 2;
937   double    flangeLength   = fmd3->GetFlangeLength() / 2;
938   double    flangeWidth    = fmd3->GetFlangeWidth() / 2;
939   new TGeoBBox("FMD3_flange_base", flangeLength, flangeWidth, flangeDepth);
940
941   // Fiducial holes 
942   const TObjArray& fiducialHoles  = fmd3->FiducialHoles();
943   double           fiducialRadius = fmd3->GetFiducialRadius();
944   TGeoTube*        fiducialShape  = new TGeoTube("FMD3_fiducial_hole", 0, 
945                                                  fiducialRadius, 
946                                                  flangeDepth+.1);
947   Int_t            nFiducialHoles = fiducialHoles.GetEntriesFast();
948   double           flangeAngle    = TMath::Pi() / 4;
949   double           flangeX        = r5->Y()+flangeLength;
950   TVector2         flangeC(flangeX * TMath::Cos(flangeAngle), 
951                            flangeX * TMath::Sin(flangeAngle));
952   TString          flangeComb("FMD3_flange_base-(");
953 #if 0// For debugging geometry 
954   TGeoVolume* fiducialVolume = new TGeoVolume("FMD3_fiducial", fiducialShape);
955   fiducialVolume->SetLineColor(kGreen);
956 #else
957   (void*)fiducialShape;
958 #endif
959   for (Int_t i = 0; i < nFiducialHoles; i++) { 
960     TVector2&        v  =  *(static_cast<TVector2*>(fiducialHoles.At(i)));
961     v                   -= flangeC;
962     TVector2         r  =  v.Rotate(-flangeAngle);
963     TGeoTranslation* t1 =  new TGeoTranslation(r.X(),  r.Y(), 0);
964     TGeoTranslation* t2 =  new TGeoTranslation(r.X(), -r.Y(), 0);
965     t1->SetName(Form("FMD3_fiducial_hole_rot%d", 2*i+0));
966     t2->SetName(Form("FMD3_fiducial_hole_rot%d", 2*i+1));
967     t1->RegisterYourself();
968     t2->RegisterYourself();
969     flangeComb.Append(Form("FMD3_fiducial_hole:FMD3_fiducial_hole_rot%d+"
970                            "FMD3_fiducial_hole:FMD3_fiducial_hole_rot%d%c",
971                            2*i+0, 2*i+1, (i == nFiducialHoles-1 ? ')' : '+')));
972 #if 1 // For debugging geometry 
973     // support->AddNode(fiducialVolume, 2*i+0, t1);
974     // support->AddNode(fiducialVolume, 2*i+1, t2);
975 #endif
976   }
977   
978   // Final flange shape, and at to full shape 
979   TGeoCompositeShape* flangeShape = new TGeoCompositeShape(flangeComb.Data());
980   flangeShape->SetName("FMD3_flange");
981   for (Int_t i = 0; i < 2; i++) { 
982     TGeoRotation* rot = new TGeoRotation();
983     rot->RotateZ((i+.5)*90);
984     TVector2 v(flangeX, 0);
985     TVector2 w = v.Rotate((i+.5) * 2 * flangeAngle);
986     TGeoCombiTrans* trans = new TGeoCombiTrans(w.X(),w.Y(),
987                                                r4->X()+flangeDepth, rot);
988     trans->SetName(Form("FMD3_flange_matrix%d", i));
989     trans->RegisterYourself();
990     coneComb.Append(Form("+FMD3_flange:FMD3_flange_matrix%d", i));
991   }
992   coneComb.Append(")-(");
993   
994   //__________________________________________________________________
995   // Holes 
996   Double_t holeL  = (fmd3->GetHoleLength()-1)/2;
997   Double_t holeD  = fmd3->GetHoleDepth()/2;
998   Double_t holeLW = fmd3->GetHoleLowWidth()/2;
999   Double_t holeHW = fmd3->GetHoleHighWidth()/2;
1000   Double_t holeA  = fmd3->GetConeOuterAngle() - 1 * TMath::Pi() / 180;
1001   Double_t holeZ  = (fmd3->GetHoleOffset() 
1002                      + holeL * TMath::Cos(holeA)
1003                      - holeD * TMath::Sin(holeA));
1004   Double_t holeX  = (fmd3->ConeR(-holeZ + fmd3->GetInnerZ() + fmd3->GetNoseZ()) 
1005                      - holeD * TMath::Cos(holeA));
1006   Double_t plateZ  = (fmd3->GetHoleOffset() 
1007                      + holeL * TMath::Cos(holeA)
1008                      - 0.033 * TMath::Sin(holeA));
1009   Double_t plateX = (fmd3->ConeR(-plateZ + fmd3->GetInnerZ()+fmd3->GetNoseZ()) 
1010                      - 0.033 * TMath::Cos(holeA));
1011   TGeoTrd1* holeShape = new TGeoTrd1("FMD3_cone_hole", 
1012                                      holeLW, holeHW, holeD, holeL);
1013   TGeoTrd1* plateShape = new TGeoTrd1("FMD3_cooling_plate", 
1014                                       holeLW, holeHW, .033, holeL);
1015   TGeoRotation* holeRot = new TGeoRotation();
1016   holeRot->SetName("FMD3_cone_hole_rotation");
1017   holeRot->RotateZ(90);
1018   holeRot->RotateY(holeA*180/TMath::Pi());
1019   TGeoCombiTrans* holeBaseTrans = new TGeoCombiTrans(holeX, 0, holeZ, holeRot);
1020   holeBaseTrans->SetName("FMD3_cone_hole_base_matrix");
1021   TGeoCombiTrans* plateBaseTrans = new TGeoCombiTrans(plateX, 0,plateZ,holeRot);
1022   (void*)holeShape;
1023   TGeoVolume* plateVolume = new TGeoVolume("F3CO", plateShape, fAl);
1024   plateShape->SetTitle("FMD3 cooling plate");
1025   plateVolume->SetTitle("FMD3 cooling plate");
1026   for (Int_t i = 0; i < 4; i++) { 
1027     Double_t        ang   = 360. / 8 * (i + .5);
1028     TGeoCombiTrans* trans = new TGeoCombiTrans(*holeBaseTrans);
1029     trans->RotateZ(ang);
1030     trans->SetName(Form("FMD3_cone_hole_matrix%d", i));
1031     trans->RegisterYourself();
1032     trans = new TGeoCombiTrans(*plateBaseTrans);
1033     trans->RotateZ(ang);
1034     trans->SetName(Form("FMD3_cooling_plate_matrix%d", i));
1035     coneComb.Append(Form("FMD3_cone_hole:FMD3_cone_hole_matrix%d+", i));
1036     support->AddNode(plateVolume, i, trans);
1037   }
1038   
1039   //__________________________________________________________________
1040   // Bolts
1041   Double_t boltRadius = fmd3->GetBoltRadius();
1042   Double_t boltLength = fmd3->GetBoltLength() / 2;
1043   Double_t boltZ1     = fmd3->GetInnerZ()+fmd3->GetNoseZ()-10;
1044   Double_t boltZ2     = fmd3->GetInnerZ()+fmd3->GetNoseZ()-20;
1045   Double_t boltXE     = 2*boltLength*TMath::Cos(fmd3->GetConeOuterAngle());
1046   Double_t boltX1     = (fmd3->ConeR(boltZ1) - boltXE);
1047   Double_t boltX2     = (fmd3->ConeR(boltZ2) - boltXE);
1048   
1049   new TGeoTube("FMD3_bolt_hole", 0, boltRadius, boltLength+.2);
1050   TGeoTube* boltShape = new TGeoTube("FMD3_bolt", 0, boltRadius, boltLength);
1051   TGeoRotation* boltRot = new TGeoRotation();
1052   boltRot->RotateY(-fmd3->GetConeOuterAngle()*180/TMath::Pi());
1053   TGeoCombiTrans* boltTrans1 = new TGeoCombiTrans(boltX1, 0, 10, boltRot);
1054   TGeoCombiTrans* boltTrans2 = new TGeoCombiTrans(boltX2, 0, 20, boltRot);
1055   TGeoCombiTrans* boltTrans3 = new TGeoCombiTrans(*boltTrans1);
1056   TGeoCombiTrans* boltTrans4 = new TGeoCombiTrans(*boltTrans2);
1057   boltTrans3->RotateZ(180);
1058   boltTrans4->RotateZ(180);
1059   boltTrans1->SetName("FMD3_bolt_matrix1");
1060   boltTrans2->SetName("FMD3_bolt_matrix2");
1061   boltTrans3->SetName("FMD3_bolt_matrix3");
1062   boltTrans4->SetName("FMD3_bolt_matrix4");
1063   boltTrans1->RegisterYourself();
1064   boltTrans2->RegisterYourself();
1065   boltTrans3->RegisterYourself();
1066   boltTrans4->RegisterYourself();
1067   coneComb.Append("FMD3_bolt_hole:FMD3_bolt_matrix1"
1068                   "+FMD3_bolt_hole:FMD3_bolt_matrix2"
1069                   "+FMD3_bolt_hole:FMD3_bolt_matrix3"
1070                   "+FMD3_bolt_hole:FMD3_bolt_matrix4)");
1071   TGeoVolume*     boltVolume = new TGeoVolume("F3SB", boltShape, fSteel);
1072   support->AddNode(boltVolume, 1, boltTrans1);
1073   support->AddNode(boltVolume, 2, boltTrans2);
1074   boltShape->SetTitle("FMD3 steering bolt");
1075   boltVolume->SetTitle("FMD3 steering bolt");
1076   
1077   //__________________________________________________________________
1078   // Final cone
1079   TGeoCompositeShape* coneShape = new TGeoCompositeShape(coneComb.Data());
1080   coneShape->SetName("FMD3_cone");
1081   coneShape->SetTitle("FMD3 cone");
1082   TGeoVolume*  coneVolume = new TGeoVolume("F3SC", coneShape, fC);
1083   coneVolume->SetLineColor(kRed);
1084   support->AddNode(coneVolume, 0, new TGeoTranslation(0, 0, 0));
1085
1086   //__________________________________________________________________
1087   // Tension boxes. 
1088   new TGeoBBox("FMD3_tension_outer", .5, 3, 5);
1089   new TGeoBBox("FMD3_tension_inner", .51, 2.7, 4.3);
1090   TString tensionExpr("FMD3_tension_outer-FMD3_tension_inner");
1091   TGeoCompositeShape* tensionShape = new TGeoCompositeShape(tensionExpr.Data());
1092   tensionShape->SetName("FMD3_tension_box");
1093   tensionShape->SetTitle("FMD3 tension box");
1094   TGeoVolume*      tensionFrame = new TGeoVolume("FMD3_tension_frame",
1095                                                   tensionShape, fAl);
1096   TGeoTube*        springShape  = new TGeoTube("FMD3_tension_spring", 
1097                                                0, .3, 4.3/2);
1098   TGeoVolume*      springVolume = new TGeoVolume("FMD3_tension_spring", 
1099                                                  springShape, fSteel);
1100   TGeoVolume*      tensionBox   = new TGeoVolume("FMD3_tension_box", 
1101                                                  tensionShape, fAir);
1102   tensionBox->AddNode(tensionFrame, 0);
1103   tensionBox->AddNode(springVolume, 0, new TGeoTranslation(0,0,4.3/2));
1104   
1105   Double_t         tensionD     = 5*TMath::Cos(fmd3->GetConeOuterAngle());
1106   Double_t         tensionZ     = r4->Z()-2*tensionD;
1107   Double_t         tensionX     = (fmd3->ConeR(fmd3->GetInnerZ()
1108                                                +fmd3->GetNoseZ() 
1109                                                -tensionZ) + 
1110                                    .5 * TMath::Sin(fmd3->GetConeOuterAngle())); 
1111   TGeoRotation*    tensionRot   = new TGeoRotation();
1112   tensionRot->RotateY(180/TMath::Pi()*fmd3->GetConeOuterAngle());
1113   TGeoCombiTrans*  tensionBase  = new TGeoCombiTrans(tensionX, 0, tensionZ, 
1114                                                      tensionRot);
1115   
1116   Double_t         wireR1       = fmd3->ConeR(fmd3->GetInnerZ()
1117                                               +fmd3->GetNoseZ());
1118   Double_t         wireR2       = fmd3->ConeR(fmd3->GetInnerZ()
1119                                               +fmd3->GetNoseZ()-
1120                                               tensionZ+tensionD);
1121   Double_t         wireL        = TMath::Sqrt(TMath::Power(wireR1-wireR2,2)+ 
1122                                               TMath::Power(tensionZ-
1123                                                            tensionD,2));
1124   Double_t         wireAngle    = TMath::ATan2(wireR2-wireR1,tensionZ-tensionD);
1125   TGeoTube*        wireShape    = new TGeoTube("FMD3_wire", 0, .1, wireL/2);
1126   TGeoVolume*      wireVolume   = new TGeoVolume("FMD3_wire", wireShape,fSteel);
1127   TGeoRotation*    wireRot      = new TGeoRotation();
1128   wireRot->RotateY(180/TMath::Pi()*wireAngle);
1129   TGeoCombiTrans*  wireBase     = new TGeoCombiTrans((wireR2-wireR1)/2+wireR1,
1130                                                      0, (tensionZ-tensionD)/2,
1131                                                      wireRot);
1132   for (Int_t i = 0; i < 2; i++) { 
1133     Double_t        thisAngle = (i+.5) * 90;
1134     TGeoCombiTrans* thisTrans = new TGeoCombiTrans(*tensionBase);
1135     thisTrans->RotateZ(thisAngle);
1136     support->AddNode(tensionBox, i, thisTrans);
1137     thisTrans = new TGeoCombiTrans(*wireBase);
1138     thisTrans->RotateZ(thisAngle);
1139     support->AddNode(wireVolume, i, thisTrans);
1140   }
1141
1142   //__________________________________________________________________
1143   // Place support volumes in half-detector volumes 
1144   Double_t         z  = fmd3->GetInnerZ();
1145   TGeoTranslation* t1 = new TGeoTranslation(0, 0, -fmd3->GetNoseZ());
1146   fmd3TopVolume->AddNode(support, 1, t1);
1147   TGeoCombiTrans*  t2 = new TGeoCombiTrans(*t1);
1148   t2->RotateZ(180);
1149   fmd3BotVolume->AddNode(support, 2, t2);
1150
1151   TGeoRotation*   rot        = new TGeoRotation("FMD3 rotatation");
1152   rot->RotateY(180);
1153   TGeoVolume*     top        = gGeoManager->GetVolume("ALIC");
1154   TGeoMatrix* mmatrix        = new TGeoCombiTrans("FMD3 trans", 0, 0, z, rot);
1155   AliFMDDebug(5, ("Placing volumes %s and %s in ALIC at z=%f", 
1156                    fmd3TopVolume->GetName(), fmd3BotVolume->GetName(), z));
1157   top->AddNode(fmd3TopVolume, fmd3->GetId(), mmatrix);
1158   top->AddNode(fmd3BotVolume, fmd3->GetId(), mmatrix);
1159
1160   return 0;
1161 }
1162
1163 #else
1164 //____________________________________________________________________
1165 TGeoVolume*
1166 AliFMDGeometryBuilder::FMD3Geometry(AliFMD3* fmd3, 
1167                                     TGeoVolume* innerTop, 
1168                                     TGeoVolume* innerBot, 
1169                                     TGeoVolume* outerTop,
1170                                     TGeoVolume* outerBot) 
1171 {
1172   // Setup the FMD3 geometry.  The FMD2 has a rather elaborate support
1173   // structure, as the support will also support the vacuum
1174   // beam-pipe. 
1175   // 
1176   // See also AliFMDGeometryBuilder::DetectorGeometry 
1177   // 
1178   if (!fmd3 || !innerTop || !innerBot || !outerTop || !outerBot) return 0;
1179   Double_t nlen    = fmd3->GetNoseLength();
1180   Double_t nz      = fmd3->GetNoseZ();
1181   Double_t noser1  = fmd3->GetNoseLowR();
1182   Double_t noser2  = fmd3->GetNoseHighR();
1183   Double_t conet   = fmd3->GetBeamThickness();
1184   Double_t conel   = fmd3->GetConeLength();
1185   Double_t backl   = fmd3->GetBackLength();
1186   // Double_t backr1  = fmd3->GetBackLowR();
1187   Double_t backr2  = fmd3->GetBackHighR();
1188   Double_t zdist   = conel -  backl - nlen;
1189   Double_t tdist   = backr2 - noser2;
1190   // Double_t beaml   = TMath::Sqrt(zdist * zdist + tdist * tdist);
1191   Double_t theta   = -180. * TMath::ATan2(tdist, zdist) / TMath::Pi();
1192   Double_t flanger = fmd3->GetFlangeR();
1193   Double_t z       = fmd3->GetInnerZ(); // fmd3->GetZ();
1194
1195   TGeoVolume* fmd3TopVolume = new TGeoVolumeAssembly(Form(fgkFMDName, 
1196                                                           fmd3->GetId(), 'T'));
1197   TGeoVolume* fmd3BotVolume = new TGeoVolumeAssembly(Form(fgkFMDName, 
1198                                                           fmd3->GetId(), 'B'));
1199   fmd3TopVolume->SetTitle("FMD3 top half");
1200   fmd3BotVolume->SetTitle("FMD3 bottom half");
1201   
1202   
1203   DetectorGeometry(fmd3, fmd3TopVolume, fmd3BotVolume, z, 
1204                    innerTop, innerBot, outerTop, outerBot);
1205
1206   
1207   TGeoVolumeAssembly* support = new TGeoVolumeAssembly("F3SU");
1208   support->SetTitle("FMD3 support");
1209   
1210   // Cone shape 
1211   TGeoPcon*        coneBase = new TGeoPcon("FMD3 cone base", 0, 180, 6);
1212   const TObjArray& radii    = fmd3.ConeRadii();
1213   TVector3*        v1       = 0;
1214   TVector3*        v4       = 0;
1215   for (Int_t i = 0; i < radii.GetEntriesFast(); i++) { 
1216     TVector3* v = static_cast<TVector3*>(radii.At(i));
1217     coneBase->DefineSection(i,  v->X(), v->Y(), v->Z());
1218     if (i == 1) v1 = v;
1219     if (i == 4) v4 = v;
1220     
1221   }
1222   Double_t  holeL = TMath::Sqrt(TMath::Power(v4->Z()-v1->Z(),2) + 
1223                                 TMath::Power(v4->X()-v1->X(),2));
1224   
1225   TGeoTrd1*       coneHole  = new TGeoTrd1("F3SC_hole",2,8,conet*3,
1226                                            (conel-2-2)/2);
1227   
1228
1229
1230   // Nose volume 
1231   TGeoTubeSeg* noseShape  = new TGeoTubeSeg(noser1, noser2, nlen / 2, 0, 180);
1232   TGeoVolume*  noseVolume = new TGeoVolume(fgkNoseName, noseShape, fC);
1233   support->AddNode(noseVolume, 0, new TGeoTranslation(0, 0, nlen/2));
1234   noseShape->SetName(fgkNoseName);
1235   noseShape->SetTitle("FMD3 nose");
1236   noseVolume->SetTitle("FMD3 nose");
1237   
1238   // Steel bolts 
1239   TGeoTube*       boltShape  = new TGeoTube("F3SB", 0, 0.3, conet / 2);
1240   TGeoVolume*     boltVolume = new TGeoVolume("F3SB", boltShape, fSteel);
1241   Double_t        z1         = -10;
1242   Double_t        x1         = (fmd3->ConeR(nz+z1));
1243   TGeoRotation*   r1         = new TGeoRotation();
1244   r1->RotateY(theta);
1245   TGeoCombiTrans* t          = new TGeoCombiTrans("F3SB1",x1,0,-z1,r1);
1246   support->AddNode(boltVolume, 1, t);
1247   z1                         = -20;
1248   x1                         = (fmd3->ConeR(nz+z1));
1249   t                          = new TGeoCombiTrans("F3SB2",x1,0,-z1,r1);
1250   support->AddNode(boltVolume, 2, t);
1251   boltShape->SetTitle("FMD3 steering bolt");
1252   boltVolume->SetTitle("FMD3 steering bolt");
1253
1254   // Cooling plates
1255   TGeoTrd1*   plateShape  = new TGeoTrd1(2, 8, 0.1, (conel-2-2)/2-.1);
1256   TGeoVolume* plateVolume = new TGeoVolume("F3CO", plateShape, fAl);
1257   plateShape->SetName("F3C0");
1258   plateShape->SetTitle("FMD3 cooling plate");
1259   plateVolume->SetTitle("FMD3 cooling plate");
1260
1261   // Shape for carbon half-cone
1262   TGeoConeSeg*    innerCone = new TGeoConeSeg("F3SC_inner", conel/2,
1263                                               noser2-conet, noser2, 
1264                                               backr2-conet, backr2, 0., 180.);
1265   innerCone->SetTitle("FMD3 cone inner");
1266   TGeoTrd1*       coneHole  = new TGeoTrd1("F3SC_hole",2,8,conet*3,
1267                                            (conel-2-2)/2);
1268   coneHole->SetTitle("FMD3 cone hole");
1269   Double_t        holeAng   = TMath::ATan2(backr2 - noser2, conel);
1270   Double_t        holeX     = ((conel-2) / 2 * TMath::Sin(holeAng) +
1271                                conet     * TMath::Cos(holeAng) +
1272                                noser2);
1273   TGeoRotation*   holeRot   = new TGeoRotation();
1274   holeRot->SetName("FMD3 cone hole rotation");
1275   holeRot->RotateZ(90);
1276   holeRot->RotateY(holeAng*180./TMath::Pi());
1277   TGeoCombiTrans* holeTrans = new TGeoCombiTrans(holeX, 0, -2, holeRot);
1278   holeRot->SetName("FMD3 cone hole");
1279
1280   // Build-up the composite shape for the cone, and add cooling plates
1281   // at the same time. 
1282   TString coneExp("F3SC_inner-(");
1283   for (int i = 0; i < 4; i++) { 
1284     Double_t        thisAng   = 360. / 8 * (i + .5);
1285     TGeoCombiTrans* thisTrans = new TGeoCombiTrans(*holeTrans);
1286     thisTrans->RotateZ(thisAng);
1287     thisTrans->SetName(Form("F3SC_rot%d", i));
1288     thisTrans->RegisterYourself();
1289     coneExp.Append(Form("F3SC_hole:F3SC_rot%d+", i));
1290
1291     const Double_t* tt         = thisTrans->GetTranslation();
1292     Double_t        x          = tt[0]+1*TMath::Cos(thisAng*TMath::Pi()/180);
1293     Double_t        y          = tt[1]+1*TMath::Sin(thisAng*TMath::Pi()/180);
1294     TGeoCombiTrans* plateTrans = new TGeoCombiTrans(x,y,tt[2]-1+nlen+conel/2,
1295                                                     thisTrans->GetRotation());
1296     support->AddNode(plateVolume, i, plateTrans);
1297   }
1298   // Remove bolt holes 
1299   coneExp.Append("F3SB:F3SB1+F3SB:F3SB2)");
1300
1301   // Finalize the half-cone shape and add volume
1302   TGeoCompositeShape* coneShape  = new TGeoCompositeShape(coneExp.Data());
1303   TGeoVolume*         coneVolume = new TGeoVolume("F3SC", coneShape, fC);
1304   coneShape->SetName("F3SC");
1305   coneShape->SetTitle("FMD3 cone");
1306   coneVolume->SetTitle("FMD3 cone");
1307   support->AddNode(coneVolume,1,new TGeoTranslation(0,0,nlen+conel/2));
1308   
1309   // The flanges 
1310   TGeoBBox* flangeShape    = new TGeoBBox((flanger - backr2) / 2, 
1311                                           fmd3->GetBeamWidth() / 2,
1312                                           backl / 2);
1313   TGeoVolume* flangeVolume = new TGeoVolume(Form(fgkFlangeName, fmd3->GetId()),
1314                                             flangeShape, fC);
1315   flangeShape->SetName(Form(fgkFlangeName, fmd3->GetId()));
1316   flangeShape->SetTitle("FMD3 flange");
1317   flangeVolume->SetTitle("FMD3 flange");
1318   
1319   Int_t    n               = fmd3->GetNFlange();
1320   Double_t r               = backr2 + (flanger - backr2) / 2;
1321   for (Int_t i = 0; i  < n/2; i++) {
1322     Double_t phi       = 360. / n * i + 180. / n;
1323     Double_t x         = r * TMath::Cos(TMath::Pi() / 180 * phi);
1324     Double_t y         = r * TMath::Sin(TMath::Pi() / 180 * phi);
1325     TGeoRotation* rot  = new TGeoRotation;
1326     rot->RotateZ(phi);
1327     TGeoMatrix* matrix = new TGeoCombiTrans(x, y, nlen+conel-backl/2, rot);
1328     matrix->SetName(Form("FMD3_flange_%02d", i));
1329     matrix->SetTitle(Form("FMD3_flange_%2d", i));
1330     support->AddNode(flangeVolume, i, matrix);
1331   }
1332
1333   // Place support volumes in half-detector volumes 
1334   z                          = fmd3->GetInnerZ();
1335   z1                         = z-nz;
1336   fmd3TopVolume->AddNode(support, 1, new TGeoTranslation(0,0,z1));
1337   r1                         = new TGeoRotation();
1338   r1->RotateZ(180);
1339   t                          = new TGeoCombiTrans(0,0,z1,r1);
1340   fmd3BotVolume->AddNode(support, 2, t);
1341
1342   TGeoRotation*   rot        = new TGeoRotation("FMD3 rotatation");
1343   rot->RotateY(180);
1344   TGeoVolume*     top        = gGeoManager->GetVolume("ALIC");
1345   TGeoMatrix* mmatrix        = new TGeoCombiTrans("FMD3 trans", 0, 0, z, rot);
1346   AliFMDDebug(5, ("Placing volumes %s and %s in ALIC at z=%f", 
1347                    fmd3TopVolume->GetName(), fmd3BotVolume->GetName(), z));
1348   top->AddNode(fmd3TopVolume, fmd3->GetId(), mmatrix);
1349   top->AddNode(fmd3BotVolume, fmd3->GetId(), mmatrix);
1350
1351   return 0;
1352 }
1353 #endif
1354
1355 //____________________________________________________________________
1356 void
1357 AliFMDGeometryBuilder::Exec(Option_t*) 
1358 {
1359   // Setup up the FMD geometry. 
1360   AliFMDDebug(1, ("\tGeometry options: %s",
1361                     (fDetailed  ? "divided into strips" : "one volume")));
1362   if (!gGeoManager) {
1363     AliFatal("No TGeoManager defined");
1364     return;
1365   }
1366
1367   fSi      = gGeoManager->GetMedium("FMD_Si$");
1368   fC       = gGeoManager->GetMedium("FMD_Carbon$");
1369   fAl      = gGeoManager->GetMedium("FMD_Aluminum$");
1370   fChip    = gGeoManager->GetMedium("FMD_Si Chip$");
1371   fAir     = gGeoManager->GetMedium("FMD_Air$");
1372   fPCB     = gGeoManager->GetMedium("FMD_PCB$");
1373   fPlastic = gGeoManager->GetMedium("FMD_Plastic$");
1374   fCopper  = gGeoManager->GetMedium("FMD_Copper$");
1375   fSteel   = gGeoManager->GetMedium("FMD_Steel$");
1376
1377   if (!fSi||!fC||!fAl||!fChip||!fAir||!fPCB||!fPlastic||!fCopper||!fSteel) {
1378     AliError("Failed to get some or all tracking mediums");
1379     return;
1380   }    
1381   AliFMDGeometry* fmd = AliFMDGeometry::Instance();
1382   AliFMDRing* inner = fmd->GetInner();
1383   AliFMDRing* outer = fmd->GetOuter();
1384   RingGeometry(inner);
1385   RingGeometry(outer);
1386   TGeoVolume* innerTop = gGeoManager->GetVolume(Form(fgkRingTopName, 
1387                                                      inner->GetId()));
1388   TGeoVolume* innerBot = gGeoManager->GetVolume(Form(fgkRingBotName, 
1389                                                      inner->GetId()));
1390   TGeoVolume* outerTop = gGeoManager->GetVolume(Form(fgkRingTopName, 
1391                                                      outer->GetId()));
1392   TGeoVolume* outerBot = gGeoManager->GetVolume(Form(fgkRingBotName, 
1393                                                      outer->GetId()));
1394   
1395   FMD1Geometry(fmd->GetFMD1(), innerTop, innerBot);
1396   FMD2Geometry(fmd->GetFMD2(), innerTop, innerBot, outerTop, outerBot);
1397   FMD3Geometry(fmd->GetFMD3(), innerTop, innerBot, outerTop, outerBot);
1398 #ifndef USE_PRE_MOVE
1399   fmd->SetSectorOff(fSectorOff);
1400   fmd->SetModuleOff(fModuleOff);
1401   fmd->SetRingOff(fRingOff);
1402   fmd->SetDetectorOff(fDetectorOff);
1403   fmd->SetActive(fActiveId.fArray, fActiveId.fN);
1404 #endif
1405   // fmd->ExtractGeomInfo();
1406   
1407 }
1408
1409
1410 //____________________________________________________________________
1411 //
1412 // EOF
1413 //