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