]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSGeometrySSDCone.cxx
Removing obsolete macros
[u/mrichter/AliRoot.git] / ITS / AliITSGeometrySSDCone.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
16 /*
17 $Log$
18 Revision 1.3  2003/04/17 22:29:23  nilsen
19 Geometry bug fixes and the like. Work still progressing.
20
21 Revision 1.2  2003/03/25 23:27:19  nilsen
22 ITS new Geometry files. Not yet ready for uses, committed to allow additional
23 development.
24
25 Revision 1.1  2003/02/10 17:03:52  nilsen
26 New version and structure of ITS V11 geometry. Work still in progress.
27
28 $Id$
29 */
30
31 #include <Riostream.h>
32 #include <stdio.h>
33 #include <stdlib.h>
34 #include <TMath.h>
35 #include <TGeometry.h>
36 #include <TNode.h>
37 #include <TTUBE.h>
38 #include <TTUBS.h>
39 #include <TPCON.h>
40 #include <TFile.h>    // only required for Tracking function?
41 #include <TCanvas.h>
42 #include <TObjArray.h>
43 #include <TLorentzVector.h>
44 #include <TObjString.h>
45 #include <TClonesArray.h>
46 #include <TBRIK.h>
47 #include <TSystem.h>
48 #include <TVector3.h>
49 #include <AliRun.h>
50 #include <AliITS.h>
51
52 #include "AliITSGeometrySSDCone.h"
53
54 ClassImp(AliITSGeometrySSDCone)
55
56 //______________________________________________________________________
57 AliITSGeometrySSDCone::AliITSGeometrySSDCone() : AliITSBaseGeometry(){
58     //Default Constructor for SSD Cone geometry
59
60     SetScalemm();
61 }
62 //______________________________________________________________________
63 AliITSGeometrySSDCone::AliITSGeometrySSDCone(AliITS *its,TVector3 &tran,
64                                              const char *moth,Int_t mat0):
65     AliITSBaseGeometry(its,0){
66     //Standard Constructor for SSD Cone geometry
67     // Inputs:
68     //   Double_t z0  Z-axis shift of this volume
69     // Outputs:
70     //   none.
71     // Return:
72     //   none.
73     Double_t t; // some general angle and coordinates [degrees].
74     Double_t Z,Rmin,Rmax; // additional point not needed in call to pcons.
75
76     fThickness = 13.0; //mm, Thickness of Rohacell+carbon fiber
77     fCthick=1.5; //mm, Carbon finber thickness
78     fRcurv=15.0; // mm, Radius of curvature.
79     fTc=51.0; // angle of SSD cone [degrees].
80     fSintc=Sind(fTc);fCostc=Cosd(fTc);fTantc=Tand(fTc);
81     fZ0=0.0;fZouterMilled=13.5-5.0;fZcylinder=170.0;fZposts=196.0;
82     fNspoaks=12;fNposts=4;fNmounts=4;
83     fRoutMax=0.5*985.0;fRoutHole=0.5*965.0;fRoutMin=0.5*945.0;
84     fRholeMax=0.5*890.0;fRholeMin=0.5*740.0;
85     fRpostMin=316.0;fdRpost=23.0;fZpostMax=196.0;fPhi0Post=30.0;
86     fRinMax=0.5*590.0;fRinCylinder=0.5*597.0;fRinHole=0.5*575.0;
87     fRinMin=0.5*562.0;fdZin=15.0;
88     // SSD-SDD Thermal/Mechanical cylinder mounts
89     fNinScrews=40;
90     fPhi0Screws=0.5*360.0/((Double_t)fNinScrews);fRcylinderScrews= 0.5*570.0;
91     fDscrewHead=8.0;fDscrewShaft=4.6;fThScrewHeadHole=8.5;
92     // SDD mounting bracket, SSD part
93     fNssdSupports=3;fPhi0SDDsupports=90.0;
94     fRsddSupportPlate = 0.5*585.0;fThSDDsupportPlate=4.0;
95     fWsddSupportPlate = 70.0;
96     fSSDcf=26; // SSD support cone Carbon Fiber materal number.
97     fSSDfs=25; // SSD support cone inserto stesalite 4411w.
98     fSSDfo=68; // SSD support cone foam, Rohacell 50A.
99     fSSDsw=14; // SSD support cone screw material,Stainless steal
100     fNcD=0; // number of screw ends (copy number)
101     fNcE=0; // number of pin end (copy number)
102
103     SetScalemm();
104     // Poly-cone Volume A. Top part of SSD cone Carbon Fiber.
105     fA.Size(7,"SSD Suport cone Carbon Fiber Surface outer left");
106      // Poly-cone Volume B. Stesalite inside volume A.
107     fB.Size(6,"SSD Suport cone Inserto Stesalite left edge");
108     // Poly-cone Volume C. Foam inside volume A.
109     fC.Size(4,"SSD Suport cone Rohacell foam left edge");
110     fD.SetName("Screw+stud used to mount things to the SSD support cone");
111     fE.SetName("pin used to mount things to the SSD support cone");
112     // Poly-cone Volume F. Foam in spoak reagion, inside volume A.
113     fF.Size(4,"SSD Top Suport cone Rohacell foam Spoak");
114     fG.Size(4,"SSD spoak carbon fiber surfaces"); // Poly-cone Volume G.
115     fH.Size(4,"SSD support cone Rohacell foam Spoak"); // Poly-cone Volume H.
116     fI.Size(9,"SSD lower/inner right part of SSD cone"); //Poly-cone Volume I.
117     fJ.Size(4,"SSD inner most foam core"); // Poly-cone Volume J.
118     fK.Size(6,"SSD inner most inserto material"); // Poly-cone Volume K.
119     fL.Size(4,"SSD Bottom cone Rohacell foam Spoak"); // Poly-cone Volume L.
120     fM.Size(4,"SSD mounting post foam substitute, Inserto");//Poly-cone Vol. M
121     fN.Size(4,"SSD mounting post CF subsititute, Inserto");//Poly-cone Vol. N
122     fO.Size(3,"SSD mounting post, carbon fiber"); // Poly-cone Volume O.
123     fP.Size(3,"SSD mounting post, Inserto"); // Poly-cone Volume P.
124     fQ.Size(4,"SSD Thermal sheal stainless steel bolts");//Poly-cone Volume Q.
125     fR.SetName("Air in front of bolt (in stasolit)");
126     fS.SetName("Air in front of Stainless Steal Screw end, N6");
127     fT.Size(2,"SSD-SDD mounting bracket Inserto-> Al."); //Poly-cone Volume T.
128     fU.Size(4,"SSD-SDD mounting bracket CF->Al."); // Poly-cone Volume U.
129     // Lets start with the upper left outer carbon fiber surface.
130     // Between za[2],rmaxa[2] and za[4],rmaxa[4] there is a curved section
131     // given by rmaxa = rmaxa[2]-r*Sind(t) for 0<=t<=fTc and 
132     // za = za[2] + r*Cosd(t) for 0<=t<=fTc. Simularly between za[1],rmina[1
133     // and za[3],rmina[3] there is a curve section given by
134     // rmina = rmina[1]-r*Sind(t) for 0<=t<=fTc and za = za[1]+r&Sind(t)
135     // for t<=0<=fTc. These curves have been replaced by straight lines
136     // between the equivelent points for simplicity.
137     Double_t dza = fThickness/fSintc-(fRoutMax-fRoutMin)/fTantc;
138     if(dza<=0){ // The number or order of the points are in error for a proper
139         // call to pcons!
140         Error("SSDcone","The definition of the points for a call to PCONS is"
141               " in error. abort.");
142         return;
143     } // end if
144     fA.P0() = 0.0;
145     fA.dP() = 360.0;
146     fA.Z(0) = fZ0;
147     fA.Rn(0) = fRoutMin;
148     fA.Rx(0) = fRoutMax;
149     fA.Z(1)  = fA.ZAt(0)+fZouterMilled - dza; // za[2] - dza.
150     fA.Rn(1) = fA.Rmin(0);
151     fA.Rx(1) = fA.Rmax(0);
152     fA.Z(2)  = fA.ZAt(0)+fZouterMilled; //From Drawing ALR-0767 and ALR-0767/3
153     fA.Rx(2) = fA.Rmax(0);
154     RadiusOfCurvature(fRcurv,0.0,fA.ZAt(1),fA.Rmin(1),fTc,fA.Z(3),fA.Rn(3));
155     fA.Rn(2) = RminFrom2Points(fA,3,1,fA.ZAt(2));
156     RadiusOfCurvature(fRcurv,0.0,fA.ZAt(2),fA.Rmax(2),fTc,fA.Z(4),fA.Rx(4));
157     fA.Rn(4) = RminFromZSSDcone(fA.ZAt(4));
158     fA.Rx(3) = RmaxFrom2Points(fA,4,2,fA.ZAt(3));
159     fA.Rn(5) = fRholeMax;
160     fA.Z(5)  = Zfrom2MinPoints(fA,4,3,fA.Rmin(5));
161     fA.Rx(5) = RmaxFromZSSDcone(fA.ZAt(5));
162     fA.Rn(6) = fRholeMax;
163     fA.Rx(6) = fA.Rmin(6);
164     fA.Z(6)  = ZFromRmaxSSDcone(fA.Rmax(6));
165     //
166     // Now lets define the Inserto Stesalite 4411w material volume.
167     fB.P0() = 0.0;
168     fB.dP() = 360.0;
169     fB.Z(0) = fA.ZAt(0);
170     fB.Rn(0) = fA.Rmin(0)+fCthick;
171     fB.Rx(0) = fA.Rmax(0)-fCthick;
172     fB.Z(1)  = fA.ZAt(1);
173     fB.Rn(1) = fB.Rmin(0);
174     fB.Rx(1) = fB.Rmax(0);
175     fB.Z(2)  = fA.ZAt(2);
176     fB.Rx(2) = fB.Rmax(1);
177     RadiusOfCurvature(fRcurv-fCthick,0.,fB.ZAt(2),fB.Rmax(2),
178                                     fTc,fB.Z(3),fB.Rx(3));
179     RadiusOfCurvature(fRcurv+fCthick,0.,fB.ZAt(1),fB.Rmin(1),
180                                     fTc,fB.Z(4),fB.Rn(4));
181     fB.Rn(2) = RminFrom2Points(fB,4,1,fB.ZAt(2));
182     fB.Rn(3) = RminFrom2Points(fB,4,1,fB.ZAt(3));
183     fB.Z(5)  = fB.ZAt(4)+(fThickness-2.0*fCthick)/fSintc;
184     fB.Rn(5) = RmaxFromZSSDcone(fB.ZAt(5),-fCthick);
185     fB.Rx(5) = fB.Rmin(5);
186     fB.Rx(4) = RmaxFrom2Points(fB,5,3,fB.ZAt(4));
187     //
188     // Now lets define the Rohacell foam material volume.
189     fC.P0() = 0.0;
190     fC.dP() = 360.0;
191     fC.Z(0) = fB.ZAt(4);
192     fC.Rn(0) = fB.Rmin(4);
193     fC.Rx(0) = fC.Rmin(0);
194     fC.Z(1)  = fB.ZAt(5);
195     fC.Rx(1) = fB.Rmin(5);
196     fC.Rn(2) = fA.Rmin(5)+fCthick;//leave space for carbon fiber covering hole
197     fC.Z(2)  = ZFromRminSSDcone(fC.Rmin(2),+fCthick);
198     fC.Rn(1) = RminFrom2Points(fC,2,0,fC.ZAt(1));
199     fC.Rx(3) = fA.Rmin(6)+fCthick;
200     fC.Rn(3) = fC.Rmax(3);
201     fC.Z(3)  = ZFromRmaxSSDcone(fC.Rx(3),-fCthick);
202     fC.Rx(2) = RmaxFrom2Points(fC,3,1,fC.ZAt(2));
203     //
204     // In volume SCB, th Inserto Stesalite 4411w material volume, there
205     // are a number of Stainless steel screw and pin studs which will be
206     // filled with screws/studs.
207     fD.Rn()=0.0,fD.Rx()=6.0,fD.Z()=0.5*10.0; // mm
208     fE.Rn()=0.0;fE.Rx()=6.0;fE.Z()=0.5*12.0; // mm
209     //
210     // There is no carbon fiber between this upper left section and the
211     // SSD spoaks. We remove it by replacing it with Rohacell foam.
212     t = fCthick/(0.5*(fRholeMax+fRholeMin));// It is not posible to get the
213     // carbon fiber thickness uniform in this phi direction. We can only
214     // make it a fixed angular thickness.
215     t *= 180.0/TMath::Pi();
216     fF.P0()  = 12.5+t; // degrees see drawing ALR-0767.
217     fF.dP()  = 5.0 - 2.0*t; // degrees
218     fF.Z(0)  = fC.ZAt(2);
219     fF.Rn(0) = fC.Rmin(3);
220     fF.Rx(0) = fF.Rmin(0);
221     fF.Rn(1) = fA.Rmin(5);
222     fF.Rx(1) = fF.Rmin(0);
223     fF.Z(1)  = ZFromRminSSDcone(fF.Rmin(1),+fCthick);
224     fF.Z(2)  = fC.ZAt(3);
225     fF.Rn(2) = fF.Rmin(1);
226     fF.Rx(2) = fF.Rmax(1);
227     fF.Rn(3) = fA.Rmin(6);
228     fF.Rx(3) = fF.Rmin(3);
229     fF.Z(3)  = ZFromRmaxSSDcone(fF.Rmax(3),-fCthick);
230     //=================================================================
231     // Now for the spoak part of the SSD cone.
232     // It is not posible to inclue the radius of curvature between
233     // the spoak part and the upper left part of the SSD cone or lowwer right
234     // part. This would be discribed by the following curves.
235     // R = Rmax - (5mm)*Sin(t) phi = phi0+(5mm*180/(Pi*fRoutHole))*Sin(t) 
236     // where 0<=t<=90 For the inner curve a simular equiation holds.
237     fG.P0()  = 12.5; // degrees see drawing ALR-0767.
238     fG.dP()  = 5.0; // degrees
239     fG.Z(0)  = fA.ZAt(5);
240     fG.Rn(0) = fA.Rmin(5);
241     fG.Rx(0) = fG.Rn(0);
242     fG.Z(1)  = fA.ZAt(6);
243     fG.Rn(1) = RminFromZSSDcone(fG.ZAt(1));
244     fG.Rx(1) = fG.Rmax(0);
245     fG.Rn(2) = fRholeMin;
246     fG.Z(2)  = ZFromRminSSDcone(fG.Rmin(2));
247     fG.Rx(2) = RmaxFromZSSDcone(fG.ZAt(2));
248     fG.Rn(3) = fG.Rmin(2);
249     fG.Rx(3) = fG.Rmin(3);
250     fG.Z(3)  = ZFromRmaxSSDcone(fG.Rmax(3));
251     // For the foam core.
252     t = fCthick/(0.5*(fRholeMax+fRholeMin));// It is not posible to get the
253     // carbon fiber thickness uniform in this phi direction. We can only
254     // make it a fixed angular thickness.
255     t *= 180.0/TMath::Pi();
256     fH.P0()  = 12.5+t; // degrees
257     fH.dP()  = 5.0 - 2.0*t; // degrees see drawing ALR-0767.
258     fH.Z(0)  = fF.ZAt(1);
259     fH.Rn(0) = fG.Rmin(0);
260     fH.Rx(0) = fH.Rmin(0);
261     fH.Z(1)  = fF.ZAt(3);
262     fH.Rn(1) = RminFromZSSDcone(fH.Z(1),+fCthick);
263     fH.Rx(1) = fH.Rmax(0);
264     fH.Z(2)  = ZFromRminSSDcone(fG.Rmin(2),+fCthick);
265     fH.Rn(2) = fG.Rmin(2);
266     fH.Rx(2) = RmaxFromZSSDcone(fH.Z(2),-fCthick);
267     fH.Z(3)  = ZFromRmaxSSDcone(fG.Rmin(3),-fCthick);
268     fH.Rn(3) = fG.Rmin(3);
269     fH.Rx(3) = fH.Rn(3);
270     //
271     //==================================================================
272     // Now for the Inner most part of the SSD cone.
273     fI.P0()  = 0.0;
274     fI.dP()  = 360.0;
275     fI.Z(0)  = fG.ZAt(2);
276     fI.Rn(0) = fG.Rmin(2);
277     fI.Rx(0) = fI.Rmin(0);
278     fI.Z(1)  = fG.ZAt(3);
279     fI.Rn(1) = RminFromZSSDcone(fI.ZAt(1));
280     fI.Rx(1) = fI.Rmax(0);
281     fI.Rn(4) = fRinMin;
282     fI.Rn(5) = fRinMin;
283     RadiusOfCurvature(fRcurv,90.0,0.0,fRinMax,90.0-fTc,Z,fI.Rx(5)); // z dummy
284     fI.Z(5)  = ZFromRmaxSSDcone(fI.Rx(5));
285     fI.Z(6)  = fZcylinder;
286     fI.Rn(6) = fRinMin;
287     fI.Z(7)  = fI.Z(6);
288     fI.Rn(7) = fRinCylinder;
289     fI.Rn(8) = fRinCylinder;
290     fI.Rx(8) = fI.Rmin(8);
291     Rmin = fI.Rmin(5);
292     RadiusOfCurvature(fRcurv,90.0-fTc,fI.Z(5),fI.Rmax(5),90.0,Z,Rmax);
293     Rmax = fRinMax;
294     fI.Z(8)  = Z+(fI.ZAt(5)-Z)*(fI.Rmax(8)-Rmax)/(fI.Rmax(5)-Rmax);
295     fI.Rx(6) = RmaxFrom2Points(fI,8,5,fI.ZAt(6));
296     fI.Rx(7) = fI.Rmax(6);
297     fI.Z(3)  = Z-fdZin;
298     fI.Z(4)  = fI.ZAt(3);
299     fI.Rx(3) = RmaxFromZSSDcone(fI.ZAt(3));
300     fI.Rx(4) = fI.Rx(3);
301     //rmin dummy
302     RadiusOfCurvature(fRcurv,90.,fI.ZAt(3),0.,90.-fTc,fI.Z(2),Rmin);
303     fI.Rn(2) = RminFromZSSDcone(fI.ZAt(2));
304     fI.Rx(2) = RmaxFromZSSDcone(fI.ZAt(2));
305     // z dummy
306     RadiusOfCurvature(fRcurv,90.-fTc,0.0,fI.Rmin(2),90.0,Z,fI.Rn(3)); 
307     // Now for Inserto volume at the inner most radius.
308     fK.P0()  = 0.0;
309     fK.dP()  = 360.0;
310     fK.Z(1)  = fI.ZAt(3)+fCthick;
311     fK.Rn(1) = fI.Rmin(3);
312     fK.Z(2)  = fK.ZAt(1);
313     fK.Rn(2) = fI.Rmin(4);
314     fK.Rn(3) = fK.Rmin(2);
315     fK.Rn(4) = fK.Rmin(2);
316     fK.Rx(4) = fI.Rmax(5)-fCthick*fSintc;
317     RadiusOfCurvature(fRcurv+fCthick,90.0,fK.ZAt(1),fK.Rmin(1),
318                                 90.0-fTc,fK.Z(0),fK.Rn(0));
319     fK.Rx(0) = fK.Rmin(0);
320     fK.Z(3)  = fK.ZAt(0)+(fThickness-2.0*fCthick)*fCostc;;
321     fK.Rx(3) = fK.Rmax(0)+(fThickness-2.0*fCthick)*fSintc;
322     fK.Rx(1) = RmaxFrom2Points(fK,3,0,fK.ZAt(1));
323     fK.Rx(2) = fK.Rmax(1);
324     fK.Z(4)  = ZFromRmaxSSDcone(fK.Rmax(4),-fCthick);
325     fK.Rn(5) = fK.Rmin(2);
326     fK.Z(5)  = fI.ZAt(6);
327     fK.Rx(5) = (fI.Rmax(5)-fI.Rmax(8))/(fI.ZAt(5)-fI.ZAt(8))*
328                (fK.ZAt(5)-fK.ZAt(4)) + fK.Rmax(4);
329     // Now for foam core at the inner most radius.
330     fJ.P0() = 0.0;
331     fJ.dP() = 360.0;
332     fJ.Rn(0) = fI.Rmin(0)-fCthick;
333     fJ.Z(0)  = ZFromRminSSDcone(fJ.Rmin(0),+fCthick);
334     fJ.Rx(0) = fJ.Rmin(0);
335     fJ.Rx(1) = fJ.Rmax(0);
336     fJ.Z(1)  = ZFromRmaxSSDcone(fJ.Rmax(1),-fCthick);
337     fJ.Rn(1) = RminFromZSSDcone(fJ.ZAt(1),+fCthick);
338     fJ.Z(2)  = fK.ZAt(0);
339     fJ.Rn(2) = fK.Rmin(0);
340     fJ.Rx(2) = RmaxFromZSSDcone(fJ.ZAt(2),-fCthick);
341     fJ.Z(3)  = fK.ZAt(3);
342     fJ.Rn(3) = fK.Rmax(3);
343     fJ.Rx(3) = fJ.Rmin(3);
344     // Now for foam core at the top of the inner most radius where 
345     // the spoaks are.
346     t = fCthick/(0.5*(fRholeMax+fRholeMin));// It is not posible to get the
347     // carbon fiber thickness uniform in this phi direction. We can only
348     // make it a fixed angular thickness.
349     t *= 180.0/TMath::Pi();
350     fL.P0() = 12.5+t; // degrees
351     fL.dP() = 5.0 - 2.0*t; // degrees see drawing ALR-0767.
352     fL.Z(0) = fH.ZAt(2);
353     fL.Rn(0) = fH.Rmin(2);
354     fL.Rx(0) = fL.Rmin(0);
355     fL.Z(1)  = fJ.ZAt(0);
356     fL.Rn(1) = fJ.Rmin(0);
357     fL.Rx(1) = fI.Rmax(1);
358     fL.Z(2)  = fH.ZAt(3);
359     fL.Rn(2) = fL.Rmin(1);
360     fL.Rx(2) = fL.Rmax(1);
361     fL.Z(3)  = fJ.ZAt(1);
362     fL.Rn(3) = fL.Rmin(2);
363     fL.Rx(3) = fL.Rmin(3);
364     // Now for the SSD mounting posts
365     fO.P0()  = fPhi0Post; // degrees
366     fO.dP()  = 180.0*fdRpost/(fRpostMin+0.5*fdRpost)/TMath::Pi(); //
367     fO.Rn(0) = fRpostMin+fdRpost;
368     fO.Rx(0) = fO.Rmin(0);
369     fO.Z(0)  = ZFromRmaxSSDcone(fO.Rmax(0));
370     fO.Rn(1) = fRpostMin;
371     fO.Z(1)  = ZFromRmaxSSDcone(fO.Rmin(1));
372     fO.Rx(1) = fO.Rmax(0);
373     fO.Z(2)  = fZ0+fZpostMax;
374     fO.Rn(2) = fRpostMin;
375     fO.Rx(2) = fO.Rmin(2)+fdRpost;
376     // Now for the SSD mounting posts
377     t = 180.0*fCthick/(fRpostMin+0.5*fdRpost)/TMath::Pi();
378     fP.dP()  = fO.DPhi()-2.0*t; // degrees
379     fP.P0()  = fO.Phi0()+t; //
380     fP.Rn(0) = fO.Rmin(0)-fCthick;
381     fP.Rx(0) = fP.Rmin(0);
382     fP.Z(0)  = ZFromRmaxSSDcone(fP.Rmax(0));
383     fP.Rn(1) = fO.Rmin(1)+fCthick;
384     fP.Rx(1) = fO.Rmin(0)-fCthick;
385     fP.Z(1)  = ZFromRmaxSSDcone(fP.Rmin(1));
386     fP.Rn(2) = fP.Rmin(1);
387     fP.Rx(2) = fP.Rmax(1);
388     fP.Z(2)  = fZ0+fZpostMax;
389     // This insrto continues into the SSD cone displacing the foam
390     // and the carbon fiber surface at those points where the posts are.
391     fM.P0()  = fP.Phi0();
392     fM.dP()  = fP.DPhi();
393     fM.Rn(0) = fRpostMin+fdRpost-fCthick;
394     fM.Rx(0) = fM.Rmin(0);
395     fM.Z(0)  = ZFromRminSSDcone(fM.Rmin(0),+fCthick);
396     fM.Rx(1) = fM.Rmax(0);
397     fM.Z(1)  = ZFromRmaxSSDcone(fM.Rmax(1),-fCthick);
398     fM.Rn(1) = RminFromZSSDcone(fM.ZAt(1),+fCthick);
399     fM.Rn(2) = fRpostMin+fCthick;
400     fM.Z(2)  = ZFromRminSSDcone(fM.Rmin(2),+fCthick);
401     fM.Rx(2) = RmaxFromZSSDcone(fM.ZAt(2),-fCthick);
402     fM.Rn(3) = fM.Rmin(2);
403     fM.Rx(3) = fM.Rmin(3);
404     fM.Z(3)  = ZFromRmaxSSDcone(fM.Rmax(3),-fCthick);
405     //
406     fN.P0()  = fP.Phi0();
407     fN.dP()  = fP.DPhi();
408     fN.Z(0)  = fM.ZAt(1);
409     fN.Rn(0) = fM.Rmax(1);
410     fN.Rx(0) = fN.Rmin(0);
411     fN.Rx(1) = fN.Rmax(0);
412     fN.Z(1)  = ZFromRmaxSSDcone(fN.Rmax(1));
413     fN.Rn(1) = RmaxFromZSSDcone(fN.ZAt(1),-fCthick);
414     fN.Z(2)  = fM.ZAt(3);
415     fN.Rn(2) = fM.Rmin(3);
416     fN.Rx(2) = RmaxFromZSSDcone(fN.ZAt(2));
417     fN.Rn(3) = fN.Rmin(2);
418     fN.Rx(3) = fN.Rmin(3);
419     fN.Z(3)  = ZFromRmaxSSDcone(fN.Rmax(3));
420     // Bolt heads holding the SSD-SDD tube to the SSD cone.
421     // Bolt -- PolyCone
422     fQ.P0()  = 0.0;
423     fQ.dP()  = 360.0;
424     fQ.Z(0)  = fI.ZAt(4)-fThSDDsupportPlate;
425     fQ.Rn(0) = 0.0;
426     fQ.Rx(0) = 0.5*fDscrewHead;
427     fQ.Z(1)  = fI.ZAt(4)-fThScrewHeadHole;
428     fQ.Rn(1) = 0.0;
429     fQ.Rx(1) = 0.5*fDscrewHead;
430     fQ.Z(2)  = fQ.ZAt(1);
431     fQ.Rn(2) = 0.0;
432     fQ.Rx(2) = 0.5*fDscrewShaft;
433     fQ.Z(3)  = fQ.ZAt(2);
434     fQ.Rn(3) = 0.0;
435     fQ.Rx(3) = fQ.Rmax(2);
436     // air infront of bolt (stasolit Volume K) -- Tube
437     fR.Z() = 0.5*(fThickness-fThScrewHeadHole);
438     fR.Rn() = 0.0;
439     fR.Rx() = 0.5*fDscrewHead;
440     // air infront of bolt (carbon fiber volume I) -- Tube
441     fS.Z() = 0.5*fThickness;
442     fS.Rn() = 0.0;
443     fS.Rx() = fR.Rmax();
444     // SDD support plate, SSD side.
445     fT.dP()  = 180.0*fWsddSupportPlate/(fRsddSupportPlate*TMath::Pi());
446     fT.P0()  = fPhi0SDDsupports=-0.5*fT.DPhi();
447     fT.Z(0)  = fK.ZAt(2);
448     fT.Rn(0) = fI.Rmin(4);
449     fT.Rx(0) = fRsddSupportPlate;
450     fT.Z(1)  = fI.ZAt(4) - fThSDDsupportPlate;
451     fT.Rn(1) = fT.Rmin(0);
452     fT.Rx(1) = fT.Rmax(0);
453     //
454     fU.dP() = fT.DPhi();
455     fU.P0() = fT.Phi0();
456     fU.Z(2) = fI.ZAt(4);
457     fU.Rn(2) = fT.Rmin(0);
458     fU.Rx(2) = fT.Rmax(0);
459     fU.Z(3)  = fT.ZAt(0);
460     fU.Rn(3) = fU.Rmin(2);
461     fU.Rx(3) = fU.Rmax(2);
462     fU.Z(1)  = fU.ZAt(2);
463     fU.Rn(1) = fI.Rmin(3);
464     fU.Rx(1) = fU.Rmax(3);
465     fU.Rx(0) = fT.Rmax(0);
466     fU.Rn(0) = fU.Rmax(0);
467     fU.Z(0)  = Zfrom2MinPoints(fI,2,3,fU.Rmax(0));
468     // Debuging
469     Print(&cout);
470 }
471 //______________________________________________________________________
472 void AliITSGeometrySSDCone::CreateG3Geometry(const char *moth,
473                                              TVector3 &trans){
474     // Calls Geant 3 geometry inilization routines with the information
475     // stored in this class.
476     // Inputs:
477     //    none.
478     // Outputs:
479     //    none.
480     // Return:
481     //    none.
482
483     PolyCone(fA,fSSDcf);
484     PolyCone(fB,fSSDfs);
485     PolyCone(fC,fSSDfo);
486     Tube(fD,fSSDsw);
487     Tube(fE,fSSDsw);
488     PolyCone(fF,fSSDfo);
489     PolyCone(fG,fSSDcf);
490     PolyCone(fH,fSSDfo);
491     PolyCone(fI,fSSDcf);
492     PolyCone(fJ,fSSDfo);
493     PolyCone(fK,fSSDfs);
494     PolyCone(fL,fSSDfo);
495     PolyCone(fM,fSSDfo);
496     PolyCone(fN,fSSDfo);
497     PolyCone(fO,fSSDcf);
498     PolyCone(fP,fSSDfs);
499     PolyCone(fQ,fSSDfo);
500     Tube(fR,fSSDsw);
501     Tube(fS,fSSDsw);
502     PolyCone(fT,fSSDfo);
503     PolyCone(fU,fSSDfo);
504     return;
505 }
506 //______________________________________________________________________
507 void AliITSGeometrySSDCone::PositionG3Geometry(AliITSBaseVolParams &moth,
508                                                Int_t cn,TVector3 &trans,
509                                                Int_t irot){
510     // Positions ths whole object at t with rotatin irot coply number cn
511     // into volume moth.
512     // Inputs:
513     //   const AliITSBaseVolParams *moth   Mother volume where this object 
514     //                                     is to be placed.
515     //   Int_t      cn      Copy number.
516     //   TVector3   &t      Translation vector for this whole volume
517     //   Int_t      irot    rotation matrix number to be applyed to this
518     //                      volume.
519     // Output:
520     //   none.
521     // Return:
522     //   none.
523     Int_t i,j,k,l,irotSpoaks,irotPost;
524     Double_t t;
525     Bool_t init=kFALSE;
526     TVector3 zero(0.0,0.0,0.0),v(0.0,0.0,0.0);
527
528     if(cn<=0) return;
529     if(cn==1) init=kTRUE;
530     Pos(fA,cn,moth,trans,0);
531     Pos(fI,cn,moth,trans,0);
532     Pos(fG,fNspoaks*(cn-1)+1,moth,trans,0);
533     irotSpoaks = irot;
534     j = 1;
535     for(i=fNspoaks*(cn-1)+2;i<fNspoaks*cn+1;i++){
536         ZMatrix(++irot,((Double_t)j)*360./((Double_t)fNspoaks));
537         Pos(fG,i,moth,trans,irot);
538         j++;
539     } // end for i
540     Pos(fO,fNposts*(cn-1)+1,moth,trans,0);
541     irotPost = irot;
542     j = 1;
543     for(i=fNposts*(cn-1)+2;i<fNposts*cn+1;i++){
544         ZMatrix(++irot,((Double_t)j)*360./((Double_t)fNposts));
545         Pos(fO,i,moth,trans,irot);
546         j++;
547     } // end for
548     if(!init) return;
549     // Inside volume A.
550     Pos(fB,1,fA,zero,0);
551     Pos(fC,1,fA,zero,0);
552     // Inside Volume B
553     k=l=0;
554     for(i=0;i<2;i++){ // position for ITS-TPC mounting brackets
555         for(j=0;j<2;j++){ // 2 screws per bracket
556             fNcD++;
557             t = -5.0+10.0*((Double_t)j)+180.*((Double_t)i);
558             v.SetX(fRoutHole*Sind(t));
559             v.SetY(fRoutHole*Cosd(t));
560             v.SetZ(fD.DzAt());
561             Pos(fD,fNcD,fB,v,0);
562         } // end for j
563         for(j=0;j<3;j++){ // 3 pins per bracket
564             fNcE++;
565             t = -3.0+3.0*((Double_t)j)+180.*((Double_t)i);
566             v.SetX(fRoutHole*Sind(t));
567             v.SetY(fRoutHole*Cosd(t));
568             v.SetZ(fE.DzAt());
569             Pos(fE,fNcE,fB,v,0);
570         } // end for j
571     } // end for i
572     for(i=0;i<2;i++){ // position for ITS-rail mounting brackets
573         for(j=0;j<4;j++){ // 4 screws per bracket
574             Double_t a[4]={0.0,2.0,5.0,7.0}; // Relative angles.
575             fNcD++;
576             t = 90.0-a[j]+187.*((Double_t)i);
577             v.SetX(fRoutHole*Sind(t));
578             v.SetY(fRoutHole*Cosd(t));
579             v.SetZ(fD.DzAt());
580             Pos(fD,fNcD,fB,v,0);
581         } // end for j
582         for(j=0;j<2;j++){ // 2 pins per bracket
583             fNcE++;
584             t = 88+7.0*((Double_t)j)+184.*((Double_t)i);
585             v.SetX(fRoutHole*Sind(t));
586             v.SetY(fRoutHole*Cosd(t));
587             v.SetZ(fE.DzAt());
588             Pos(fE,fNcE,fB,v,0);
589         } // end for j
590     } // end for i
591     for(i=0;i<fNmounts;i++){ // mounting holes/screws for beam 
592         // pipe support and SPD cone support (dump side,non-dump 
593         // side has them to).
594         for(j=0;j<2;j++){ // 2 screws per bracket
595             fNcD++;
596             t = 180.*20./(fRoutHole*TMath::Pi());
597             t = 45.0+((Double_t)(j-1))*t+90.*((Double_t)i);
598             v.SetX(fRoutHole*Sind(t));
599             v.SetY(fRoutHole*Cosd(t));
600             v.SetZ(fD.DzAt());
601             Pos(fD,fNcD,fB,v,0);
602         } // end for j
603         for(j=0;j<1;j++){ // 1 pins per bracket
604             fNcE++;
605             t = 45.0+90.*((Double_t)i);
606             v.SetX(fRoutHole*Sind(t));
607             v.SetY(fRoutHole*Cosd(t));
608             v.SetZ(fE.DzAt());
609             Pos(fE,fNcE,fB,v,0);
610         } // end for j
611     } // end for i
612     Pos(fF,1,fA,zero,0);
613     Pos(fL,1,fI,zero,0);
614     for(i=1;i<fNspoaks;i++){
615         Pos(fF,i+1,fA,zero,irotSpoaks+i);
616         Pos(fL,i+1,fA,zero,irotSpoaks+i);
617     } // end for i
618     Pos(fH,1,fG,zero,0);
619     Pos(fK,1,fI,zero,0);
620     Pos(fJ,1,fI,zero,0);
621     Pos(fP,1,fO,zero,0);
622     Pos(fM,1,fJ,zero,0);
623     Pos(fN,1,fI,zero,0);
624     for(i=1;i<fNposts;i++){
625         Pos(fN,i+1,fI,zero,irotPost+i);
626         Pos(fM,i+1,fJ,zero,irotPost+i);
627     } // end for i
628     // Add SDD thermal sheald screws, and SDD suppport bracket.
629     t = fPhi0Screws;
630     v.SetX(fRcylinderScrews*Sind(t));
631     v.SetY(fRcylinderScrews*Cosd(t));
632     v.SetZ(0.0);
633     Pos(fQ,1,fK,v,0);
634     v.SetZ(fK.ZAt(2)-fR.DzAt());
635     Pos(fR,1,fK,v,0);
636     v.SetZ(fI.ZAt(4)-fS.DzAt());
637     Pos(fS,1,fI,v,0);
638     Pos(fT,1,fK,zero,0);
639     Pos(fU,1,fI,zero,0);
640     Int_t ir=2,it=2;
641     Double_t d = 360./((Double_t)fNssdSupports);
642     cout << fPhi0SDDsupports<<" " << fNssdSupports<< " " << fPhi0Screws<< endl;
643     for(i=1;i<fNinScrews;i++){
644         t = fPhi0Screws + ((Double_t)i)*360./((Double_t)fNinScrews);
645         v.SetX(fRcylinderScrews*Sind(t));
646         v.SetY(fRcylinderScrews*Cosd(t));
647         v.SetZ(0.0);
648         Pos(fQ,i+1,fK,v,0);
649         cout << " t="<<t<<endl;
650         if((t>=-0.5*fT.DPhi()&&t<=+0.5*fT.DPhi())||
651            (t>=d-0.5*fT.DPhi()&&t<=d+0.5*fT.DPhi())||
652            (t>=2.*d-0.5*fT.DPhi()&&t<=2.*d+0.5*fT.DPhi())){
653             // SDD support here
654             if((t>=-0.5*fT.DPhi()&&t<=+0.5*fT.DPhi())) continue;
655             t = ((Double_t)(it-1))*d;
656             cout << " Zmatrix t=" << t << endl;
657             ZMatrix(++irot,t);
658             Pos(fT,it,fK,zero,irot);
659             Pos(fU,it,fI,zero,irot);
660             it++;
661         }else{ // Air here
662             v.SetZ(fK.ZAt(2)-fR.DzAt());
663             Pos(fR,ir,fK,v,0);
664             v.SetZ(fI.ZAt(4)-fS.DzAt());
665             Pos(fS,ir,fI,v,0);
666             ir++;
667         } // end if
668     } // end for i
669     // Air holes and bolts in SSD mounting studs, volumes I,M,N,O.
670     return;
671 }
672 //______________________________________________________________________
673 void AliITSGeometrySSDCone::CreateG3Materials(){
674     // Fills the Geant 3 banks with Material and Medium definisions.
675     // Inputs:
676     //   none.
677     // Outputs:
678     //   none.
679     // Returns:
680     //   none.
681     Int_t i;
682     Int_t Z[5],N[5];
683     Double_t W[5],dens;
684
685     // epoxy
686     dens = 10.*GetA(1)+13.*GetA(6)+3.*GetA(8);
687     Z[0] = 1; W[0] = 10.*GetA(Z[0])/dens; // Hydrogen Content
688     Z[1] = 6; W[1] = 13.*GetA(Z[1])/dens; // Carbon Content
689     Z[2] = 8; W[2] =  3.*GetA(Z[2])/dens; // Oxegen
690     // Carbon fiber is about 64% carbon fiber and 36% epoxy by volume.
691     // Now need to add in the carbon fiber
692     W[0] *= 0.36*GetA(Z[0]);
693     W[1]  = 0.36*W[1]*dens*GetA(Z[1]) + 0.64*GetA(Z[1]);
694     W[2] *= 0.36*GetA(Z[2]);
695     // Renormilize the weights
696     dens = 0.0;
697     for(i=0;i<3;i++){dens += W[i];}
698     for(i=0;i<3;i++){W[i] /= dens;}
699     dens = 1.7; // grams/cm^3 taken as density of G10 PDG book.
700     MixtureByWeight(fSSDcf,"Carbon Fiber for SSD support cone",Z,W,dens,3,0);
701     // epoxy 
702     dens = 10.*GetA(1)+13.*GetA(6)+3.*GetA(8);
703     Z[0] = 1; W[0] = 10.*GetA(Z[0])/dens; // Hydrogen Content
704     Z[1] = 6; W[1] = 13.*GetA(Z[1])/dens; // Carbon Content
705     Z[2] = 8; W[2] =  3.*GetA(Z[2])/dens; // Oxegen
706     Z[3] = 14;W[3] = 0.0; // no Silicon in epoxy.
707     // glass fiber is about 64% carbon fiber and 36% epoxy by volume.
708     // Now need to add in the glass fiber
709     W[0] *= 0.36*GetA(Z[0]);
710     W[1] *= 0.36*GetA(Z[1]);
711     W[2]  = 0.36*W[2]*dens*GetA(Z[2]) + 0.64*2.0*GetA(Z[2]);
712     W[3]  = 0.64*GetA(Z[3]); // Si
713     // Renormilize the weights
714     dens = 0.0;
715     for(i=0;i<4;i++){dens += W[i];}
716     for(i=0;i<4;i++){W[i] /= dens;}
717     dens = 1.7; // grams/cm^3 taken as density of G10 PDG book.
718     MixtureByWeight(fSSDfs,"Inserto stealite 4411w for SSD support cone",
719                     Z,W,dens,4,0);
720     // Rohacell 51 C14 H10 N2 O6 from Flavio Tosello
721     // http://cesweb.grantadesign.com/demo/index.do
722     Z[0] = 1; N[0] = 10; // Hydrogen Content
723     Z[1] = 6; N[1] = 14; // Carbon Content
724     Z[2] = 7; N[2] =  2; // Nitrogen Content
725     Z[3] = 8; N[3] =  6; // Oxigen Content
726     dens = 0.0513; // grams/cm^3 From Flavio Tosello 
727     //  http://www.emkayplatics.co.uk/roh51.html
728     MixtureByNumber(fSSDfo,"Foam core (Rohacell 51) for SSD support cone",
729                     Z,N,dens,4,0);
730     // Stainless steel. Temperary values.
731     Z[0] =  6; W[0] = 0.5; // Carbon Content
732     Z[1] = 25; W[1] = 0.5; // Iron Content
733     dens = 7.87; // Grams/cm^3  density of iron used.
734     MixtureByWeight(fSSDsw,"Stainless steal screw, pin, and stud material",
735                     Z,W,dens,2,0);
736 }
737 //______________________________________________________________________
738 void AliITSGeometrySSDCone::BuildDisplayGeometry(){
739     // Fill Root geometry banks for fast simple ITS simulation event
740     // display. See Display.C, and related code, for more details.
741     // Inputs:
742     //    none.
743     // Outputs:
744     //   none.
745     // Return:
746     //  none.
747
748     // No need to display ITS cones.
749 }
750 //______________________________________________________________________
751 void AliITSGeometrySSDCone::Print(ostream *os){
752     // Prints out the data kept in this class
753     // Inputs:
754     //   ostream *os pointer to the output stream
755     // Outputs:
756     //   none.
757     // Return:
758     //   none.
759
760     *os << "Object AliITSGeometrySSDCone" << endl;
761     *os << " Object fA" << endl << fA << endl;
762     *os << " Object fB" << endl << fB << endl;
763     *os << " Object fC" << endl << fC << endl;
764     *os << " Object fD" << endl << fD << endl;
765     *os << " Object fE" << endl << fE << endl;
766     *os << " Object fF" << endl << fF << endl;
767     *os << " Object fG" << endl << fG << endl;
768     *os << " Object fH" << endl << fH << endl;
769     *os << " Object fI" << endl << fI << endl;
770     *os << " Object fJ" << endl << fJ << endl;
771     *os << " Object fK" << endl << fK << endl;
772     *os << " Object fL" << endl << fL << endl;
773     *os << " Object fM" << endl << fM << endl;
774     *os << " Object fN" << endl << fN << endl;
775     *os << " Object fO" << endl << fO << endl;
776     *os << " Object fP" << endl << fP << endl;
777     *os << " Object fQ" << endl << fQ << endl;
778     *os << " Object fR" << endl << fR << endl;
779     *os << " Object fS" << endl << fS << endl;
780     *os << " Object fT" << endl << fT << endl;
781     *os << " Object fU" << endl << fU << endl;
782     return;
783 }
784 //______________________________________________________________________
785 void AliITSGeometrySSDCone::Read(istream *is){
786     // Read in data written with Print above.
787     // Inputs:
788     //   istream *is Input stream pointer
789     // Output:
790     //   none.
791     // Return:
792     //   none.
793     char s[50];
794
795     is->getline(s,49);
796     is->getline(s,49);
797     *is >> fA;
798     is->getline(s,49);
799     *is >> fB;
800     is->getline(s,49);
801     *is >> fC;
802     is->getline(s,49);
803     *is >> fD;
804     is->getline(s,49);
805     *is >> fE;
806     is->getline(s,49);
807     *is >> fF;
808     is->getline(s,49);
809     *is >> fG;
810     is->getline(s,49);
811     *is >> fH;
812     is->getline(s,49);
813     *is >> fI;
814     is->getline(s,49);
815     *is >> fJ;
816     is->getline(s,49);
817     *is >> fK;
818     is->getline(s,49);
819     *is >> fL;
820     is->getline(s,49);
821     *is >> fM;
822     is->getline(s,49);
823     *is >> fN;
824     is->getline(s,49);
825     *is >> fO;
826     is->getline(s,49);
827     *is >> fP;
828     is->getline(s,49);
829     *is >> fQ;
830     is->getline(s,49);
831     *is >> fR;
832     is->getline(s,49);
833     *is >> fS;
834     is->getline(s,49);
835     *is >> fT;
836     is->getline(s,49);
837     *is >> fU;
838     return;
839 }
840 //______________________________________________________________________
841 ostream &operator<<(ostream &os,AliITSGeometrySSDCone &s){
842     // Operator << for C++ like output of AliITSGeometrySSDCone class.
843     // Inputs:
844     //   ostream &os              The output stream
845     //   AliITSGeometrySSDCone &s The class to be outputed
846     // Outputs:
847     //   none.
848     // Return:
849     //  ostream &os  The address of the output stream
850
851     s.Print(&os);
852     return os;
853 }
854 //______________________________________________________________________
855 istream &operator>>(istream &is,AliITSGeometrySSDCone &s){
856     // Operator >> for C++ like input of AliITSGeometrySSDCone class.
857     // Inputs:
858     //   istream &is              The input stream
859     //   AliITSGeometrySSDCone &s The class to be inputed
860     // Outputs:
861     //   none.
862     // Return:
863     //  istream &is  The address of the input stream
864
865     s.Read(&is);
866     return is;
867 }