]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMD.cxx
Several changes:
[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 /* $Id$ */
16 /** @file    AliFMD.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Sun Mar 26 17:59:18 2006
19     @brief   Implementation of AliFMD base class 
20 */
21 //____________________________________________________________________
22 //                                                                          
23 // Forward Multiplicity Detector based on Silicon wafers. This class
24 // is the driver for especially simulation. 
25 //
26 // The Forward Multiplicity Detector consists of 3 sub-detectors FMD1,
27 // FMD2, and FMD3, each of which has 1 or 2 rings of silicon sensors. 
28 //                                                       
29 // This is the base class for all FMD manager classes. 
30 //                    
31 // The actual code is done by various separate classes.   Below is
32 // diagram showing the relationship between the various FMD classes
33 // that handles the simulation
34 //
35 //
36 //       +----------+   +----------+   
37 //       | AliFMDv1 |   | AliFMDv0 |   
38 //       +----------+   +----------+   
39 //            |              |                    +-----------------+
40 //       +----+--------------+                 +--| AliFMDDigitizer |
41 //       |                                     |  +-----------------+
42 //       |           +---------------------+   |
43 //       |        +--| AliFMDBaseDigitizer |<--+
44 //       V     1  |  +---------------------+   |
45 //  +--------+<>--+                            |  +------------------+
46 //  | AliFMD |                                 +--| AliFMDSDigitizer |    
47 //  +--------+<>--+                               +------------------+       
48 //             1  |  +---------------------+
49 //                +--| AliFMDReconstructor |
50 //                   +---------------------+
51 //
52 // *  AliFMD 
53 //    This defines the interface for the various parts of AliROOT that
54 //    uses the FMD, like AliFMDSimulator, AliFMDDigitizer, 
55 //    AliFMDReconstructor, and so on. 
56 //
57 // *  AliFMDv0
58 //    This is a concrete implementation of the AliFMD interface. 
59 //    It is the responsibility of this class to create the FMD
60 //    geometry.
61 //
62 // *  AliFMDv1 
63 //    This is a concrete implementation of the AliFMD interface. 
64 //    It is the responsibility of this class to create the FMD
65 //    geometry, process hits in the FMD, and serve hits and digits to
66 //    the various clients. 
67 //  
68 // *  AliFMDSimulator
69 //    This is the base class for the FMD simulation tasks.   The
70 //    simulator tasks are responsible to implment the geoemtry, and
71 //    process hits. 
72 //                                                                          
73 // *  AliFMDReconstructor
74 //    This is a concrete implementation of the AliReconstructor that
75 //    reconstructs pseudo-inclusive-multiplicities from digits (raw or
76 //    from simulation)
77 //
78 // Calibration and geometry parameters are managed by separate
79 // singleton managers.  These are AliFMDGeometry and
80 // AliFMDParameters.  Please refer to these classes for more
81 // information on these.
82 //
83
84 // These files are not in the same directory, so there's no reason to
85 // ask the preprocessor to search in the current directory for these
86 // files by including them with `#include "..."' 
87 #include <cmath>                // __CMATH__
88 #include <TClonesArray.h>       // ROOT_TClonesArray
89 #include <TGeometry.h>          // ROOT_TGeomtry
90 #include <TNode.h>              // ROOT_TNode
91 #include <TXTRU.h>              // ROOT_TXTRU
92 #include <TRotMatrix.h>         // ROOT_TRotMatrix
93 #include <TTUBE.h>              // ROOT_TTUBE
94 #include <TTree.h>              // ROOT_TTree
95 #include <TBrowser.h>           // ROOT_TBrowser
96 #include <TVirtualMC.h>         // ROOT_TVirtualMC
97 #include <TVector2.h>           // ROOT_TVector2 
98 #include <TGeoManager.h>        // ROOT_TGeoManager
99
100 #include <AliRunDigitizer.h>    // ALIRUNDIGITIZER_H
101 #include <AliLoader.h>          // ALILOADER_H
102 #include <AliRun.h>             // ALIRUN_H
103 #include <AliMC.h>              // ALIMC_H
104 #include <AliMagF.h>            // ALIMAGF_H
105 // #include <AliLog.h>          // ALILOG_H
106 #include "AliFMDDebug.h" // Better debug macros
107 #include "AliFMD.h"             // ALIFMD_H
108 #include "AliFMDDigit.h"        // ALIFMDDIGIT_H
109 #include "AliFMDSDigit.h"       // ALIFMDSDIGIT_H
110 #include "AliFMDHit.h"          // ALIFMDHIT_H
111 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
112 #include "AliFMDDetector.h"     // ALIFMDDETECTOR_H
113 #include "AliFMDRing.h"         // ALIFMDRING_H
114 #include "AliFMDDigitizer.h"    // ALIFMDDIGITIZER_H
115 #include "AliFMDHitDigitizer.h" // ALIFMDSDIGITIZER_H
116 // #define USE_SSDIGITIZER 
117 //#ifdef USE_SSDIGITIZER
118 //# include "AliFMDSSDigitizer.h"       // ALIFMDSDIGITIZER_H
119 //#endif
120 // #include "AliFMDGeometryBuilder.h"
121 #include "AliFMDRawWriter.h"    // ALIFMDRAWWRITER_H
122
123 //____________________________________________________________________
124 ClassImp(AliFMD)
125 #if 0
126   ; // This is to keep Emacs from indenting the next line 
127 #endif 
128
129 //____________________________________________________________________
130 AliFMD::AliFMD()
131   : AliDetector(),
132     fSDigits(0), 
133     fNsdigits(0),
134     fDetailed(kTRUE),
135     fUseOld(kFALSE),
136     fUseAssembly(kTRUE),
137     fBad(0) 
138 {
139   //
140   // Default constructor for class AliFMD
141   //
142   AliFMDDebug(10, ("\tDefault CTOR"));
143   fHits        = 0;
144   fDigits      = 0;
145   fIshunt      = 0;
146   fBad         = new TClonesArray("AliFMDHit");
147 }
148
149 //____________________________________________________________________
150 AliFMD::AliFMD(const char *name, const char *title)
151   : AliDetector (name, title),
152     fSDigits(0),
153     fNsdigits(0),
154     fDetailed(kTRUE),
155     fUseOld(kFALSE),
156     fUseAssembly(kFALSE),
157     fBad(0)
158 {
159   //
160   // Standard constructor for Forward Multiplicity Detector
161   //
162   AliFMDDebug(10, ("\tStandard CTOR"));
163   fBad         = new TClonesArray("AliFMDHit");
164   
165   // Initialise Hit array
166   HitsArray();
167   gAlice->GetMCApp()->AddHitList(fHits);
168
169   // (S)Digits for the detectors disk
170   DigitsArray();
171   SDigitsArray();
172   
173   // CHC: What is this?
174   fIshunt = 0;
175   //PH  SetMarkerColor(kRed);
176   //PH  SetLineColor(kYellow);
177 }
178
179 //____________________________________________________________________
180 AliFMD::~AliFMD ()
181 {
182   // Destructor for base class AliFMD
183   if (fHits) {
184     fHits->Delete();
185     delete fHits;
186     fHits = 0;
187   }
188   if (fDigits) {
189     fDigits->Delete();
190     delete fDigits;
191     fDigits = 0;
192   }
193   if (fSDigits) {
194     fSDigits->Delete();
195     delete fSDigits;
196     fSDigits = 0;
197   }
198   if (fBad) {
199     fBad->Delete();
200     delete fBad;
201     fBad = 0;
202   }
203 }
204
205
206 //====================================================================
207 //
208 // GEometry ANd Traking
209 //
210 //____________________________________________________________________
211 void 
212 AliFMD::CreateGeometry()
213 {
214   //
215   // Create the geometry of Forward Multiplicity Detector.  The actual
216   // construction of the geometry is delegated to the class
217   // AliFMDGeometryBuilder, invoked by the singleton manager
218   // AliFMDGeometry. 
219   //
220   AliFMDGeometry*  fmd = AliFMDGeometry::Instance();
221   fmd->SetDetailed(fDetailed);
222   fmd->UseAssembly(fUseAssembly);
223   fmd->Build();
224 }    
225
226 //____________________________________________________________________
227 void AliFMD::CreateMaterials() 
228 {
229   // Define the materials and tracking mediums needed by the FMD
230   // simulation.   These mediums are made by sending the messages
231   // AliMaterial, AliMixture, and AliMedium to the passed AliModule
232   // object module.   The defined mediums are 
233   // 
234   //    FMD Si$         Silicon (active medium in sensors)
235   //    FMD C$          Carbon fibre (support cone for FMD3 and vacuum pipe)
236   //    FMD Al$         Aluminium (honeycomb support plates)
237   //    FMD PCB$        Printed Circuit Board (FEE board with VA1_3)
238   //    FMD Chip$       Electronics chips (currently not used)
239   //    FMD Air$        Air (Air in the FMD)
240   //    FMD Plastic$    Plastic (Support legs for the hybrid cards)
241   //
242   // The geometry builder should really be the one that creates the
243   // materials, but the architecture of AliROOT makes that design
244   // akward.  What should happen, was that the AliFMDGeometryBuilder
245   // made the mediums, and that this class retrives pointers from the
246   // TGeoManager, and registers the mediums here.  Alas, it's not
247   // really that easy. 
248   //
249   AliFMDDebug(10, ("\tCreating materials"));
250   // Get pointer to geometry singleton object. 
251   AliFMDGeometry* geometry = AliFMDGeometry::Instance();
252   geometry->Init();
253 #if 0
254   if (gGeoManager && gGeoManager->GetMedium("FMD Si$")) {
255     // We need to figure out the some stuff about the geometry
256     fmd->ExtractGeomInfo();
257     return;
258   }
259 #endif  
260   Int_t    id;
261   Double_t a                = 0;
262   Double_t z                = 0;
263   Double_t density          = 0;
264   Double_t radiationLength  = 0;
265   Double_t absorbtionLength = 999;
266   Int_t    fieldType        = gAlice->Field()->Integ();     // Field type 
267   Double_t maxField         = gAlice->Field()->Max();     // Field max.
268   Double_t maxBending       = 0;     // Max Angle
269   Double_t maxStepSize      = 0.001; // Max step size 
270   Double_t maxEnergyLoss    = 1;     // Max Delta E
271   Double_t precision        = 0.001; // Precision
272   Double_t minStepSize      = 0.001; // Minimum step size 
273  
274   // Silicon 
275   a                = 28.0855;
276   z                = 14.;
277   density          = geometry->GetSiDensity();
278   radiationLength  = 9.36;
279   maxBending       = 1;
280   maxStepSize      = .001;
281   precision        = .001;
282   minStepSize      = .001;
283   id               = kSiId;
284   AliMaterial(id, "Si$", a, z, density, radiationLength, absorbtionLength);
285   AliMedium(kSiId, "Si$", id,1,fieldType,maxField,maxBending,
286             maxStepSize,maxEnergyLoss,precision,minStepSize);
287   
288
289   // Carbon 
290   a                = 12.011;
291   z                = 6.;
292   density          = 2.265;
293   radiationLength  = 18.8;
294   maxBending       = 10;
295   maxStepSize      = .01;
296   precision        = .003;
297   minStepSize      = .003;
298   id               = kCarbonId;
299   AliMaterial(id, "Carbon$", a, z, density, radiationLength, absorbtionLength);
300   AliMedium(kCarbonId, "Carbon$", id,0,fieldType,maxField,maxBending,
301                     maxStepSize,maxEnergyLoss,precision,minStepSize);
302
303   // Aluminum
304   a                = 26.981539;
305   z                = 13.;
306   density          = 2.7;
307   radiationLength  = 8.9;
308   id               = kAlId;
309   AliMaterial(id, "Aluminum$",a,z, density, radiationLength, absorbtionLength);
310   AliMedium(kAlId, "Aluminum$", id, 0, fieldType, maxField, maxBending,
311             maxStepSize, maxEnergyLoss, precision, minStepSize);
312   
313   
314   // Copper 
315   a                = 63.546;
316   z                = 29;
317   density          =  8.96;
318   radiationLength  =  1.43;
319   id               = kCopperId;
320   AliMaterial(id, "Copper$", 
321                       a, z, density, radiationLength, absorbtionLength);
322   AliMedium(kCopperId, "Copper$", id, 0, fieldType, maxField, maxBending,
323             maxStepSize, maxEnergyLoss, precision, minStepSize);
324   
325
326   // Silicon chip 
327   {
328     Float_t as[] = { 12.0107,      14.0067,      15.9994,
329                       1.00794,     28.0855,     107.8682 };
330     Float_t zs[] = {  6.,           7.,           8.,
331                       1.,          14.,          47. };
332     Float_t ws[] = {  0.039730642,  0.001396798,  0.01169634,
333                       0.004367771,  0.844665,     0.09814344903 };
334     density          = 2.36436;
335     maxBending       = 10;
336     maxStepSize      = .01;
337     precision        = .003;
338     minStepSize      = .003;
339     id               = kSiChipId;
340     AliMixture(id, "Si Chip$", as, zs, density, 6, ws);
341     AliMedium(kSiChipId, "Si Chip$",  id, 0, fieldType, maxField, maxBending, 
342               maxStepSize, maxEnergyLoss, precision, minStepSize);
343   }
344   
345   // Kaption
346   {
347     Float_t as[] = { 1.00794,  12.0107,  14.010,   15.9994};
348     Float_t zs[] = { 1.,        6.,       7.,       8.};
349     Float_t ws[] = { 0.026362,  0.69113,  0.07327,  0.209235};
350     density          = 1.42;
351     maxBending       = 1;
352     maxStepSize      = .001;
353     precision        = .001;
354     minStepSize      = .001;
355     id               = kKaptonId;
356     AliMixture(id, "Kaption$", as, zs, density, 4, ws);
357     AliMedium(kKaptonId, "Kaption$", id,0,fieldType,maxField,maxBending,
358               maxStepSize,maxEnergyLoss,precision,minStepSize);
359   }
360
361   // Air
362   {
363     Float_t as[] = { 12.0107, 14.0067,   15.9994,  39.948 };
364     Float_t zs[] = {  6.,      7.,       8.,       18. };
365     Float_t ws[] = { 0.000124, 0.755267, 0.231781, 0.012827 }; 
366     density      = .00120479;
367     maxBending   = 1;
368     maxStepSize  = .001;
369     precision    = .001;
370     minStepSize  = .001;
371     id           = kAirId;
372     AliMixture(id, "Air$", as, zs, density, 4, ws);
373     AliMedium(kAirId, "Air$", id,0,fieldType,maxField,maxBending,
374               maxStepSize,maxEnergyLoss,precision,minStepSize);
375   }
376   
377   // PCB
378   {
379     Float_t zs[] = { 14.,         20.,         13.,         12.,
380                       5.,         22.,         11.,         19.,
381                      26.,          9.,          8.,          6.,
382                       7.,          1.};
383     Float_t as[] = { 28.0855,     40.078,      26.981538,   24.305, 
384                      10.811,      47.867,      22.98977,    39.0983,
385                      55.845,      18.9984,     15.9994,     12.0107,
386                      14.0067,      1.00794};
387     Float_t ws[] = {  0.15144894,  0.08147477,  0.04128158,  0.00904554, 
388                       0.01397570,  0.00287685,  0.00445114,  0.00498089,
389                       0.00209828,  0.00420000,  0.36043788,  0.27529426,
390                       0.01415852,  0.03427566};
391     density      = 1.8;
392     maxBending   = 1;
393     maxStepSize  = .001;
394     precision    = .001;
395     minStepSize  = .001;
396     id           = kPcbId;
397     AliMixture(id, "PCB$", as, zs, density, 14, ws);
398     AliMedium(kPcbId, "PCB$", id,0,fieldType,maxField,maxBending,
399               maxStepSize,maxEnergyLoss,precision,minStepSize);
400   }
401   
402   // Stainless steel
403   {
404     Float_t as[] = { 55.847, 51.9961, 58.6934, 28.0855 };
405     Float_t zs[] = { 26.,    24.,     28.,     14.     };
406     Float_t ws[] = { .715,   .18,     .1,      .005    };
407     density      = 7.88;
408     id           = kSteelId;
409     AliMixture(id, "Steel$", as, zs, density, 4, ws);
410     AliMedium(kSteelId, "Steel$", id, 0, fieldType, maxField, maxBending, 
411               maxStepSize, maxEnergyLoss, precision, minStepSize);
412   }
413   // Plastic 
414   {
415     Float_t as[] = { 1.01, 12.01 };
416     Float_t zs[] = { 1.,   6.    };
417     Float_t ws[] = { 1.,   1.    };
418     density      = 1.03;
419     maxBending   = 10;
420     maxStepSize  = .01;
421     precision    = .003;
422     minStepSize  = .003;
423     id           = kPlasticId;
424     AliMixture(id, "Plastic$", as, zs, density, -2, ws);
425     AliMedium(kPlasticId, "Plastic$", id,0,fieldType,maxField,maxBending,
426               maxStepSize,maxEnergyLoss,precision,minStepSize);
427   }
428
429 }
430
431 //____________________________________________________________________
432 void  
433 AliFMD::Init()
434 {
435   // Initialize the detector 
436   // 
437   AliFMDDebug(1, ("Initialising FMD detector object"));
438   TVirtualMC*      mc     = TVirtualMC::GetMC();
439   AliFMDGeometry*  fmd    = AliFMDGeometry::Instance();
440   const TArrayI&   actGeo = fmd->ActiveIds();
441   TArrayI          actVmc(actGeo.fN);
442   for (Int_t i = 0; i < actGeo.fN; i++) {
443     TGeoVolume *sens = gGeoManager->GetVolume(actGeo[i]);
444     if (!sens) {
445       AliError(Form("No TGeo volume for sensitive volume ID=%d",actGeo[i]));
446       continue;
447     }   
448     actVmc[i] = mc->VolId(sens->GetName());
449     AliFMDDebug(1, ("Active vol id # %d: %d changed to %d", 
450                     i, actGeo[i], actVmc[i]));
451   }
452   fmd->SetActive(actVmc.fArray, actVmc.fN);
453   // fmd->InitTransformations();
454 }
455
456 //____________________________________________________________________
457 void
458 AliFMD::FinishEvent()
459 {
460   // Called at the end of the an event in simulations.  If the debug
461   // level is high enough, then the `bad' hits are printed.
462   // 
463   if (AliLog::GetDebugLevel("FMD", "AliFMD") < 10) return;
464   if (fBad && fBad->GetEntries() > 0) {
465     AliWarning((Form("EndEvent", "got %d 'bad' hits", fBad->GetEntries())));
466     TIter next(fBad);
467     AliFMDHit* hit;
468     while ((hit = static_cast<AliFMDHit*>(next()))) hit->Print("D");
469     fBad->Clear();
470   }
471 }
472
473
474
475 //====================================================================
476 //
477 // Hit and Digit managment 
478 //
479 //____________________________________________________________________
480 void 
481 AliFMD::MakeBranch(Option_t * option)
482 {
483   // Create Tree branches for the FMD.
484   //
485   // Options:
486   //
487   //    H          Make a branch of TClonesArray of AliFMDHit's
488   //    D          Make a branch of TClonesArray of AliFMDDigit's
489   //    S          Make a branch of TClonesArray of AliFMDSDigit's
490   // 
491   const Int_t kBufferSize = 16000;
492   TString branchname(GetName());
493   TString opt(option);
494   
495   if (opt.Contains("H", TString::kIgnoreCase)) {
496     HitsArray();
497     AliDetector::MakeBranch(option); 
498   }
499   if (opt.Contains("D", TString::kIgnoreCase)) { 
500     DigitsArray();
501     MakeBranchInTree(fLoader->TreeD(), branchname.Data(),
502                      &fDigits, kBufferSize, 0);
503   }
504   if (opt.Contains("S", TString::kIgnoreCase)) { 
505     SDigitsArray();
506     MakeBranchInTree(fLoader->TreeS(), branchname.Data(),
507                      &fSDigits, kBufferSize, 0);
508   }
509 }
510
511 //____________________________________________________________________
512 void 
513 AliFMD::SetTreeAddress()
514 {
515   // Set branch address for the Hits, Digits, and SDigits Tree.
516   if (fLoader->TreeH()) HitsArray();
517   AliDetector::SetTreeAddress();
518
519   TTree *treeD = fLoader->TreeD();
520   if (treeD) {
521     DigitsArray();
522     TBranch* branch = treeD->GetBranch ("FMD");
523     if (branch) branch->SetAddress(&fDigits);
524   }
525
526   TTree *treeS = fLoader->TreeS();
527   if (treeS) {
528     SDigitsArray();
529     TBranch* branch = treeS->GetBranch ("FMD");
530     if (branch) branch->SetAddress(&fSDigits);
531   }
532 }
533
534 //____________________________________________________________________
535 void 
536 AliFMD::SetHitsAddressBranch(TBranch *b)
537 {
538   // Set the TClonesArray to read hits into. 
539   b->SetAddress(&fHits);
540 }
541 //____________________________________________________________________
542 void 
543 AliFMD::SetSDigitsAddressBranch(TBranch *b)
544 {
545   // Set the TClonesArray to read hits into. 
546   b->SetAddress(&fSDigits);
547 }
548
549 //____________________________________________________________________
550 void 
551 AliFMD::AddHit(Int_t track, Int_t *vol, Float_t *hits) 
552 {
553   // Add a hit to the hits tree 
554   // 
555   // The information of the two arrays are decoded as 
556   // 
557   // Parameters
558   //    track                Track #
559   //    ivol[0]  [UShort_t ] Detector # 
560   //    ivol[1]  [Char_t   ] Ring ID 
561   //    ivol[2]  [UShort_t ] Sector #
562   //    ivol[3]  [UShort_t ] Strip # 
563   //    hits[0]  [Float_t  ] Track's X-coordinate at hit 
564   //    hits[1]  [Float_t  ] Track's Y-coordinate at hit
565   //    hits[3]  [Float_t  ] Track's Z-coordinate at hit
566   //    hits[4]  [Float_t  ] X-component of track's momentum             
567   //    hits[5]  [Float_t  ] Y-component of track's momentum             
568   //    hits[6]  [Float_t  ] Z-component of track's momentum            
569   //    hits[7]  [Float_t  ] Energy deposited by track                  
570   //    hits[8]  [Int_t    ] Track's particle Id # 
571   //    hits[9]  [Float_t  ] Time when the track hit
572   // 
573   // 
574   AddHitByFields(track, 
575                  UShort_t(vol[0]),  // Detector # 
576                  Char_t(vol[1]),    // Ring ID
577                  UShort_t(vol[2]),  // Sector # 
578                  UShort_t(vol[3]),  // Strip # 
579                  hits[0],           // X
580                  hits[1],           // Y
581                  hits[2],           // Z
582                  hits[3],           // Px
583                  hits[4],           // Py
584                  hits[5],           // Pz
585                  hits[6],           // Energy loss 
586                  Int_t(hits[7]),    // PDG 
587                  hits[8]);          // Time
588 }
589
590 //____________________________________________________________________
591 AliFMDHit*
592 AliFMD::AddHitByFields(Int_t    track, 
593                        UShort_t detector, 
594                        Char_t   ring, 
595                        UShort_t sector, 
596                        UShort_t strip, 
597                        Float_t  x, 
598                        Float_t  y, 
599                        Float_t  z,
600                        Float_t  px, 
601                        Float_t  py, 
602                        Float_t  pz,
603                        Float_t  edep,
604                        Int_t    pdg,
605                        Float_t  t, 
606                        Float_t  l, 
607                        Bool_t   stop)
608 {
609   // Add a hit to the list
610   //
611   // Parameters:
612   // 
613   //    track     Track #
614   //    detector  Detector # (1, 2, or 3)                      
615   //    ring      Ring ID ('I' or 'O')
616   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
617   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
618   //    x         Track's X-coordinate at hit
619   //    y         Track's Y-coordinate at hit
620   //    z         Track's Z-coordinate at hit
621   //    px        X-component of track's momentum 
622   //    py        Y-component of track's momentum
623   //    pz        Z-component of track's momentum
624   //    edep      Energy deposited by track
625   //    pdg       Track's particle Id #
626   //    t         Time when the track hit 
627   //    l         Track length through the material. 
628   //    stop      Whether track was stopped or disappeared
629   // 
630   TClonesArray& a = *(HitsArray());
631   // Search through the list of already registered hits, and see if we
632   // find a hit with the same parameters.  If we do, then don't create
633   // a new hit, but rather update the energy deposited in the hit.
634   // This is done, so that a FLUKA based simulation will get the
635   // number of hits right, not just the enerrgy deposition. 
636   AliFMDHit* hit = 0;
637   for (Int_t i = 0; i < fNhits; i++) {
638     if (!a.At(i)) continue;
639     hit = static_cast<AliFMDHit*>(a.At(i));
640     if (hit->Detector() == detector 
641         && hit->Ring() == ring
642         && hit->Sector() == sector 
643         && hit->Strip() == strip
644         && hit->Track() == track) {
645       AliFMDDebug(1, ("already had a hit in FMD%d%c[%2d,%3d] for track # %d,"
646                        " adding energy (%f) to that hit (%f) -> %f", 
647                        detector, ring, sector, strip, track, edep, hit->Edep(),
648                        hit->Edep() + edep));
649       hit->SetEdep(hit->Edep() + edep);
650       return hit;
651     }
652   }
653   // If hit wasn't already registered, do so know. 
654   hit = new (a[fNhits]) AliFMDHit(fIshunt, track, detector, ring, sector, 
655                                   strip, x, y, z, px, py, pz, edep, pdg, t, 
656                                   l, stop);
657   fNhits++;
658   return hit;
659 }
660
661 //____________________________________________________________________
662 void 
663 AliFMD::AddDigit(Int_t* digits, Int_t*)
664 {
665   // Add a digit to the Digit tree 
666   // 
667   // Paramters 
668   //
669   //    digits[0]  [UShort_t] Detector #
670   //    digits[1]  [Char_t]   Ring ID
671   //    digits[2]  [UShort_t] Sector #
672   //    digits[3]  [UShort_t] Strip #
673   //    digits[4]  [UShort_t] ADC Count 
674   //    digits[5]  [Short_t]  ADC Count, -1 if not used
675   //    digits[6]  [Short_t]  ADC Count, -1 if not used 
676   // 
677   AddDigitByFields(UShort_t(digits[0]),  // Detector #
678                    Char_t(digits[1]),    // Ring ID
679                    UShort_t(digits[2]),  // Sector #
680                    UShort_t(digits[3]),  // Strip #
681                    UShort_t(digits[4]),  // ADC Count1 
682                    Short_t(digits[5]),   // ADC Count2 
683                    Short_t(digits[6]),   // ADC Count3 
684                    Short_t(digits[7])); 
685 }
686
687 //____________________________________________________________________
688 void 
689 AliFMD::AddDigitByFields(UShort_t detector, 
690                          Char_t   ring, 
691                          UShort_t sector, 
692                          UShort_t strip, 
693                          UShort_t count1, 
694                          Short_t  count2,
695                          Short_t  count3, 
696                          Short_t  count4)
697 {
698   // add a real digit - as coming from data
699   // 
700   // Parameters 
701   //
702   //    detector  Detector # (1, 2, or 3)                      
703   //    ring      Ring ID ('I' or 'O')
704   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
705   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
706   //    count1    ADC count (a 10-bit word)
707   //    count2    ADC count (a 10-bit word), or -1 if not used
708   //    count3    ADC count (a 10-bit word), or -1 if not used
709   TClonesArray& a = *(DigitsArray());
710   
711   new (a[fNdigits++]) 
712     AliFMDDigit(detector, ring, sector, strip, count1, count2, count3, count4);
713   AliFMDDebug(15, ("Adding digit # %5d/%5d for FMD%d%c[%2d,%3d]=(%d,%d,%d,%d)",
714                    fNdigits-1, a.GetEntriesFast(),
715                    detector, ring, sector, strip, 
716                    count1, count2, count3, count4));
717   
718 }
719
720 //____________________________________________________________________
721 void 
722 AliFMD::AddSDigit(Int_t* digits)
723 {
724   // Add a digit to the SDigit tree 
725   // 
726   // Paramters 
727   //
728   //    digits[0]  [UShort_t] Detector #
729   //    digits[1]  [Char_t]   Ring ID
730   //    digits[2]  [UShort_t] Sector #
731   //    digits[3]  [UShort_t] Strip #
732   //    digits[4]  [Float_t]  Total energy deposited 
733   //    digits[5]  [UShort_t] ADC Count 
734   //    digits[6]  [Short_t]  ADC Count, -1 if not used
735   //    digits[7]  [Short_t]  ADC Count, -1 if not used 
736   // 
737   AddSDigitByFields(UShort_t(digits[0]),   // Detector #
738                     Char_t(digits[1]),     // Ring ID
739                     UShort_t(digits[2]),   // Sector #
740                     UShort_t(digits[3]),   // Strip #
741                     Float_t(digits[4]),    // Edep
742                     UShort_t(digits[5]),   // ADC Count1 
743                     Short_t(digits[6]),    // ADC Count2 
744                     Short_t(digits[7]),    // ADC Count3 
745                     Short_t(digits[8]),    // ADC Count4
746                     UShort_t(digits[9]),   // N particles
747                     UShort_t(digits[10])); // N primaries
748    
749 }
750
751 //____________________________________________________________________
752 void 
753 AliFMD::AddSDigitByFields(UShort_t detector, 
754                           Char_t   ring, 
755                           UShort_t sector, 
756                           UShort_t strip, 
757                           Float_t  edep,
758                           UShort_t count1, 
759                           Short_t  count2,
760                           Short_t  count3, 
761                           Short_t  count4, 
762                           UShort_t ntot, 
763                           UShort_t nprim)
764 {
765   // add a summable digit
766   // 
767   // Parameters 
768   //
769   //    detector  Detector # (1, 2, or 3)                      
770   //    ring      Ring ID ('I' or 'O')
771   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
772   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
773   //    edep      Total energy deposited
774   //    count1    ADC count (a 10-bit word)
775   //    count2    ADC count (a 10-bit word), or -1 if not used
776   //    count3    ADC count (a 10-bit word), or -1 if not used
777   //
778   TClonesArray& a = *(SDigitsArray());
779   // AliFMDDebug(0, ("Adding sdigit # %d", fNsdigits));
780   
781   new (a[fNsdigits++]) 
782     AliFMDSDigit(detector, ring, sector, strip, edep, 
783                  count1, count2, count3, count4, ntot, nprim);
784 }
785
786 //____________________________________________________________________
787 void 
788 AliFMD::ResetSDigits()
789 {
790   // Reset number of digits and the digits array for this detector. 
791   //
792   fNsdigits   = 0;
793   if (fSDigits) fSDigits->Clear();
794 }
795
796
797 //____________________________________________________________________
798 TClonesArray*
799 AliFMD::HitsArray() 
800 {
801   // Initialize hit array if not already, and return pointer to it. 
802   if (!fHits) { 
803     fHits = new TClonesArray("AliFMDHit", 1000);
804     fNhits = 0;
805   }
806   return fHits;
807 }
808
809 //____________________________________________________________________
810 TClonesArray*
811 AliFMD::DigitsArray() 
812 {
813   // Initialize digit array if not already, and return pointer to it. 
814   if (!fDigits) { 
815     fDigits = new TClonesArray("AliFMDDigit", 1000);
816     fNdigits = 0;
817   }
818   return fDigits;
819 }
820
821 //____________________________________________________________________
822 TClonesArray*
823 AliFMD::SDigitsArray() 
824 {
825   // Initialize digit array if not already, and return pointer to it. 
826   if (!fSDigits) { 
827     fSDigits = new TClonesArray("AliFMDSDigit", 1000);
828     fNsdigits = 0;
829   }
830   return fSDigits;
831 }
832
833 //====================================================================
834 //
835 // Digitization 
836 //
837 //____________________________________________________________________
838 void 
839 AliFMD::Hits2Digits() 
840 {
841   // Create AliFMDDigit's from AliFMDHit's.  This is done by making a
842   // AliFMDDigitizer, and executing that code.
843   // 
844   AliFMDHitDigitizer digitizer(this, AliFMDHitDigitizer::kDigits);
845   digitizer.Init();
846   digitizer.Exec("");
847 }
848
849 //____________________________________________________________________
850 void 
851 AliFMD::Hits2SDigits() 
852 {
853   // Create AliFMDSDigit's from AliFMDHit's.  This is done by creating
854   // an AliFMDSDigitizer object, and executing it. 
855   // 
856   AliFMDHitDigitizer digitizer(this, AliFMDHitDigitizer::kSDigits);
857   digitizer.Init();
858   digitizer.Exec("");
859 }
860
861   
862 //____________________________________________________________________
863 AliDigitizer* 
864 AliFMD::CreateDigitizer(AliRunDigitizer* manager) const
865 {
866   // Create a digitizer object 
867   
868   /* This is what we probably _should_ do */
869   AliFMDBaseDigitizer* digitizer = 0;
870   
871 #ifdef USE_SSDIGITIZER
872   digitizer = new AliFMDSSDigitizer(manager);
873 #else 
874   /* This is what we actually do, and will work */
875 #if 0
876   AliInfo("SDigit->Digit conversion not really supported, "
877           "doing Hit->Digit conversion instead");
878 #endif
879   digitizer = new AliFMDDigitizer(manager);
880 #endif
881   return digitizer;
882 }
883
884 //====================================================================
885 //
886 // Raw data simulation 
887 //
888 //__________________________________________________________________
889 void 
890 AliFMD::Digits2Raw() 
891 {
892   // Turn digits into raw data. 
893   // 
894   // This uses the class AliFMDRawWriter to do the job.   Please refer
895   // to that class for more information. 
896   AliFMDRawWriter writer(this);
897   writer.Exec();
898 }
899
900
901 //====================================================================
902 //
903 // Utility 
904 //
905 //__________________________________________________________________
906 void 
907 AliFMD::Browse(TBrowser* b) 
908 {
909   // Browse this object. 
910   //
911   AliFMDDebug(30, ("\tBrowsing the FMD"));
912   AliDetector::Browse(b);
913   b->Add(AliFMDGeometry::Instance());
914 }
915
916 //____________________________________________________________________  
917 void
918 AliFMD::AddAlignableVolumes() const
919 {
920   //
921   // Create entries for alignable volumes associating the symbolic volume
922   // name with the corresponding volume path. Needs to be syncronized with
923   // eventual changes in the geometry.
924   // 
925   // This code was made by Raffaele Grosso <rgrosso@mail.cern.ch>.  I
926   // (cholm) will probably want to change it.   For one, I think it
927   // should be the job of the geometry manager to deal with this. 
928   AliInfo("Add FMD alignable volumes");
929   AliFMDGeometry::Instance()->SetAlignableVolumes();
930 #if 0  
931   for(size_t f = 1; f <= 3; f++){ // Detector 1,2,3
932     for(size_t tb =  0; tb <2 ; tb++){ // Top/Bottom 
933       char     stb = tb == 0 ? 'T' : 'B';
934       unsigned min = tb == 0 ? 0   : 5;
935
936       TString halfVol(Form("/ALIC_1/F%dM%c_%d", f, stb, f));
937       TString halfSym(halfVol);
938       if(!gGeoManager->SetAlignableEntry(halfSym.Data(),halfVol.Data()))
939         AliFatal(Form("Alignable entry %s not created. "
940                       "Volume path %s not valid", 
941                       halfSym.Data(),halfVol.Data()));
942       for(size_t io = 0; io < 2; io++){ // inner, outer 
943         if (f==1 && io==1) continue; // Only one ring in FMD1 
944         if(tb == 1 && io==1) min=10;
945         char     sio = (io == 0 ? 'I' : 'O');
946         unsigned nio = (io == 0 ? 3   : 9);
947         unsigned max = (io == 0 ? 5   : 10) + min;
948         
949         for(size_t i = min; i < max; i++) { // Modules
950           TString modVol(Form("%s/F%c%cV_7%d/F%cSE_%d", halfVol.Data(), 
951                               sio, stb, nio, sio, i));
952           TString modSym(modVol);
953           if(!gGeoManager->SetAlignableEntry(modSym.Data(),modVol.Data()))
954             AliFatal(Form("Alignable entry %s not created. "
955                           "Volume path %s not valid", 
956                           modSym.Data(), modVol.Data()));
957         }
958       }
959     }
960   }
961 #endif
962 }
963 //___________________________________________________________________
964 //
965 // EOF
966 //