]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONv0.cxx
new class to take into account ITS material distribution in tracking v1
[u/mrichter/AliRoot.git] / MUON / AliMUONv0.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 $Log$
18 Revision 1.15  2001/01/17 20:57:45  hristov
19 Unused variable removed
20
21 Revision 1.14  2000/12/21 22:42:55  morsch
22 Constructor contains default set-up for segmentation.
23 Record charged particles only.
24
25 Revision 1.13  2000/10/06 10:03:38  morsch
26 Call to gMC->VolId() moved to Init()
27
28 Revision 1.12  2000/10/02 21:28:09  fca
29 Removal of useless dependecies via forward declarations
30
31 Revision 1.11  2000/06/27 07:31:07  morsch
32 fChambers = 0; deleted from constructor.
33
34 Revision 1.10  2000/06/26 14:02:38  morsch
35 Add class AliMUONConstants with MUON specific constants using static memeber data and access methods.
36
37 Revision 1.9  2000/06/15 07:58:49  morsch
38 Code from MUON-dev joined
39
40 Revision 1.8.4.9  2000/06/12 19:20:49  morsch
41 Constructor sets default geometry, segmentation and response parameters.
42
43 Revision 1.8.4.8  2000/06/09 21:55:28  morsch
44 Most coding rule violations corrected.
45
46 Revision 1.8.4.7  2000/05/02 13:15:18  morsch
47 Coding rule violations RS3, RN13 corected
48
49 Revision 1.8.4.6  2000/05/02 10:24:26  morsch
50 Public access to fdGas and fdAlu of AliMUONChamber replaced by getters.
51
52 Revision 1.8.4.5  2000/04/26 19:58:47  morsch
53 Obsolete reference to trig_ removed.
54
55 Revision 1.8.4.4  2000/04/19 19:42:47  morsch
56 change NCH to kNCH
57
58 Revision 1.8.4.3  2000/02/17 08:17:43  morsch
59 Gammas and neutrons are also scored in the stepmanager
60 */
61
62 /////////////////////////////////////////////////////////
63 //  Manager and hits classes for set:MUON version 0    //
64 /////////////////////////////////////////////////////////
65
66 #include <TTUBE.h>
67 #include <TNode.h> 
68 #include <TRandom.h> 
69 #include <TLorentzVector.h> 
70 #include <iostream.h>
71
72 #include "AliMUONv0.h"
73 #include "AliMUONChamber.h"
74 #include "AliRun.h"
75 #include "AliMC.h"
76 #include "AliMagF.h"
77 #include "AliMUONHit.h"
78 #include "AliMUONPadHit.h"
79 #include "AliCallf77.h"
80 #include "AliConst.h" 
81 #include "AliMUONResponseV0.h"
82 #include "AliMUONResponseTrigger.h"
83 #include "AliMUONSegmentationV0.h"
84 #include "AliMUONSegmentationV01.h"
85 #include "AliMUONSegmentationV02.h"
86 #include "AliMUONSegmentationV04.h"
87 #include "AliMUONSegmentationV05.h"
88 #include "AliMUONSegmentationSlat.h"
89 #include "AliMUONSegmentationSlatN.h"
90 #include "AliMUONSegmentationTrigger.h"
91 #include "AliMUONSegmentationTriggerX.h"
92 #include "AliMUONSegmentationTriggerY.h"
93 #include "AliMUONConstants.h"
94
95 ClassImp(AliMUONv0)
96  
97 //___________________________________________
98 AliMUONv0::AliMUONv0() : AliMUON()
99 {
100 // Constructor
101     fChambers = 0;
102 }
103  
104 //___________________________________________
105 AliMUONv0::AliMUONv0(const char *name, const char *title)
106        : AliMUON(name,title)
107 {
108 // Constructor
109     SetIshunt(0);
110     SetMaxStepGas(0.1);
111     SetMaxStepAlu(0.1);
112     
113 //
114 // First define the number of planes that are segmented (1 or 2) by a call
115 // to SetNsec. 
116 // Then chose for each chamber (chamber plane) the segmentation 
117 // and response model.
118 // They should be equal for the two chambers of each station. In a future
119 // version this will be enforced.
120 //
121 //  
122     Int_t chamber;
123     // Default response: 5 mm of gas
124     AliMUONResponseV0* response0 = new AliMUONResponseV0;
125     response0->SetSqrtKx3AndDeriveKx2Kx4(0.7131); // sqrt(0.5085)
126     response0->SetSqrtKy3AndDeriveKy2Ky4(0.7642); // sqrt(0.5840)
127     response0->SetPitch(0.25); // anode-cathode distance
128     response0->SetSigmaIntegration(10.);
129     response0->SetChargeSlope(50);
130     response0->SetChargeSpread(0.18, 0.18);
131     response0->SetMaxAdc(4096);
132     response0->SetZeroSuppression(6);
133     
134     // Response for 4 mm of gas (station 1)
135     // automatic consistency with width of sensitive medium in CreateGeometry ????
136     AliMUONResponseV0* responseSt1 = new AliMUONResponseV0;
137     // Mathieson parameters from L.Kharmandarian's thesis, page 190
138     responseSt1->SetSqrtKx3AndDeriveKx2Kx4(0.7000); // sqrt(0.4900)
139     responseSt1->SetSqrtKy3AndDeriveKy2Ky4(0.7550); // sqrt(0.5700)
140     responseSt1->SetPitch(0.20); // anode-cathode distance
141     responseSt1->SetSigmaIntegration(10.);
142     // ChargeSlope larger to compensate for the smaller anode-cathode distance
143     // and keep the same most probable ADC channel for mip's
144     responseSt1->SetChargeSlope(62.5); 
145     // assumed proportionality to anode-cathode distance for ChargeSpread
146     responseSt1->SetChargeSpread(0.144, 0.144);
147     responseSt1->SetMaxAdc(4096);
148     responseSt1->SetZeroSuppression(6);
149     
150 //--------------------------------------------------------
151 // Configuration for Chamber TC1/2  (Station 1) ----------           
152 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
153     Float_t rseg1[4]={17.5, 55.2, 71.3, 95.5};
154     Int_t   nseg1[4]={4, 4, 2, 1};
155 //
156     chamber=1;
157 //^^^^^^^^^
158     SetNsec(chamber-1,2);
159 //
160     AliMUONSegmentationV01 *seg11=new AliMUONSegmentationV01(4);
161     
162     seg11->SetSegRadii(rseg1);
163 //  seg11->SetPadSize(3, 0.5);
164     seg11->SetPadSize(2.4, 0.4); // smaller pad size
165 //  seg11->SetDAnod(3.0/3./4);
166     seg11->SetDAnod(0.20); // smaller distance between anode wires
167     seg11->SetPadDivision(nseg1);
168     
169     SetSegmentationModel(chamber-1, 1, seg11);
170 //
171     AliMUONSegmentationV02 *seg12=new AliMUONSegmentationV02(4);
172     seg12->SetSegRadii(rseg1); 
173 //  seg12->SetPadSize(0.75, 2.0);
174     seg12->SetPadSize(0.6, 1.6); // smaller pad size
175 //  seg12->SetDAnod(3.0/3./4);
176     seg12->SetDAnod(0.20); // smaller distance between anode wires
177     seg12->SetPadDivision(nseg1);
178     
179     SetSegmentationModel(chamber-1, 2, seg12);
180     
181 //  SetResponseModel(chamber-1, response0);         
182     SetResponseModel(chamber-1, responseSt1); // special response           
183     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
184     
185     chamber=2;
186 //^^^^^^^^^
187 //
188     SetNsec(chamber-1,2);
189 //
190     AliMUONSegmentationV01 *seg21=new AliMUONSegmentationV01(4);
191     seg21->SetSegRadii(rseg1);
192 //  seg21->SetPadSize(3, 0.5);
193     seg21->SetPadSize(2.4, 0.4); // smaller pad size
194 //  seg21->SetDAnod(3.0/3./4);
195     seg21->SetDAnod(0.20); // smaller distance between anode wires
196     seg21->SetPadDivision(nseg1);
197     SetSegmentationModel(chamber-1, 1, seg21);
198 //
199     AliMUONSegmentationV02 *seg22=new AliMUONSegmentationV02(4);
200     seg22->SetSegRadii(rseg1); 
201 //  seg22->SetPadSize(0.75, 2.);
202     seg22->SetPadSize(0.6, 1.6); // smaller pad size
203 //  seg22->SetDAnod(3.0/3./4);
204     seg22->SetDAnod(0.20); // smaller distance between anode wires
205     seg22->SetPadDivision(nseg1);
206     SetSegmentationModel(chamber-1, 2, seg22);
207     
208 //  SetResponseModel(chamber-1, response0);         
209     SetResponseModel(chamber-1, responseSt1); // special response
210     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
211     
212 //
213 //--------------------------------------------------------
214 // Configuration for Chamber TC3/4 (Station 2) -----------
215 ///^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
216 // Float_t rseg2[4]={23.5, 87.7, 122.4, 122.5};
217     Float_t rseg2[4]={23.5, 47.1, 87.7, 122.5};
218     Int_t   nseg2[4]={4, 4, 2, 1};
219 //
220     chamber=3;
221 //^^^^^^^^^
222     SetNsec(chamber-1,2);
223 //
224     AliMUONSegmentationV01 *seg31=new AliMUONSegmentationV01(4);
225     seg31->SetSegRadii(rseg2);
226     seg31->SetPadSize(3.0, 0.5);
227     seg31->SetDAnod(3.0/3./4);
228     seg31->SetPadDivision(nseg2);
229     SetSegmentationModel(chamber-1, 1, seg31);
230 //
231     AliMUONSegmentationV02 *seg32=new AliMUONSegmentationV02(4);
232     seg32->SetSegRadii(rseg2); 
233     seg32->SetPadSize(0.75, 2.0);
234     seg32->SetPadDivision(nseg2);
235     seg32->SetDAnod(3.0/3./4);
236     
237     SetSegmentationModel(chamber-1, 2, seg32);
238     
239     SetResponseModel(chamber-1, response0);         
240     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
241     
242     chamber=4;
243 //^^^^^^^^^
244 //
245     SetNsec(chamber-1,2);
246 //
247     AliMUONSegmentationV01 *seg41=new AliMUONSegmentationV01(4);
248     seg41->SetSegRadii(rseg2);
249     seg41->SetPadSize(3.0, 0.5);
250     seg41->SetDAnod(3.0/3./4);
251     seg41->SetPadDivision(nseg2);
252     SetSegmentationModel(chamber-1, 1, seg41);
253 //
254     AliMUONSegmentationV02 *seg42=new AliMUONSegmentationV02(4);
255     seg42->SetSegRadii(rseg2); 
256     seg42->SetPadSize(0.75, 2.0);
257     seg42->SetPadDivision(nseg2);
258     seg42->SetDAnod(3.0/3./4);
259     
260     SetSegmentationModel(chamber-1, 2, seg42);
261     
262     SetResponseModel(chamber-1, response0);         
263     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
264     
265     
266 //--------------------------------------------------------
267 // Configuration for Chamber TC5/6  (Station 3) ----------          
268 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
269     Int_t   nseg3[4]={4, 4, 2, 1};
270     Int_t   npcb5[36] = {0,0,2,0,
271                          0,0,3,0,
272                          0,1,3,0,
273                          0,2,2,0,
274                          0,1,2,0, 
275                          0,2,2,0, 
276                          0,1,3,0, 
277                          0,0,3,0,
278                          0,0,2,0};
279
280     Float_t shift = 1.5/2.;
281     // Float_t xpos5[8]    = {2., 2., 2., 42., 42., 2., 2., 2.};
282     Float_t xpos5[9]    = {2., 2., 2., 2.,32., 2., 2., 2., 2.};
283     Float_t ypos5       = -(20.+4.*(40.-2.*shift));
284     
285     chamber=5;
286     SetNsec(chamber-1,2);
287     AliMUONSegmentationSlat *seg51=new AliMUONSegmentationSlat(4);
288     seg51->SetNSlats(9); 
289     seg51->SetShift(shift);  
290     seg51->SetNPCBperSector(npcb5); 
291     seg51->SetSlatXPositions(xpos5);
292     seg51->SetSlatYPosition(ypos5);
293     seg51->SetPadSize(10.,0.5);
294     seg51->SetDAnod(0.25);
295     seg51->SetPadDivision(nseg3);
296     SetSegmentationModel(chamber-1, 1, seg51);
297     
298     AliMUONSegmentationSlatN *seg52=new AliMUONSegmentationSlatN(4);
299     seg52->SetNSlats(9); 
300     seg52->SetShift(shift);  
301     seg52->SetNPCBperSector(npcb5); 
302     seg52->SetSlatXPositions(xpos5);
303     seg52->SetSlatYPosition(ypos5);
304     seg52->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
305     seg52->SetDAnod(0.25);
306     seg52->SetPadDivision(nseg3);
307     SetSegmentationModel(chamber-1, 2, seg52);
308     SetResponseModel(chamber-1, response0);         
309     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
310     
311     chamber=6;
312     SetNsec(chamber-1,2);
313     AliMUONSegmentationSlat *seg61=new AliMUONSegmentationSlat(4);
314     seg61->SetNSlats(9); 
315     seg61->SetShift(shift);  
316     seg61->SetNPCBperSector(npcb5); 
317     seg61->SetSlatXPositions(xpos5);
318     seg61->SetSlatYPosition(ypos5);
319     seg61->SetPadSize(10.,0.5);
320     seg61->SetDAnod(0.25);
321     seg61->SetPadDivision(nseg3);
322     SetSegmentationModel(chamber-1, 1, seg61);
323     
324     AliMUONSegmentationSlatN *seg62=new AliMUONSegmentationSlatN(4);
325     seg62->SetNSlats(9); 
326     seg62->SetShift(shift);  
327     seg62->SetNPCBperSector(npcb5); 
328     seg62->SetSlatXPositions(xpos5);
329     seg62->SetSlatYPosition(ypos5);
330     seg62->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
331     seg62->SetDAnod(0.25);
332     seg62->SetPadDivision(nseg3);
333     SetSegmentationModel(chamber-1, 2, seg62);
334     SetResponseModel(chamber-1, response0);         
335     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
336     
337 //--------------------------------------------------------
338 // Configuration for Chamber TC7/8  (Station 4) ----------           
339 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
340
341     Int_t   nseg4[4]={4, 4, 2, 1};
342
343     chamber=7;
344 //^^^^^^^^^
345
346     SetNsec(chamber-1,2);
347 //
348     AliMUONSegmentationSlat *seg71=new AliMUONSegmentationSlat(4);
349     Int_t npcb7[44] = {0,0,0,3,
350                        0,0,2,2,
351                        0,0,3,2,
352                        0,2,2,1,
353                        0,2,2,1,
354                        0,1,2,1, 
355                        0,2,2,1, 
356                        0,2,2,1, 
357                        0,0,3,2, 
358                        0,0,2,2, 
359                        0,0,0,3};
360     Float_t xpos7[11]   = {2., 2., 2., 2., 2., 39.5, 2., 2., 2., 2., 2.};
361     Float_t ypos7       = -(20.+5.*(40.-2.*shift));
362     
363     seg71->SetNSlats(11);  
364     seg71->SetShift(shift);  
365     seg71->SetNPCBperSector(npcb7); 
366     seg71->SetSlatXPositions(xpos7);
367     seg71->SetSlatYPosition(ypos7);
368     
369     seg71->SetPadSize(10.,0.5);
370     seg71->SetDAnod(0.25);
371     seg71->SetPadDivision(nseg4);
372     SetSegmentationModel(chamber-1, 1, seg71);
373     
374     AliMUONSegmentationSlatN *seg72=new AliMUONSegmentationSlatN(4);
375     
376     SetSegmentationModel(chamber-1, 2, seg72);
377     seg72->SetNSlats(11);  
378     seg72->SetShift(shift);   
379     seg72->SetNPCBperSector(npcb7); 
380     seg72->SetSlatXPositions(xpos7);
381     seg72->SetSlatYPosition(ypos7);
382     seg72->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
383     seg72->SetDAnod(0.25);
384     seg72->SetPadDivision(nseg4);
385     
386     SetResponseModel(chamber-1, response0);         
387     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
388     
389     chamber=8;
390 //^^^^^^^^^
391     SetNsec(chamber-1,2);
392 //
393     AliMUONSegmentationSlat *seg81=new AliMUONSegmentationSlat(4);
394     
395     seg81->SetNSlats(11);  
396     seg81->SetShift(shift);  
397     seg81->SetNPCBperSector(npcb7); 
398     seg81->SetSlatXPositions(xpos7);
399     seg81->SetSlatYPosition(ypos7);
400     seg81->SetPadSize(10.,0.5);
401     seg81->SetDAnod(0.25);
402     seg81->SetPadDivision(nseg4);
403     SetSegmentationModel(chamber-1, 1, seg81);
404     
405     AliMUONSegmentationSlat *seg82=new AliMUONSegmentationSlatN(4);
406
407     SetSegmentationModel(chamber-1, 2, seg82);
408     seg82->SetNSlats(11);  
409     seg82->SetShift(shift);  
410     seg82->SetNPCBperSector(npcb7); 
411     seg82->SetSlatXPositions(xpos7);
412     seg82->SetSlatYPosition(ypos7);
413     seg82->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
414     seg82->SetDAnod(0.25);
415     seg82->SetPadDivision(nseg4);
416     
417     SetResponseModel(chamber-1, response0);         
418     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
419     
420
421 //--------------------------------------------------------
422 // Configuration for Chamber TC9/10  (Station 5) ---------           
423 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
424     chamber=9;
425 //^^^^^^^^^
426     
427     SetNsec(chamber-1,2);
428 //
429     AliMUONSegmentationSlat *seg91=new AliMUONSegmentationSlat(4);
430     Int_t   npcb9[52] = {0,0,0,3,
431                          0,0,0,4,
432                          0,0,2,3,
433                          0,0,3,3,
434                          0,2,2,2,
435                          0,2,2,2,
436                          0,1,2,2, 
437                          0,2,2,2, 
438                          0,2,2,2, 
439                          0,0,3,3, 
440                          0,0,2,3, 
441                          0,0,0,4, 
442                          0,0,0,3};   
443     
444     // Float_t xpos9[13]   = {2., 2., 2., 2., 2., 2., 39.5 , 2., 2., 2., 2., 2., 2.};
445     Float_t xpos9[13]   = {2., 2., 2., 2., 2., 2., 39.5, 2., 2., 2., 2., 2., 2.};
446     Float_t ypos9       = -(20.+6.*(40.-2.*shift));
447     
448     seg91->SetNSlats(13);  
449     seg91->SetShift(shift);  
450     seg91->SetNPCBperSector(npcb9); 
451     seg91->SetSlatXPositions(xpos9);
452     seg91->SetSlatYPosition(ypos9);
453     seg91->SetPadSize(10.,0.5);
454     seg91->SetDAnod(0.25);
455     seg91->SetPadDivision(nseg4);
456     SetSegmentationModel(chamber-1, 1, seg91);
457     
458     AliMUONSegmentationSlatN *seg92=new AliMUONSegmentationSlatN(4);
459     
460     SetSegmentationModel(chamber-1, 2, seg92);
461     seg92->SetNSlats(13);  
462     seg92->SetShift(shift);   
463     seg92->SetNPCBperSector(npcb9); 
464     seg92->SetSlatXPositions(xpos9);
465     seg92->SetSlatYPosition(ypos9);
466     seg92->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
467     seg92->SetDAnod(0.25);
468     seg92->SetPadDivision(nseg4);
469     
470     SetResponseModel(chamber-1, response0);         
471     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
472     
473     chamber=10;
474 //^^^^^^^^^
475     SetNsec(chamber-1,2);
476 //
477     AliMUONSegmentationSlat *seg101=new AliMUONSegmentationSlat(4);
478     
479     seg101->SetNSlats(13);  
480     seg101->SetShift(shift);  
481     seg101->SetNPCBperSector(npcb9); 
482     seg101->SetSlatXPositions(xpos9);
483     seg101->SetSlatYPosition(ypos9);
484     seg101->SetPadSize(10.,0.5);
485     seg101->SetDAnod(0.25);
486     seg101->SetPadDivision(nseg4);
487     SetSegmentationModel(chamber-1, 1, seg101);
488     
489     AliMUONSegmentationSlatN *seg102=new AliMUONSegmentationSlatN(4);
490     
491     SetSegmentationModel(chamber-1, 2, seg102);
492     seg102->SetNSlats(13);  
493     seg102->SetShift(shift);   
494     seg102->SetNPCBperSector(npcb9); 
495     seg102->SetSlatXPositions(xpos9);
496     seg102->SetSlatYPosition(ypos9);
497     seg102->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
498     seg102->SetDAnod(0.25);
499     seg102->SetPadDivision(nseg4);
500     
501     SetResponseModel(chamber-1, response0);         
502     Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
503     
504 //--------------------------------------------------------
505 // Configuration for Trigger Stations -------------------- 
506 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
507 // Cluster-size off
508     AliMUONResponseTrigger* responseTrigger0 =  new AliMUONResponseTrigger;
509 // Cluster-size on  
510 // AliMUONResponseTriggerV1* responseTrigger0 =  new AliMUONResponseTriggerV1;
511  
512     chamber=11;
513     SetNsec(chamber-1,2);
514     AliMUONSegmentationTriggerX *seg111=new AliMUONSegmentationTriggerX;
515     SetSegmentationModel(chamber-1, 1, seg111);
516     AliMUONSegmentationTriggerY *seg112=new AliMUONSegmentationTriggerY;
517     SetSegmentationModel(chamber-1, 2, seg112);
518     
519     SetResponseModel(chamber-1, responseTrigger0);      
520     Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
521     
522  
523     chamber=12;
524     SetNsec(chamber-1,2);
525     AliMUONSegmentationTriggerX *seg121=new AliMUONSegmentationTriggerX;
526     SetSegmentationModel(chamber-1, 1, seg121);
527     AliMUONSegmentationTriggerY *seg122=new AliMUONSegmentationTriggerY;
528     SetSegmentationModel(chamber-1, 2, seg122);
529     
530     SetResponseModel(chamber-1, responseTrigger0);
531     Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
532     
533     chamber=13;
534     SetNsec(chamber-1,2);
535     AliMUONSegmentationTriggerX *seg131=new AliMUONSegmentationTriggerX;
536     SetSegmentationModel(chamber-1, 1, seg131);
537     AliMUONSegmentationTriggerY *seg132=new AliMUONSegmentationTriggerY;
538     SetSegmentationModel(chamber-1, 2, seg132);
539     SetResponseModel(chamber-1, responseTrigger0);      
540     Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
541     
542     chamber=14;
543     SetNsec(chamber-1,2);
544     AliMUONSegmentationTriggerX *seg141=new AliMUONSegmentationTriggerX;
545     SetSegmentationModel(chamber-1, 1, seg141);
546     AliMUONSegmentationTriggerY *seg142=new AliMUONSegmentationTriggerY;
547     SetSegmentationModel(chamber-1, 2, seg142);
548     
549     SetResponseModel(chamber-1, responseTrigger0); 
550     Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
551 }
552
553 void AliMUONv0::CreateGeometry()
554 {
555 // Creates coarse geometry for hit density simulations
556     Int_t *idtmed = fIdtmed->GetArray()-1099;
557 //
558      Float_t zpos, dAlu, tpar[3];
559      Int_t idAir=idtmed[1100];
560      Int_t idAlu=idtmed[1103];     
561
562      AliMUONChamber *iChamber;
563      // Loop over all chambers (tracking and trigger)
564      for (Int_t ch = 0; ch < AliMUONConstants::NCh(); ch++) {
565          char alu[8];
566          char gas[8];
567      
568          iChamber=(AliMUONChamber*) (*fChambers)[ch];
569          // Z of the chamber
570          zpos=iChamber->Z(); 
571          dAlu=iChamber->DAlu();
572          if (ch < AliMUONConstants::NTrackingCh()) {
573            // tracking chambers
574              sprintf(alu,"CA0%1d",ch);
575              sprintf(gas,"CG0%1d",ch);   
576          } else {
577            // trigger chambers
578              sprintf(alu,"CA%2d",ch);
579              sprintf(gas,"CG%2d",ch);    
580          }
581 //
582          tpar[0] = iChamber->RInner(); 
583          tpar[1] = iChamber->ROuter();
584          tpar[2] = (dAlu+0.2)/2.;
585          if (ch !=4 && ch !=5) {
586              gMC->Gsvolu(alu, "TUBE", idAlu, tpar, 3);
587              tpar[2] = 0.1;
588              gMC->Gsvolu(gas, "TUBE", idAir, tpar, 3);
589          } else {
590              gMC->Gsvolu(alu, "TUBE", idAlu, tpar, 3);
591              tpar[2] = 0.1;
592              gMC->Gsvolu(gas, "TUBE", idAir, tpar, 3);
593          }
594          gMC->Gspos(gas, 1, alu,  0., 0., 0., 0, "ONLY");
595          gMC->Gspos(alu, 1, "ALIC", 0., 0., zpos, 0, "ONLY");
596      }
597 }
598
599 //___________________________________________
600 void AliMUONv0::CreateMaterials()
601 {
602 // Creates materials for coarse geometry
603     AliMaterial(15, "AIR$      ", 14.61, 7.3, .001205, 30423.24, 67500);
604     AliMaterial(9, "ALUMINIUM$", 26.98, 13., 2.7, 8.9, 37.2);
605
606     Float_t epsil  = .001; // Tracking precision, 
607     Float_t stemax = -1.;  // Maximum displacement for multiple scat 
608     Float_t tmaxfd = -20.; // Maximum angle due to field deflection 
609     Float_t deemax = -.3;  // Maximum fractional energy loss, DLS 
610     Float_t stmin  = -.8;
611     Int_t isxfld   = gAlice->Field()->Integ();
612     Float_t sxmgmx = gAlice->Field()->Max();
613
614     //
615     //    Air 
616     AliMedium(1, "AIR_CH_US         ", 15, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
617     AliMedium(4, "ALU_CH_US          ", 9, 0, isxfld, sxmgmx, tmaxfd, fMaxStepAlu, 
618             fMaxDestepAlu, epsil, stmin);
619
620 }
621
622 void AliMUONv0::Init()
623 {
624    // 
625    // Initialize Tracking Chambers
626    //
627     char vName[8];
628     printf("\n\n\n Start Init for version 0 - CPC chamber type\n\n\n");
629     for (Int_t i=0; i<AliMUONConstants::NCh(); i++) {
630 // Initialise chamber
631         ((AliMUONChamber*) (*fChambers)[i])->Init();
632 // Set sensitive volume Id
633         if (i < AliMUONConstants::NTrackingCh()) {
634             // tracking chambers
635             sprintf(vName,"CG0%1d",i);   
636         } else {
637             // trigger chambers
638             sprintf(vName,"CG%2d",i);    
639         }
640         ((AliMUONChamber*) (*fChambers)[i])->SetGid(gMC->VolId(vName));
641     }
642 }
643
644 void AliMUONv0::StepManager()
645 {
646 //
647 // Step manager for hit density simulations
648   Int_t          copy, id;
649   static Int_t   idvol;
650   static Int_t   vol[2];
651   Int_t          ipart;
652   TLorentzVector pos;
653   TLorentzVector mom;
654   Float_t        theta,phi;
655   
656   //  modifs perso
657   static Float_t hits[15];
658
659   TClonesArray &lhits = *fHits;
660   //
661   // Only gas gap inside chamber
662   // Tag chambers and record hits when track enters 
663   idvol=-1;
664   id=gMC->CurrentVolID(copy);
665   
666     for (Int_t i=1; i<=AliMUONConstants::NCh(); i++) {
667       if(id==((AliMUONChamber*)(*fChambers)[i-1])->GetGid()){ 
668           vol[0]=i; 
669           idvol=i-1;
670       }
671     }
672     if (idvol == -1) return;
673   //
674   // Get current particle id (ipart), track position (pos)  and momentum (mom) 
675   gMC->TrackPosition(pos);
676   gMC->TrackMomentum(mom);
677
678   ipart  = gMC->TrackPid();
679   //
680   // record hits when track enters ...
681   if( !(gMC->TrackCharge()) ) return; 
682   if( gMC->IsTrackEntering()) {
683       Double_t tc = mom[0]*mom[0]+mom[1]*mom[1];
684       Double_t rt = TMath::Sqrt(tc);
685       theta   = Float_t(TMath::ATan2(rt,Double_t(mom[2])))*kRaddeg;
686       phi     = Float_t(TMath::ATan2(Double_t(mom[1]),Double_t(mom[0])))*kRaddeg;
687       hits[0] = Float_t(ipart);             // Geant3 particle type
688       hits[1] = pos[0];                     // X-position for hit
689       hits[2] = pos[1];                     // Y-position for hit
690       hits[3] = pos[2];                     // Z-position for hit
691       hits[4] = theta;                      // theta angle of incidence
692       hits[5] = phi;                        // phi angle of incidence 
693       hits[8] = -1;                         // first padhit
694       hits[9] = -1;                         // last pad hit
695
696       // modifs personel
697       hits[10] = mom[3]; // hit Energy
698       hits[11] = mom[0]; // Px
699       hits[12] = mom[1]; // Py
700       hits[13] = mom[2]; // Pz
701       hits[14] = gMC->TrackTime();
702       
703       // fin modifs perso
704       new(lhits[fNhits++]) 
705           AliMUONHit(fIshunt,gAlice->CurrentTrack(),vol,hits);
706
707   }
708 }
709
710
711
712
713
714
715
716
717