]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHv3.cxx
Code cleaning, all wranings removed with new Makefile options
[u/mrichter/AliRoot.git] / RICH / AliRICHv3.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 /* $Id$ */
17
18 #include <Riostream.h>
19
20 #include <TBRIK.h>
21 #include <TGeometry.h>
22 #include <TLorentzVector.h>
23 #include <TNode.h>
24 #include <TParticle.h>
25 #include <TVector3.h>
26 #include <TVirtualMC.h>
27 #include <TPDGCode.h> //for kNuetron
28
29 #include "AliConst.h"
30 #include "AliMagF.h"
31 #include "AliPDG.h"
32 #include "AliRICHGeometry.h"
33 #include "AliRICHResponseV0.h"
34 #include "AliRICHSegmentationV1.h"
35 #include "AliRICHv3.h"
36 #include "AliRun.h"
37
38 ClassImp(AliRICHv3)
39
40 //______________________________________________________________
41 //    Implementation of the RICH version 3 with azimuthal rotation
42
43
44 AliRICHv3::AliRICHv3(const char *sName, const char *sTitle)
45           :AliRICH(sName,sTitle)
46 {
47 // The named ctor currently creates a single copy of 
48 // AliRICHGeometry AliRICHSegmentationV1 AliRICHResponseV0
49 // and initialises the corresponding models of all 7 chambers with these stuctures.
50 // Note: all chambers share the single copy of models. MUST be changed later (???).
51   if(GetDebug())Info("named ctor","Start.");
52
53    fCkovNumber=fFreonProd=0;
54    
55    AliRICHGeometry       *pRICHGeometry    =new AliRICHGeometry;           // ??? to be moved to AlRICHChamber::named ctor
56    AliRICHSegmentationV1 *pRICHSegmentation=new AliRICHSegmentationV1;     // ??? to be moved to AlRICHChamber::named ctor
57    AliRICHResponseV0     *pRICHResponse    =new AliRICHResponseV0;         // ??? to be moved to AlRICHChamber::named ctor
58      
59    fChambers = new TObjArray(kNCH);
60    for (Int_t i=0; i<kNCH; i++){    
61       fChambers->AddAt(new AliRICHChamber,i); // ??? to be changed to named ctor of AliRICHChamber
62       SetGeometryModel(i,pRICHGeometry);
63       SetSegmentationModel(i,pRICHSegmentation);
64       SetResponseModel(i,pRICHResponse);
65       ((AliRICHChamber*)fChambers->At(i))->Init(i); // ??? to be removed     
66    }
67   if(GetDebug())Info("named ctor","Stop.");
68 }//AliRICHv3::ctor(const char *pcName, const char *pcTitle)
69
70 AliRICHv3::~AliRICHv3()
71 {
72 // Dtor deletes RICH models. In future (???) AliRICHChamber will be responsible for that.
73    if(GetDebug()) cout<<ClassName()<<"::dtor()>\n";
74       
75    if(fChambers) {
76      AliRICHChamber *ch = GetChamber(0); 
77      if(ch) {
78        delete ch->GetGeometryModel();
79        delete ch->GetResponseModel();
80        delete ch->GetSegmentationModel();
81      }
82    }
83 }//AliRICHv3::dtor()
84
85
86 void AliRICHv3::CreateGeometry()
87 {
88 // Provides geometry structure for simulation (currently GEANT volumes tree)         
89    if(GetDebug()) cout<<ClassName()<<"::CreateGeometry()>\n";
90
91   AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH"); 
92   AliRICHSegmentationV0*  segmentation;
93   AliRICHGeometry*  geometry;
94   AliRICHChamber*       iChamber;
95
96   iChamber = &(pRICH->Chamber(0));
97   segmentation=(AliRICHSegmentationV0*) iChamber->GetSegmentationModel();
98   geometry=iChamber->GetGeometryModel();
99
100   Float_t distance;
101   distance = geometry->GetFreonThickness()/2 + geometry->GetQuartzThickness() + geometry->GetGapThickness();
102   geometry->SetRadiatorToPads(distance);
103     
104   //Opaque quartz thickness
105   Float_t oqua_thickness = .5;
106   //CsI dimensions
107
108
109   Float_t csi_width = segmentation->Npx()*segmentation->Dpx() + segmentation->DeadZone();
110   Float_t csi_length = segmentation->Npy()*segmentation->Dpy() + 2*segmentation->DeadZone();
111   
112   
113   Int_t *idtmed = fIdtmed->GetArray()-999;
114     
115     Int_t i;
116     Float_t zs;
117     Int_t idrotm[1099];
118     Float_t par[3];
119     
120     // --- Define the RICH detector 
121     //     External aluminium box 
122     par[0] = 68.8;
123     par[1] = 13;                 //Original Settings
124     par[2] = 70.86;
125     gMC->Gsvolu("RICH", "BOX ", idtmed[1009], par, 3);
126     
127     //     Air 
128     par[0] = 66.3;
129     par[1] = 13;                 //Original Settings
130     par[2] = 68.35;
131     gMC->Gsvolu("SRIC", "BOX ", idtmed[1000], par, 3);
132     
133     //    Air 2 (cutting the lower part of the box)
134     
135     par[0] = 1.25;
136     par[1] = 3;                 //Original Settings
137     par[2] = 70.86;
138     gMC->Gsvolu("AIR2", "BOX ", idtmed[1000], par, 3);
139
140     //    Air 3 (cutting the lower part of the box)
141     
142     par[0] = 66.3;
143     par[1] = 3;                 //Original Settings
144     par[2] = 1.2505;
145     gMC->Gsvolu("AIR3", "BOX ", idtmed[1000], par, 3);
146     
147     //     Honeycomb 
148     par[0] = 66.3;
149     par[1] = .188;                 //Original Settings
150     par[2] = 68.35;
151     gMC->Gsvolu("HONE", "BOX ", idtmed[1001], par, 3);
152     
153     //     Aluminium sheet 
154     par[0] = 66.3;
155     par[1] = .025;                 //Original Settings
156     par[2] = 68.35;
157     /*par[0] = 66.5;
158     par[1] = .025;
159     par[2] = 63.1;*/
160     gMC->Gsvolu("ALUM", "BOX ", idtmed[1009], par, 3);
161     
162     //     Quartz 
163     par[0] = geometry->GetQuartzWidth()/2;
164     par[1] = geometry->GetQuartzThickness()/2;
165     par[2] = geometry->GetQuartzLength()/2;
166     gMC->Gsvolu("QUAR", "BOX ", idtmed[1002], par, 3);
167     
168     //     Spacers (cylinders) 
169     par[0] = 0.;
170     par[1] = .5;
171     par[2] = geometry->GetFreonThickness()/2;
172     gMC->Gsvolu("SPAC", "TUBE", idtmed[1002], par, 3);
173     
174     //     Feet (freon slabs supports)
175
176     par[0] = .7;
177     par[1] = .3;
178     par[2] = 1.9;
179     gMC->Gsvolu("FOOT", "BOX", idtmed[1009], par, 3);
180
181     //     Opaque quartz 
182     par[0] = geometry->GetQuartzWidth()/2;
183     par[1] = .2;
184     par[2] = geometry->GetQuartzLength()/2;
185     gMC->Gsvolu("OQUA", "BOX ", idtmed[1007], par, 3);
186   
187     //     Frame of opaque quartz
188     par[0] = geometry->GetOuterFreonWidth()/2;
189     par[1] = geometry->GetFreonThickness()/2;
190     par[2] = geometry->GetOuterFreonLength()/2; 
191     gMC->Gsvolu("OQF1", "BOX ", idtmed[1007], par, 3);
192
193     par[0] = geometry->GetInnerFreonWidth()/2;
194     par[1] = geometry->GetFreonThickness()/2;
195     par[2] = geometry->GetInnerFreonLength()/2; 
196     gMC->Gsvolu("OQF2", "BOX ", idtmed[1007], par, 3);
197     
198     
199     //     Freon 
200     par[0] = geometry->GetOuterFreonWidth()/2 - oqua_thickness;
201     par[1] = geometry->GetFreonThickness()/2;
202     par[2] = geometry->GetOuterFreonLength()/2 - 2*oqua_thickness; 
203     gMC->Gsvolu("FRE1", "BOX ", idtmed[1003], par, 3);
204
205     par[0] = geometry->GetInnerFreonWidth()/2 - oqua_thickness;
206     par[1] = geometry->GetFreonThickness()/2;
207     par[2] = geometry->GetInnerFreonLength()/2 - 2*oqua_thickness; 
208     gMC->Gsvolu("FRE2", "BOX ", idtmed[1003], par, 3);
209     
210     //     Methane 
211     par[0] = csi_width/2;
212     par[1] = geometry->GetGapThickness()/2;
213     par[2] = csi_length/2;
214     gMC->Gsvolu("META", "BOX ", idtmed[1004], par, 3);
215     
216     //     Methane gap 
217     par[0] = csi_width/2;
218     par[1] = geometry->GetProximityGapThickness()/2;
219     par[2] = csi_length/2;
220     gMC->Gsvolu("GAP ", "BOX ", idtmed[1008], par, 3);
221     
222     //     CsI photocathode 
223     par[0] = csi_width/2;
224     par[1] = .25;
225     par[2] = csi_length/2;
226     gMC->Gsvolu("CSI ", "BOX ", idtmed[1005], par, 3);
227     
228     //     Anode grid 
229     par[0] = 0.;
230     par[1] = .001;
231     par[2] = 20.;
232     gMC->Gsvolu("GRID", "TUBE", idtmed[1006], par, 3);
233
234     // Wire supports
235     // Bar of metal
236     
237     par[0] = csi_width/2;
238     par[1] = 1.05;
239     par[2] = 1.05;
240     gMC->Gsvolu("WSMe", "BOX ", idtmed[1009], par, 3);
241
242     // Ceramic pick up (base)
243     
244     par[0] =  csi_width/2;
245     par[1] = .25;
246     par[2] = 1.05;
247     gMC->Gsvolu("WSG1", "BOX ", idtmed[1010], par, 3);
248
249     // Ceramic pick up (head)
250
251     par[0] = csi_width/2;
252     par[1] = .1;
253     par[2] = .1;
254     gMC->Gsvolu("WSG2", "BOX ", idtmed[1010], par, 3);
255
256     // Aluminium supports for methane and CsI
257     // Short bar
258
259     par[0] = csi_width/2;
260     par[1] = geometry->GetGapThickness()/2 + .25;
261     par[2] = (68.35 - csi_length/2)/2;
262     gMC->Gsvolu("SMSH", "BOX", idtmed[1009], par, 3);
263     
264     // Long bar
265
266     par[0] = (66.3 - csi_width/2)/2;
267     par[1] = geometry->GetGapThickness()/2 + .25;
268     par[2] = csi_length/2 + 68.35 - csi_length/2;
269     gMC->Gsvolu("SMLG", "BOX", idtmed[1009], par, 3);
270     
271     // Aluminium supports for freon
272     // Short bar
273
274     par[0] = geometry->GetQuartzWidth()/2;
275     par[1] = .3;
276     par[2] = (68.35 - geometry->GetQuartzLength()/2)/2;
277     gMC->Gsvolu("SFSH", "BOX", idtmed[1009], par, 3);
278     
279     // Long bar
280
281     par[0] = (66.3 - geometry->GetQuartzWidth()/2)/2;
282     par[1] = .3;
283     par[2] = geometry->GetQuartzLength()/2 + 68.35 - geometry->GetQuartzLength()/2;
284     gMC->Gsvolu("SFLG", "BOX", idtmed[1009], par, 3);
285     
286     // PCB backplane
287     
288     par[0] = csi_width/2;
289     par[1] = .25;
290     par[2] = csi_length/4 -.5025;
291     gMC->Gsvolu("PCB ", "BOX", idtmed[1011], par, 3);
292
293     
294     // Backplane supports
295
296     // Aluminium slab
297     
298     par[0] = 33.15;
299     par[1] = 2;
300     par[2] = 21.65;
301     gMC->Gsvolu("BACK", "BOX", idtmed[1009], par, 3);
302     
303     // Big hole
304     
305     par[0] = 9.05;
306     par[1] = 2;
307     par[2] = 4.4625;
308     gMC->Gsvolu("BKHL", "BOX", idtmed[1000], par, 3);
309
310     // Small hole
311     
312     par[0] = 5.7;
313     par[1] = 2;
314     par[2] = 4.4625;
315     gMC->Gsvolu("BKHS", "BOX", idtmed[1000], par, 3);
316
317     // Place holes inside backplane support
318
319     gMC->Gspos("BKHS", 1, "BACK", .8 + 5.7,0., .6 + 4.4625, 0, "ONLY");
320     gMC->Gspos("BKHS", 2, "BACK", -.8 - 5.7,0., .6 + 4.4625, 0, "ONLY");
321     gMC->Gspos("BKHS", 3, "BACK", .8 + 5.7,0., -.6 - 4.4625, 0, "ONLY");
322     gMC->Gspos("BKHS", 4, "BACK", -.8 - 5.7,0., -.6 - 4.4625, 0, "ONLY");
323     gMC->Gspos("BKHS", 5, "BACK", .8 + 5.7,0., .6 + 8.925 + 1.2 + 4.4625, 0, "ONLY");
324     gMC->Gspos("BKHS", 6, "BACK", -.8 - 5.7,0., .6 + 8.925 + 1.2 + 4.4625, 0, "ONLY");
325     gMC->Gspos("BKHS", 7, "BACK", .8 + 5.7,0., -.6 - 8.925 - 1.2 - 4.4625, 0, "ONLY");
326     gMC->Gspos("BKHS", 8, "BACK", -.8 - 5.7,0., -.6 - 8.925 - 1.2 - 4.4625, 0, "ONLY");
327     gMC->Gspos("BKHL", 1, "BACK", .8 + 11.4 + 1.6 + 9.05, 0., .6 + 4.4625, 0, "ONLY");
328     gMC->Gspos("BKHL", 2, "BACK", -.8 - 11.4 - 1.6 - 9.05, 0., .6 + 4.4625, 0, "ONLY");
329     gMC->Gspos("BKHL", 3, "BACK", .8 + 11.4 + 1.6 + 9.05, 0., -.6 - 4.4625, 0, "ONLY");
330     gMC->Gspos("BKHL", 4, "BACK", -.8 - 11.4 - 1.6 - 9.05, 0., -.6 - 4.4625, 0, "ONLY");
331     gMC->Gspos("BKHL", 5, "BACK", .8 + 11.4+ 1.6 + 9.05, 0., .6 + 8.925 + 1.2 + 4.4625, 0, "ONLY");
332     gMC->Gspos("BKHL", 6, "BACK", -.8 - 11.4 - 1.6 - 9.05, 0., .6 + 8.925 + 1.2 + 4.4625, 0, "ONLY");
333     gMC->Gspos("BKHL", 7, "BACK", .8 + 11.4 + 1.6 + 9.05, 0., -.6 - 8.925 - 1.2 - 4.4625, 0, "ONLY");
334     gMC->Gspos("BKHL", 8, "BACK", -.8 - 11.4 - 1.6 - 9.05, 0., -.6 - 8.925 - 1.2 - 4.4625, 0, "ONLY");
335
336     
337   
338     // --- Places the detectors defined with GSVOLU 
339     //     Place material inside RICH 
340     gMC->Gspos("SRIC", 1, "RICH", 0.,0., 0., 0, "ONLY");
341     gMC->Gspos("AIR2", 1, "RICH", 66.3 + 1.2505, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .6 - .05 - .376 -.5 - 3.35, 0., 0, "ONLY");
342     gMC->Gspos("AIR2", 2, "RICH", -66.3 - 1.2505, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .6 - .05 - .376 -.5 - 3.35, 0., 0, "ONLY");
343     gMC->Gspos("AIR3", 1, "RICH", 0.,  1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .6 - .05 - .376 -.5 - 3.35, -68.35 - 1.25, 0, "ONLY");
344     gMC->Gspos("AIR3", 2, "RICH", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .6 - .05 - .376 -.5 - 3.35,  68.35 + 1.25, 0, "ONLY");
345     
346       
347     gMC->Gspos("ALUM", 1, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .6 - .05 - .376 -.025, 0., 0, "ONLY");
348     gMC->Gspos("HONE", 1, "SRIC", 0., 1.276- geometry->GetGapThickness()/2  - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .6 - .05 - .188, 0., 0, "ONLY");
349     gMC->Gspos("ALUM", 2, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .6 - .025, 0., 0, "ONLY");
350     gMC->Gspos("FOOT", 1, "SRIC", 64.95, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .3, 36.9, 0, "ONLY");
351     gMC->Gspos("FOOT", 2, "SRIC", 21.65, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .3 , 36.9, 0, "ONLY");
352     gMC->Gspos("FOOT", 3, "SRIC", -21.65, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .3, 36.9, 0, "ONLY");
353     gMC->Gspos("FOOT", 4, "SRIC", -64.95, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .3, 36.9, 0, "ONLY");
354     gMC->Gspos("FOOT", 5, "SRIC", 64.95, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .3, -36.9, 0, "ONLY");
355     gMC->Gspos("FOOT", 6, "SRIC", 21.65, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .3, -36.9, 0, "ONLY");
356     gMC->Gspos("FOOT", 7, "SRIC", -21.65, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .3, -36.9, 0, "ONLY");
357     gMC->Gspos("FOOT", 8, "SRIC", -64.95, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .4 - .3, -36.9, 0, "ONLY");
358     gMC->Gspos("OQUA", 1, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()- .2, 0., 0, "ONLY");
359     
360     // Supports placing
361
362     // Methane supports
363     gMC->Gspos("SMLG", 1, "SRIC", csi_width/2 + (66.3 - csi_width/2)/2, 1.276 + .25, 0., 0, "ONLY");
364     gMC->Gspos("SMLG", 2, "SRIC", - csi_width/2 - (66.3 - csi_width/2)/2, 1.276 + .25, 0., 0, "ONLY");
365     gMC->Gspos("SMSH", 1, "SRIC", 0., 1.276 + .25, csi_length/2 + (68.35 - csi_length/2)/2, 0, "ONLY");
366     gMC->Gspos("SMSH", 2, "SRIC", 0., 1.276 + .25, - csi_length/2 - (68.35 - csi_length/2)/2, 0, "ONLY");
367
368     //Freon supports
369
370     Float_t supp_y = 1.276 - geometry->GetGapThickness()/2- geometry->GetQuartzThickness() -geometry->GetFreonThickness() - .2 + .3; //y position of freon supports
371
372     gMC->Gspos("SFLG", 1, "SRIC", geometry->GetQuartzWidth()/2 + (66.3 - geometry->GetQuartzWidth()/2)/2, supp_y, 0., 0, "ONLY");
373     gMC->Gspos("SFLG", 2, "SRIC", - geometry->GetQuartzWidth()/2 - (66.3 - geometry->GetQuartzWidth()/2)/2, supp_y, 0., 0, "ONLY");
374     gMC->Gspos("SFSH", 1, "SRIC", 0., supp_y, geometry->GetQuartzLength()/2 + (68.35 - geometry->GetQuartzLength()/2)/2, 0, "ONLY");
375     gMC->Gspos("SFSH", 2, "SRIC", 0., supp_y, - geometry->GetQuartzLength()/2 - (68.35 - geometry->GetQuartzLength()/2)/2, 0, "ONLY");
376     
377     AliMatrix(idrotm[1019], 0., 0., 90., 0., 90., 90.);
378     
379
380     Int_t nspacers = 30;
381     
382     for (i = 0; i < nspacers/3; i++) {
383         zs = -11.6/2 + (TMath::Abs(nspacers/6) - i) * 12.2;
384         gMC->Gspos("SPAC", i, "FRE1", 10.5, 0., zs, idrotm[1019], "ONLY");  //Original settings 
385     }
386     
387     for (i = nspacers/3; i < (nspacers*2)/3; i++) {
388         zs = -11.6/2 + (nspacers/3 + TMath::Abs(nspacers/6) - i) * 12.2;
389         gMC->Gspos("SPAC", i, "FRE1", 0, 0., zs, idrotm[1019], "ONLY");  //Original settings 
390     }
391     
392     for (i = (nspacers*2)/3; i < nspacers; ++i) {
393         zs = -11.6/2 + ((nspacers*2)/3 + TMath::Abs(nspacers/6) - i) * 12.2;
394         gMC->Gspos("SPAC", i, "FRE1", -10.5, 0., zs, idrotm[1019], "ONLY"); //Original settings  
395     }
396
397     for (i = 0; i < nspacers/3; i++) {
398         zs = -11.6/2 + (TMath::Abs(nspacers/6) - i) * 12.2;
399         gMC->Gspos("SPAC", i, "FRE2", 10.5, 0., zs, idrotm[1019], "ONLY");  //Original settings 
400     }
401     
402     for (i = nspacers/3; i < (nspacers*2)/3; i++) {
403         zs = -11.6/2 + (nspacers/3 + TMath::Abs(nspacers/6) - i) * 12.2;
404         gMC->Gspos("SPAC", i, "FRE2", 0, 0., zs, idrotm[1019], "ONLY");  //Original settings 
405     }
406     
407     for (i = (nspacers*2)/3; i < nspacers; ++i) {
408         zs = -11.6/2 + ((nspacers*2)/3 + TMath::Abs(nspacers/6) - i) * 12.2;
409         gMC->Gspos("SPAC", i, "FRE2", -10.5, 0., zs, idrotm[1019], "ONLY"); //Original settings  
410     }
411
412     
413     gMC->Gspos("FRE1", 1, "OQF1", 0., 0., 0., 0, "ONLY");
414     gMC->Gspos("FRE2", 1, "OQF2", 0., 0., 0., 0, "ONLY");
415     gMC->Gspos("OQF1", 1, "SRIC", geometry->GetOuterFreonWidth()/2 + geometry->GetInnerFreonWidth()/2 + 2, 1.276 - geometry->GetGapThickness()/2- geometry->GetQuartzThickness() -geometry->GetFreonThickness()/2, 0., 0, "ONLY"); //Original settings (31.3)
416     gMC->Gspos("OQF2", 2, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()/2, 0., 0, "ONLY");          //Original settings 
417     gMC->Gspos("OQF1", 3, "SRIC", - (geometry->GetOuterFreonWidth()/2 + geometry->GetInnerFreonWidth()/2) - 2, 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness() - geometry->GetFreonThickness()/2, 0., 0, "ONLY");       //Original settings (-31.3)
418     gMC->Gspos("QUAR", 1, "SRIC", 0., 1.276 - geometry->GetGapThickness()/2 - geometry->GetQuartzThickness()/2, 0., 0, "ONLY");
419     gMC->Gspos("GAP ", 1, "META", 0., geometry->GetGapThickness()/2 - geometry->GetProximityGapThickness()/2 - 0.0001, 0., 0, "ONLY");
420     gMC->Gspos("META", 1, "SRIC", 0., 1.276, 0., 0, "ONLY");
421     gMC->Gspos("CSI ", 1, "SRIC", 0., 1.276 + geometry->GetGapThickness()/2 + .25, 0., 0, "ONLY");
422     printf("CSI pos: %f\n",1.276 + geometry->GetGapThickness()/2 + .25);
423    
424     // Wire support placing
425
426     gMC->Gspos("WSG2", 1, "GAP ", 0., geometry->GetProximityGapThickness()/2 - .1, 0., 0, "ONLY");
427     gMC->Gspos("WSG1", 1, "CSI ", 0., 0., 0., 0, "ONLY");
428     gMC->Gspos("WSMe", 1, "SRIC ", 0., 1.276 + geometry->GetGapThickness()/2 + .5 + 1.05, 0., 0, "ONLY");
429
430     // Backplane placing
431     
432     gMC->Gspos("BACK", 1, "SRIC ", -33.15, 1.276 + geometry->GetGapThickness()/2 + .5 + 2.1 + 2, 43.3, 0, "ONLY");
433     gMC->Gspos("BACK", 2, "SRIC ", 33.15, 1.276 + geometry->GetGapThickness()/2 + .5 + 2.1 + 2 , 43.3, 0, "ONLY");
434     gMC->Gspos("BACK", 3, "SRIC ", -33.15, 1.276 + geometry->GetGapThickness()/2 + .5 + 2.1 + 2, 0., 0, "ONLY");
435     gMC->Gspos("BACK", 4, "SRIC ", 33.15, 1.276 + geometry->GetGapThickness()/2 + .5 + 2.1 + 2, 0., 0, "ONLY");
436     gMC->Gspos("BACK", 5, "SRIC ", 33.15, 1.276 + geometry->GetGapThickness()/2 + .5 + 2.1 + 2, -43.3, 0, "ONLY");
437     gMC->Gspos("BACK", 6, "SRIC ", -33.15, 1.276 + geometry->GetGapThickness()/2 + .5 + 2.1 + 2, -43.3, 0, "ONLY");
438
439     // PCB placing
440     
441     gMC->Gspos("PCB ", 1, "SRIC ", 0.,  1.276 + geometry->GetGapThickness()/2 + .5 + 1.05, csi_width/4 + .5025 + 2.5, 0, "ONLY");
442     gMC->Gspos("PCB ", 2, "SRIC ", 0.,  1.276 + geometry->GetGapThickness()/2 + .5 + 1.05, -csi_width/4 - .5025 - 2.5, 0, "ONLY");
443
444 // Place chambers into mother volume ALIC
445            
446    Double_t dOffset        = geometry->GetOffset() - geometry->GetGapThickness()/2;  // distance from center of mother volume ALIC to methane
447    
448    Double_t dAlpha         = geometry->GetAlphaAngle(); // angle between centers of chambers - y-z plane
449    Double_t dAlphaRad      = dAlpha*kDegrad;
450    
451    Double_t dBeta          = geometry->GetBetaAngle();   // angle between center of chambers - y-x plane
452    Double_t dBetaRad       = dBeta*kDegrad;
453    
454    Double_t dRotAngle      = geometry->GetRotationAngle();     // the whole RICH is to be rotated in x-y plane + means clockwise rotation 
455    Double_t dRotAngleRad   = dRotAngle*kDegrad;
456     
457    
458    TRotMatrix *pRotMatrix; // tmp pointer
459    
460    TVector3 vector(0,dOffset,0); // Position of chamber 2 without rotation
461     
462 // Chamber 0  standalone (no other chambers in this row) 
463    pRotMatrix = new TRotMatrix("rot993","rot993", 0., 0., 0.,0.,0.,0.);
464    const Double_t* r   = pRotMatrix->SetAngles(90., 0., 90.-dAlpha , 90.,  dAlpha, -90.);
465    Double_t* rr  = RotateXY(r, -dRotAngleRad);
466    AliMatrix(idrotm[1000], rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
467    pRotMatrix->SetAngles(rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
468
469    vector.SetXYZ(0,dOffset,0);  vector.RotateX(dAlphaRad); 
470    vector.RotateZ(-dRotAngleRad);
471    
472    gMC->Gspos("RICH",1,"ALIC",vector.X(),vector.Y(),vector.Z(),idrotm[1000], "ONLY");           
473    Chamber(0).SetChamberTransform(vector.X(),vector.Y(),vector.Z(),pRotMatrix);
474   if(GetDebug()) Info("CreateGeometry 0","%8.3f %8.3f %8.3f %8.3f %8.3f %8.3f",rr[0],rr[1],rr[2],rr[3],rr[4],rr[5]);
475   if(GetDebug()) Info("CreateGeometry 0","x=%8.3f y=%8.3f z=%8.3f",vector.X(),vector.Y(),vector.Z());   
476 // Chamber 1   
477    pRotMatrix = new TRotMatrix("rot994","rot994", 0., 0., 0.,0.,0.,0.);
478    r   = pRotMatrix->SetAngles(90., -dBeta, 90., 90.-dBeta,  0., 0.);
479    rr  = RotateXY(r, -dRotAngleRad);
480    AliMatrix(idrotm[1001], rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
481    pRotMatrix->SetAngles(rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
482    vector.SetXYZ(0,dOffset,0);  vector.RotateZ(-dBetaRad); 
483    vector.RotateZ(-dRotAngleRad);
484    
485    gMC->Gspos("RICH",2,"ALIC",vector.X(),vector.Y(),vector.Z(),idrotm[1001], "ONLY");           
486    Chamber(1).SetChamberTransform(vector.X(),vector.Y(),vector.Z(),pRotMatrix);
487   if(GetDebug()) Info("CreateGeometry 1","%8.3f %8.3f %8.3f %8.3f %8.3f %8.3f",rr[0],rr[1],rr[2],rr[3],rr[4],rr[5]);
488   if(GetDebug()) Info("CreateGeometry 1","x=%8.3f y=%8.3f z=%8.3f",vector.X(),vector.Y(),vector.Z());
489 // Chamber 2   the top one with no Alpha-Beta rotation
490    pRotMatrix = new TRotMatrix("rot995","rot995", 0., 0., 0.,0.,0.,0.);
491    r   = pRotMatrix->SetAngles(90., 0., 90., 90.,  0., 0.);
492    rr  = RotateXY(r, -dRotAngleRad);
493    AliMatrix(idrotm[1002], rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
494    pRotMatrix->SetAngles(rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
495    vector.SetXYZ(0,dOffset,0);
496    vector.RotateZ(-dRotAngleRad);
497    gMC->Gspos("RICH",3,"ALIC",vector.X(),vector.Y(),vector.Z(),idrotm[1002], "ONLY");           
498    Chamber(2).SetChamberTransform(vector.X(),vector.Y(),vector.Z(),pRotMatrix);
499   if(GetDebug()) Info("CreateGeometry 2","%8.3f %8.3f %8.3f %8.3f %8.3f %8.3f",rr[0],rr[1],rr[2],rr[3],rr[4],rr[5]);
500   if(GetDebug()) Info("CreateGeometry 2","x=%8.3f y=%8.3f z=%8.3f",vector.X(),vector.Y(),vector.Z());   
501 // Chamber 3
502    pRotMatrix = new TRotMatrix("rot996","rot996", 0., 0., 0.,0.,0.,0.);
503    r   = pRotMatrix->SetAngles(90., dBeta, 90., 90.+dBeta,  0., 0.);
504    rr  = RotateXY(r, -dRotAngleRad);
505    AliMatrix(idrotm[1003], rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
506    pRotMatrix->SetAngles(rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
507    vector.SetXYZ(0,dOffset,0);  vector.RotateZ(dBetaRad); 
508    vector.RotateZ(-dRotAngleRad);
509    
510    gMC->Gspos("RICH",4,"ALIC",vector.X(),vector.Y(),vector.Z(),idrotm[1003], "ONLY");           
511    Chamber(3).SetChamberTransform(vector.X(),vector.Y(),vector.Z(),pRotMatrix);
512   if(GetDebug()) Info("CreateGeometry 3","%8.3f %8.3f %8.3f %8.3f %8.3f %8.3f",rr[0],rr[1],rr[2],rr[3],rr[4],rr[5]);
513   if(GetDebug()) Info("CreateGeometry 3","x=%8.3f y=%8.3f z=%8.3f",vector.X(),vector.Y(),vector.Z());
514 // Chamber 4   
515    pRotMatrix = new TRotMatrix("rot997","rot997", 0., 0., 0.,0.,0.,0.);
516    r   = pRotMatrix->SetAngles(90., 360.-dBeta, 108.2, 90.-dBeta,  18.2, 90.-dBeta);
517    rr  = RotateXY(r, -dRotAngleRad);
518    AliMatrix(idrotm[1004], rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
519    pRotMatrix->SetAngles(rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
520    vector.SetXYZ(0,dOffset,0);  vector.RotateZ(-dBetaRad); vector.RotateX(-dAlphaRad); 
521    vector.RotateZ(-dRotAngleRad);
522    
523    gMC->Gspos("RICH",5,"ALIC",vector.X(),vector.Y(),vector.Z(),idrotm[1004], "ONLY");
524    Chamber(4).SetChamberTransform(vector.X(),vector.Y(),vector.Z(),pRotMatrix);
525   if(GetDebug()) Info("CreateGeometry 4","%8.3f %8.3f %8.3f %8.3f %8.3f %8.3f",rr[0],rr[1],rr[2],rr[3],rr[4],rr[5]);
526   if(GetDebug()) Info("CreateGeometry 4","x=%8.3f y=%8.3f z=%8.3f",vector.X(),vector.Y(),vector.Z());
527 // Chamber 5   
528    pRotMatrix = new TRotMatrix("rot998","rot998", 0., 0., 0.,0.,0.,0.);
529    r   = pRotMatrix->SetAngles(90., 0., 90.+dAlpha, 90.,  dAlpha, 90.);
530    rr  = RotateXY(r, -dRotAngleRad);
531    AliMatrix(idrotm[1005], rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
532    pRotMatrix->SetAngles(rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);   
533    vector.SetXYZ(0,dOffset,0); vector.RotateX(-dAlphaRad); 
534    vector.RotateZ(-dRotAngleRad);
535       
536    gMC->Gspos("RICH",6,"ALIC",vector.X(),vector.Y(),vector.Z(),idrotm[1005], "ONLY");           
537    Chamber(5).SetChamberTransform(vector.X(),vector.Y(),vector.Z(),pRotMatrix);
538   if(GetDebug()) Info("CreateGeometry 5","%8.3f %8.3f %8.3f %8.3f %8.3f %8.3f",rr[0],rr[1],rr[2],rr[3],rr[4],rr[5]);
539   if(GetDebug()) Info("CreateGeometry 5","x=%8.3f y=%8.3f z=%8.3f",vector.X(),vector.Y(),vector.Z());
540 // Chamber 6          
541    pRotMatrix = new TRotMatrix("rot999","rot999", 0., 0., 0.,0.,0.,0.);
542    r   = pRotMatrix->SetAngles(90., dBeta, 108.2, 90.+dBeta,  18.2, 90.+dBeta);
543    rr  = RotateXY(r, -dRotAngleRad);
544    AliMatrix(idrotm[1006], rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
545    pRotMatrix->SetAngles(rr[0], rr[1], rr[2], rr[3], rr[4], rr[5]);
546    vector.SetXYZ(0,dOffset,0);  vector.RotateZ(dBetaRad); vector.RotateX(-dAlphaRad); 
547    vector.RotateZ(-dRotAngleRad);
548       
549    gMC->Gspos("RICH",7,"ALIC",vector.X(),vector.Y(),vector.Z(),idrotm[1006], "ONLY");
550    Chamber(6).SetChamberTransform(vector.X(),vector.Y(),vector.Z(),pRotMatrix);
551   if(GetDebug()) Info("CreateGeometry 6","%8.3f %8.3f %8.3f %8.3f %8.3f %8.3f",rr[0],rr[1],rr[2],rr[3],rr[4],rr[5]);
552   if(GetDebug()) Info("CreateGeometry 6","x=%8.3f y=%8.3f z=%8.3f",vector.X(),vector.Y(),vector.Z());
553       
554 }//void AliRICHv3::CreateGeometry()
555 //______________________________________________________________________________
556 void AliRICHv3::Init()
557 {//Makes nothing for a while   
558   if(GetDebug())Info("Init","Start.");
559   if(GetDebug())Info("Init","Stop.");    
560 }
561 //______________________________________________________________________________
562 void AliRICHv3::BuildGeometry()    
563 {//Provides geometry structure for event display (ROOT TNode tree)
564   if(GetDebug())Info("BuildGeometry","Start.");
565   
566     TNode *node, *subnode, *top;
567     
568     const int kColorRICH = kRed;
569     //
570     top=gAlice->GetGeometry()->GetNode("alice");
571
572     AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH"); 
573     AliRICHChamber*       iChamber;
574     AliRICHGeometry*  geometry;
575  
576     iChamber = &(pRICH->Chamber(0));
577     AliRICHSegmentationV1* segmentation=(AliRICHSegmentationV1*) iChamber->GetSegmentationModel();
578     geometry=iChamber->GetGeometryModel();
579     
580     new TBRIK("S_RICH","S_RICH","void",71.09999,11.5,73.15);
581
582     Float_t padplane_width = segmentation->GetPadPlaneWidth();
583     Float_t padplane_length = segmentation->GetPadPlaneLength();
584
585
586     new TBRIK("PHOTO","PHOTO","void", padplane_width/2,.1,padplane_length/2);
587
588 // Chamber 0             
589     top->cd();
590     node = new TNode("RICH1","RICH1","S_RICH",Chamber(0).GetX(),Chamber(0).GetY(),Chamber(0).GetZ(),"rot993");
591     node->SetLineColor(kColorRICH);
592     node->cd();
593     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
594     subnode->SetLineColor(kGreen);
595     fNodes->Add(subnode);
596     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,padplane_length/2 + segmentation->DeadZone()/2,"");
597     subnode->SetLineColor(kGreen);
598     fNodes->Add(subnode);
599     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
600     subnode->SetLineColor(kGreen);
601     fNodes->Add(subnode);
602     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
603     subnode->SetLineColor(kGreen);
604     fNodes->Add(subnode);
605     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,-padplane_length/2 - segmentation->DeadZone()/2,"");
606     subnode->SetLineColor(kGreen);
607     fNodes->Add(subnode);
608     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
609     subnode->SetLineColor(kGreen);
610     fNodes->Add(subnode);
611     fNodes->Add(node);
612
613 // Chamber 1
614     top->cd(); 
615     node = new TNode("RICH2","RICH2","S_RICH",Chamber(1).GetX(),Chamber(1).GetY(),Chamber(1).GetZ(),"rot994");
616     node->SetLineColor(kColorRICH);
617     node->cd();
618     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
619     subnode->SetLineColor(kGreen);
620     fNodes->Add(subnode);
621     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,padplane_length/2 + segmentation->DeadZone()/2,"");
622     subnode->SetLineColor(kGreen);
623     fNodes->Add(subnode);
624     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
625     subnode->SetLineColor(kGreen);
626     fNodes->Add(subnode);
627     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
628     subnode->SetLineColor(kGreen);
629     fNodes->Add(subnode);
630     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,-padplane_length/2 - segmentation->DeadZone()/2,"");
631     subnode->SetLineColor(kGreen);
632     fNodes->Add(subnode);
633     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
634     subnode->SetLineColor(kGreen);
635     fNodes->Add(subnode);
636     fNodes->Add(node);
637
638 // Chamber 2
639     top->cd();
640     node = new TNode("RICH3","RICH3","S_RICH",Chamber(2).GetX(),Chamber(2).GetY(),Chamber(2).GetZ(),"rot995");
641     node->SetLineColor(kColorRICH);
642     node->cd();
643     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
644     subnode->SetLineColor(kGreen);
645     fNodes->Add(subnode);
646     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,padplane_length/2 + segmentation->DeadZone()/2,"");
647     subnode->SetLineColor(kGreen);
648     fNodes->Add(subnode);
649     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
650     subnode->SetLineColor(kGreen);
651     fNodes->Add(subnode);
652     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
653     subnode->SetLineColor(kGreen);
654     fNodes->Add(subnode);
655     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,-padplane_length/2 - segmentation->DeadZone()/2,"");
656     subnode->SetLineColor(kGreen);
657     fNodes->Add(subnode);
658     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
659     subnode->SetLineColor(kGreen);
660     fNodes->Add(subnode);
661     fNodes->Add(node);
662     
663 // Chamber 3
664     top->cd();
665     node = new TNode("RICH4","RICH4","S_RICH",Chamber(3).GetX(),Chamber(3).GetY(),Chamber(3).GetZ(),"rot996");
666     node->SetLineColor(kColorRICH);
667     node->cd();
668     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
669     subnode->SetLineColor(kGreen);
670     fNodes->Add(subnode);
671     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,padplane_length/2 + segmentation->DeadZone()/2,"");
672     subnode->SetLineColor(kGreen);
673     fNodes->Add(subnode);
674     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
675     subnode->SetLineColor(kGreen);
676     fNodes->Add(subnode);
677     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
678     subnode->SetLineColor(kGreen);
679     fNodes->Add(subnode);
680     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,-padplane_length/2 - segmentation->DeadZone()/2,"");
681     subnode->SetLineColor(kGreen);
682     fNodes->Add(subnode);
683     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
684     subnode->SetLineColor(kGreen);
685     fNodes->Add(subnode);
686     fNodes->Add(node);
687
688 // Chamber 4
689     top->cd();
690     node = new TNode("RICH5","RICH5","S_RICH",Chamber(4).GetX(),Chamber(4).GetY(),Chamber(4).GetZ(),"rot997");
691     node->SetLineColor(kColorRICH);
692     node->cd();
693     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
694     subnode->SetLineColor(kGreen);
695     fNodes->Add(subnode);
696     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,padplane_length/2 + segmentation->DeadZone()/2,"");
697     subnode->SetLineColor(kGreen);
698     fNodes->Add(subnode);
699     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
700     subnode->SetLineColor(kGreen);
701     fNodes->Add(subnode);
702     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
703     subnode->SetLineColor(kGreen);
704     fNodes->Add(subnode);
705     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,-padplane_length/2 - segmentation->DeadZone()/2,"");
706     subnode->SetLineColor(kGreen);
707     fNodes->Add(subnode);
708     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
709     subnode->SetLineColor(kGreen);
710     fNodes->Add(subnode);
711     fNodes->Add(node);
712
713 // Chamber 5
714     top->cd();
715     node = new TNode("RICH6","RICH6","S_RICH",Chamber(5).GetX(),Chamber(5).GetY(),Chamber(5).GetZ(),"rot998");
716     node->SetLineColor(kColorRICH);
717     fNodes->Add(node);node->cd();
718     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
719     subnode->SetLineColor(kGreen);
720     fNodes->Add(subnode);
721     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,padplane_length/2 + segmentation->DeadZone()/2,"");
722     subnode->SetLineColor(kGreen);
723     fNodes->Add(subnode);
724     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
725     subnode->SetLineColor(kGreen);
726     fNodes->Add(subnode);
727     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
728     subnode->SetLineColor(kGreen);
729     fNodes->Add(subnode);
730     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,-padplane_length/2 - segmentation->DeadZone()/2,"");
731     subnode->SetLineColor(kGreen);
732     fNodes->Add(subnode);
733     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
734     subnode->SetLineColor(kGreen);
735     fNodes->Add(subnode);
736
737 // Chamber 6
738     top->cd();
739     node = new TNode("RICH7","RICH7","S_RICH",Chamber(6).GetX(),Chamber(6).GetY(),Chamber(6).GetZ(),"rot999");
740     node->SetLineColor(kColorRICH);
741     node->cd();
742     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
743     subnode->SetLineColor(kGreen);
744     fNodes->Add(subnode);
745     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,padplane_length/2 + segmentation->DeadZone()/2,"");
746     subnode->SetLineColor(kGreen);
747     fNodes->Add(subnode);
748     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,padplane_length/2 + segmentation->DeadZone()/2,"");
749     subnode->SetLineColor(kGreen);
750     fNodes->Add(subnode);
751     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",padplane_width + segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
752     subnode->SetLineColor(kGreen);
753     fNodes->Add(subnode);
754     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,-padplane_length/2 - segmentation->DeadZone()/2,"");
755     subnode->SetLineColor(kGreen);
756     fNodes->Add(subnode);
757     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-padplane_width - segmentation->DeadZone(),5,-padplane_length/2 - segmentation->DeadZone()/2,"");
758     subnode->SetLineColor(kGreen);
759     fNodes->Add(subnode);
760     fNodes->Add(node); 
761   if(GetDebug())Info("BuildGeometry","Stop.");    
762 }//AliRICHv3::BuildGeometry()
763 //______________________________________________________________________________
764 Double_t* AliRICHv3::RotateXY(const Double_t* r, Double_t a)
765 {
766     // Rotatation in xy-plane
767     // by angle a
768     // The resulting rotation matrix is given back in the G3 notation. 
769     Double_t* rr = new Double_t[6];
770     Double_t m[9];
771     Int_t i,j,k;
772     
773     for (i = 0; i < 3; i++) {
774         j = 3*i;
775         m[j]   = r[j] * TMath::Cos(a) - r[j+1] * TMath::Sin(a);
776         m[j+1] = r[j] * TMath::Sin(a) + r[j+1] * TMath::Cos(a);
777         m[j+2] = r[j+2];
778     }
779     
780     for (i = 0; i < 3; i++) {
781             j = 3*i;
782             k = 2*i;
783             rr[k]    = TMath::ACos(m[j+2])        * kRaddeg;
784             rr[k+1]  = TMath::ATan2(m[j+1], m[j]) * kRaddeg;
785     }
786     return rr;
787 }//Double_t* AliRICHv3::RotateXY(const Double_t* r, Double_t a)
788 //______________________________________________________________________________
789 void AliRICHv3::StepManager()
790 {//Full Step Manager
791
792     Int_t          copy, id;
793     static Int_t   idvol;
794     static Int_t   vol[2];
795     Int_t          ipart;
796     static Float_t hits[22];
797     static Float_t ckovData[19];
798     TLorentzVector position;
799     TLorentzVector momentum;
800     Float_t        pos[3];
801     Float_t        mom[4];
802     Float_t        localPos[3];
803     Float_t        localMom[4];
804     Float_t        localTheta,localPhi;
805     Float_t        theta,phi;
806     Float_t        destep, step;
807     Double_t        ranf[2];
808     Int_t          nPads=1;
809     Float_t        coscerenkov;
810     static Float_t eloss, xhit, yhit, tlength;
811     const  Float_t kBig=1.e10;
812        
813     TClonesArray &lhits = *fHits;
814     TParticle *current = (TParticle*)(*gAlice->Particles())[gAlice->GetCurrentTrackNumber()];
815
816  //if (current->Energy()>1)
817    //{
818         
819     // Only gas gap inside chamber
820     // Tag chambers and record hits when track enters 
821     
822  
823     id=gMC->CurrentVolID(copy);
824     idvol = copy-1;
825     Float_t cherenkovLoss=0;
826     //gAlice->KeepTrack(gAlice->GetCurrentTrackNumber());
827     
828     gMC->TrackPosition(position);
829     pos[0]=position(0);
830     pos[1]=position(1);
831     pos[2]=position(2);
832     //bzero((char *)ckovData,sizeof(ckovData)*19);
833     ckovData[1] = pos[0];                 // X-position for hit
834     ckovData[2] = pos[1];                 // Y-position for hit
835     ckovData[3] = pos[2];                 // Z-position for hit
836     ckovData[6] = 0;                      // dummy track length
837     //ckovData[11] = gAlice->GetCurrentTrackNumber();
838     
839     //printf("\n+++++++++++\nTrack: %d\n++++++++++++\n",gAlice->GetCurrentTrackNumber());
840
841     //AliRICH *RICH = (AliRICH *) gAlice->GetDetector("RICH"); 
842     
843     /********************Store production parameters for Cerenkov photons************************/ 
844 //is it a Cerenkov photon? 
845     if (gMC->TrackPid() == 50000050) { 
846
847       //if (gMC->VolId("GAP ")==gMC->CurrentVolID(copy))
848         //{                    
849           Float_t ckovEnergy = current->Energy();
850           //energy interval for tracking
851           if  (ckovEnergy > 5.6e-09 && ckovEnergy < 7.8e-09 )       
852             //if (ckovEnergy > 0)
853             {
854               if (gMC->IsTrackEntering()){        //is track entering?
855                 //printf("Track entered (1)\n");
856                 if (gMC->VolId("FRE1")==gMC->CurrentVolID(copy) || gMC->VolId("FRE2")==gMC->CurrentVolID(copy))
857                   {                                                          //is it in freo?
858                     if (gMC->IsNewTrack()){                          //is it the first step?
859                       //printf("I'm in!\n");
860                       Int_t mother = current->GetFirstMother(); 
861                       
862                       //printf("Second Mother:%d\n",current->GetSecondMother());
863                       
864                       ckovData[10] = mother;
865                       ckovData[11] = gAlice->GetCurrentTrackNumber();
866                       ckovData[12] = 1;             //Media where photon was produced 1->Freon, 2->Quarz
867                       //printf("Produced in FREO\n");
868                       fCkovNumber++;
869                       fFreonProd=1;
870                       //printf("Index: %d\n",fCkovNumber);
871                     }    //first step question
872                   }        //freo question
873                 
874                 if (gMC->IsNewTrack()){                                  //is it first step?
875                   if (gMC->VolId("QUAR")==gMC->CurrentVolID(copy))             //is it in quarz?
876                     {
877                       ckovData[12] = 2;
878                       //printf("Produced in QUAR\n");
879                     }    //quarz question
880                 }        //first step question
881                 
882                 //printf("Before %d\n",fFreonProd);
883               }   //track entering question
884               
885               if (ckovData[12] == 1)                                        //was it produced in Freon?
886                 //if (fFreonProd == 1)
887                 {
888                   if (gMC->IsTrackEntering()){                                     //is track entering?
889                     //printf("Track entered (2)\n");
890                     //printf("Current volume (should be META): %s\n",gMC->CurrentVolName());
891                     //printf("VolId: %d, CurrentVolID: %d\n",gMC->VolId("META"),gMC->CurrentVolID(copy));
892                     if (gMC->VolId("META")==gMC->CurrentVolID(copy))                //is it in gap?      
893                       {
894                         //printf("Got in META\n");
895                         gMC->TrackMomentum(momentum);
896                         mom[0]=momentum(0);
897                         mom[1]=momentum(1);
898                         mom[2]=momentum(2);
899                         mom[3]=momentum(3);
900                         
901                         gMC->Gmtod(mom,localMom,2);
902                         Float_t cophi = TMath::Cos(TMath::ATan2(localMom[0], localMom[1]));
903                         Float_t t = (1. - .025 / cophi) * (1. - .05 /  cophi);
904                         /**************** Photons lost in second grid have to be calculated by hand************/ 
905                         gMC->GetRandom()->RndmArray(1,ranf);
906                         if (ranf[0] > t) {
907                           gMC->StopTrack();
908                           ckovData[13] = 5;
909                           AddCerenkov(gAlice->GetCurrentTrackNumber(),vol,ckovData);
910                           //printf("Added One (1)!\n");
911                           //printf("Lost one in grid\n");
912                         }
913                         /**********************************************************************************/
914                       }    //gap
915                     
916                     //printf("Current volume (should be CSI) (1): %s\n",gMC->CurrentVolName());
917                     //printf("VolId: %d, CurrentVolID: %d\n",gMC->VolId("CSI "),gMC->CurrentVolID(copy));
918                     if (gMC->VolId("CSI ")==gMC->CurrentVolID(copy))             //is it in csi?      
919                       {
920                         //printf("Got in CSI\n");
921                         gMC->TrackMomentum(momentum);
922                         mom[0]=momentum(0);
923                         mom[1]=momentum(1);
924                         mom[2]=momentum(2);
925                         mom[3]=momentum(3);
926
927                         gMC->Gmtod(mom,localMom,2);
928                         /********* Photons lost by Fresnel reflection have to be calculated by hand********/ 
929                         /***********************Cerenkov phtons (always polarised)*************************/
930                         Double_t localTc = localMom[0]*localMom[0]+localMom[2]*localMom[2];
931                         Double_t localRt = TMath::Sqrt(localTc);
932                         localTheta   = Float_t(TMath::ATan2(localRt,Double_t(localMom[1])));
933                         Double_t cotheta = TMath::Abs(cos(localTheta));
934                         Float_t t = Fresnel(ckovEnergy*1e9,cotheta,1);
935                             gMC->GetRandom()->RndmArray(1,ranf);
936                             if (ranf[0] < t) {
937                               gMC->StopTrack();
938                               ckovData[13] = 6;
939                               AddCerenkov(gAlice->GetCurrentTrackNumber(),vol,ckovData);
940                                 
941                               //printf("Added One (2)!\n");
942                               //printf("Lost by Fresnel\n");
943                             }
944                             /**********************************************************************************/
945                       }
946                   } //track entering?
947                   
948                   
949                   /********************Evaluation of losses************************/
950                   /******************still in the old fashion**********************/
951                   
952                   TArrayI procs;
953                   Int_t i1 = gMC->StepProcesses(procs);            //number of physics mechanisms acting on the particle
954                   for (Int_t i = 0; i < i1; ++i) {
955                     //        Reflection loss 
956                     if (procs[i] == kPLightReflection) {        //was it reflected
957                       ckovData[13]=10;
958                       if (gMC->VolId("FRE1")==gMC->CurrentVolID(copy) || gMC->VolId("FRE2")==gMC->CurrentVolID(copy)) 
959                         ckovData[13]=1;
960                       if (gMC->CurrentVolID(copy) == gMC->VolId("QUAR")) 
961                         ckovData[13]=2;
962                       //gMC->StopTrack();
963                       //AddCerenkov(gAlice->GetCurrentTrackNumber(),vol,ckovData);
964                     } //reflection question
965                      
966                     //        Absorption loss 
967                     else if (procs[i] == kPLightAbsorption) {              //was it absorbed?
968                       //printf("Got in absorption\n");
969                       ckovData[13]=20;
970                       if (gMC->VolId("FRE1")==gMC->CurrentVolID(copy) || gMC->VolId("FRE2")==gMC->CurrentVolID(copy)) 
971                         ckovData[13]=11;
972                       if (gMC->CurrentVolID(copy) == gMC->VolId("QUAR")) 
973                         ckovData[13]=12;
974                       if (gMC->CurrentVolID(copy) == gMC->VolId("META")) 
975                         ckovData[13]=13;
976                       if (gMC->CurrentVolID(copy) == gMC->VolId("GAP ")) 
977                         ckovData[13]=13;
978                       
979                       if (gMC->CurrentVolID(copy) == gMC->VolId("SRIC")) 
980                         ckovData[13]=15;
981                       
982                       //        CsI inefficiency 
983                       if (gMC->CurrentVolID(copy) == gMC->VolId("CSI ")) {
984                         ckovData[13]=16;
985                       }
986                       gMC->StopTrack();
987                       AddCerenkov(gAlice->GetCurrentTrackNumber(),vol,ckovData);
988                       //printf("Added One (3)!\n");
989                       //printf("Added cerenkov %d\n",fCkovNumber);
990                     } //absorption question 
991                     
992                     
993                     //        Photon goes out of tracking scope 
994                     else if (procs[i] == kPStop) {                 //is it below energy treshold?
995                       ckovData[13]=21;
996                       gMC->StopTrack();
997                       AddCerenkov(gAlice->GetCurrentTrackNumber(),vol,ckovData);
998                       //printf("Added One (4)!\n");
999                     }   // energy treshold question         
1000                   }  //number of mechanisms cycle
1001                   /**********************End of evaluation************************/
1002                 } //freon production question
1003             } //energy interval question
1004         //}//inside the proximity gap question
1005     } //cerenkov photon question
1006       
1007     /**************************************End of Production Parameters Storing*********************/ 
1008     
1009     
1010     /*******************************Treat photons that hit the CsI (Ckovs and Feedbacks)************/ 
1011     
1012     if (gMC->TrackPid() == 50000050 || gMC->TrackPid() == 50000051) {
1013       //printf("Cerenkov\n");
1014       
1015       //if (gMC->TrackPid() == 50000051)
1016         //printf("Tracking a feedback\n");
1017       
1018       if (gMC->VolId("CSI ")==gMC->CurrentVolID(copy))
1019         {
1020           //printf("Current volume (should be CSI) (2): %s\n",gMC->CurrentVolName());
1021           //printf("VolId: %d, CurrentVolID: %d\n",gMC->VolId("CSI "),gMC->CurrentVolID(copy));
1022           //printf("Got in CSI\n");
1023           //printf("Tracking a %d\n",gMC->TrackPid());
1024           if (gMC->Edep() > 0.){
1025                 gMC->TrackPosition(position);
1026                 gMC->TrackMomentum(momentum);
1027                 pos[0]=position(0);
1028                 pos[1]=position(1);
1029                 pos[2]=position(2);
1030                 mom[0]=momentum(0);
1031                 mom[1]=momentum(1);
1032                 mom[2]=momentum(2);
1033                 mom[3]=momentum(3);
1034                 Double_t tc = mom[0]*mom[0]+mom[1]*mom[1];
1035                 Double_t rt = TMath::Sqrt(tc);
1036                 theta   = Float_t(TMath::ATan2(rt,Double_t(mom[2])))*kRaddeg;
1037                 phi     = Float_t(TMath::ATan2(Double_t(mom[1]),Double_t(mom[0])))*kRaddeg;
1038                 
1039                 gMC->CurrentVolOffID(2,copy);
1040                 vol[0]=copy;
1041                 idvol=vol[0]-1;
1042                 
1043
1044                 gMC->Gmtod(pos,localPos,1);
1045
1046                 //Chamber(idvol).GlobaltoLocal(pos,localPos);
1047                                                                     
1048                 gMC->Gmtod(mom,localMom,2);
1049
1050                 //Chamber(idvol).GlobaltoLocal(mom,localMom);
1051                 
1052                 gMC->CurrentVolOffID(2,copy);
1053                 vol[0]=copy;
1054                 idvol=vol[0]-1;
1055
1056                 //Int_t sector=((AliRICHChamber*) (*fChambers)[idvol])
1057                         //->Sector(localPos[0], localPos[2]);
1058                 //printf("Sector:%d\n",sector);
1059
1060                 /*if (gMC->TrackPid() == 50000051){
1061                   fFeedbacks++;
1062                   printf("Feedbacks:%d\n",fFeedbacks);
1063                 }*/     
1064                 
1065         //PH            ((AliRICHChamber*) (*fChambers)[idvol])
1066                 ((AliRICHChamber*)fChambers->At(idvol))
1067                     ->SigGenInit(localPos[0], localPos[2], localPos[1]);
1068                 if(idvol<kNCH) {        
1069                     ckovData[0] = gMC->TrackPid();        // particle type
1070                     ckovData[1] = pos[0];                 // X-position for hit
1071                     ckovData[2] = pos[1];                 // Y-position for hit
1072                     ckovData[3] = pos[2];                 // Z-position for hit
1073                     ckovData[4] = theta;                      // theta angle of incidence
1074                     ckovData[5] = phi;                      // phi angle of incidence 
1075                     ckovData[8] = (Float_t) fNsdigits;      // first sdigit
1076                     ckovData[9] = -1;                       // last pad hit
1077                     ckovData[13] = 4;                       // photon was detected
1078                     ckovData[14] = mom[0];
1079                     ckovData[15] = mom[1];
1080                     ckovData[16] = mom[2];
1081                     
1082                     destep = gMC->Edep();
1083                     gMC->SetMaxStep(kBig);
1084                     cherenkovLoss  += destep;
1085                     ckovData[7]=cherenkovLoss;
1086                     
1087                     //nPads = Hits2SDigits(localPos[0],localPos[2],cherenkovLoss,idvol,kCerenkov);//for photons in CsI kir
1088                                     
1089                     if (fNsdigits > (Int_t)ckovData[8]) {
1090                         ckovData[8]= ckovData[8]+1;
1091                         ckovData[9]= (Float_t) fNsdigits;
1092                     }
1093
1094                     //printf("Cerenkov loss: %f\n", cherenkovLoss);
1095
1096                     ckovData[17] = nPads;
1097                     //printf("nPads:%d",nPads);
1098                     
1099                     //TClonesArray *Hits = RICH->Hits();
1100                     AliRICHhit *mipHit =  (AliRICHhit*) (fHits->UncheckedAt(0));
1101                     if (mipHit)
1102                       {
1103                         mom[0] = current->Px();
1104                         mom[1] = current->Py();
1105                         mom[2] = current->Pz();
1106                         Float_t mipPx = mipHit->MomX();
1107                         Float_t mipPy = mipHit->MomY();
1108                         Float_t mipPz = mipHit->MomZ();
1109                         
1110                         Float_t r = mom[0]*mom[0] + mom[1]*mom[1] + mom[2]*mom[2];
1111                         Float_t rt = TMath::Sqrt(r);
1112                         Float_t mipR = mipPx*mipPx + mipPy*mipPy + mipPz*mipPz; 
1113                         Float_t mipRt = TMath::Sqrt(mipR);
1114                         if ((rt*mipRt) > 0)
1115                           {
1116                             coscerenkov = (mom[0]*mipPx + mom[1]*mipPy + mom[2]*mipPz)/(rt*mipRt);
1117                           }
1118                         else
1119                           {
1120                             coscerenkov = 0;
1121                           }
1122                         Float_t cherenkov = TMath::ACos(coscerenkov);
1123                         ckovData[18]=cherenkov;
1124                       }
1125                     //if (sector != -1)
1126                     //{
1127                     AddHit(gAlice->GetCurrentTrackNumber(),vol,ckovData);
1128                     AddCerenkov(gAlice->GetCurrentTrackNumber(),vol,ckovData);
1129                     //printf("Added One (5)!\n");
1130                     //}
1131                 }
1132             }
1133         }
1134     }
1135     
1136     /***********************************************End of photon hits*********************************************/
1137     
1138
1139     /**********************************************Charged particles treatment*************************************/
1140
1141     else if (gMC->TrackCharge())
1142     //else if (1 == 1)
1143       {
1144 //If MIP
1145         /*if (gMC->IsTrackEntering())
1146           {                
1147             hits[13]=20;//is track entering?
1148           }*/
1149         if (gMC->VolId("FRE1")==gMC->CurrentVolID(copy) || gMC->VolId("FRE2")==gMC->CurrentVolID(copy))
1150           {
1151             gMC->TrackMomentum(momentum);
1152             mom[0]=momentum(0);
1153             mom[1]=momentum(1);
1154             mom[2]=momentum(2);
1155             mom[3]=momentum(3);
1156             hits [19] = mom[0];
1157             hits [20] = mom[1];
1158             hits [21] = mom[2];
1159             fFreonProd=1;
1160           }
1161
1162         if (gMC->VolId("GAP ")== gMC->CurrentVolID(copy)) {
1163 // Get current particle id (ipart), track position (pos)  and momentum (mom)
1164             
1165             gMC->CurrentVolOffID(3,copy);
1166             vol[0]=copy;
1167             idvol=vol[0]-1;
1168
1169             //Int_t sector=((AliRICHChamber*) (*fChambers)[idvol])
1170                         //->Sector(localPos[0], localPos[2]);
1171             //printf("Sector:%d\n",sector);
1172             
1173             gMC->TrackPosition(position);
1174             gMC->TrackMomentum(momentum);
1175             pos[0]=position(0);
1176             pos[1]=position(1);
1177             pos[2]=position(2);
1178             mom[0]=momentum(0);
1179             mom[1]=momentum(1);
1180             mom[2]=momentum(2);
1181             mom[3]=momentum(3);
1182
1183             gMC->Gmtod(pos,localPos,1);
1184             
1185             //Chamber(idvol).GlobaltoLocal(pos,localPos);
1186                                                                     
1187             gMC->Gmtod(mom,localMom,2);
1188
1189             //Chamber(idvol).GlobaltoLocal(mom,localMom);
1190             
1191             ipart  = gMC->TrackPid();
1192             //
1193             // momentum loss and steplength in last step
1194             destep = gMC->Edep();
1195             step   = gMC->TrackStep();
1196   
1197             //
1198             // record hits when track enters ...
1199             if( gMC->IsTrackEntering()) {
1200 //              gMC->SetMaxStep(fMaxStepGas);
1201                 Double_t tc = mom[0]*mom[0]+mom[1]*mom[1];
1202                 Double_t rt = TMath::Sqrt(tc);
1203                 theta   = Float_t(TMath::ATan2(rt,Double_t(mom[2])))*kRaddeg;
1204                 phi     = Float_t(TMath::ATan2(Double_t(mom[1]),Double_t(mom[0])))*kRaddeg;
1205                 
1206
1207                 Double_t localTc = localMom[0]*localMom[0]+localMom[2]*localMom[2];
1208                 Double_t localRt = TMath::Sqrt(localTc);
1209                 localTheta   = Float_t(TMath::ATan2(localRt,Double_t(localMom[1])))*kRaddeg;                       
1210                 localPhi     = Float_t(TMath::ATan2(Double_t(localMom[2]),Double_t(localMom[0])))*kRaddeg;    
1211                 
1212                 hits[0] = Float_t(ipart);         // particle type
1213                 hits[1] = localPos[0];                 // X-position for hit
1214                 hits[2] = localPos[1];                 // Y-position for hit
1215                 hits[3] = localPos[2];                 // Z-position for hit
1216                 hits[4] = localTheta;                  // theta angle of incidence
1217                 hits[5] = localPhi;                    // phi angle of incidence 
1218                 hits[8] = (Float_t) fNsdigits;    // first sdigit
1219                 hits[9] = -1;                     // last pad hit
1220                 hits[13] = fFreonProd;           // did id hit the freon?
1221                 hits[14] = mom[0];
1222                 hits[15] = mom[1];
1223                 hits[16] = mom[2];
1224                 hits[18] = 0;               // dummy cerenkov angle
1225
1226                 tlength = 0;
1227                 eloss   = 0;
1228                 fFreonProd = 0;
1229         
1230                 Chamber(idvol).LocaltoGlobal(localPos,hits+1);
1231            
1232                 
1233                 //To make chamber coordinates x-y had to pass localPos[0], localPos[2]
1234                 xhit    = localPos[0];
1235                 yhit    = localPos[2];
1236                 // Only if not trigger chamber
1237                 if(idvol<kNCH) {
1238                     //
1239                     //  Initialize hit position (cursor) in the segmentation model 
1240           //PH              ((AliRICHChamber*) (*fChambers)[idvol])
1241                     ((AliRICHChamber*)fChambers->At(idvol))
1242                         ->SigGenInit(localPos[0], localPos[2], localPos[1]);
1243                 }
1244             }
1245             
1246             // 
1247             // Calculate the charge induced on a pad (disintegration) in case 
1248             //
1249             // Mip left chamber ...
1250             if( gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared()){
1251                 gMC->SetMaxStep(kBig);
1252                 eloss   += destep;
1253                 tlength += step;
1254                 
1255                                 
1256                 // Only if not trigger chamber
1257                 if(idvol<kNCH) {
1258                   if (eloss > 0) 
1259                     {
1260                       if(gMC->TrackPid() == kNeutron)
1261                         printf("\n\n\n\n\n Neutron Making Pad Hit!!! \n\n\n\n");
1262                       //nPads = Hits2SDigits(xhit,yhit,eloss,idvol,kMip); //for MIP kir
1263                       hits[17] = nPads;
1264                       //printf("nPads:%d",nPads);
1265                     }
1266                 }
1267                 
1268                 hits[6]=tlength;
1269                 hits[7]=eloss;
1270                 if (fNsdigits > (Int_t)hits[8]) {
1271                     hits[8]= hits[8]+1;
1272                     hits[9]= (Float_t) fNsdigits;
1273                 }
1274                 
1275                 //if(sector !=-1)
1276                 new(lhits[fNhits++]) AliRICHhit(fIshunt,gAlice->GetCurrentTrackNumber(),vol,hits);
1277                 eloss = 0; 
1278                 //
1279                 // Check additional signal generation conditions 
1280                 // defined by the segmentation
1281                 // model (boundary crossing conditions) 
1282             } else if 
1283           //PH          (((AliRICHChamber*) (*fChambers)[idvol])
1284                 (((AliRICHChamber*)fChambers->At(idvol))
1285                  ->SigGenCond(localPos[0], localPos[2], localPos[1]))
1286             {
1287           //PH          ((AliRICHChamber*) (*fChambers)[idvol])
1288                 ((AliRICHChamber*)fChambers->At(idvol))
1289                     ->SigGenInit(localPos[0], localPos[2], localPos[1]);
1290                 if (eloss > 0) 
1291                   {
1292                     if(gMC->TrackPid() == kNeutron)
1293                       printf("\n\n\n\n\n Neutron Making Pad Hit!!! \n\n\n\n");
1294                     //nPads = Hits2SDigits(xhit,yhit,eloss,idvol,kMip);//for N kir
1295                     hits[17] = nPads;
1296                     //printf("Npads:%d",NPads);
1297                   }
1298                 xhit     = localPos[0];
1299                 yhit     = localPos[2]; 
1300                 eloss    = destep;
1301                 tlength += step ;
1302                 //
1303                 // nothing special  happened, add up energy loss
1304             } else {        
1305                 eloss   += destep;
1306                 tlength += step ;
1307             }
1308         }
1309       }
1310     /*************************************************End of MIP treatment**************************************/
1311 }//void AliRICHv3::StepManager()