]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMD.cxx
Added functionalities in display the RICH event and pattern recognition. Minor bugs...
[u/mrichter/AliRoot.git] / FMD / AliFMD.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 //////////////////////////////////////////////////////////////////////////////
19 //                                                                          
20 //  Forward Multiplicity Detector based on Silicon wafers This class
21 //  contains the base procedures for the Forward Multiplicity detector
22 //  Detector consists of 5 Si volumes covered pseudorapidity interval
23 //  from 1.7 to 5.1.
24 //                                                                           
25 // The actual code is done by various separate classes.   Below is
26 // diagram showing the relationship between the various FMD classes
27 // that handles the geometry 
28 //
29 //
30 //       +----------+   +----------+   
31 //       | AliFMDv1 |   | AliFMDv1 |   
32 //       +----------+   +----------+   
33 //            |              |
34 //       +----+--------------+
35 //       |
36 //       |           +------------+ 1  +---------------+
37 //       |        +- | AliFMDRing |<>--| AliFMDPolygon | 
38 //       V     2  |  +------------+    +---------------+   
39 //  +--------+<>--+        |
40 //  | AliFMD |             ^                       
41 //  +--------+<>--+        V 1..2                     
42 //             3  | +-------------------+ 
43 //                +-| AliFMDSubDetector | 
44 //                  +-------------------+
45 //                           ^              
46 //                           |
47 //             +-------------+-------------+
48 //             |             |             |          
49 //        +---------+   +---------+   +---------+
50 //        | AliFMD1 |   | AliFMD2 |   | AliFMD3 |
51 //        +---------+   +---------+   +---------+
52 //      
53 //
54 // *  AliFMD 
55 //    This defines the interface for the various parts of AliROOT that
56 //    uses the FMD, like AliFMDDigitizer, AliFMDReconstructor, and so
57 //    on. 
58 //
59 // *  AliFMDv1 
60 //    This is a concrete implementation of the AliFMD interface. 
61 //    It is the responsibility of this class to create the FMD
62 //    geometry, process hits in the FMD, and serve hits and digits to
63 //    the various clients. 
64 //  
65 //    It uses the objects of class AliFMDSubDetector to do the various
66 //    stuff for FMD1, 2, and 3 
67 //
68 // *  AliFMDRing 
69 //    This class contains all stuff needed to do with a ring.  It's
70 //    used by the AliFMDSubDetector objects to instantise inner and
71 //    outer rings.  The AliFMDRing objects are shared by the
72 //    AliFMDSubDetector objects, and owned by the AliFMDv1 object. 
73 //
74 // *  AliFMDPolygon 
75 //    The code I lifted from TGeoPolygon to help with the geometry of
76 //    the modules, as well as to decide wether a hit is actually with
77 //    in the real module shape.  The point is, that the shape of the
78 //    various ring modules are really polygons (much like the lid of a
79 //    coffin), but it's segmented at constant radius.  That is very
80 //    hard to implement using GEANT 3.21 shapes, so instead the
81 //    modules are implemented as TUBS (tube sections), and in the step
82 //    procedure we do the test whether the track was inside the real
83 //    shape of the module.  
84 //
85 // *  AliFMD1, AliFMD2, and AliFMD3 
86 //    These are specialisation of AliFMDSubDetector, that contains the
87 //    particularities of each of the sub-detector system.  It is
88 //    envisioned that the classes should also define the support
89 //    volumes and material for each of the detectors.                          
90 //                                                                          
91 //Begin_Html
92 /*
93    </pre>
94    <br clear=left>
95    <p>
96      The responsible person for this module is
97      <a href="mailto:Alla.Maevskaia@cern.ch">Alla Maevskaia</a>.
98    </p>
99    <p>
100      Many modifications by <a href="mailto:cholm@nbi.dk">Christian
101        Holm Christensen</a>.
102    </p>
103    <pre>
104 */
105 //End_Html
106
107 #ifndef ROOT_TClonesArray
108 #include <TClonesArray.h>
109 #endif
110 #ifndef ROOT_TGeomtry
111 # include <TGeometry.h>
112 #endif
113 #ifndef ROOT_TNode
114 # include <TNode.h>
115 #endif
116 #ifndef ROOT_TTUBE
117 # include <TTUBE.h>
118 #endif
119 #ifndef ROOT_TTree
120 # include <TTree.h>
121 #endif
122 #ifndef ROOT_TVirtualMC
123 # include <TVirtualMC.h>
124 #endif
125 #ifndef ROOT_TBrowser
126 # include <TBrowser.h>
127 #endif
128 #ifndef ROOT_TMath
129 # include <TMath.h>
130 #endif
131
132 #ifndef ALIRUNDIGITIZER_H
133 # include "AliRunDigitizer.h"
134 #endif
135 #ifndef ALILOADER_H
136 # include "AliLoader.h"
137 #endif
138 #ifndef ALIRUN_H
139 # include "AliRun.h"
140 #endif
141 #ifndef ALIMC_H
142 # include "AliMC.h"
143 #endif
144 #ifndef ALILog_H
145 # include "AliLog.h"
146 #endif
147 #ifndef ALIMAGF_H
148 # include "AliMagF.h"
149 #endif
150 #ifndef ALIFMD_H
151 # include "AliFMD.h"
152 #endif
153 #ifndef ALIFMDDIGIG_H
154 # include "AliFMDDigit.h"
155 #endif
156 #ifndef ALIFMDHIT_H
157 # include "AliFMDHit.h"
158 #endif
159 #ifndef ALIFMDDIGITIZER_H
160 # include "AliFMDDigitizer.h"
161 #endif
162 #ifndef ALIFMD1_H
163 # include "AliFMD1.h"
164 #endif
165 #ifndef ALIFMD2_H
166 # include "AliFMD2.h"
167 #endif
168 #ifndef ALIFMD3_H
169 # include "AliFMD3.h"
170 #endif
171 #ifndef ALIALTROBUFFER_H
172 # include "AliAltroBuffer.h"
173 #endif
174
175 //____________________________________________________________________
176 ClassImp(AliFMD);
177
178 //____________________________________________________________________
179 AliFMD::AliFMD()
180   : fInner(0), 
181     fOuter(0),
182     fFMD1(0),
183     fFMD2(0), 
184     fFMD3(0)
185 {
186   //
187   // Default constructor for class AliFMD
188   //
189   AliDebug(0, "Default CTOR");
190   fHits     = 0;
191   fDigits   = 0;
192   fSDigits  = 0;
193   fNsdigits = 0;
194   fIshunt   = 0;
195 }
196
197 //____________________________________________________________________
198 AliFMD::AliFMD(const char *name, const char *title, bool detailed)
199   : AliDetector (name, title),
200     fInner(0), 
201     fOuter(0),
202     fFMD1(0),
203     fFMD2(0), 
204     fFMD3(0)
205 {
206   //
207   // Standard constructor for Forward Multiplicity Detector
208   //
209   AliDebug(0, "Standard CTOR");
210
211   // Initialise Hit array
212   HitsArray();
213   gAlice->GetMCApp()->AddHitList(fHits);
214
215   // (S)Digits for the detectors disk
216   DigitsArray();
217   SDigitsArray();
218   
219   // CHC: What is this?
220   fIshunt = 0;
221   SetMarkerColor(kRed);
222   SetLineColor(kYellow);
223   SetSiDensity();
224
225   // Create sub-volume managers 
226   fInner = new AliFMDRing('I', detailed);
227   fOuter = new AliFMDRing('O', detailed);
228   fFMD1  = new AliFMD1();
229   fFMD2  = new AliFMD2();
230   fFMD3  = new AliFMD3();
231
232   // Specify parameters of sub-volume managers 
233   fFMD1->SetInner(fInner);
234   fFMD1->SetOuter(0);
235
236   fFMD2->SetInner(fInner);
237   fFMD2->SetOuter(fOuter);
238   
239   fFMD3->SetInner(fInner);
240   fFMD3->SetOuter(fOuter);
241
242   SetLegLength();
243   SetLegRadius();
244   SetLegOffset();
245   SetModuleSpacing();
246   
247   fInner->SetLowR(4.3);
248   fInner->SetHighR(17.2);
249   fInner->SetWaferRadius(13.4/2);
250   fInner->SetTheta(36/2);
251   fInner->SetNStrips(512);
252   fInner->SetSiThickness(.03);
253   fInner->SetPrintboardThickness(.11);
254   fInner->SetBondingWidth(.5);
255
256   fOuter->SetLowR(15.6);
257   fOuter->SetHighR(28.0);
258   fOuter->SetWaferRadius(13.4/2);
259   fOuter->SetTheta(18/2);
260   fOuter->SetNStrips( 256);
261   fOuter->SetSiThickness(.03);
262   fOuter->SetPrintboardThickness(.1);
263   fOuter->SetBondingWidth(.5);
264   
265   
266   fFMD1->SetHoneycombThickness(1);
267   fFMD1->SetInnerZ(340.0);
268   
269   fFMD2->SetHoneycombThickness(1);
270   fFMD2->SetInnerZ(83.4);
271   fFMD2->SetOuterZ(75.2);
272
273   fFMD3->SetHoneycombThickness(1);
274   fFMD3->SetInnerZ(-62.8);
275   fFMD3->SetOuterZ(-75.2);
276 }
277
278 //____________________________________________________________________
279 AliFMD::~AliFMD ()
280 {
281   // Destructor for base class AliFMD
282   if (fHits) {
283     fHits->Delete();
284     delete fHits;
285     fHits = 0;
286   }
287   if (fDigits) {
288     fDigits->Delete();
289     delete fDigits;
290     fDigits = 0;
291   }
292   if (fSDigits) {
293     fSDigits->Delete();
294     delete fSDigits;
295     fSDigits = 0;
296   }
297 }
298
299 //====================================================================
300 //
301 // GEometry ANd Traking
302 //
303 //____________________________________________________________________
304 void 
305 AliFMD::CreateGeometry()
306 {
307  //
308   // Create the geometry of Forward Multiplicity Detector version 0
309   //
310   // DebugGuard guard("AliFMD::CreateGeometry");
311   AliDebug(10, "Creating geometry");
312
313   fInner->Init();
314   fOuter->Init();
315
316   TString name;
317   Double_t par[3];
318
319   par[0]      =  fLegRadius - .1;
320   par[1]      =  fLegRadius;
321   par[2]      =  fLegLength / 2;
322   name        =  "SLEG";
323   fShortLegId =  gMC->Gsvolu(name.Data(),"TUBE",(*fIdtmed)[kPlasticId],par,3);
324   
325   par[2]      += fModuleSpacing / 2;
326   name        = "LLEG";
327   fLongLegId  =  gMC->Gsvolu(name.Data(),"TUBE",(*fIdtmed)[kPlasticId],par,3);
328
329   fInner->SetupGeometry((*fIdtmed)[kAirId], 
330                         (*fIdtmed)[kSiId], 
331                         (*fIdtmed)[kPcbId], 
332                         fPrintboardRotationId, 
333                         fIdentityRotationId);
334   fOuter->SetupGeometry((*fIdtmed)[kAirId], 
335                         (*fIdtmed)[kSiId], 
336                         (*fIdtmed)[kPcbId], 
337                         fPrintboardRotationId, 
338                         fIdentityRotationId);
339
340   fFMD1->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
341   fFMD2->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
342   fFMD3->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
343   
344   fFMD1->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);
345   fFMD2->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);
346   fFMD3->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);    
347 }    
348
349 //____________________________________________________________________
350 void AliFMD::CreateMaterials() 
351 {
352   // Register various materials and tracking mediums with the
353   // backend.   
354   // 
355   // Currently defined materials and mediums are 
356   // 
357   //    FMD Air         Normal air 
358   //    FMD Si          Active silicon of sensors 
359   //    FMD Carbon      Normal carbon used in support, etc. 
360   //    FMD Kapton      Carbon used in Honeycomb
361   //    FMD PCB         Printed circuit board material 
362   //    FMD Plastic     Material for support legs 
363   // 
364   // Also defined are two rotation matricies. 
365   //
366   // DebugGuard guard("AliFMD::CreateMaterials");
367   AliDebug(10, "Creating materials");
368   Int_t    id;
369   Double_t a                = 0;
370   Double_t z                = 0;
371   Double_t density          = 0;
372   Double_t radiationLength  = 0;
373   Double_t absorbtionLength = 999;
374   Int_t    fieldType        = gAlice->Field()->Integ();     // Field type 
375   Double_t maxField         = gAlice->Field()->Max();     // Field max.
376   Double_t maxBending       = 0;     // Max Angle
377   Double_t maxStepSize      = 0.001; // Max step size 
378   Double_t maxEnergyLoss    = 1;     // Max Delta E
379   Double_t precision        = 0.001; // Precision
380   Double_t minStepSize      = 0.001; // Minimum step size 
381  
382   // Silicon 
383   a                = 28.0855;
384   z                = 14.;
385   density          = fSiDensity;
386   radiationLength  = 9.36;
387   maxBending       = 1;
388   maxStepSize      = .001;
389   precision        = .001;
390   minStepSize      = .001;
391   id               = kSiId;
392   AliMaterial(id, "FMD Si$", a, z, density, radiationLength, absorbtionLength);
393   AliMedium(kSiId, "FMD Si$",id,1,fieldType,maxField,maxBending,
394             maxStepSize,maxEnergyLoss,precision,minStepSize);
395
396
397   // Carbon 
398   a                = 12.011;
399   z                = 6.;
400   density          = 2.265;
401   radiationLength  = 18.8;
402   maxBending       = 10;
403   maxStepSize      = .01;
404   precision        = .003;
405   minStepSize      = .003;
406   id               = kCarbonId;
407   AliMaterial(id, "FMD Carbon$", a, z, density, radiationLength, 
408               absorbtionLength);
409   AliMedium(kCarbonId, "FMD Carbon$",id,0,fieldType,maxField,maxBending,
410             maxStepSize,maxEnergyLoss,precision,minStepSize);
411
412   // Silicon chip 
413   {
414     Float_t as[] = { 12.0107,      14.0067,      15.9994,
415                      1.00794,      28.0855,     107.8682 };
416     Float_t zs[] = {  6.,           7.,           8.,
417                       1.,          14.,          47. };
418     Float_t ws[] = {  0.039730642,  0.001396798,  0.01169634,
419                       0.004367771,  0.844665,     0.09814344903 };
420     density = 2.36436;
421     maxBending       = 10;
422     maxStepSize      = .01;
423     precision        = .003;
424     minStepSize      = .003;
425     id = kSiChipId;
426     AliMixture(id, "FMD Si Chip$", as, zs, density, 6, ws);
427     AliMedium(kSiChipId, "FMD Si Chip$", id, 0, fieldType, maxField, 
428               maxBending, maxStepSize, maxEnergyLoss, precision, minStepSize);
429   }
430   
431
432   // Kaption 
433   {
434     Float_t as[] = { 1.00794,  12.0107,  14.010,   15.9994};
435     Float_t zs[] = { 1.,        6.,       7.,       8.};
436     Float_t ws[] = { 0.026362,  0.69113,  0.07327,  0.209235};
437     density          = 1.42;
438     maxBending       = 1;
439     maxStepSize      = .001;
440     precision        = .001;
441     minStepSize      = .001;
442     id               = kKaptionId;
443     AliMixture(id, "FMD Kaption$", as, zs, density, 4, ws);
444     AliMedium(kKaptionId, "FMD Kaption$",id,0,fieldType,maxField,maxBending,
445               maxStepSize,maxEnergyLoss,precision,minStepSize);
446   }
447   
448   // Air
449   {
450     Float_t as[] = { 12.0107, 14.0067,   15.9994,  39.948 };
451     Float_t zs[] = {  6.,      7.,       8.,       18. };
452     Float_t ws[] = { 0.000124, 0.755267, 0.231781, 0.012827 }; 
453     density      = .00120479;
454     maxBending   = 1;
455     maxStepSize  = .001;
456     precision    = .001;
457     minStepSize  = .001;
458     id           = kAirId;
459     AliMixture(id, "FMD Air$", as, zs, density, 4, ws);
460     AliMedium(kAirId, "FMD Air$", id,0,fieldType,maxField,maxBending,
461               maxStepSize,maxEnergyLoss,precision,minStepSize);
462   }
463   
464   // PCB
465   {
466     Float_t zs[] = { 14.,         20.,         13.,         12.,
467                       5.,         22.,         11.,         19.,
468                      26.,          9.,          8.,          6.,
469                       7.,          1.};
470     Float_t as[] = { 28.0855,     40.078,      26.981538,   24.305, 
471                      10.811,      47.867,      22.98977,    39.0983,
472                      55.845,      18.9984,     15.9994,     12.0107,
473                      14.0067,      1.00794};
474     Float_t ws[] = {  0.15144894,  0.08147477,  0.04128158,  0.00904554, 
475                       0.01397570,  0.00287685,  0.00445114,  0.00498089,
476                       0.00209828,  0.00420000,  0.36043788,  0.27529426,
477                       0.01415852,  0.03427566};
478     density      = 1.8;
479     maxBending   = 1;
480     maxStepSize  = .001;
481     precision    = .001;
482     minStepSize  = .001;
483     id           = kPcbId;
484     AliMixture(id, "FMD PCB$", as, zs, density, 14, ws);
485     AliMedium(kPcbId, "FMD PCB$", id,1,fieldType,maxField,maxBending,
486               maxStepSize,maxEnergyLoss,precision,minStepSize);
487   }
488   
489   // Plastic 
490   {
491     Float_t as[] = { 1.01, 12.01 };
492     Float_t zs[] = { 1.,   6.    };
493     Float_t ws[] = { 1.,   1.    };
494     density      = 1.03;
495     maxBending   = 10;
496     maxStepSize  = .01;
497     precision    = .003;
498     minStepSize  = .003;
499     id           = kPlasticId;
500     AliMixture(id, "FMD Plastic$", as, zs, density, -2, ws);
501     AliMedium(kPlasticId, "FMD Plastic$", id,0,fieldType,maxField,maxBending,
502                 maxStepSize,maxEnergyLoss,precision,minStepSize);
503   }
504   AliMatrix(fPrintboardRotationId, 90, 90, 0, 90, 90, 0);
505   AliMatrix(fIdentityRotationId, 90, 0, 90, 90, 0, 0);
506 }
507
508 //____________________________________________________________________
509 void  
510 AliFMD::Init()
511 {
512   //
513   // Initialis the FMD after it has been built
514   Int_t i;
515   //
516   if (fDebug) {
517     std::cout << "\n" << ClassName() << ": " << std::flush;
518     for (i = 0; i < 35; i++) std::cout << "*";
519     std::cout << " FMD_INIT ";
520     for (i = 0; i < 35; i++) std::cout << "*";
521     std::cout << "\n" << ClassName() << ": " << std::flush;
522     //
523     // Here the FMD initialisation code (if any!)
524     for (i = 0; i < 80; i++) std::cout << "*";
525     std::cout << std::endl;
526   }
527   //
528   //
529 }
530
531 //====================================================================
532 //
533 // Graphics and event display
534 //
535 //____________________________________________________________________
536 void 
537 AliFMD::BuildGeometry()
538 {
539   //
540   // Build simple ROOT TNode geometry for event display
541   //
542   // Build a simplified geometry of the FMD used for event display  
543   // 
544   AliDebug(10, "Creating a simplified geometry");
545
546   TNode* top = gAlice->GetGeometry()->GetNode("alice");
547   
548   fFMD1->SimpleGeometry(fNodes, top, GetLineColor(), 0);
549   fFMD2->SimpleGeometry(fNodes, top, GetLineColor(), 0);
550   fFMD3->SimpleGeometry(fNodes, top, GetLineColor(), 0);
551 }
552
553 //____________________________________________________________________
554 void 
555 AliFMD::DrawDetector()
556 {
557   //
558   // Draw a shaded view of the Forward multiplicity detector version 0
559   //
560   // DebugGuard guard("AliFMD::DrawDetector");
561   AliDebug(10, "Draw detector");
562   
563   //Set ALIC mother transparent
564   gMC->Gsatt("ALIC","SEEN",0);
565
566   //Set volumes visible
567   fFMD1->Gsatt();
568   fFMD2->Gsatt();
569   fFMD3->Gsatt();
570   fInner->Gsatt();
571   fOuter->Gsatt();
572
573   //
574   gMC->Gdopt("hide", "on");
575   gMC->Gdopt("shad", "on");
576   gMC->Gsatt("*", "fill", 7);
577   gMC->SetClipBox(".");
578   gMC->SetClipBox("*", 0, 1000, -1000, 1000, -1000, 1000);
579   gMC->DefaultRange();
580   gMC->Gdraw("alic", 40, 30, 0, 12, 12, .055, .055);
581   gMC->Gdhead(1111, "Forward Multiplicity Detector");
582   gMC->Gdman(16, 10, "MAN");
583   gMC->Gdopt("hide", "off");
584 }
585
586 //____________________________________________________________________
587 const Int_t 
588 AliFMD::DistanceToPrimitive(Int_t, Int_t)
589 {
590   //
591   // Calculate the distance from the mouse to the FMD on the screen
592   // Dummy routine
593   //
594   return 9999;
595 }
596
597 //====================================================================
598 //
599 // Hit and Digit managment 
600 //
601 //____________________________________________________________________
602 void 
603 AliFMD::MakeBranch(Option_t * option)
604 {
605   // Create Tree branches for the FMD.
606   const Int_t kBufferSize = 16000;
607   TString branchname(GetName());
608   TString opt(option);
609   
610   if (opt.Contains("H", TString::kIgnoreCase)) {
611     HitsArray();
612     AliDetector::MakeBranch(option); 
613   }
614   if (opt.Contains("D", TString::kIgnoreCase)) { 
615     DigitsArray();
616     MakeBranchInTree(fLoader->TreeD(), branchname.Data(),
617                      &fDigits, kBufferSize, 0);
618   }
619   if (opt.Contains("S", TString::kIgnoreCase)) { 
620     SDigitsArray();
621     MakeBranchInTree(fLoader->TreeS(), branchname.Data(),
622                      &fSDigits, kBufferSize, 0);
623   }
624 }
625
626 //____________________________________________________________________
627 void 
628 AliFMD::SetTreeAddress()
629 {
630   // Set branch address for the Hits and Digits Tree.
631
632   if (fLoader->TreeH()) HitsArray();
633   AliDetector::SetTreeAddress();
634
635   TTree *treeD = fLoader->TreeD();
636   if (treeD) {
637     DigitsArray();
638     TBranch* branch = treeD->GetBranch ("FMD");
639     if (branch) branch->SetAddress(&fDigits);
640   }
641
642   TTree *treeS = fLoader->TreeS();
643   if (treeS) {
644     SDigitsArray();
645     TBranch* branch = treeS->GetBranch ("FMD");
646     if (branch) branch->SetAddress(&fSDigits);
647   }
648 }
649
650
651
652 //____________________________________________________________________
653 void 
654 AliFMD::SetHitsAddressBranch(TBranch *b)
655 {
656   b->SetAddress(&fHits);
657 }
658
659 //____________________________________________________________________
660 void 
661 AliFMD::AddHit(Int_t track, Int_t *vol, Float_t *hits) 
662 {
663   // Add a hit to the hits tree 
664   // 
665   // The information of the two arrays are decoded as 
666   // 
667   // Parameters
668   //    track                Track #
669   //    ivol[0]  [UShort_t ] Detector # 
670   //    ivol[1]  [Char_t   ] Ring ID 
671   //    ivol[2]  [UShort_t ] Sector #
672   //    ivol[3]  [UShort_t ] Strip # 
673   //    hits[0]  [Float_t  ] Track's X-coordinate at hit 
674   //    hits[1]  [Float_t  ] Track's Y-coordinate at hit
675   //    hits[3]  [Float_t  ] Track's Z-coordinate at hit
676   //    hits[4]  [Float_t  ] X-component of track's momentum             
677   //    hits[5]  [Float_t  ] Y-component of track's momentum             
678   //    hits[6]  [Float_t  ] Z-component of track's momentum            
679   //    hits[7]  [Float_t  ] Energy deposited by track                  
680   //    hits[8]  [Int_t    ] Track's particle Id # 
681   //    hits[9]  [Float_t  ] Time when the track hit                    
682   AddHit(track, 
683          UShort_t(vol[0]),  // Detector # 
684          Char_t(vol[1]),    // Ring ID
685          UShort_t(vol[2]),  // Sector # 
686          UShort_t(vol[3]),  // Strip # 
687          hits[0],           // X
688          hits[1],           // Y
689          hits[2],           // Z
690          hits[3],           // Px
691          hits[4],           // Py
692          hits[5],           // Pz
693          hits[6],           // Energy loss 
694          Int_t(hits[7]),    // PDG 
695          hits[8]);          // Time
696 }
697
698 //____________________________________________________________________
699 void 
700 AliFMD::AddHit(Int_t    track, 
701                UShort_t detector, 
702                Char_t   ring, 
703                UShort_t sector, 
704                UShort_t strip, 
705                Float_t  x, 
706                Float_t  y, 
707                Float_t  z,
708                Float_t  px, 
709                Float_t  py, 
710                Float_t  pz,
711                Float_t  edep,
712                Int_t    pdg,
713                Float_t  t)
714 {
715   //
716   // Add a hit to the list
717   //
718   // Parameters:
719   // 
720   //    track     Track #
721   //    detector  Detector # (1, 2, or 3)                      
722   //    ring      Ring ID ('I' or 'O')
723   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
724   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
725   //    x         Track's X-coordinate at hit
726   //    y         Track's Y-coordinate at hit
727   //    z         Track's Z-coordinate at hit
728   //    px        X-component of track's momentum 
729   //    py        Y-component of track's momentum
730   //    pz        Z-component of track's momentum
731   //    edep      Energy deposited by track
732   //    pdg       Track's particle Id #
733   //    t         Time when the track hit 
734   // 
735   TClonesArray& a = *(HitsArray());
736   // Search through the list of already registered hits, and see if we
737   // find a hit with the same parameters.  If we do, then don't create
738   // a new hit, but rather update the energy deposited in the hit.
739   // This is done, so that a FLUKA based simulation will get the
740   // number of hits right, not just the enerrgy deposition. 
741   for (Int_t i = 0; i < fNhits; i++) {
742     if (!a.At(i)) continue;
743     AliFMDHit* hit = static_cast<AliFMDHit*>(a.At(i));
744     if (hit->Detector() == detector 
745         && hit->Ring() == ring
746         && hit->Sector() == sector 
747         && hit->Strip() == strip
748         && hit->Track() == track) {
749       hit->SetEdep(hit->Edep() + edep);
750       return;
751     }
752   }
753   // If hit wasn't already registered, do so know. 
754   new (a[fNhits]) AliFMDHit(fIshunt, track, detector, ring, sector, strip, 
755                             x, y, z, px, py, pz, edep, pdg, t);
756   fNhits++;
757 }
758
759 //____________________________________________________________________
760 void 
761 AliFMD::AddDigit(Int_t* digits)
762 {
763   // Add a digit to the Digit tree 
764   // 
765   // Paramters 
766   //
767   //    digits[0]  [UShort_t] Detector #
768   //    digits[1]  [Char_t]   Ring ID
769   //    digits[2]  [UShort_t] Sector #
770   //    digits[3]  [UShort_t] Strip #
771   //    digits[4]  [UShort_t] ADC Count 
772   //    digits[5]  [Short_t]  ADC Count, -1 if not used
773   //    digits[6]  [Short_t]  ADC Count, -1 if not used 
774   // 
775   AddDigit(UShort_t(digits[0]),  // Detector #
776            Char_t(digits[1]),    // Ring ID
777            UShort_t(digits[2]),  // Sector #
778            UShort_t(digits[3]),  // Strip #
779            UShort_t(digits[4]),  // ADC Count1 
780            Short_t(digits[5]),   // ADC Count2 
781            Short_t(digits[6]));  // ADC Count3 
782 }
783
784 //____________________________________________________________________
785 void 
786 AliFMD::AddDigit(UShort_t detector, 
787                  Char_t   ring, 
788                  UShort_t sector, 
789                  UShort_t strip, 
790                  UShort_t count1, 
791                  Short_t  count2,
792                  Short_t  count3)
793 {
794   // add a real digit - as coming from data
795   // 
796   // Parameters 
797   //
798   //    detector  Detector # (1, 2, or 3)                      
799   //    ring      Ring ID ('I' or 'O')
800   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
801   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
802   //    count1    ADC count (a 10-bit word)
803   //    count2    ADC count (a 10-bit word), or -1 if not used
804   //    count3    ADC count (a 10-bit word), or -1 if not used
805   TClonesArray& a = *(DigitsArray());
806   
807   new (a[fNdigits++]) 
808     AliFMDDigit(detector, ring, sector, strip, count1, count2, count3);
809 }
810
811 //____________________________________________________________________
812 void 
813 AliFMD::AddSDigit(Int_t* digits)
814 {
815   // Add a digit to the SDigit tree 
816   // 
817   // Paramters 
818   //
819   //    digits[0]  [UShort_t] Detector #
820   //    digits[1]  [Char_t]   Ring ID
821   //    digits[2]  [UShort_t] Sector #
822   //    digits[3]  [UShort_t] Strip #
823   //    digits[4]  [Float_t]  Total energy deposited 
824   //    digits[5]  [UShort_t] ADC Count 
825   //    digits[6]  [Short_t]  ADC Count, -1 if not used
826   //    digits[7]  [Short_t]  ADC Count, -1 if not used 
827   // 
828   AddSDigit(UShort_t(digits[0]),  // Detector #
829             Char_t(digits[1]),    // Ring ID
830             UShort_t(digits[2]),  // Sector #
831             UShort_t(digits[3]),  // Strip #
832             Float_t(digits[4]),   // Edep
833             UShort_t(digits[5]),  // ADC Count1 
834             Short_t(digits[6]),   // ADC Count2 
835             Short_t(digits[7]));  // ADC Count3 
836 }
837
838 //____________________________________________________________________
839 void 
840 AliFMD::AddSDigit(UShort_t detector, 
841                   Char_t   ring, 
842                   UShort_t sector, 
843                   UShort_t strip, 
844                   Float_t  edep,
845                   UShort_t count1, 
846                   Short_t  count2,
847                   Short_t  count3)
848 {
849   // add a summable digit
850   // 
851   // Parameters 
852   //
853   //    detector  Detector # (1, 2, or 3)                      
854   //    ring      Ring ID ('I' or 'O')
855   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
856   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
857   //    edep      Total energy deposited
858   //    count1    ADC count (a 10-bit word)
859   //    count2    ADC count (a 10-bit word), or -1 if not used
860   //    count3    ADC count (a 10-bit word), or -1 if not used
861   TClonesArray& a = *(SDigitsArray());
862   
863   new (a[fNsdigits++]) 
864     AliFMDSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
865 }
866
867 //____________________________________________________________________
868 void 
869 AliFMD::ResetSDigits()
870 {
871   //
872   // Reset number of digits and the digits array for this detector
873   //
874   fNsdigits   = 0;
875   if (fSDigits) fSDigits->Clear();
876 }
877
878
879 //____________________________________________________________________
880 TClonesArray*
881 AliFMD::HitsArray() 
882 {
883   // Initialize hit array if not already, and return pointer to it. 
884   if (!fHits) { 
885     fHits = new TClonesArray("AliFMDHit", 1000);
886     fNhits = 0;
887   }
888   return fHits;
889 }
890
891 //____________________________________________________________________
892 TClonesArray*
893 AliFMD::DigitsArray() 
894 {
895   // Initialize digit array if not already, and return pointer to it. 
896   if (!fDigits) { 
897     fDigits = new TClonesArray("AliFMDDigit", 1000);
898     fNdigits = 0;
899   }
900   return fDigits;
901 }
902
903 //____________________________________________________________________
904 TClonesArray*
905 AliFMD::SDigitsArray() 
906 {
907   // Initialize digit array if not already, and return pointer to it. 
908   if (!fSDigits) { 
909     fSDigits = new TClonesArray("AliFMDSDigit", 1000);
910     fNsdigits = 0;
911   }
912   return fSDigits;
913 }
914
915 //====================================================================
916 //
917 // Digitization 
918 //
919 //____________________________________________________________________
920 void 
921 AliFMD::Hits2Digits() 
922 {
923   AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
924   manager->SetInputStream(0, "galice.root");
925   manager->SetOutputFile("H2Dfile");
926   
927   /* AliDigitizer* dig =*/ CreateDigitizer(manager);
928   manager->Exec("");
929 }
930
931 //____________________________________________________________________
932 void 
933 AliFMD::Hits2SDigits() 
934 {
935   AliDigitizer* sdig = new AliFMDSDigitizer("galice.root");
936   sdig->Exec("");
937 }
938
939   
940 //____________________________________________________________________
941 AliDigitizer* 
942 AliFMD::CreateDigitizer(AliRunDigitizer* manager) const
943 {
944   // Create a digitizer object 
945   return new AliFMDDigitizer(manager);
946 }
947
948 //====================================================================
949 //
950 // Raw data simulation 
951 //
952 //__________________________________________________________________
953 void 
954 AliFMD::Digits2Raw() 
955 {
956   AliFMD* fmd = static_cast<AliFMD*>(gAlice->GetDetector(GetName()));
957   fLoader->LoadDigits();
958   TTree* digitTree = fLoader->TreeD();
959   if (!digitTree) {
960     Error("Digits2Raw", "no digit tree");
961     return;
962   }
963   
964   TClonesArray* digits = new TClonesArray("AliFMDDigit", 1000);
965   fmd->SetTreeAddress();
966   TBranch* digitBranch = digitTree->GetBranch(GetName());
967   if (!digitBranch) {
968     Error("Digits2Raw", "no branch for %s", GetName());
969     return;
970   }
971   digitBranch->SetAddress(&digits);
972   
973   Int_t nEvents = Int_t(digitTree->GetEntries());
974   for (Int_t event = 0; event < nEvents; event++) {
975     fmd->ResetDigits();
976     digitTree->GetEvent(event);
977     
978     Int_t nDigits = digits->GetEntries();
979     if (nDigits < 1) continue;
980
981
982     UShort_t prevDetector = 0;
983     Char_t   prevRing     = '\0';
984     UShort_t prevSector   = 0;
985     // UShort_t prevStrip    = 0;
986
987     // The first seen strip number for a channel 
988     UShort_t startStrip   = 0;
989     
990     // Which channel number in the ALTRO channel we're at 
991     UShort_t offset       = 0;
992
993     // How many times the ALTRO Samples one VA1_ALICE channel 
994     Int_t sampleRate = 1;
995
996     // A buffer to hold 1 ALTRO channel - Normally, one ALTRO channel
997     // holds 128 VA1_ALICE channels, sampled at a rate of `sampleRate' 
998     TArrayI channel(128 * sampleRate);
999     
1000     // The Altro buffer 
1001     AliAltroBuffer* altro = 0;
1002     
1003     // Loop over the digits in the event.  Note, that we assume the
1004     // the digits are in order in the branch.   If they were not, we'd
1005     // have to cache all channels before we could write the data to
1006     // the ALTRO buffer, or we'd have to set up a map of the digits. 
1007     for (Int_t i = 0; i < nDigits; i++) {
1008       // Get the digit
1009       AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
1010
1011       UShort_t det    = digit->Detector();
1012       Char_t   ring   = digit->Ring();
1013       UShort_t sector = digit->Sector();
1014       UShort_t strip  = digit->Strip();
1015       if (det != prevDetector) {
1016         AliDebug(10, Form("FMD: New DDL, was %d, now %d",
1017                           kBaseDDL + prevDetector - 1,
1018                           kBaseDDL + det - 1));
1019         // If an altro exists, delete the object, flushing the data to
1020         // disk, and closing the file. 
1021         if (altro) { 
1022           // When the first argument is false, we write the real
1023           // header. 
1024           AliDebug(10, Form("New altro: Write channel at %d Strip: %d "
1025                             "Sector: %d  Ring: %d", 
1026                             i, startStrip, prevSector, prevRing));
1027           // TPC to FMD translations 
1028           // 
1029           //    TPC                FMD
1030           //    ----------+-----------
1031           //    pad       |      strip
1032           //    row       |     sector
1033           //    sector    |       ring
1034           // 
1035           altro->WriteChannel(Int_t(startStrip), 
1036                               Int_t(prevSector), 
1037                               Int_t((prevRing == 'I' ? 0 : 1)), 
1038                               channel.fN, channel.fArray, 0);
1039           altro->Flush();
1040           altro->WriteDataHeader(kFALSE, kFALSE);
1041           delete altro;
1042           altro = 0;
1043         }
1044
1045         prevDetector = det;
1046         // Need to open a new DDL! 
1047         Int_t ddlId = kBaseDDL + det - 1;
1048         TString filename(Form("%s_%d.ddl", GetName(),  ddlId));
1049
1050         AliDebug(10, Form("New altro buffer with DDL file %s", 
1051                           filename.Data()));
1052         AliDebug(10, Form("New altro at %d", i));
1053         // Create a new altro buffer - a `1' as the second argument
1054         // means `write mode' 
1055         altro = new AliAltroBuffer(filename.Data(), 1);
1056         
1057         // Write a dummy (first argument is true) header to the DDL
1058         // file - later on, when we close the file, we write the real
1059         // header
1060         altro->WriteDataHeader(kTRUE, kFALSE);
1061
1062         // Figure out the sample rate 
1063         if (digit->Count2() > 0) sampleRate = 2;
1064         if (digit->Count3() > 0) sampleRate = 3;
1065
1066         channel.Set(128 * sampleRate);
1067         offset     = 0;
1068         prevRing   = ring;
1069         prevSector = sector;
1070         startStrip = strip;
1071       }
1072       else if (offset == 128                        
1073                || digit->Ring() != prevRing 
1074                || digit->Sector() != prevSector) {
1075         // Force a new Altro channel
1076         AliDebug(10, Form("Flushing channel to disk because %s",
1077                           (offset == 128 ? "channel is full" :
1078                            (ring != prevRing ? "new ring up" :
1079                             "new sector up"))));
1080         AliDebug(10, Form("New Channel: Write channel at %d Strip: %d "
1081                           "Sector: %d  Ring: %d", 
1082                           i, startStrip, prevSector, prevRing));
1083         altro->WriteChannel(Int_t(startStrip), 
1084                             Int_t(prevSector), 
1085                             Int_t((prevRing == 'I' ? 0 : 1)), 
1086                             channel.fN, channel.fArray, 0);
1087         // Reset and update channel variables 
1088         channel.Reset(0);
1089         offset     = 0; 
1090         startStrip = strip;
1091         prevRing   = ring;
1092         prevSector = sector;
1093       }
1094
1095       // Store the counts of the ADC in the channel buffer 
1096       channel[offset * sampleRate] = digit->Count1();
1097       if (sampleRate > 1) 
1098         channel[offset * sampleRate + 1] = digit->Count2();
1099       if (sampleRate > 2) 
1100         channel[offset * sampleRate + 2] = digit->Count3();
1101       offset++;
1102     }
1103     // Finally, we need to close the final ALTRO buffer if it wasn't
1104     // already 
1105     if (altro) {
1106       altro->Flush();
1107       altro->WriteDataHeader(kFALSE, kFALSE);
1108       delete altro;
1109     }
1110   }
1111   fLoader->UnloadDigits();
1112 }
1113
1114 //==================================================================
1115 //
1116 // Various setter functions for the common paramters 
1117 //
1118
1119 //__________________________________________________________________
1120 void 
1121 AliFMD::SetLegLength(Double_t length) 
1122 {
1123  // Set lenght of plastic legs that hold the hybrid (print board and
1124   // silicon sensor) onto the honeycomp support
1125   //
1126   // DebugGuard guard("AliFMD::SetLegLength");
1127    AliDebug(10, "AliFMD::SetLegLength");
1128   fLegLength = length;
1129   fInner->SetLegLength(fLegLength);
1130   fOuter->SetLegLength(fLegLength);
1131 }
1132
1133 //__________________________________________________________________
1134 void 
1135 AliFMD::SetLegOffset(Double_t offset) 
1136 {
1137   // Set offset from edge of hybrid to plastic legs that hold the
1138   // hybrid (print board and silicon sensor) onto the honeycomp
1139   // support 
1140   //
1141   // DebugGuard guard("AliFMD::SetLegOffset");
1142   AliDebug(10, "AliFMD::SetLegOffset");
1143   fInner->SetLegOffset(offset);
1144   fOuter->SetLegOffset(offset);
1145 }
1146
1147 //__________________________________________________________________
1148 void 
1149 AliFMD::SetLegRadius(Double_t radius) 
1150 {
1151   // Set the diameter of the plastic legs that hold the hybrid (print
1152   // board and silicon sensor) onto the honeycomp support
1153   //
1154   // DebugGuard guard("AliFMD::SetLegRadius");
1155   AliDebug(10, "AliFMD::SetLegRadius");
1156   fLegRadius = radius;
1157   fInner->SetLegRadius(fLegRadius);
1158   fOuter->SetLegRadius(fLegRadius);
1159 }
1160
1161 //__________________________________________________________________
1162 void 
1163 AliFMD::SetModuleSpacing(Double_t spacing) 
1164 {
1165   // Set the distance between the front and back sensor modules
1166   // (module staggering). 
1167   //
1168   // DebugGuard guard("AliFMD::SetModuleSpacing");
1169   AliDebug(10, "AliFMD::SetModuleSpacing");  
1170   fModuleSpacing = spacing;
1171   fInner->SetModuleSpacing(fModuleSpacing);
1172   fOuter->SetModuleSpacing(fModuleSpacing);
1173 }
1174
1175 //====================================================================
1176 //
1177 // Utility 
1178 //
1179 //__________________________________________________________________
1180 void 
1181 AliFMD::Browse(TBrowser* b) 
1182 {
1183   AliDebug(10, "AliFMD::Browse");
1184   AliDetector::Browse(b);
1185   if (fInner) b->Add(fInner, "Inner Ring");
1186   if (fOuter) b->Add(fOuter, "Outer Ring");
1187   if (fFMD1)  b->Add(fFMD1,  "FMD1 SubDetector");
1188   if (fFMD2)  b->Add(fFMD2,  "FMD2 SubDetector");
1189   if (fFMD3)  b->Add(fFMD3,  "FMD3 SubDetector");
1190 }
1191
1192
1193 //********************************************************************
1194 //
1195 // AliFMDv0
1196 //
1197 //__________________________________________________________________
1198
1199 ClassImp(AliFMDv0);
1200
1201 //********************************************************************
1202 //
1203 // AliFMDv1
1204 //
1205 //__________________________________________________________________
1206
1207 ClassImp(AliFMDv1);
1208
1209
1210 //_//____________________________________________________________________
1211 void 
1212 AliFMDv1::StepManager()
1213 {
1214   //
1215   // Called for every step in the Forward Multiplicity Detector
1216   //
1217   // The procedure is as follows: 
1218   // 
1219   //   - IF NOT track is alive THEN RETURN ENDIF
1220   //   - IF NOT particle is charged THEN RETURN ENDIF
1221   //   - IF NOT volume name is "STRI" or "STRO" THEN RETURN ENDIF 
1222   //   - Get strip number (volume copy # minus 1)
1223   //   - Get phi division number (mother volume copy #)
1224   //   - Get module number (grand-mother volume copy #)
1225   //   - section # = 2 * module # + phi division # - 1
1226   //   - Get ring Id from volume name 
1227   //   - Get detector # from grand-grand-grand-mother volume name 
1228   //   - Get pointer to sub-detector object. 
1229   //   - Get track position 
1230   //   - IF track is entering volume AND track is inside real shape THEN
1231   //   -   Reset energy deposited 
1232   //   -   Get track momentum 
1233   //   -   Get particle ID # 
1234   ///  - ENDIF
1235   //   - IF track is inside volume AND inside real shape THEN 
1236   ///  -   Update energy deposited 
1237   //   - ENDIF 
1238   //   - IF track is inside real shape AND (track is leaving volume,
1239   //         or it died, or it is stopped  THEN
1240   //   -   Create a hit 
1241   //   - ENDIF
1242   //     
1243   //
1244   // DebugGuard guard("AliFMDv1::StepManager");
1245   AliDebug(10, "AliFMDv1::StepManager");
1246   // return;
1247
1248   // If the track is gone, return
1249   if (!gMC->IsTrackAlive()) return;
1250   
1251   // Only process charged particles 
1252   if(TMath::Abs(gMC->TrackCharge()) <= 0) return; 
1253
1254   // Only do stuff is the track is in one of the strips. 
1255   TString vol(gMC->CurrentVolName());
1256   if (!vol.Contains("STR")) return;
1257
1258
1259   // Get the strip number.  Note, that GEANT numbers divisions from 1,
1260   // so we subtract one 
1261   Int_t strip;             
1262   gMC->CurrentVolID(strip);
1263   strip--;                 
1264
1265   // Get the phi division of the module 
1266   Int_t phiDiv;                         // * The phi division number (1 or 2)
1267   gMC->CurrentVolOffID(1, phiDiv);      //   in the module  
1268
1269   // Active volume number - not used. 
1270   // Int_t active;                         
1271   // gMC->CurrentVolOffID(2, active);      
1272
1273   // Get the module number in the ring. 
1274   Int_t module;                    
1275   gMC->CurrentVolOffID(3, module); 
1276   
1277   // Ring copy number - the same as the detector number - not used
1278   // Int_t ringCopy;                       // * Ring copy number
1279   // gMC->CurrentVolOffID(4, ringCopy);    //   Same as detector number 
1280   
1281   // Get the detector number from the path name 
1282   Int_t detector = Int_t((gMC->CurrentVolOffName(5)[3]) - 48);
1283
1284   // The sector number, calculated from module and phi division # 
1285   Int_t  sector =  2 * module + phiDiv - 1;
1286
1287   // The ring ID is encoded in the volume name 
1288   Char_t ring = vol[3];
1289
1290   // Get a pointer to the sub detector structure 
1291   AliFMDSubDetector* det = 0;
1292   switch (detector) {
1293   case 1: det = fFMD1; break;
1294   case 2: det = fFMD2; break;
1295   case 3: det = fFMD3; break;
1296   }
1297   if (!det) return;
1298
1299   // Get the current track position 
1300   TLorentzVector v;
1301   gMC->TrackPosition(v);
1302   // Check that the track is actually within the active area 
1303   Bool_t isWithin = det->CheckHit(ring, module, v.X(), v.Y());
1304   Bool_t entering = gMC->IsTrackEntering() && isWithin;
1305   Bool_t inside   = gMC->IsTrackInside()   && isWithin;
1306   Bool_t out      = (gMC->IsTrackExiting() 
1307                      || gMC->IsTrackDisappeared() 
1308                      || gMC->IsTrackStop() 
1309                      || !isWithin);
1310 // Reset the energy deposition for this track, and update some of
1311   // our parameters.
1312   if (entering) {
1313     fCurrentDeltaE = 0;
1314
1315     // Get production vertex and momentum of the track 
1316     fCurrentV = v;
1317     gMC->TrackMomentum(fCurrentP);
1318     fCurrentPdg = gMC->IdFromPDG(gMC->TrackPid());
1319
1320     // if (fAnalyser) 
1321     //   fAnalyser->Update(detector, ring, isWithin, v.X(), v.Y());
1322   }
1323   
1324   // If the track is inside, then update the energy deposition
1325   if (inside && fCurrentDeltaE >= 0) 
1326     fCurrentDeltaE += 1000 * gMC->Edep();
1327
1328   // The track exits the volume, or it disappeared in the volume, or
1329   // the track is stopped because it no longer fulfills the cuts
1330   // defined, then we create a hit. 
1331   if (out && fCurrentDeltaE >= 0) {
1332     fCurrentDeltaE += 1000 * gMC->Edep();
1333
1334     AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber(),
1335            detector, ring,  sector, strip,
1336            fCurrentV.X(), fCurrentV.Y(), fCurrentV.Z(),
1337            fCurrentP.X(), fCurrentP.Y(), fCurrentP.Z(), 
1338            fCurrentDeltaE, fCurrentPdg, fCurrentV.T());
1339     fCurrentDeltaE = -1;
1340   }
1341 }
1342 //___________________________________________________________________
1343 //
1344 // EOF
1345 //