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