]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDv2.cxx
removing temporary file
[u/mrichter/AliRoot.git] / FMD / AliFMDv2.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 // Forward Multiplicity detector based on Silicon version 0        //
18 //
19 //Begin Html       
20 /*
21 <img src="gif/AliFMDv0Class.gif">
22 */
23 //End Html
24 //                                                                  //
25 //                                                                  //
26 //////////////////////////////////////////////////////////////////////
27 #include <TMath.h>
28 #include <TGeometry.h>
29 #include <TTUBE.h>
30 #include <TFile.h>
31 #include <TTree.h>
32 #include <TNode.h>
33 #include <TClonesArray.h>
34 #include <TLorentzVector.h>
35 #include <TDirectory.h>
36 #include "AliFMDv2.h"
37 #include "AliFMDv0.h"
38 #include "AliRun.h"
39 #include <Riostream.h>
40 #include "AliMagF.h"
41 #include "AliFMDhit.h"
42 #include "AliFMDdigit.h"
43 #include <stdlib.h>
44
45 ClassImp(AliFMDv2)
46
47 //--------------------------------------------------------------------
48 AliFMDv2::AliFMDv2(const char *name, const char *title):
49  AliFMD(name,title)
50 {
51   //
52   // Standart constructor for Forward Multiplicity Detector version 0
53   //
54   fIdSens1=0;
55   fIdSens2=0;
56   fIdSens3=0;
57   fIdSens4=0;
58   fIdSens5=0;
59 //  setBufferSize(128000);
60 }
61 //-------------------------------------------------------------------------
62 void AliFMDv2::CreateGeometry()
63 {
64  //
65   // Create the geometry of Forward Multiplicity Detector version 0
66   //
67   //Detector consists of 6 volumes: 
68   // 1st covered pseudorapidity interval from 3.3 to 2.0
69   // and placed on 65cm in Z-direction;
70   // 2nd - from 2.0 to 1.6 and Z=85 cm;
71   // 3d  - the same pseudorapidity interval as the 1st 
72   // but on the other side from the interaction point z=-65cm;
73   // 4th - simmetricaly with the 2nd : 
74   // pseudorapidity from 2.0 to 1.6, Z=-85cm   
75   // 5th - from 3.6 to 4.7, Z=-270cm
76   // 6th - from 4.5 to 5.5 , Z=-630cm.
77   // Each part has 400mkm Si (sensetive area, detector itself),
78   // 0.75cm of plastic simulated electronics material,
79   // Al support ring 2cm thickness and 1cm width placed on 
80   // the outer radius of each Si disk;
81   //    
82   // begin Html
83   /*
84    <img src="gif/AliFMDv0.gif">
85    */
86   //
87
88   Int_t *idtmed = fIdtmed->GetArray();
89    
90   Int_t ifmd;
91   Int_t idrotm[999];
92   Float_t zFMD,par[3],ppcon[5];
93     //ppcon1[15],ppcon2[15],ppcon3[15],ppcon4[15],ppcon5[15];
94   Float_t z[5]={62.8, 75.2, -83.4, -75.2, -340.};
95   Float_t NylonTube[3]={0.2,0.6,0.45};
96   Float_t zPCB=0.12; Float_t zHoneyComb=0.5; 
97   Float_t zSi=0.03;
98  
99 //-------------------------------------------------------------------
100  //  FMD 
101  //------------------------------------------------------------------
102   //  cout<<" !!!!!!!!!!!New FMD geometry !!!!!!!!!"<<endl;
103   
104
105   char nameFMD[5], nameSi[5], nameSector[5], nameRing[5];
106   Char_t nameHoney[5], nameHoneyIn[5], nameHoneyOut[5];
107   Char_t namePCB[5], nameCopper[5], nameChips[5], nameG10[5];
108   Char_t nameLPCB[5], nameLCopper[5], nameLChips[5], nameGL10[5];;
109   Float_t rin[5]={4.2,15.4,4.2,15.4,4.2};
110   Float_t rout[5]={17.4,28.4,17.4,28.4,17.4};
111   Float_t RinHoneyComb[5] ={ 5.15,16.4,  5.15,16.4,  5.15};
112   Float_t RoutHoneyComb[5]={20.63,34.92,22.3, 32.02,20.63};
113   Float_t zInside;
114   Float_t zCooper=0.01; Float_t zChips=0.01;
115   Float_t yNylonTube[5]={10,20,10,20,10};
116
117
118   AliMatrix(idrotm[901], 90, 0, 90, 90, 180, 0);
119   
120   
121   // Nylon tubes
122    gMC->Gsvolu("GNYL","TUBE", idtmed[1], NylonTube, 3);  //support nylon tube
123    Float_t wideSupport=zSi+3*zPCB+2*NylonTube[2]+zHoneyComb;
124    //   cout<<" wideSupport "<<wideSupport<<endl;
125
126   for (ifmd=0; ifmd<5; ifmd++)
127     {
128       sprintf(nameFMD,"FMD%d",ifmd+1);
129       ppcon[0]=0;
130       ppcon[1]=360;
131       ppcon[2]=4;
132       
133       //      ppcon[3]=z[ifmd];
134       ppcon[3]=-wideSupport/2.;
135       ppcon[4]=rin[ifmd];
136       ppcon[5]=rout[ifmd];
137       
138       ppcon[6]=ppcon[3]+zSi+2*zPCB+2*NylonTube[2];
139       ppcon[7]=rin[ifmd];
140       ppcon[8]=rout[ifmd];
141       
142       ppcon[9]=ppcon[6];
143       ppcon[10]=RinHoneyComb[ifmd];
144       ppcon[11]=RoutHoneyComb[ifmd];
145
146       ppcon[12]=ppcon[9]+zHoneyComb+zPCB;
147       ppcon[13]=RinHoneyComb[ifmd];
148       ppcon[14]=RoutHoneyComb[ifmd];
149       zFMD=z[ifmd]+wideSupport/2.;
150       //      cout<<" Si "<<ifmd+1<<" "<<zFMD<<" "<<ppcon[12]<<endl;
151       gMC->Gsvolu(nameFMD,"PCON",idtmed[0],ppcon,15);
152       if (z[ifmd] >0){  
153         gMC->Gspos(nameFMD,1,"ALIC",0,0,z[ifmd],0, "ONLY");}
154       else { 
155         gMC->Gspos(nameFMD,1,"ALIC",0,0,z[ifmd],idrotm[901], "ONLY");}
156      //silicon
157       sprintf(nameSi,"GSI%d",ifmd+1);
158       sprintf(nameSector,"GSC%d",ifmd+1);
159       sprintf(nameRing,"GRN%d",ifmd+1);
160       //honeycomb support
161       sprintf(nameHoney,"GSU%d",ifmd+1);
162       gMC->Gsvolu(nameHoney,"TUBE", idtmed[0], par, 0);  //honeycomb 
163       sprintf(nameHoneyIn,"GHI%d",ifmd+1);
164       gMC->Gsvolu(nameHoneyIn,"TUBE", idtmed[7], par, 0);  //honey comb inside 
165       sprintf(nameHoneyOut,"GHO%d",ifmd+1);
166       gMC->Gsvolu(nameHoneyOut,"TUBE", idtmed[6], par, 0);  //honey comb skin
167       //PCB
168       sprintf(namePCB,"GPC%d",ifmd+1);
169       gMC->Gsvolu(namePCB,"TUBE", idtmed[0], par, 0); //PCB
170       sprintf(nameCopper,"GCO%d",ifmd+1);
171       gMC->Gsvolu(nameCopper,"TUBE", idtmed[3], par, 0);  // Cooper
172       sprintf(nameChips,"GCH%d",ifmd+1);
173       gMC->Gsvolu(nameChips,"TUBE", idtmed[5], par, 0); // Si chips
174       sprintf(nameG10,"G10%d",ifmd+1);
175       gMC->Gsvolu(nameG10,"TUBE", idtmed[2], par, 0);  //G10 plate
176       //last PCB
177       sprintf(nameLPCB,"GPL%d",ifmd+1);
178       gMC->Gsvolu(nameLPCB,"TUBE", idtmed[0], par, 0); //PCB
179       sprintf(nameLCopper,"GCL%d",ifmd+1);
180       gMC->Gsvolu(nameLCopper,"TUBE", idtmed[3], par, 0);  // Cooper
181       sprintf(nameLChips,"GHL%d",ifmd+1);
182       gMC->Gsvolu(nameLChips,"TUBE", idtmed[5], par, 0); // Si chips
183       sprintf(nameGL10,"G1L%d",ifmd+1);
184       gMC->Gsvolu(nameGL10,"TUBE", idtmed[2], par, 0);  //G10 plate
185      
186      
187       par[0]=rin[ifmd]; // pipe size
188       par[1]=rout[ifmd];
189       par[2]=zSi/2;
190       gMC->Gsvolu(nameSi,"TUBE", idtmed[4], par, 3);
191       printf ("rin %f rout %f ZFMD %f\n",par[0],par[1],z[ifmd]);
192       zInside=ppcon[3]+par[2];
193       //      cout<<" zInside "<<zInside<<endl;
194       gMC->Gspos(nameSi,ifmd+1,nameFMD,0,0,zInside,0, "ONLY");
195       //PCB 1
196       zInside += par[2]+zPCB/2;
197       par[2]=zPCB/2;
198       //   cout<<"FMD "<<ifmd<<" zInside PCB1 = "<<zInside<<endl;
199       gMC->Gsposp(namePCB,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
200       zInside += zPCB;
201       //  cout<<"FMD "<<ifmd<<" zInside PCB2 = "<<zInside<<endl;
202       gMC->Gsposp(namePCB,2,nameFMD,0,0,zInside,0, "ONLY",par,3);
203       Float_t NulonTubeBegin=zInside+2*zPCB;
204       par[2]=zPCB/2-0.02;
205       Float_t zInPCB = -zPCB/2+par[2];
206       gMC->Gsposp(nameG10,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
207       zInPCB+=par[2]+zCooper/2 ;
208       par[2]=zCooper/2;
209       gMC->Gsposp(nameCopper,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
210       zInPCB += zCooper/2 + zChips/2;
211       par[2]=zChips/2;
212       gMC->Gsposp(nameChips,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
213       //HoneyComb
214       zHoneyComb=0.8;   
215       par[0] = RinHoneyComb[ifmd];
216       par[1] = RoutHoneyComb[ifmd];
217       par[2] = zHoneyComb/2;
218       //  cout<<"FMD "<<ifmd<<" zInside before HoneyComb = "<<zInside<<" par[2] "<<par[2]<<endl;
219       zInside += 2*NylonTube[2]+par[2];
220       // cout<<" zInside HoneyComb"<<zInside<<endl;
221       gMC->Gsposp(nameHoney,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
222       par[2]=0.1/2;
223       Float_t zHoney=-zHoneyComb/2+par[2];
224       //  cout<<" zHoney shurka 1 "<<zHoney <<endl;
225       gMC->Gsposp(nameHoneyOut,1,nameHoney,0,0,zHoney,0,
226                   "ONLY",par,3); //shkurki
227       zHoney=zHoneyComb/2-par[2];
228       //  cout<<" zHoney shurka 2 "<<zHoney <<endl;
229       gMC->Gsposp(nameHoneyOut,2,nameHoney,0,0,zHoney,0, "ONLY",par,3);
230       par[2]=(zHoneyComb-2.*0.1)/2; //soty vnutri
231       gMC->Gsposp(nameHoneyIn,1,nameHoney,0,0,0,0, "ONLY",par,3);
232       
233       gMC->Gspos("GNYL",1,nameFMD,0,yNylonTube[ifmd],
234                  NulonTubeBegin+NylonTube[2]/2.,0, "ONLY");
235       gMC->Gspos("GNYL",2,nameFMD,0,-yNylonTube[ifmd],
236                  NulonTubeBegin+NylonTube[2]/2.,0, "ONLY");
237          
238       //last PCB
239       par[0]=RoutHoneyComb[ifmd]-9;
240       par[1]=RoutHoneyComb[ifmd];
241       par[2]=zPCB/2;
242       //  cout<<" rin "<<par[0]<<" rout "<<par[1]<<endl;
243       zInside += zHoneyComb/2+par[2];
244       gMC->Gsposp(nameLPCB,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
245       
246        par[2]=zPCB/2-0.02;
247        zInPCB = -zPCB/2+par[2];
248        gMC->Gsposp(nameGL10,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
249        zInPCB+=par[2]+zCooper/2 ;
250        par[2]=zCooper/2;
251        gMC->Gsposp(nameLCopper,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
252        zInPCB += zCooper/2 + zChips/2;
253        par[2]=zChips/2;
254        gMC->Gsposp(nameLChips,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
255       
256        /*
257       if (ifmd==3) {
258         par[0]=0.1;
259         par[1]=7.;
260         par[2]=0.1;
261         gMC->Gsvolu("GST1","BOX ", idtmed[6], par, 3);  //support carbon stick
262         gMC->Gspos("GST1",1,nameFMD,0.,49.2-par[1],z[3]-0.7,0,"ONLY");
263         gMC->Gspos("GST1",2,nameFMD,0.,-49.2+par[1],z[3]-0.7,0,"ONLY");
264         
265         Float_t pSupportCone[5];
266         pSupportCone[0]=1.5;
267         pSupportCone[1]=34.4;
268         pSupportCone[2]=34.5;
269         pSupportCone[3]=30.;
270         pSupportCone[4]=29.9;
271         gMC->Gsvolu("GCY1","CONE", idtmed[6],pSupportCone , 5);  //support carbon stick
272         gMC->Gspos("GCY1",1,nameFMD,0.,0.,z[3]+1.5,0,"ONLY");
273         
274       }
275      if (ifmd==2) {
276         par[0]=0.1;
277         par[1]=1.3;
278         par[2]=0.1;
279         gMC->Gsvolu("GST3","BOX ", idtmed[6], par, 3);  //support carbon stick
280         gMC->Gspos("GST3",1,nameFMD,0.,23.-par[1],z[2]+2.1,0,"ONLY");
281         gMC->Gspos("GST3",2,nameFMD,0.,-23.+par[1],z[2]+2.1,0,"ONLY");
282
283         Float_t pSupportCone[5];
284         pSupportCone[0]=1.5;
285         pSupportCone[1]=30.;
286         pSupportCone[2]=29.9;
287         pSupportCone[3]=21.;
288         pSupportCone[4]=20.9;
289         gMC->Gsvolu("GCY2","CONE", idtmed[6],pSupportCone , 5);  //support carbon stick
290         gMC->Gspos("GCY2",1,nameFMD,0.,0.,z[2]+1.5,0,"ONLY");
291         
292      }
293        */
294        
295      //Granularity
296     fSectorsSi1=20;
297     fRingsSi1=256*3;
298     //fRingsSi1=3; // for drawing only
299     fSectorsSi2=40;
300     fRingsSi2=128*3;
301     // fRingsSi2=3; //for  drawing onl
302     if(ifmd==1||ifmd==3)
303       { 
304         gMC->Gsdvn(nameSector, nameSi , fSectorsSi2, 2);
305         gMC->Gsdvn(nameRing, nameSector, fRingsSi2, 1);
306       }
307     else
308       {
309         gMC->Gsdvn(nameSector, nameSi , fSectorsSi1, 2);
310         gMC->Gsdvn(nameRing, nameSector , fRingsSi1, 1);
311       }
312     
313     }
314 }    
315
316
317 //------------------------------------------------------------------------
318
319
320 void AliFMDv2::CreateMaterials()
321 {
322  Int_t isxfld   = gAlice->Field()->Integ();
323  Float_t sxmgmx = gAlice->Field()->Max();
324
325  // Plastic CH
326  Float_t aPlastic[2]={1.01,12.01};
327  Float_t zPlastic[2]={1,6};
328  Float_t wPlastic[2]={1,1};
329  Float_t denPlastic=1.03;
330  //     60% SiO2 , 40% G10FR4 
331  // PC board
332  Float_t apcb[3]  = { 28.0855,15.9994,17.749 };
333  Float_t zpcb[3]  = { 14.,8.,8.875 };
334  Float_t wpcb[3]  = { .28,.32,.4 };
335  Float_t denspcb  = 1.8;
336    //
337  
338  //*** Definition Of avaible FMD materials ***
339  AliMaterial(0, "FMD Air$", 14.61, 7.3, .001205, 30423.,999); 
340  AliMixture(1, "Plastic$",aPlastic,zPlastic,denPlastic,-2,wPlastic);
341  AliMixture(2, "SSD PCB$",   apcb, zpcb, denspcb, 3, wpcb);
342  AliMaterial(3, "SSD Copper$", 63.546, 29., 8.96, 1.43, 999.);
343  AliMaterial(4, "SSD Si$",      28.0855, 14., 2.33, 9.36, 999.);
344  AliMaterial(5, "SSD Si chip$", 28.0855, 14., 2.33, 9.36, 999.);
345  AliMaterial(6, "SSD C$",       12.011,   6., 2.265,18.8, 999.);
346  AliMaterial(7, "SSD Kapton$", 12.011, 6., 0.01, 31.27, 999.);//honeycomb
347   AliMaterial(8, "SSD G10FR4$", 17.749, 8.875, 1.8, 21.822, 999.);
348    
349
350 //**
351  AliMedium(0, "FMD air$", 0, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
352  AliMedium(1, "Plastic$", 1, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
353  AliMedium(2, "SSD PCB$", 2, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
354  AliMedium(3, "SSD Copper$", 3, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
355  AliMedium(4, "SSD Si$", 4, 1, isxfld, sxmgmx, 1., .001, 1., .001, .001);
356  AliMedium(5, "SSD Si chip$", 5, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
357  AliMedium(6, "SSD C$", 6, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
358  AliMedium(7, "SSD Kapton$", 7, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
359  AliMedium(8, "SSD G10FR4$", 8, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
360  
361
362
363 }
364
365 //---------------------------------------------------------------------
366 void AliFMDv2::DrawDetector()
367 {
368 //
369 // Draw a shaded view of the Forward multiplicity detector version 0
370 //
371
372
373 //Set ALIC mother transparent
374 gMC->Gsatt("ALIC","SEEN",0);
375 //
376 //Set volumes visible
377 gMC->Gsatt("FMD1","SEEN",1);
378 gMC->Gsatt("FMD2","SEEN",1);
379 gMC->Gsatt("FMD3","SEEN",1);
380 gMC->Gsatt("FMD4","SEEN",1);
381 gMC->Gsatt("FMD5","SEEN",1);
382
383 //
384 gMC->Gdopt("hide","on");
385 gMC->Gdopt("shad","on");
386 gMC->SetClipBox(".");
387 gMC->SetClipBox("*",0,1000,-1000,1000,-1000,1000);
388 gMC->DefaultRange();
389 gMC->Gdraw("alic",40,30,0,12,9.5,.2,0.2);
390 gMC->Gdhead(1111,"Forward multiplicity detector");
391 gMC->Gdopt("hide","off");
392 }
393 //-------------------------------------------------------------------
394 void AliFMDv2::Init()
395 {
396 // Initialises version 0 of the Forward Multiplicity Detector
397 //
398 AliFMD::Init();
399 fIdSens1=gMC->VolId("GRN1");
400 fIdSens2=gMC->VolId("GRN2");
401 fIdSens3=gMC->VolId("GRN3");
402 fIdSens4=gMC->VolId("GRN4");
403 fIdSens5=gMC->VolId("GRN5");
404 printf("*** FMD version 2 initialized ***\n");
405 }
406
407 //-------------------------------------------------------------------
408
409 void AliFMDv2::StepManager()
410 {
411   //
412   // Called for every step in the Forward Multiplicity Detector
413   //
414   Int_t id,copy,copy1,copy2;
415   static Float_t hits[9];
416   static Int_t vol[3];
417   static Float_t de;
418   TLorentzVector pos;
419   TLorentzVector mom;
420
421
422   TClonesArray &lhits = *fHits;
423   if(!gMC->IsTrackAlive()) return; // particle has disappeared
424
425   Float_t charge = gMC->TrackCharge();
426   if(TMath::Abs(charge)<=0.) return; //take only charged particles
427
428   //  printf(" in StepManeger \n");
429   id=gMC->CurrentVolID(copy);
430   
431 // Check the sensetive volume
432    if(id==fIdSens1||id==fIdSens2||id==fIdSens3||id==fIdSens4||id==fIdSens5)
433      {
434        if(gMC->IsTrackEntering())
435          {
436            vol[2]=copy;
437            gMC->CurrentVolOffID(1,copy1);
438            vol[1]=copy1;
439            gMC->CurrentVolOffID(2,copy2);
440            vol[0]=copy2;
441
442            gMC->TrackPosition(pos);
443            hits[0]=pos[0];
444            hits[1]=pos[1];
445            hits[2]=pos[2];
446
447            gMC->TrackMomentum(mom);
448            hits[3]=mom[0];
449            hits[4]=mom[1];
450            hits[5]=mom[2];
451
452            Int_t iPart= gMC->TrackPid();
453            Int_t partId=gMC->IdFromPDG(iPart);
454            hits[7]=partId;
455            hits[8]=1e9*gMC->TrackTime();
456            de=0.;
457          }
458        if(gMC->IsTrackInside()){
459            de=de+1000.*gMC->Edep();
460        }
461        
462        if(gMC->IsTrackExiting()
463           ||gMC->IsTrackDisappeared()||
464           gMC->IsTrackStop())
465          {
466              hits[6]=de+1000.*gMC->Edep();
467       new(lhits[fNhits++]) AliFMDhit(fIshunt,gAlice->CurrentTrack(),vol,hits);
468          } // IsTrackExiting()
469      }
470   }
471
472 //--------------------------------------------------------------------------
473
474 void AliFMDv2::Response( Float_t Edep)
475 {
476   Float_t I=1.664*0.04*2.33/22400; // = 0.69e-6;
477   Float_t chargeOnly=Edep/I;
478   //Add noise ~500electrons
479   Int_t charge=500;
480   if (Edep>0)
481      charge=Int_t(gRandom->Gaus(chargeOnly,500));       
482  }   
483
484
485
486
487
488
489