]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDv1.cxx
new FMD geometry
[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 <iostream.h>
42 #include <fstream.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     printf(name,nameSi);
129     
130     zfmd=TMath::Abs(z[ifmd]);
131     printf("zfmd %f z[ifmd] %f",zfmd,z[ifmd]);
132     AliFMD::Eta2Radius(etain[ifmd],zfmd,&rin[ifmd]);
133     AliFMD::Eta2Radius(etaout[ifmd],zfmd,&rout[ifmd]);
134     
135     par[0]=rin[ifmd]; // pipe size
136     par[1]=rout[ifmd];
137     par[2]=zFMD/2;
138     gMC->Gsvolu(name,"TUBE", idtmed[3], par, 3);
139     //    par[2]=zDet/2;
140     gMC->Gsvolu(nameSi,"TUBE", idtmed[1], par, 0);
141     
142     printf ("rin %f rout %f ZFMD %f\n",par[0],par[1],z[ifmd]);
143     if (z[ifmd] < 0){  
144       gMC->Gspos(name,1,"ALIC",0,0,z[ifmd],0, "ONLY");}
145     else { 
146       gMC->Gspos(name,1,"ALIC",0,0,z[ifmd],idrotm[901], "ONLY");}
147   //Silicon detector
148     par[2]=zDet/2;
149     zpos=zFMD/2 -par[2];
150     gMC->Gsposp(nameSi,ifmd+1,name,0,0,zpos,0, "ONLY",par,3);
151     cout<<" Si "<<nameSi<<" ifmd "<<ifmd<<" rin "<<par[0]<<" rout "<<par[1]<<
152       " zDet "<<par[2]<<endl;
153     //Granularity
154     cout<<"fSectorsSi1 "<<fSectorsSi1<<
155       " fRingsSi1 "<<fRingsSi1<<
156       " fSectorsSi2 "<<fSectorsSi2<<
157       " fRingsSi2 "<<fRingsSi2<<endl;
158    if(ifmd==1||ifmd==3)
159       { 
160         gMC->Gsdvn(nameSector, nameSi , fSectorsSi2, 2);
161         gMC->Gsdvn(nameRing, nameSector, fRingsSi2, 1);
162       }
163     else
164       {
165         gMC->Gsdvn(nameSector, nameSi , fSectorsSi1, 2);
166         gMC->Gsdvn(nameRing, nameSector , fRingsSi1, 1);
167       }
168
169     //Plastic slice for electronics
170     par[2]=zElectronic/2;
171     zpos=zpos-zDet/2-par[2];
172     gMC->Gsposp("GEL ",ifmd+1,name,0,0,zpos,0, "ONLY",par,3);
173
174    //Simple Al support
175    par[1]=rout[ifmd];
176    par[0]=rout[ifmd]-2;
177    par[2]=zSupport/2;
178    zpos=zpos-zElectronic/2-par[2];
179    //   gMC->Gsposp("GSUP",ifmd+1,name,0,0,zpos,0, "ONLY",par,3);
180    
181
182   }  
183
184 }    
185
186 //------------------------------------------------------------------------
187 void AliFMDv1::CreateMaterials()
188 {
189  Int_t isxfld   = gAlice->Field()->Integ();
190  Float_t sxmgmx = gAlice->Field()->Max();
191
192  // Plastic CH
193  Float_t aPlastic[2]={1.01,12.01};
194  Float_t zPlastic[2]={1,6};
195  Float_t wPlastic[2]={1,1};
196  Float_t denPlastic=1.03;
197    //
198  
199  //*** Definition Of avaible FMD materials ***
200  AliMaterial(0, "Si chip$", 28.0855,14.,2.33,9.36,999);
201  AliMaterial(1, "Al supprt$", 26.980,13.,2.70,8.9,999);
202  AliMaterial(2, "FMD Air$", 14.61, 7.3, .001205, 30423.,999); 
203  AliMixture( 5, "Plastic$",aPlastic,zPlastic,denPlastic,-2,wPlastic);
204  
205
206 //**
207  AliMedium(1, "Si chip$", 0, 1, isxfld, sxmgmx, 1., .001, 1., .001, .001);
208  AliMedium(2, "Al support$", 1, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
209  AliMedium(3, "FMD air$", 2, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
210  AliMedium(4, "Plastic$", 5, 0,isxfld, sxmgmx,  10., .01, 1., .003, .003);
211  
212
213
214 }
215 //---------------------------------------------------------------------
216 void AliFMDv1::DrawDetector()
217 {
218 //
219 // Draw a shaded view of the Forward multiplicity detector version 0
220 //
221
222 AliMC* pMC = AliMC::GetMC();
223
224 //Set ALIC mother transparent
225 pMC->Gsatt("ALIC","SEEN",0);
226 //
227 //Set volumes visible
228 gMC->Gsatt("FMD1","SEEN",1);
229 gMC->Gsatt("FMD2","SEEN",1);
230 gMC->Gsatt("FMD3","SEEN",1);
231 gMC->Gsatt("FMD4","SEEN",1);
232 gMC->Gsatt("FMD5","SEEN",1);
233
234 //
235 gMC->Gdopt("hide","on");
236 gMC->Gdopt("shad","on");
237 gMC->SetClipBox(".");
238 gMC->SetClipBox("*",0,1000,-1000,1000,-1000,1000);
239 gMC->DefaultRange();
240 gMC->Gdraw("alic",40,30,0,12,9.5,.2,0.2);
241 gMC->Gdhead(1111,"Forward multiplicity detector");
242 gMC->Gdopt("hide","off");
243 }
244 //-------------------------------------------------------------------
245 void AliFMDv1::Init()
246 {
247 // Initialises version 0 of the Forward Multiplicity Detector
248 //
249 AliMC* gMC=AliMC::GetMC();
250 AliFMD::Init();
251 fIdSens1=gMC->VolId("GRN1");
252 fIdSens2=gMC->VolId("GRN2");
253 fIdSens3=gMC->VolId("GRN3");
254 fIdSens4=gMC->VolId("GRN4");
255 fIdSens5=gMC->VolId("GRN5");
256 printf("*** FMD version 1 initialized ***\n");
257 }
258
259 //-------------------------------------------------------------------
260
261 void AliFMDv1::StepManager()
262 {
263   //
264   // Called for every step in the Forward Multiplicity Detector
265   //
266   Int_t id,copy,copy1,copy2;
267   static Float_t hits[9];
268   static Int_t vol[3];
269   static Float_t de;
270   TLorentzVector pos;
271   TLorentzVector mom;
272
273
274   TClonesArray &lhits = *fHits;
275   AliMC* gMC=AliMC::GetMC();
276   if(!gMC->IsTrackAlive()) return; // particle has disappeared
277
278   Float_t charge = gMC->TrackCharge();
279   if(TMath::Abs(charge)<=0.) return; //take only charged particles
280
281   //  printf(" in StepManeger \n");
282   id=gMC->CurrentVolID(copy);
283   //((TGeant3*)gMC)->Gpcxyz();
284   
285 // Check the sensetive volume
286    if(id==fIdSens1||id==fIdSens2||id==fIdSens3||id==fIdSens4||id==fIdSens5)
287      {
288        if(gMC->IsTrackEntering())
289          {
290            vol[2]=copy;
291            gMC->CurrentVolOffID(1,copy1);
292            vol[1]=copy1;
293            gMC->CurrentVolOffID(2,copy2);
294            vol[0]=copy2;
295            //  printf("vol0 %d vol1 %d vol2 %d\n",vol[0],vol[1],vol[2]); 
296            gMC->TrackPosition(pos);
297            hits[0]=pos[0];
298            hits[1]=pos[1];
299            hits[2]=pos[2];
300            //      printf(" Zpos %f \n",hits[2]);
301            gMC->TrackMomentum(mom);
302            hits[3]=mom[0];
303            hits[4]=mom[1];
304            hits[5]=mom[2];
305
306            Int_t iPart= gMC->TrackPid();
307            Int_t partId=gMC->IdFromPDG(iPart);
308            hits[7]=partId;
309            hits[8]=1e9*gMC->TrackTime();
310            de=0.;
311          }
312        if(gMC->IsTrackInside()){
313            de=de+1000.*gMC->Edep();
314            //      cout<<" de "<<de<<endl;
315            // cout<<" step "<<gMC->TrackStep()<<endl;
316        }
317        
318        if(gMC->IsTrackExiting()
319           ||gMC->IsTrackDisappeared()||
320           gMC->IsTrackStop())
321          {
322              hits[6]=de+1000.*gMC->Edep();
323              //     cout<<" idSens "<<id<<endl;
324              //cout<<" hits "<<hits[6]<<endl;
325              //     for(Int_t i=0; i<9; i++) printf(" i %d hits %f \n",i,hits[i]);
326       new(lhits[fNhits++]) AliFMDhit(fIshunt,gAlice->CurrentTrack(),vol,hits);
327          } // IsTrackExiting()
328      }
329   }
330 //--------------------------------------------------------------------------
331
332 void AliFMDv1::Response( Float_t Edep)
333 {
334   Float_t I=1.664*0.04*2.33/22400; // = 0.69e-6;
335   Float_t chargeOnly=Edep/I;
336   //Add noise ~500electrons
337   Int_t charge=500;
338   if (Edep>0)
339      charge=Int_t(gRandom->Gaus(chargeOnly,500));       
340  }   
341
342
343
344
345
346