Corrections to obey the coding conventions
[u/mrichter/AliRoot.git] / FMD / AliFMDv2.cxx
CommitLineData
0d5ea2e8 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#include <TMath.h>
28#include <TGeometry.h>
29#include <TTUBE.h>
30#include <TFile.h>
31#include <TTree.h>
32#include <TNode.h>
33#include <TClonesArray.h>
34#include <TLorentzVector.h>
35#include <TDirectory.h>
36#include "AliFMDv2.h"
37#include "AliFMDv0.h"
38#include "AliRun.h"
39#include <Riostream.h>
40#include "AliMagF.h"
41#include "AliFMDhit.h"
42#include "AliFMDdigit.h"
43#include <stdlib.h>
44
45ClassImp(AliFMDv2)
46
47//--------------------------------------------------------------------
48AliFMDv2::AliFMDv2(const char *name, const char *title):
49 AliFMD(name,title)
50{
51 //
52 // Standart constructor for Forward Multiplicity Detector version 0
53 //
54 fIdSens1=0;
55 fIdSens2=0;
56 fIdSens3=0;
57 fIdSens4=0;
58 fIdSens5=0;
59// setBufferSize(128000);
60}
61//-------------------------------------------------------------------------
62void AliFMDv2::CreateGeometry()
63{
64 //
65 // Create the geometry of Forward Multiplicity Detector version 0
66 //
67 //Detector consists of 6 volumes:
68 // 1st covered pseudorapidity interval from 3.3 to 2.0
69 // and placed on 65cm in Z-direction;
70 // 2nd - from 2.0 to 1.6 and Z=85 cm;
71 // 3d - the same pseudorapidity interval as the 1st
72 // but on the other side from the interaction point z=-65cm;
73 // 4th - simmetricaly with the 2nd :
74 // pseudorapidity from 2.0 to 1.6, Z=-85cm
75 // 5th - from 3.6 to 4.7, Z=-270cm
76 // 6th - from 4.5 to 5.5 , Z=-630cm.
77 // Each part has 400mkm Si (sensetive area, detector itself),
78 // 0.75cm of plastic simulated electronics material,
79 // Al support ring 2cm thickness and 1cm width placed on
80 // the outer radius of each Si disk;
81 //
82 // begin Html
83 /*
84 <img src="gif/AliFMDv0.gif">
85 */
86 //
87
88 Int_t *idtmed = fIdtmed->GetArray();
89
90 Int_t ifmd;
91 Int_t idrotm[999];
92 Float_t zFMD,par[3],ppcon[5];
93 //ppcon1[15],ppcon2[15],ppcon3[15],ppcon4[15],ppcon5[15];
94 Float_t z[5]={62.8, 75.2, -83.4, -75.2, -340.};
95 Float_t NylonTube[3]={0.2,0.6,0.45};
96 Float_t zPCB=0.12; Float_t zHoneyComb=0.5;
97 Float_t zSi=0.03;
98
99//-------------------------------------------------------------------
100 // FMD
101 //------------------------------------------------------------------
102 // cout<<" !!!!!!!!!!!New FMD geometry !!!!!!!!!"<<endl;
103
104
105 char nameFMD[5], nameSi[5], nameSector[5], nameRing[5];
106 Char_t nameHoney[5], nameHoneyIn[5], nameHoneyOut[5];
107 Char_t namePCB[5], nameCopper[5], nameChips[5], nameG10[5];
108 Char_t nameLPCB[5], nameLCopper[5], nameLChips[5], nameGL10[5];;
109 Float_t rin[5]={4.2,15.4,4.2,15.4,4.2};
110 Float_t rout[5]={17.4,28.4,17.4,28.4,17.4};
111 Float_t RinHoneyComb[5] ={ 5.15,16.4, 5.15,16.4, 5.15};
112 Float_t RoutHoneyComb[5]={20.63,34.92,22.3, 32.02,20.63};
113 Float_t zInside;
114 Float_t zCooper=0.01; Float_t zChips=0.01;
115 Float_t yNylonTube[5]={10,20,10,20,10};
116
117
118 AliMatrix(idrotm[901], 90, 0, 90, 90, 180, 0);
119
120
121 // Nylon tubes
122 gMC->Gsvolu("GNYL","TUBE", idtmed[1], NylonTube, 3); //support nylon tube
123 Float_t wideSupport=zSi+3*zPCB+2*NylonTube[2]+zHoneyComb;
124 // cout<<" wideSupport "<<wideSupport<<endl;
125
126 for (ifmd=0; ifmd<5; ifmd++)
127 {
128 sprintf(nameFMD,"FMD%d",ifmd+1);
129 ppcon[0]=0;
130 ppcon[1]=360;
131 ppcon[2]=4;
132
23ab7963 133 ppcon[3]=-wideSupport;
134 ppcon[4]=rin[ifmd]+0.1;
135 ppcon[5]=rout[ifmd]+0.1;
0d5ea2e8 136
137 ppcon[6]=ppcon[3]+zSi+2*zPCB+2*NylonTube[2];
23ab7963 138 ppcon[7]=rin[ifmd]+0.1;
139 ppcon[8]=rout[ifmd]+0.1;
0d5ea2e8 140
141 ppcon[9]=ppcon[6];
23ab7963 142 ppcon[10]=RinHoneyComb[ifmd]+0.1;
143 ppcon[11]=RoutHoneyComb[ifmd]+0.1;
0d5ea2e8 144
23ab7963 145 ppcon[12]=ppcon[9]+2*zHoneyComb+zPCB;
146 ppcon[13]=RinHoneyComb[ifmd]+0.1;
147 ppcon[14]=RoutHoneyComb[ifmd]+0.1;
0d5ea2e8 148 zFMD=z[ifmd]+wideSupport/2.;
149 // cout<<" Si "<<ifmd+1<<" "<<zFMD<<" "<<ppcon[12]<<endl;
150 gMC->Gsvolu(nameFMD,"PCON",idtmed[0],ppcon,15);
151 if (z[ifmd] >0){
152 gMC->Gspos(nameFMD,1,"ALIC",0,0,z[ifmd],0, "ONLY");}
153 else {
154 gMC->Gspos(nameFMD,1,"ALIC",0,0,z[ifmd],idrotm[901], "ONLY");}
155 //silicon
156 sprintf(nameSi,"GSI%d",ifmd+1);
157 sprintf(nameSector,"GSC%d",ifmd+1);
158 sprintf(nameRing,"GRN%d",ifmd+1);
23ab7963 159
0d5ea2e8 160 //honeycomb support
161 sprintf(nameHoney,"GSU%d",ifmd+1);
162 gMC->Gsvolu(nameHoney,"TUBE", idtmed[0], par, 0); //honeycomb
163 sprintf(nameHoneyIn,"GHI%d",ifmd+1);
164 gMC->Gsvolu(nameHoneyIn,"TUBE", idtmed[7], par, 0); //honey comb inside
165 sprintf(nameHoneyOut,"GHO%d",ifmd+1);
166 gMC->Gsvolu(nameHoneyOut,"TUBE", idtmed[6], par, 0); //honey comb skin
167 //PCB
168 sprintf(namePCB,"GPC%d",ifmd+1);
169 gMC->Gsvolu(namePCB,"TUBE", idtmed[0], par, 0); //PCB
170 sprintf(nameCopper,"GCO%d",ifmd+1);
171 gMC->Gsvolu(nameCopper,"TUBE", idtmed[3], par, 0); // Cooper
172 sprintf(nameChips,"GCH%d",ifmd+1);
173 gMC->Gsvolu(nameChips,"TUBE", idtmed[5], par, 0); // Si chips
174 sprintf(nameG10,"G10%d",ifmd+1);
175 gMC->Gsvolu(nameG10,"TUBE", idtmed[2], par, 0); //G10 plate
176 //last PCB
177 sprintf(nameLPCB,"GPL%d",ifmd+1);
178 gMC->Gsvolu(nameLPCB,"TUBE", idtmed[0], par, 0); //PCB
179 sprintf(nameLCopper,"GCL%d",ifmd+1);
180 gMC->Gsvolu(nameLCopper,"TUBE", idtmed[3], par, 0); // Cooper
181 sprintf(nameLChips,"GHL%d",ifmd+1);
182 gMC->Gsvolu(nameLChips,"TUBE", idtmed[5], par, 0); // Si chips
183 sprintf(nameGL10,"G1L%d",ifmd+1);
184 gMC->Gsvolu(nameGL10,"TUBE", idtmed[2], par, 0); //G10 plate
185
186
187 par[0]=rin[ifmd]; // pipe size
188 par[1]=rout[ifmd];
189 par[2]=zSi/2;
190 gMC->Gsvolu(nameSi,"TUBE", idtmed[4], par, 3);
0d5ea2e8 191 zInside=ppcon[3]+par[2];
23ab7963 192 gMC->Gspos(nameSi,ifmd+1,nameFMD,0,0,zInside,0, "ONLY");
0d5ea2e8 193 //PCB 1
194 zInside += par[2]+zPCB/2;
195 par[2]=zPCB/2;
0d5ea2e8 196 gMC->Gsposp(namePCB,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
197 zInside += zPCB;
0d5ea2e8 198 gMC->Gsposp(namePCB,2,nameFMD,0,0,zInside,0, "ONLY",par,3);
199 Float_t NulonTubeBegin=zInside+2*zPCB;
200 par[2]=zPCB/2-0.02;
201 Float_t zInPCB = -zPCB/2+par[2];
202 gMC->Gsposp(nameG10,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
203 zInPCB+=par[2]+zCooper/2 ;
204 par[2]=zCooper/2;
205 gMC->Gsposp(nameCopper,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
206 zInPCB += zCooper/2 + zChips/2;
207 par[2]=zChips/2;
208 gMC->Gsposp(nameChips,1,namePCB,0,0,zInPCB,0, "ONLY",par,3);
209 //HoneyComb
210 zHoneyComb=0.8;
211 par[0] = RinHoneyComb[ifmd];
212 par[1] = RoutHoneyComb[ifmd];
213 par[2] = zHoneyComb/2;
0d5ea2e8 214 zInside += 2*NylonTube[2]+par[2];
0d5ea2e8 215 gMC->Gsposp(nameHoney,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
216 par[2]=0.1/2;
217 Float_t zHoney=-zHoneyComb/2+par[2];
0d5ea2e8 218 gMC->Gsposp(nameHoneyOut,1,nameHoney,0,0,zHoney,0,
219 "ONLY",par,3); //shkurki
220 zHoney=zHoneyComb/2-par[2];
0d5ea2e8 221 gMC->Gsposp(nameHoneyOut,2,nameHoney,0,0,zHoney,0, "ONLY",par,3);
222 par[2]=(zHoneyComb-2.*0.1)/2; //soty vnutri
223 gMC->Gsposp(nameHoneyIn,1,nameHoney,0,0,0,0, "ONLY",par,3);
224
225 gMC->Gspos("GNYL",1,nameFMD,0,yNylonTube[ifmd],
226 NulonTubeBegin+NylonTube[2]/2.,0, "ONLY");
227 gMC->Gspos("GNYL",2,nameFMD,0,-yNylonTube[ifmd],
228 NulonTubeBegin+NylonTube[2]/2.,0, "ONLY");
229
230 //last PCB
231 par[0]=RoutHoneyComb[ifmd]-9;
232 par[1]=RoutHoneyComb[ifmd];
233 par[2]=zPCB/2;
0d5ea2e8 234 zInside += zHoneyComb/2+par[2];
235 gMC->Gsposp(nameLPCB,1,nameFMD,0,0,zInside,0, "ONLY",par,3);
236
237 par[2]=zPCB/2-0.02;
238 zInPCB = -zPCB/2+par[2];
239 gMC->Gsposp(nameGL10,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
240 zInPCB+=par[2]+zCooper/2 ;
241 par[2]=zCooper/2;
242 gMC->Gsposp(nameLCopper,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
243 zInPCB += zCooper/2 + zChips/2;
244 par[2]=zChips/2;
245 gMC->Gsposp(nameLChips,1,nameLPCB,0,0,zInPCB,0, "ONLY",par,3);
246
23ab7963 247
0d5ea2e8 248 //Granularity
249 fSectorsSi1=20;
250 fRingsSi1=256*3;
23ab7963 251 // fRingsSi1=3; // for drawing only
0d5ea2e8 252 fSectorsSi2=40;
253 fRingsSi2=128*3;
254 // fRingsSi2=3; //for drawing onl
255 if(ifmd==1||ifmd==3)
256 {
257 gMC->Gsdvn(nameSector, nameSi , fSectorsSi2, 2);
258 gMC->Gsdvn(nameRing, nameSector, fRingsSi2, 1);
259 }
260 else
261 {
262 gMC->Gsdvn(nameSector, nameSi , fSectorsSi1, 2);
263 gMC->Gsdvn(nameRing, nameSector , fRingsSi1, 1);
264 }
265
266 }
267}
268
269
270//------------------------------------------------------------------------
271
272
273void AliFMDv2::CreateMaterials()
274{
275 Int_t isxfld = gAlice->Field()->Integ();
276 Float_t sxmgmx = gAlice->Field()->Max();
277
278 // Plastic CH
279 Float_t aPlastic[2]={1.01,12.01};
280 Float_t zPlastic[2]={1,6};
281 Float_t wPlastic[2]={1,1};
282 Float_t denPlastic=1.03;
283 // 60% SiO2 , 40% G10FR4
284 // PC board
285 Float_t apcb[3] = { 28.0855,15.9994,17.749 };
286 Float_t zpcb[3] = { 14.,8.,8.875 };
287 Float_t wpcb[3] = { .28,.32,.4 };
288 Float_t denspcb = 1.8;
289 //
290
291 //*** Definition Of avaible FMD materials ***
292 AliMaterial(0, "FMD Air$", 14.61, 7.3, .001205, 30423.,999);
293 AliMixture(1, "Plastic$",aPlastic,zPlastic,denPlastic,-2,wPlastic);
294 AliMixture(2, "SSD PCB$", apcb, zpcb, denspcb, 3, wpcb);
295 AliMaterial(3, "SSD Copper$", 63.546, 29., 8.96, 1.43, 999.);
296 AliMaterial(4, "SSD Si$", 28.0855, 14., 2.33, 9.36, 999.);
297 AliMaterial(5, "SSD Si chip$", 28.0855, 14., 2.33, 9.36, 999.);
298 AliMaterial(6, "SSD C$", 12.011, 6., 2.265,18.8, 999.);
299 AliMaterial(7, "SSD Kapton$", 12.011, 6., 0.01, 31.27, 999.);//honeycomb
300 AliMaterial(8, "SSD G10FR4$", 17.749, 8.875, 1.8, 21.822, 999.);
301
302
303//**
304 AliMedium(0, "FMD air$", 0, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
305 AliMedium(1, "Plastic$", 1, 0,isxfld, sxmgmx, 10., .01, 1., .003, .003);
306 AliMedium(2, "SSD PCB$", 2, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
307 AliMedium(3, "SSD Copper$", 3, 0,isxfld, sxmgmx, 10., .01, 1., .003, .003);
308 AliMedium(4, "SSD Si$", 4, 1, isxfld, sxmgmx, 1., .001, 1., .001, .001);
309 AliMedium(5, "SSD Si chip$", 5, 0,isxfld, sxmgmx, 10., .01, 1., .003, .003);
310 AliMedium(6, "SSD C$", 6, 0,isxfld, sxmgmx, 10., .01, 1., .003, .003);
311 AliMedium(7, "SSD Kapton$", 7, 0, isxfld, sxmgmx, 1., .001, 1., .001, .001);
312 AliMedium(8, "SSD G10FR4$", 8, 0,isxfld, sxmgmx, 10., .01, 1., .003, .003);
313
314
315
316}
317
318//---------------------------------------------------------------------
319void AliFMDv2::DrawDetector()
320{
321//
322// Draw a shaded view of the Forward multiplicity detector version 0
323//
324
325
326//Set ALIC mother transparent
327gMC->Gsatt("ALIC","SEEN",0);
328//
329//Set volumes visible
330gMC->Gsatt("FMD1","SEEN",1);
331gMC->Gsatt("FMD2","SEEN",1);
332gMC->Gsatt("FMD3","SEEN",1);
333gMC->Gsatt("FMD4","SEEN",1);
334gMC->Gsatt("FMD5","SEEN",1);
335
336//
337gMC->Gdopt("hide","on");
338gMC->Gdopt("shad","on");
339gMC->SetClipBox(".");
340gMC->SetClipBox("*",0,1000,-1000,1000,-1000,1000);
341gMC->DefaultRange();
342gMC->Gdraw("alic",40,30,0,12,9.5,.2,0.2);
343gMC->Gdhead(1111,"Forward multiplicity detector");
344gMC->Gdopt("hide","off");
345}
346//-------------------------------------------------------------------
347void AliFMDv2::Init()
348{
349// Initialises version 0 of the Forward Multiplicity Detector
350//
351AliFMD::Init();
352fIdSens1=gMC->VolId("GRN1");
353fIdSens2=gMC->VolId("GRN2");
354fIdSens3=gMC->VolId("GRN3");
355fIdSens4=gMC->VolId("GRN4");
356fIdSens5=gMC->VolId("GRN5");
357printf("*** FMD version 2 initialized ***\n");
358}
359
360//-------------------------------------------------------------------
361
362void AliFMDv2::StepManager()
363{
364 //
365 // Called for every step in the Forward Multiplicity Detector
366 //
367 Int_t id,copy,copy1,copy2;
368 static Float_t hits[9];
369 static Int_t vol[3];
370 static Float_t de;
371 TLorentzVector pos;
372 TLorentzVector mom;
373
374
375 TClonesArray &lhits = *fHits;
376 if(!gMC->IsTrackAlive()) return; // particle has disappeared
377
378 Float_t charge = gMC->TrackCharge();
379 if(TMath::Abs(charge)<=0.) return; //take only charged particles
380
381 // printf(" in StepManeger \n");
382 id=gMC->CurrentVolID(copy);
383
384// Check the sensetive volume
385 if(id==fIdSens1||id==fIdSens2||id==fIdSens3||id==fIdSens4||id==fIdSens5)
386 {
387 if(gMC->IsTrackEntering())
388 {
389 vol[2]=copy;
390 gMC->CurrentVolOffID(1,copy1);
391 vol[1]=copy1;
392 gMC->CurrentVolOffID(2,copy2);
393 vol[0]=copy2;
394
395 gMC->TrackPosition(pos);
396 hits[0]=pos[0];
397 hits[1]=pos[1];
398 hits[2]=pos[2];
399
400 gMC->TrackMomentum(mom);
401 hits[3]=mom[0];
402 hits[4]=mom[1];
403 hits[5]=mom[2];
404
405 Int_t iPart= gMC->TrackPid();
406 Int_t partId=gMC->IdFromPDG(iPart);
407 hits[7]=partId;
408 hits[8]=1e9*gMC->TrackTime();
409 de=0.;
410 }
411 if(gMC->IsTrackInside()){
412 de=de+1000.*gMC->Edep();
413 }
414
415 if(gMC->IsTrackExiting()
416 ||gMC->IsTrackDisappeared()||
417 gMC->IsTrackStop())
418 {
419 hits[6]=de+1000.*gMC->Edep();
420 new(lhits[fNhits++]) AliFMDhit(fIshunt,gAlice->CurrentTrack(),vol,hits);
421 } // IsTrackExiting()
422 }
423 }
424
425//--------------------------------------------------------------------------
426
427void AliFMDv2::Response( Float_t Edep)
428{
429 Float_t I=1.664*0.04*2.33/22400; // = 0.69e-6;
430 Float_t chargeOnly=Edep/I;
431 //Add noise ~500electrons
432 Int_t charge=500;
433 if (Edep>0)
434 charge=Int_t(gRandom->Gaus(chargeOnly,500));
435 }
436
437
438
439
440
441
442