Transition to NewIO
[u/mrichter/AliRoot.git] / FMD / AliFMDv1.cxx
... / ...
CommitLineData
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/* $Id$ */
17
18 /////////////////////////////////////////////////////////////////////
19// //
20// Forward Multiplicity detector based on Silicon version 0 //
21//
22//Begin Html
23/*
24<img src="gif/AliFMDv0Class.gif">
25*/
26//End Html
27// //
28// //
29//////////////////////////////////////////////////////////////////////
30
31#include <Riostream.h>
32#include <stdlib.h>
33
34#include <TClonesArray.h>
35#include <TDirectory.h>
36#include <TFile.h>
37#include <TGeometry.h>
38#include <TLorentzVector.h>
39#include <TMath.h>
40#include <TNode.h>
41#include <TTUBE.h>
42#include <TTree.h>
43#include <TVirtualMC.h>
44
45#include "AliFMDdigit.h"
46#include "AliFMDhit.h"
47#include "AliFMDv0.h"
48#include "AliFMDv1.h"
49#include "AliMagF.h"
50#include "AliRun.h"
51
52ClassImp(AliFMDv1)
53
54//--------------------------------------------------------------------
55AliFMDv1::AliFMDv1(const char *name, const char *title):
56 AliFMD(name,title)
57{
58 //
59 // Standart constructor for Forward Multiplicity Detector version 0
60 //
61 fIdSens1=0;
62 fIdSens2=0;
63 fIdSens3=0;
64 fIdSens4=0;
65 fIdSens5=0;
66// setBufferSize(128000);
67}
68//-------------------------------------------------------------------------
69void AliFMDv1::CreateGeometry()
70{
71 //
72 // Create the geometry of Forward Multiplicity Detector version 0
73 //
74 //Detector consists of 6 volumes:
75 // 1st covered pseudorapidity interval from 3.3 to 2.0
76 // and placed on 65cm in Z-direction;
77 // 2nd - from 2.0 to 1.6 and Z=85 cm;
78 // 3d - the same pseudorapidity interval as the 1st
79 // but on the other side from the interaction point z=-65cm;
80 // 4th - simmetricaly with the 2nd :
81 // pseudorapidity from 2.0 to 1.6, Z=-85cm
82 // 5th - from 3.6 to 4.7, Z=-270cm
83 // 6th - from 4.5 to 5.5 , Z=-630cm.
84 // Each part has 400mkm Si (sensetive area, detector itself),
85 // 0.75cm of plastic simulated electronics material,
86 // Al support ring 2cm thickness and 1cm width placed on
87 // the outer radius of each Si disk;
88 //
89 // begin Html
90 /*
91 <img src="gif/AliFMDv0.gif">
92 */
93 //
94
95
96 Int_t *idtmed = fIdtmed->GetArray();
97
98 Int_t ifmd;
99 Int_t idrotm[999];
100 Float_t zfmd,par[3];
101 char name[5], nameSi[5], nameSector[5], nameRing[5];
102
103 Float_t rin[6], rout[6],zpos;
104
105 Float_t etain[5]= {3.40, 2.29, 3.68, 2.29, 5.09};
106 Float_t etaout[6]={2.01, 1.70, 2.28, 1.70, 3.68};
107 // Float_t z[6]={64., 85., -64., -85., -270., -630};
108 Float_t z[6]={62.8, 75.2, -83.4, -75.2, -340.};
109 Float_t zDet=0.03;
110 Float_t zElectronic=0.1;
111 Float_t zSupport=1.;
112
113 Float_t zFMD=1.;
114//-------------------------------------------------------------------
115 // FMD
116 //------------------------------------------------------------------
117 cout<<" !!!!!!!!!!!New FMD geometry !!!!!!!!!"<<endl;
118
119 AliMatrix(idrotm[901], 90, 0, 90, 90, 180, 0);
120
121 // gMC->Gsvolu("GSI","TUBE", idtmed[1], par, 0);
122 gMC->Gsvolu("GEL ","TUBE", idtmed[4], par, 0);
123 gMC->Gsvolu("GSUP","TUBE", idtmed[2], par, 0);
124
125 for (ifmd =0; ifmd < 5; ifmd++){
126
127 sprintf(name,"FMD%d",ifmd+1);
128 sprintf(nameSi,"GSI%d",ifmd+1);
129 sprintf(nameSector,"GSC%d",ifmd+1);
130 sprintf(nameRing,"GRN%d",ifmd+1);
131
132 zfmd=TMath::Abs(z[ifmd]);
133
134 AliFMD::Eta2Radius(etain[ifmd],zfmd,&rin[ifmd]);
135 AliFMD::Eta2Radius(etaout[ifmd],zfmd,&rout[ifmd]);
136
137 par[0]=rin[ifmd]; // pipe size
138 par[1]=rout[ifmd];
139 par[2]=zFMD/2;
140 gMC->Gsvolu(name,"TUBE", idtmed[3], par, 3);
141 gMC->Gsvolu(nameSi,"TUBE", idtmed[1], par, 0);
142
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
152 //Granularity
153 fSectorsSi1=20;
154 fRingsSi1=256;
155 fSectorsSi2=40;
156 fRingsSi2=128;
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//------------------------------------------------------------------------
186void 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//---------------------------------------------------------------------
215void AliFMDv1::DrawDetector()
216{
217//
218// Draw a shaded view of the Forward multiplicity detector version 0
219//
220
221//Set ALIC mother transparent
222gMC->Gsatt("ALIC","SEEN",0);
223//
224//Set volumes visible
225gMC->Gsatt("FMD1","SEEN",1);
226gMC->Gsatt("FMD2","SEEN",1);
227gMC->Gsatt("FMD3","SEEN",1);
228gMC->Gsatt("FMD4","SEEN",1);
229gMC->Gsatt("FMD5","SEEN",1);
230
231//
232gMC->Gdopt("hide","on");
233gMC->Gdopt("shad","on");
234gMC->SetClipBox(".");
235gMC->SetClipBox("*",0,1000,-1000,1000,-1000,1000);
236gMC->DefaultRange();
237gMC->Gdraw("alic",40,30,0,12,9.5,.2,0.2);
238gMC->Gdhead(1111,"Forward multiplicity detector");
239gMC->Gdopt("hide","off");
240}
241//-------------------------------------------------------------------
242void AliFMDv1::Init()
243{
244// Initialises version 0 of the Forward Multiplicity Detector
245//
246AliFMD::Init();
247fIdSens1=gMC->VolId("GRN1");
248fIdSens2=gMC->VolId("GRN2");
249fIdSens3=gMC->VolId("GRN3");
250fIdSens4=gMC->VolId("GRN4");
251fIdSens5=gMC->VolId("GRN5");
252printf("*** FMD version 1 initialized ***\n");
253}
254
255//-------------------------------------------------------------------
256
257void 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 if(!gMC->IsTrackAlive()) return; // particle has disappeared
272
273 Float_t charge = gMC->TrackCharge();
274 if(TMath::Abs(charge)<=0.) return; //take only charged particles
275
276 // printf(" in StepManeger \n");
277 id=gMC->CurrentVolID(copy);
278 //((TGeant3*)gMC)->Gpcxyz();
279
280// Check the sensetive volume
281 if(id==fIdSens1||id==fIdSens2||id==fIdSens3||id==fIdSens4||id==fIdSens5)
282 {
283 if(gMC->IsTrackEntering())
284 {
285 vol[2]=copy;
286 gMC->CurrentVolOffID(1,copy1);
287 vol[1]=copy1;
288 gMC->CurrentVolOffID(2,copy2);
289 vol[0]=copy2;
290
291 gMC->TrackPosition(pos);
292 hits[0]=pos[0];
293 hits[1]=pos[1];
294 hits[2]=pos[2];
295
296 gMC->TrackMomentum(mom);
297 hits[3]=mom[0];
298 hits[4]=mom[1];
299 hits[5]=mom[2];
300
301 Int_t iPart= gMC->TrackPid();
302 Int_t partId=gMC->IdFromPDG(iPart);
303 hits[7]=partId;
304 hits[8]=1e9*gMC->TrackTime();
305 de=0.;
306 }
307 if(gMC->IsTrackInside()){
308 de=de+1000.*gMC->Edep();
309 }
310
311 if(gMC->IsTrackExiting()
312 ||gMC->IsTrackDisappeared()||
313 gMC->IsTrackStop())
314 {
315 hits[6]=de+1000.*gMC->Edep();
316 new(lhits[fNhits++]) AliFMDhit(fIshunt,gAlice->CurrentTrack(),vol,hits);
317 } // IsTrackExiting()
318 }
319 }
320//--------------------------------------------------------------------------
321
322void AliFMDv1::Response( Float_t Edep)
323{
324 Float_t I=1.664*0.04*2.33/22400; // = 0.69e-6;
325 Float_t chargeOnly=Edep/I;
326 //Add noise ~500electrons
327 Int_t charge=500;
328 if (Edep>0)
329 charge=Int_t(gRandom->Gaus(chargeOnly,500));
330 }
331
332
333
334
335
336