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