]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDv1.cxx
volume overlap fixed
[u/mrichter/AliRoot.git] / FMD / AliFMDv1.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
28 #include <TMath.h>
29 #include <TGeometry.h>
30 #include <TTUBE.h>
31 #include <TFile.h>
32 #include <TTree.h>
33 #include <TNode.h>
34 #include <TClonesArray.h>
35 #include <TLorentzVector.h>
36 #include <TDirectory.h>
37 #include "AliFMDv1.h"
38 #include "AliFMDv0.h"
39 #include "AliRun.h"
40 #include <Riostream.h>
41 #include "AliMagF.h"
42 #include "AliFMDhit.h"
43 #include "AliFMDdigit.h"
44 #include <stdlib.h>
45 //#include "TGeant3.h"
46 //class TGeant3;
47 ClassImp(AliFMDv1)
48
49 //--------------------------------------------------------------------
50 AliFMDv1::AliFMDv1(const char *name, const char *title):
51  AliFMD(name,title)
52 {
53   //
54   // Standart constructor for Forward Multiplicity Detector version 0
55   //
56   fIdSens1=0;
57   fIdSens2=0;
58   fIdSens3=0;
59   fIdSens4=0;
60   fIdSens5=0;
61 //  setBufferSize(128000);
62 }
63 //-------------------------------------------------------------------------
64 void AliFMDv1::CreateGeometry()
65 {
66  //
67   // Create the geometry of Forward Multiplicity Detector version 0
68   //
69   //Detector consists of 6 volumes: 
70   // 1st covered pseudorapidity interval from 3.3 to 2.0
71   // and placed on 65cm in Z-direction;
72   // 2nd - from 2.0 to 1.6 and Z=85 cm;
73   // 3d  - the same pseudorapidity interval as the 1st 
74   // but on the other side from the interaction point z=-65cm;
75   // 4th - simmetricaly with the 2nd : 
76   // pseudorapidity from 2.0 to 1.6, Z=-85cm   
77   // 5th - from 3.6 to 4.7, Z=-270cm
78   // 6th - from 4.5 to 5.5 , Z=-630cm.
79   // Each part has 400mkm Si (sensetive area, detector itself),
80   // 0.75cm of plastic simulated electronics material,
81   // Al support ring 2cm thickness and 1cm width placed on 
82   // the outer radius of each Si disk;
83   //    
84   // begin Html
85   /*
86    <img src="gif/AliFMDv0.gif">
87    */
88   //
89
90   
91   Int_t *idtmed = fIdtmed->GetArray();
92    
93   Int_t ifmd;
94   Int_t idrotm[999];
95   Float_t zfmd,par[3];
96   char name[5], nameSi[5], nameSector[5], nameRing[5];
97
98   Float_t rin[6], rout[6],zpos;
99
100   Float_t etain[5]= {3.40, 2.29, 3.68, 2.29, 5.09};
101   Float_t etaout[6]={2.01, 1.70, 2.28, 1.70, 3.68};
102   //  Float_t z[6]={64., 85., -64., -85., -270., -630};
103   Float_t z[6]={62.8, 75.2, -83.4, -75.2, -340.};
104   Float_t zDet=0.03;
105   Float_t zElectronic=0.1;
106   Float_t zSupport=1.;
107
108   Float_t zFMD=1.;
109 //-------------------------------------------------------------------
110  //  FMD 
111  //------------------------------------------------------------------
112         cout<<" !!!!!!!!!!!New FMD geometry !!!!!!!!!"<<endl;
113
114   AliMatrix(idrotm[901], 90, 0, 90, 90, 180, 0);
115
116   //  gMC->Gsvolu("GSI","TUBE", idtmed[1], par, 0);
117   gMC->Gsvolu("GEL ","TUBE", idtmed[4], par, 0);
118   gMC->Gsvolu("GSUP","TUBE", idtmed[2], par, 0);
119
120   for (ifmd =0; ifmd < 5; ifmd++){
121
122     sprintf(name,"FMD%d",ifmd+1);
123     sprintf(nameSi,"GSI%d",ifmd+1);
124     sprintf(nameSector,"GSC%d",ifmd+1);
125     sprintf(nameRing,"GRN%d",ifmd+1);
126
127     zfmd=TMath::Abs(z[ifmd]);
128
129     AliFMD::Eta2Radius(etain[ifmd],zfmd,&rin[ifmd]);
130     AliFMD::Eta2Radius(etaout[ifmd],zfmd,&rout[ifmd]);
131     
132     par[0]=rin[ifmd]; // pipe size
133     par[1]=rout[ifmd];
134     par[2]=zFMD/2;
135     gMC->Gsvolu(name,"TUBE", idtmed[3], par, 3);
136     gMC->Gsvolu(nameSi,"TUBE", idtmed[1], par, 0);
137     
138     if (z[ifmd] < 0){  
139       gMC->Gspos(name,1,"ALIC",0,0,z[ifmd],0, "ONLY");}
140     else { 
141       gMC->Gspos(name,1,"ALIC",0,0,z[ifmd],idrotm[901], "ONLY");}
142   //Silicon detector
143     par[2]=zDet/2;
144     zpos=zFMD/2 -par[2];
145     gMC->Gsposp(nameSi,ifmd+1,name,0,0,zpos,0, "ONLY",par,3);
146
147     //Granularity
148     fSectorsSi1=20;
149     fRingsSi1=256;
150     fSectorsSi2=40;
151     fRingsSi2=128;
152    if(ifmd==1||ifmd==3)
153       { 
154         gMC->Gsdvn(nameSector, nameSi , fSectorsSi2, 2);
155         gMC->Gsdvn(nameRing, nameSector, fRingsSi2, 1);
156       }
157     else
158       {
159         gMC->Gsdvn(nameSector, nameSi , fSectorsSi1, 2);
160         gMC->Gsdvn(nameRing, nameSector , fRingsSi1, 1);
161       }
162
163     //Plastic slice for electronics
164     par[2]=zElectronic/2;
165     zpos=zpos-zDet/2-par[2];
166     gMC->Gsposp("GEL ",ifmd+1,name,0,0,zpos,0, "ONLY",par,3);
167
168    //Simple Al support
169    par[1]=rout[ifmd];
170    par[0]=rout[ifmd]-2;
171    par[2]=zSupport/2;
172    zpos=zpos-zElectronic/2-par[2];
173    //   gMC->Gsposp("GSUP",ifmd+1,name,0,0,zpos,0, "ONLY",par,3);
174    
175
176   }  
177
178 }    
179
180 //------------------------------------------------------------------------
181 void AliFMDv1::CreateMaterials()
182 {
183  Int_t isxfld   = gAlice->Field()->Integ();
184  Float_t sxmgmx = gAlice->Field()->Max();
185
186  // Plastic CH
187  Float_t aPlastic[2]={1.01,12.01};
188  Float_t zPlastic[2]={1,6};
189  Float_t wPlastic[2]={1,1};
190  Float_t denPlastic=1.03;
191    //
192  
193  //*** Definition Of avaible FMD materials ***
194  AliMaterial(0, "Si chip$", 28.0855,14.,2.33,9.36,999);
195  AliMaterial(1, "Al supprt$", 26.980,13.,2.70,8.9,999);
196  AliMaterial(2, "FMD Air$", 14.61, 7.3, .001205, 30423.,999); 
197  AliMixture( 5, "Plastic$",aPlastic,zPlastic,denPlastic,-2,wPlastic);
198  
199
200 //**
201  AliMedium(1, "Si chip$", 0, 1, isxfld, sxmgmx, 1., .001, 1., .001, .001);
202  AliMedium(2, "Al support$", 1, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
203  AliMedium(3, "FMD air$", 2, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
204  AliMedium(4, "Plastic$", 5, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
205  
206
207
208 }
209 //---------------------------------------------------------------------
210 void AliFMDv1::DrawDetector()
211 {
212 //
213 // Draw a shaded view of the Forward multiplicity detector version 0
214 //
215
216
217 //Set ALIC mother transparent
218 gMC->Gsatt("ALIC","SEEN",0);
219 //
220 //Set volumes visible
221 gMC->Gsatt("FMD1","SEEN",1);
222 gMC->Gsatt("FMD2","SEEN",1);
223 gMC->Gsatt("FMD3","SEEN",1);
224 gMC->Gsatt("FMD4","SEEN",1);
225 gMC->Gsatt("FMD5","SEEN",1);
226
227 //
228 gMC->Gdopt("hide","on");
229 gMC->Gdopt("shad","on");
230 gMC->SetClipBox(".");
231 gMC->SetClipBox("*",0,1000,-1000,1000,-1000,1000);
232 gMC->DefaultRange();
233 gMC->Gdraw("alic",40,30,0,12,9.5,.2,0.2);
234 gMC->Gdhead(1111,"Forward multiplicity detector");
235 gMC->Gdopt("hide","off");
236 }
237 //-------------------------------------------------------------------
238 void AliFMDv1::Init()
239 {
240 // Initialises version 0 of the Forward Multiplicity Detector
241 //
242 AliFMD::Init();
243 fIdSens1=gMC->VolId("GRN1");
244 fIdSens2=gMC->VolId("GRN2");
245 fIdSens3=gMC->VolId("GRN3");
246 fIdSens4=gMC->VolId("GRN4");
247 fIdSens5=gMC->VolId("GRN5");
248 printf("*** FMD version 1 initialized ***\n");
249 }
250
251 //-------------------------------------------------------------------
252
253 void AliFMDv1::StepManager()
254 {
255   //
256   // Called for every step in the Forward Multiplicity Detector
257   //
258   Int_t id,copy,copy1,copy2;
259   static Float_t hits[9];
260   static Int_t vol[3];
261   static Float_t de;
262   TLorentzVector pos;
263   TLorentzVector mom;
264
265
266   TClonesArray &lhits = *fHits;
267   if(!gMC->IsTrackAlive()) return; // particle has disappeared
268
269   Float_t charge = gMC->TrackCharge();
270   if(TMath::Abs(charge)<=0.) return; //take only charged particles
271
272   //  printf(" in StepManeger \n");
273   id=gMC->CurrentVolID(copy);
274   //((TGeant3*)gMC)->Gpcxyz();
275   
276 // Check the sensetive volume
277    if(id==fIdSens1||id==fIdSens2||id==fIdSens3||id==fIdSens4||id==fIdSens5)
278      {
279        if(gMC->IsTrackEntering())
280          {
281            vol[2]=copy;
282            gMC->CurrentVolOffID(1,copy1);
283            vol[1]=copy1;
284            gMC->CurrentVolOffID(2,copy2);
285            vol[0]=copy2;
286
287            gMC->TrackPosition(pos);
288            hits[0]=pos[0];
289            hits[1]=pos[1];
290            hits[2]=pos[2];
291
292            gMC->TrackMomentum(mom);
293            hits[3]=mom[0];
294            hits[4]=mom[1];
295            hits[5]=mom[2];
296
297            Int_t iPart= gMC->TrackPid();
298            Int_t partId=gMC->IdFromPDG(iPart);
299            hits[7]=partId;
300            hits[8]=1e9*gMC->TrackTime();
301            de=0.;
302          }
303        if(gMC->IsTrackInside()){
304            de=de+1000.*gMC->Edep();
305        }
306        
307        if(gMC->IsTrackExiting()
308           ||gMC->IsTrackDisappeared()||
309           gMC->IsTrackStop())
310          {
311              hits[6]=de+1000.*gMC->Edep();
312       new(lhits[fNhits++]) AliFMDhit(fIshunt,gAlice->CurrentTrack(),vol,hits);
313          } // IsTrackExiting()
314      }
315   }
316 //--------------------------------------------------------------------------
317
318 void AliFMDv1::Response( Float_t Edep)
319 {
320   Float_t I=1.664*0.04*2.33/22400; // = 0.69e-6;
321   Float_t chargeOnly=Edep/I;
322   //Add noise ~500electrons
323   Int_t charge=500;
324   if (Edep>0)
325      charge=Int_t(gRandom->Gaus(chargeOnly,500));       
326  }   
327
328
329
330
331
332