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