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