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