]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMD.cxx
Fix some small error in the building of the FMD volume paths when generating the...
[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 <math.h>               // __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 "AliFMD.h"             // ALIFMD_H
107 #include "AliFMDDigit.h"        // ALIFMDDIGIT_H
108 #include "AliFMDSDigit.h"       // ALIFMDSDIGIT_H
109 #include "AliFMDHit.h"          // ALIFMDHIT_H
110 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
111 #include "AliFMDDetector.h"     // ALIFMDDETECTOR_H
112 #include "AliFMDRing.h"         // ALIFMDRING_H
113 #include "AliFMDDigitizer.h"    // ALIFMDDIGITIZER_H
114 #include "AliFMDSDigitizer.h"   // ALIFMDSDIGITIZER_H
115 // #include "AliFMDGeometryBuilder.h"
116 #include "AliFMDRawWriter.h"    // ALIFMDRAWWRITER_H
117 #include "AliFMDPoints.h"       // ALIFMDPOINTS_H
118
119 //____________________________________________________________________
120 ClassImp(AliFMD)
121 #if 0
122   ; // This is to keep Emacs from indenting the next line 
123 #endif 
124
125 //____________________________________________________________________
126 AliFMD::AliFMD()
127   : AliDetector(),
128     fSDigits(0), 
129     fNsdigits(0),
130     fDetailed(kTRUE),
131     fUseOld(kFALSE),
132     fUseAssembly(kTRUE),
133     fBad(0) 
134 {
135   //
136   // Default constructor for class AliFMD
137   //
138   AliDebug(10, "\tDefault CTOR");
139   fHits        = 0;
140   fDigits      = 0;
141   fIshunt      = 0;
142   fBad         = new TClonesArray("AliFMDHit");
143 }
144
145 //____________________________________________________________________
146 AliFMD::AliFMD(const AliFMD& other)
147   : AliDetector(other),
148     fSDigits(other.fSDigits), 
149     fNsdigits(other.fNsdigits),
150     fDetailed(other.fDetailed),
151     fUseOld(other.fUseOld),
152     fUseAssembly(other.fUseAssembly),
153     fBad(other.fBad)
154 {
155   // Copy constructor 
156 }
157
158 //____________________________________________________________________
159 AliFMD::AliFMD(const char *name, const char *title)
160   : AliDetector (name, title),
161     fSDigits(0),
162     fNsdigits(0),
163     fDetailed(kTRUE),
164     fUseOld(kFALSE),
165     fUseAssembly(kFALSE),
166     fBad(0)
167 {
168   //
169   // Standard constructor for Forward Multiplicity Detector
170   //
171   AliDebug(10, "\tStandard CTOR");
172   fBad         = new TClonesArray("AliFMDHit");
173   
174   // Initialise Hit array
175   HitsArray();
176   gAlice->GetMCApp()->AddHitList(fHits);
177
178   // (S)Digits for the detectors disk
179   DigitsArray();
180   SDigitsArray();
181   
182   // CHC: What is this?
183   fIshunt = 0;
184   SetMarkerColor(kRed);
185   SetLineColor(kYellow);
186 }
187
188 //____________________________________________________________________
189 AliFMD::~AliFMD ()
190 {
191   // Destructor for base class AliFMD
192   if (fHits) {
193     fHits->Delete();
194     delete fHits;
195     fHits = 0;
196   }
197   if (fDigits) {
198     fDigits->Delete();
199     delete fDigits;
200     fDigits = 0;
201   }
202   if (fSDigits) {
203     fSDigits->Delete();
204     delete fSDigits;
205     fSDigits = 0;
206   }
207   if (fBad) {
208     fBad->Delete();
209     delete fBad;
210     fBad = 0;
211   }
212 }
213
214 //____________________________________________________________________
215 AliFMD&
216 AliFMD::operator=(const AliFMD& other)
217 {
218   // Assignment operator
219   AliDetector::operator=(other);
220   fSDigits              = other.fSDigits; 
221   fNsdigits             = other.fNsdigits;
222   fDetailed             = other.fDetailed;
223   fBad                  = other.fBad;
224   return *this;
225 }
226
227 //====================================================================
228 //
229 // GEometry ANd Traking
230 //
231 //____________________________________________________________________
232 void 
233 AliFMD::CreateGeometry()
234 {
235   //
236   // Create the geometry of Forward Multiplicity Detector.  The actual
237   // construction of the geometry is delegated to the class
238   // AliFMDGeometryBuilder, invoked by the singleton manager
239   // AliFMDGeometry. 
240   //
241   AliFMDGeometry*  fmd = AliFMDGeometry::Instance();
242   fmd->SetDetailed(fDetailed);
243   fmd->UseAssembly(fUseAssembly);
244   fmd->Build();
245 }    
246
247 //____________________________________________________________________
248 void AliFMD::CreateMaterials() 
249 {
250   // Define the materials and tracking mediums needed by the FMD
251   // simulation.   These mediums are made by sending the messages
252   // AliMaterial, AliMixture, and AliMedium to the passed AliModule
253   // object module.   The defined mediums are 
254   // 
255   //    FMD Si$         Silicon (active medium in sensors)
256   //    FMD C$          Carbon fibre (support cone for FMD3 and vacuum pipe)
257   //    FMD Al$         Aluminium (honeycomb support plates)
258   //    FMD PCB$        Printed Circuit Board (FEE board with VA1_3)
259   //    FMD Chip$       Electronics chips (currently not used)
260   //    FMD Air$        Air (Air in the FMD)
261   //    FMD Plastic$    Plastic (Support legs for the hybrid cards)
262   //
263   // The geometry builder should really be the one that creates the
264   // materials, but the architecture of AliROOT makes that design
265   // akward.  What should happen, was that the AliFMDGeometryBuilder
266   // made the mediums, and that this class retrives pointers from the
267   // TGeoManager, and registers the mediums here.  Alas, it's not
268   // really that easy. 
269   //
270   AliDebug(10, "\tCreating materials");
271   // Get pointer to geometry singleton object. 
272   AliFMDGeometry* geometry = AliFMDGeometry::Instance();
273   geometry->Init();
274 #if 0
275   if (gGeoManager && gGeoManager->GetMedium("FMD Si$")) {
276     // We need to figure out the some stuff about the geometry
277     fmd->ExtractGeomInfo();
278     return;
279   }
280 #endif  
281   Int_t    id;
282   Double_t a                = 0;
283   Double_t z                = 0;
284   Double_t density          = 0;
285   Double_t radiationLength  = 0;
286   Double_t absorbtionLength = 999;
287   Int_t    fieldType        = gAlice->Field()->Integ();     // Field type 
288   Double_t maxField         = gAlice->Field()->Max();     // Field max.
289   Double_t maxBending       = 0;     // Max Angle
290   Double_t maxStepSize      = 0.001; // Max step size 
291   Double_t maxEnergyLoss    = 1;     // Max Delta E
292   Double_t precision        = 0.001; // Precision
293   Double_t minStepSize      = 0.001; // Minimum step size 
294  
295   // Silicon 
296   a                = 28.0855;
297   z                = 14.;
298   density          = geometry->GetSiDensity();
299   radiationLength  = 9.36;
300   maxBending       = 1;
301   maxStepSize      = .001;
302   precision        = .001;
303   minStepSize      = .001;
304   id               = kSiId;
305   AliMaterial(id, "Si$", a, z, density, radiationLength, absorbtionLength);
306   AliMedium(kSiId, "Si$", id,1,fieldType,maxField,maxBending,
307             maxStepSize,maxEnergyLoss,precision,minStepSize);
308   
309
310   // Carbon 
311   a                = 12.011;
312   z                = 6.;
313   density          = 2.265;
314   radiationLength  = 18.8;
315   maxBending       = 10;
316   maxStepSize      = .01;
317   precision        = .003;
318   minStepSize      = .003;
319   id               = kCarbonId;
320   AliMaterial(id, "Carbon$", a, z, density, radiationLength, absorbtionLength);
321   AliMedium(kCarbonId, "Carbon$", id,0,fieldType,maxField,maxBending,
322                     maxStepSize,maxEnergyLoss,precision,minStepSize);
323
324   // Aluminum
325   a                = 26.981539;
326   z                = 13.;
327   density          = 2.7;
328   radiationLength  = 8.9;
329   id               = kAlId;
330   AliMaterial(id, "Aluminum$",a,z, density, radiationLength, absorbtionLength);
331   AliMedium(kAlId, "Aluminum$", id, 0, fieldType, maxField, maxBending,
332             maxStepSize, maxEnergyLoss, precision, minStepSize);
333   
334   
335   // Copper 
336   a                = 63.546;
337   z                = 29;
338   density          =  8.96;
339   radiationLength  =  1.43;
340   id               = kCopperId;
341   AliMaterial(id, "Copper$", 
342                       a, z, density, radiationLength, absorbtionLength);
343   AliMedium(kCopperId, "Copper$", id, 0, fieldType, maxField, maxBending,
344             maxStepSize, maxEnergyLoss, precision, minStepSize);
345   
346
347   // Silicon chip 
348   {
349     Float_t as[] = { 12.0107,      14.0067,      15.9994,
350                       1.00794,     28.0855,     107.8682 };
351     Float_t zs[] = {  6.,           7.,           8.,
352                       1.,          14.,          47. };
353     Float_t ws[] = {  0.039730642,  0.001396798,  0.01169634,
354                       0.004367771,  0.844665,     0.09814344903 };
355     density          = 2.36436;
356     maxBending       = 10;
357     maxStepSize      = .01;
358     precision        = .003;
359     minStepSize      = .003;
360     id               = kSiChipId;
361     AliMixture(id, "Si Chip$", as, zs, density, 6, ws);
362     AliMedium(kSiChipId, "Si Chip$",  id, 0, fieldType, maxField, maxBending, 
363               maxStepSize, maxEnergyLoss, precision, minStepSize);
364   }
365   
366   // Kaption
367   {
368     Float_t as[] = { 1.00794,  12.0107,  14.010,   15.9994};
369     Float_t zs[] = { 1.,        6.,       7.,       8.};
370     Float_t ws[] = { 0.026362,  0.69113,  0.07327,  0.209235};
371     density          = 1.42;
372     maxBending       = 1;
373     maxStepSize      = .001;
374     precision        = .001;
375     minStepSize      = .001;
376     id               = kKaptonId;
377     AliMixture(id, "Kaption$", as, zs, density, 4, ws);
378     AliMedium(kKaptonId, "Kaption$", id,0,fieldType,maxField,maxBending,
379               maxStepSize,maxEnergyLoss,precision,minStepSize);
380   }
381
382   // Air
383   {
384     Float_t as[] = { 12.0107, 14.0067,   15.9994,  39.948 };
385     Float_t zs[] = {  6.,      7.,       8.,       18. };
386     Float_t ws[] = { 0.000124, 0.755267, 0.231781, 0.012827 }; 
387     density      = .00120479;
388     maxBending   = 1;
389     maxStepSize  = .001;
390     precision    = .001;
391     minStepSize  = .001;
392     id           = kAirId;
393     AliMixture(id, "Air$", as, zs, density, 4, ws);
394     AliMedium(kAirId, "Air$", id,0,fieldType,maxField,maxBending,
395               maxStepSize,maxEnergyLoss,precision,minStepSize);
396   }
397   
398   // PCB
399   {
400     Float_t zs[] = { 14.,         20.,         13.,         12.,
401                       5.,         22.,         11.,         19.,
402                      26.,          9.,          8.,          6.,
403                       7.,          1.};
404     Float_t as[] = { 28.0855,     40.078,      26.981538,   24.305, 
405                      10.811,      47.867,      22.98977,    39.0983,
406                      55.845,      18.9984,     15.9994,     12.0107,
407                      14.0067,      1.00794};
408     Float_t ws[] = {  0.15144894,  0.08147477,  0.04128158,  0.00904554, 
409                       0.01397570,  0.00287685,  0.00445114,  0.00498089,
410                       0.00209828,  0.00420000,  0.36043788,  0.27529426,
411                       0.01415852,  0.03427566};
412     density      = 1.8;
413     maxBending   = 1;
414     maxStepSize  = .001;
415     precision    = .001;
416     minStepSize  = .001;
417     id           = kPcbId;
418     AliMixture(id, "PCB$", as, zs, density, 14, ws);
419     AliMedium(kPcbId, "PCB$", id,0,fieldType,maxField,maxBending,
420               maxStepSize,maxEnergyLoss,precision,minStepSize);
421   }
422   
423   // Plastic 
424   {
425     Float_t as[] = { 1.01, 12.01 };
426     Float_t zs[] = { 1.,   6.    };
427     Float_t ws[] = { 1.,   1.    };
428     density      = 1.03;
429     maxBending   = 10;
430     maxStepSize  = .01;
431     precision    = .003;
432     minStepSize  = .003;
433     id           = kPlasticId;
434     AliMixture(id, "Plastic$", as, zs, density, -2, ws);
435     AliMedium(kPlasticId, "Plastic$", id,0,fieldType,maxField,maxBending,
436               maxStepSize,maxEnergyLoss,precision,minStepSize);
437   }
438 }
439
440 //____________________________________________________________________
441 void  
442 AliFMD::Init()
443 {
444   // Initialize the detector 
445   // 
446   AliDebug(1, "Initialising FMD detector object");
447   // AliFMDGeometry*  fmd = AliFMDGeometry::Instance();
448   // fmd->InitTransformations();
449 }
450
451 //____________________________________________________________________
452 void
453 AliFMD::FinishEvent()
454 {
455   // Called at the end of the an event in simulations.  If the debug
456   // level is high enough, then the `bad' hits are printed.
457   // 
458   if (AliLog::GetDebugLevel("FMD", "AliFMD") < 10) return;
459   if (fBad && fBad->GetEntries() > 0) {
460     AliWarning((Form("EndEvent", "got %d 'bad' hits", fBad->GetEntries())));
461     TIter next(fBad);
462     AliFMDHit* hit;
463     while ((hit = static_cast<AliFMDHit*>(next()))) hit->Print("D");
464     fBad->Clear();
465   }
466 }
467
468
469 //====================================================================
470 //
471 // Graphics and event display
472 //
473 //____________________________________________________________________
474 void 
475 AliFMD::BuildGeometry()
476 {
477   //
478   // Build simple ROOT TNode geometry for event display. With the new
479   // geometry modeller, TGeoManager, this seems rather redundant. 
480   AliDebug(10, "\tCreating a simplified geometry");
481
482   AliFMDGeometry* fmd = AliFMDGeometry::Instance();
483   
484   static TXTRU*     innerShape = 0;
485   static TXTRU*     outerShape = 0;
486   static TObjArray* innerRot   = 0;
487   static TObjArray* outerRot   = 0;
488
489   if (!innerShape || !outerShape) {
490     // Make the shapes for the modules 
491     for (Int_t i = 0; i < 2; i++) {
492       AliFMDRing* r = 0;
493       switch (i) {
494       case 0: r = fmd->GetRing('I'); break;
495       case 1: r = fmd->GetRing('O'); break;
496       }
497       if (!r) {
498         AliError(Form("no ring found for i=%d", i));
499         return;
500       }
501       Double_t    siThick  = r->GetSiThickness();
502       const Int_t knv      = r->GetNVerticies();
503       Double_t    theta    = r->GetTheta();
504       Int_t       nmod     = r->GetNModules();
505       
506       TXTRU* shape = new TXTRU(r->GetName(), r->GetTitle(), "void", knv, 2);
507       for (Int_t j = 0; j < knv; j++) {
508         TVector2* vv = r->GetVertex(knv - 1 - j);
509         shape->DefineVertex(j, vv->X(), vv->Y());
510       }
511       shape->DefineSection(0, -siThick / 2, 1, 0, 0);
512       shape->DefineSection(1, +siThick / 2, 1, 0, 0);
513       shape->SetLineColor(GetLineColor());
514       
515       TObjArray* rots = new TObjArray(nmod);
516       for (Int_t j = 0; j < nmod; j++) {
517         Double_t th = (j + .5) * theta * 2;
518         TString name(Form("FMD_ring_%c_rot_%02d", r->GetId(), j));
519         TString title(Form("FMD Ring %c Rotation # %d", r->GetId(), j));
520         TRotMatrix* rot = new TRotMatrix(name.Data(), title.Data(),
521                                          90, th, 90, fmod(90+th,360), 0, 0);
522         rots->AddAt(rot, j);
523       }
524       
525       switch (r->GetId()) {
526       case 'i':
527       case 'I': innerShape = shape; innerRot = rots; break;
528       case 'o':
529       case 'O': outerShape = shape; outerRot = rots; break;
530       }
531     }
532   }
533   
534   TNode* top = gAlice->GetGeometry()->GetNode("alice");
535   
536   for (Int_t i = 1; i <= 3; i++) {
537     AliFMDDetector* det = fmd->GetDetector(i);
538     if (!det) {
539       Warning("BuildGeometry", "FMD%d seems to be disabled", i);
540       continue;
541     }
542     Double_t w  = 0;
543     Double_t rh = det->GetRing('I')->GetHighR();
544     Char_t   id = 'I';
545     if (det->GetRing('O')) {
546       w  = TMath::Abs(det->GetRingZ('O') - det->GetRingZ('I'));
547       id = (TMath::Abs(det->GetRingZ('O')) 
548             > TMath::Abs(det->GetRingZ('I')) ? 'O' : 'I');
549       rh = det->GetRing('O')->GetHighR();
550     }
551     w += (det->GetRing(id)->GetModuleSpacing() +
552           det->GetRing(id)->GetSiThickness());
553     TShape* shape = new TTUBE(det->GetName(), det->GetTitle(), "void",
554                               det->GetRing('I')->GetLowR(), rh, w / 2);
555     Double_t z = (det->GetRingZ('I') - w / 2);
556     if (z > 0) z += det->GetRing(id)->GetModuleSpacing();
557     top->cd();
558     TNode* node = new TNode(det->GetName(), det->GetTitle(), shape, 
559                             0, 0, z, 0);
560     fNodes->Add(node);
561     
562     for (Int_t j = 0; j < 2; j++) {
563       AliFMDRing* r      = 0;
564       TShape*     rshape = 0;
565       TObjArray*  rots   = 0;
566       switch (j) {
567       case 0: 
568         r = det->GetRing('I'); rshape = innerShape; rots = innerRot; break;
569       case 1: 
570         r = det->GetRing('O'); rshape = outerShape; rots = outerRot; break;
571       }
572       if (!r) continue;
573       
574       Double_t    siThick  = r->GetSiThickness();
575       Int_t       nmod     = r->GetNModules();
576       Double_t    modspace = r->GetModuleSpacing();
577       Double_t    rz       = - (z - det->GetRingZ(r->GetId()));
578       
579       for (Int_t k = 0; k < nmod; k++) {
580         node->cd();
581         Double_t    offz    = (k % 2 == 1 ? modspace : 0);
582         TRotMatrix* rot     = static_cast<TRotMatrix*>(rots->At(k));
583         TString name(Form("%s%c_module_%02d", det->GetName(), r->GetId(),k));
584         TString title(Form("%s%c Module %d", det->GetName(), r->GetId(),k));
585         TNode* mnod = new TNode(name.Data(), title.Data(), rshape, 
586                                 0, 0, rz - siThick / 2 
587                                 + TMath::Sign(offz,z), rot);
588         mnod->SetLineColor(GetLineColor());
589         fNodes->Add(mnod);
590       } // for (Int_t k = 0 ; ...)
591     } // for (Int_t j = 0 ; ...)
592   } // for (Int_t i = 1 ; ...)
593 }
594
595 //____________________________________________________________________
596 void 
597 AliFMD::LoadPoints(Int_t /* track */) 
598 {
599   // Store x, y, z of all hits in memory for display. 
600   // 
601   // Normally, the hits are drawn using TPolyMarker3D - however, that
602   // is not very useful for the FMD.  Therefor, this member function
603   // is overloaded to make TMarker3D, via the class AliFMDPoints.
604   // AliFMDPoints is a local class. 
605   //
606   if (!fHits) {
607     AliError(Form("fHits == 0. Name is %s",GetName()));
608     return;
609   }
610   Int_t nHits = fHits->GetEntriesFast();
611   if (nHits == 0) {
612     return;
613   }
614   Int_t tracks = gAlice->GetMCApp()->GetNtrack();
615   if (fPoints == 0) fPoints = new TObjArray(2 * tracks);
616
617   // Get geometry 
618   AliFMDGeometry* geom = AliFMDGeometry::Instance();
619   geom->Init();
620   geom->InitTransformations();
621
622   // Now make markers for each hit  
623   // AliInfo(Form("Drawing %d hits (have %d points) for track %d", 
624   //              nHits, fPoints->GetEntriesFast(), track));
625   for (Int_t ihit = 0; ihit < nHits; ihit++) {
626     AliFMDHit* hit = static_cast<AliFMDHit*>(fHits->At(ihit));
627     if (!hit) continue;
628     Double_t edep    = hit->Edep();
629     Double_t m       = hit->M();
630     Double_t poverm  = (m == 0 ? 0 : hit->P());
631     Double_t absQ    = TMath::Abs(hit->Q());
632     Bool_t   bad     = kFALSE;
633     // This `if' is to debug abnormal energy depositions.  We trigger on
634     // p/m approx larger than or equal to a MIP, and a large edep - more 
635     // than 1 keV - a MIP is 100 eV. 
636     if (edep > absQ * absQ && poverm > 1) bad = kTRUE;
637
638     AliFMDPoints* p1 = new AliFMDPoints(hit, GetMarkerColor());
639     // AliPoints* p1 = new AliPoints();
640     // p1->SetMarkerColor(GetMarkerColor());
641     // p1->SetMarkerSize(GetMarkerSize());
642     // p1->SetPoint(0, hit->X(), hit->Y(), hit->Z());
643     p1->SetDetector(this);
644     p1->SetParticle(hit->GetTrack());
645     fPoints->AddAt(p1, hit->GetTrack());
646     if (bad) {
647       p1->SetMarkerColor(4);
648       // p1->SetMarkerSize(2 * GetMarkerSize());
649     }
650     
651     Double_t x, y, z;
652     geom->Detector2XYZ(hit->Detector(), hit->Ring(), hit->Sector(), 
653                        hit->Strip(), x, y, z);
654     AliFMDPoints* p = new AliFMDPoints(hit, 3);
655     // AliPoints* p = new AliPoints();
656     // p->SetMarkerColor(3);
657     // p->SetMarkerSize(GetMarkerSize());
658     // p->SetPoint(0, x, y, z);
659     p->SetDetector(this);
660     p->SetParticle(hit->GetTrack());
661     p->SetXYZ(x, y, z);
662     p->SetMarkerColor(3);
663     fPoints->AddAt(p, tracks+hit->GetTrack());
664     if (bad) {
665       p->SetMarkerColor(5);
666       // p->SetMarkerSize(2 * GetMarkerSize());
667     }
668     // AliInfo(Form("Adding point at %d", tracks+hit->GetTrack()));
669   }
670 }
671
672 //____________________________________________________________________
673 void 
674 AliFMD::DrawDetector()
675 {
676   // Draw a shaded view of the Forward multiplicity detector.  This
677   // isn't really useful anymore. 
678   AliDebug(10, "\tDraw detector");
679 }
680
681 //____________________________________________________________________
682 Int_t 
683 AliFMD::DistancetoPrimitive(Int_t, Int_t)
684 {
685   // Calculate the distance from the mouse to the FMD on the screen
686   // Dummy routine.
687   //
688   return 9999;
689 }
690
691 //====================================================================
692 //
693 // Hit and Digit managment 
694 //
695 //____________________________________________________________________
696 void 
697 AliFMD::MakeBranch(Option_t * option)
698 {
699   // Create Tree branches for the FMD.
700   //
701   // Options:
702   //
703   //    H          Make a branch of TClonesArray of AliFMDHit's
704   //    D          Make a branch of TClonesArray of AliFMDDigit's
705   //    S          Make a branch of TClonesArray of AliFMDSDigit's
706   // 
707   const Int_t kBufferSize = 16000;
708   TString branchname(GetName());
709   TString opt(option);
710   
711   if (opt.Contains("H", TString::kIgnoreCase)) {
712     HitsArray();
713     AliDetector::MakeBranch(option); 
714   }
715   if (opt.Contains("D", TString::kIgnoreCase)) { 
716     DigitsArray();
717     MakeBranchInTree(fLoader->TreeD(), branchname.Data(),
718                      &fDigits, kBufferSize, 0);
719   }
720   if (opt.Contains("S", TString::kIgnoreCase)) { 
721     SDigitsArray();
722     MakeBranchInTree(fLoader->TreeS(), branchname.Data(),
723                      &fSDigits, kBufferSize, 0);
724   }
725 }
726
727 //____________________________________________________________________
728 void 
729 AliFMD::SetTreeAddress()
730 {
731   // Set branch address for the Hits, Digits, and SDigits Tree.
732   if (fLoader->TreeH()) HitsArray();
733   AliDetector::SetTreeAddress();
734
735   TTree *treeD = fLoader->TreeD();
736   if (treeD) {
737     DigitsArray();
738     TBranch* branch = treeD->GetBranch ("FMD");
739     if (branch) branch->SetAddress(&fDigits);
740   }
741
742   TTree *treeS = fLoader->TreeS();
743   if (treeS) {
744     SDigitsArray();
745     TBranch* branch = treeS->GetBranch ("FMD");
746     if (branch) branch->SetAddress(&fSDigits);
747   }
748 }
749
750 //____________________________________________________________________
751 void 
752 AliFMD::SetHitsAddressBranch(TBranch *b)
753 {
754   // Set the TClonesArray to read hits into. 
755   b->SetAddress(&fHits);
756 }
757
758 //____________________________________________________________________
759 void 
760 AliFMD::AddHit(Int_t track, Int_t *vol, Float_t *hits) 
761 {
762   // Add a hit to the hits tree 
763   // 
764   // The information of the two arrays are decoded as 
765   // 
766   // Parameters
767   //    track                Track #
768   //    ivol[0]  [UShort_t ] Detector # 
769   //    ivol[1]  [Char_t   ] Ring ID 
770   //    ivol[2]  [UShort_t ] Sector #
771   //    ivol[3]  [UShort_t ] Strip # 
772   //    hits[0]  [Float_t  ] Track's X-coordinate at hit 
773   //    hits[1]  [Float_t  ] Track's Y-coordinate at hit
774   //    hits[3]  [Float_t  ] Track's Z-coordinate at hit
775   //    hits[4]  [Float_t  ] X-component of track's momentum             
776   //    hits[5]  [Float_t  ] Y-component of track's momentum             
777   //    hits[6]  [Float_t  ] Z-component of track's momentum            
778   //    hits[7]  [Float_t  ] Energy deposited by track                  
779   //    hits[8]  [Int_t    ] Track's particle Id # 
780   //    hits[9]  [Float_t  ] Time when the track hit
781   // 
782   // 
783   AddHitByFields(track, 
784                  UShort_t(vol[0]),  // Detector # 
785                  Char_t(vol[1]),    // Ring ID
786                  UShort_t(vol[2]),  // Sector # 
787                  UShort_t(vol[3]),  // Strip # 
788                  hits[0],           // X
789                  hits[1],           // Y
790                  hits[2],           // Z
791                  hits[3],           // Px
792                  hits[4],           // Py
793                  hits[5],           // Pz
794                  hits[6],           // Energy loss 
795                  Int_t(hits[7]),    // PDG 
796                  hits[8]);          // Time
797 }
798
799 //____________________________________________________________________
800 AliFMDHit*
801 AliFMD::AddHitByFields(Int_t    track, 
802                        UShort_t detector, 
803                        Char_t   ring, 
804                        UShort_t sector, 
805                        UShort_t strip, 
806                        Float_t  x, 
807                        Float_t  y, 
808                        Float_t  z,
809                        Float_t  px, 
810                        Float_t  py, 
811                        Float_t  pz,
812                        Float_t  edep,
813                        Int_t    pdg,
814                        Float_t  t, 
815                        Float_t  l, 
816                        Bool_t   stop)
817 {
818   // Add a hit to the list
819   //
820   // Parameters:
821   // 
822   //    track     Track #
823   //    detector  Detector # (1, 2, or 3)                      
824   //    ring      Ring ID ('I' or 'O')
825   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
826   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
827   //    x         Track's X-coordinate at hit
828   //    y         Track's Y-coordinate at hit
829   //    z         Track's Z-coordinate at hit
830   //    px        X-component of track's momentum 
831   //    py        Y-component of track's momentum
832   //    pz        Z-component of track's momentum
833   //    edep      Energy deposited by track
834   //    pdg       Track's particle Id #
835   //    t         Time when the track hit 
836   //    l         Track length through the material. 
837   //    stop      Whether track was stopped or disappeared
838   // 
839   TClonesArray& a = *(HitsArray());
840   // Search through the list of already registered hits, and see if we
841   // find a hit with the same parameters.  If we do, then don't create
842   // a new hit, but rather update the energy deposited in the hit.
843   // This is done, so that a FLUKA based simulation will get the
844   // number of hits right, not just the enerrgy deposition. 
845   AliFMDHit* hit = 0;
846   for (Int_t i = 0; i < fNhits; i++) {
847     if (!a.At(i)) continue;
848     hit = static_cast<AliFMDHit*>(a.At(i));
849     if (hit->Detector() == detector 
850         && hit->Ring() == ring
851         && hit->Sector() == sector 
852         && hit->Strip() == strip
853         && hit->Track() == track) {
854       AliDebug(1, Form("already had a hit in FMD%d%c[%2d,%3d] for track # %d,"
855                        " adding energy (%f) to that hit (%f) -> %f", 
856                        detector, ring, sector, strip, track, edep, hit->Edep(),
857                        hit->Edep() + edep));
858       hit->SetEdep(hit->Edep() + edep);
859       return hit;
860     }
861   }
862   // If hit wasn't already registered, do so know. 
863   hit = new (a[fNhits]) AliFMDHit(fIshunt, track, detector, ring, sector, 
864                                   strip, x, y, z, px, py, pz, edep, pdg, t, 
865                                   l, stop);
866   fNhits++;
867   return hit;
868 }
869
870 //____________________________________________________________________
871 void 
872 AliFMD::AddDigit(Int_t* digits, Int_t*)
873 {
874   // Add a digit to the Digit tree 
875   // 
876   // Paramters 
877   //
878   //    digits[0]  [UShort_t] Detector #
879   //    digits[1]  [Char_t]   Ring ID
880   //    digits[2]  [UShort_t] Sector #
881   //    digits[3]  [UShort_t] Strip #
882   //    digits[4]  [UShort_t] ADC Count 
883   //    digits[5]  [Short_t]  ADC Count, -1 if not used
884   //    digits[6]  [Short_t]  ADC Count, -1 if not used 
885   // 
886   AddDigitByFields(UShort_t(digits[0]),  // Detector #
887                    Char_t(digits[1]),    // Ring ID
888                    UShort_t(digits[2]),  // Sector #
889                    UShort_t(digits[3]),  // Strip #
890                    UShort_t(digits[4]),  // ADC Count1 
891                    Short_t(digits[5]),   // ADC Count2 
892                    Short_t(digits[6]));  // ADC Count3 
893 }
894
895 //____________________________________________________________________
896 void 
897 AliFMD::AddDigitByFields(UShort_t detector, 
898                          Char_t   ring, 
899                          UShort_t sector, 
900                          UShort_t strip, 
901                          UShort_t count1, 
902                          Short_t  count2,
903                          Short_t  count3)
904 {
905   // add a real digit - as coming from data
906   // 
907   // Parameters 
908   //
909   //    detector  Detector # (1, 2, or 3)                      
910   //    ring      Ring ID ('I' or 'O')
911   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
912   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
913   //    count1    ADC count (a 10-bit word)
914   //    count2    ADC count (a 10-bit word), or -1 if not used
915   //    count3    ADC count (a 10-bit word), or -1 if not used
916   TClonesArray& a = *(DigitsArray());
917   
918   new (a[fNdigits++]) 
919     AliFMDDigit(detector, ring, sector, strip, count1, count2, count3);
920 }
921
922 //____________________________________________________________________
923 void 
924 AliFMD::AddSDigit(Int_t* digits)
925 {
926   // Add a digit to the SDigit tree 
927   // 
928   // Paramters 
929   //
930   //    digits[0]  [UShort_t] Detector #
931   //    digits[1]  [Char_t]   Ring ID
932   //    digits[2]  [UShort_t] Sector #
933   //    digits[3]  [UShort_t] Strip #
934   //    digits[4]  [Float_t]  Total energy deposited 
935   //    digits[5]  [UShort_t] ADC Count 
936   //    digits[6]  [Short_t]  ADC Count, -1 if not used
937   //    digits[7]  [Short_t]  ADC Count, -1 if not used 
938   // 
939   AddSDigitByFields(UShort_t(digits[0]),  // Detector #
940                     Char_t(digits[1]),    // Ring ID
941                     UShort_t(digits[2]),  // Sector #
942                     UShort_t(digits[3]),  // Strip #
943                     Float_t(digits[4]),   // Edep
944                     UShort_t(digits[5]),  // ADC Count1 
945                     Short_t(digits[6]),   // ADC Count2 
946                     Short_t(digits[7]));  // ADC Count3 
947 }
948
949 //____________________________________________________________________
950 void 
951 AliFMD::AddSDigitByFields(UShort_t detector, 
952                           Char_t   ring, 
953                           UShort_t sector, 
954                           UShort_t strip, 
955                           Float_t  edep,
956                           UShort_t count1, 
957                           Short_t  count2,
958                           Short_t  count3)
959 {
960   // add a summable digit
961   // 
962   // Parameters 
963   //
964   //    detector  Detector # (1, 2, or 3)                      
965   //    ring      Ring ID ('I' or 'O')
966   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
967   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
968   //    edep      Total energy deposited
969   //    count1    ADC count (a 10-bit word)
970   //    count2    ADC count (a 10-bit word), or -1 if not used
971   //    count3    ADC count (a 10-bit word), or -1 if not used
972   //
973   TClonesArray& a = *(SDigitsArray());
974   
975   new (a[fNsdigits++]) 
976     AliFMDSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
977 }
978
979 //____________________________________________________________________
980 void 
981 AliFMD::ResetSDigits()
982 {
983   // Reset number of digits and the digits array for this detector. 
984   //
985   fNsdigits   = 0;
986   if (fSDigits) fSDigits->Clear();
987 }
988
989
990 //____________________________________________________________________
991 TClonesArray*
992 AliFMD::HitsArray() 
993 {
994   // Initialize hit array if not already, and return pointer to it. 
995   if (!fHits) { 
996     fHits = new TClonesArray("AliFMDHit", 1000);
997     fNhits = 0;
998   }
999   return fHits;
1000 }
1001
1002 //____________________________________________________________________
1003 TClonesArray*
1004 AliFMD::DigitsArray() 
1005 {
1006   // Initialize digit array if not already, and return pointer to it. 
1007   if (!fDigits) { 
1008     fDigits = new TClonesArray("AliFMDDigit", 1000);
1009     fNdigits = 0;
1010   }
1011   return fDigits;
1012 }
1013
1014 //____________________________________________________________________
1015 TClonesArray*
1016 AliFMD::SDigitsArray() 
1017 {
1018   // Initialize digit array if not already, and return pointer to it. 
1019   if (!fSDigits) { 
1020     fSDigits = new TClonesArray("AliFMDSDigit", 1000);
1021     fNsdigits = 0;
1022   }
1023   return fSDigits;
1024 }
1025
1026 //====================================================================
1027 //
1028 // Digitization 
1029 //
1030 //____________________________________________________________________
1031 void 
1032 AliFMD::Hits2Digits() 
1033 {
1034   // Create AliFMDDigit's from AliFMDHit's.  This is done by making a
1035   // AliFMDDigitizer, and executing that code.
1036   // 
1037   Warning("Hits2Digits", "Try not to use this method.\n"
1038           "Instead, use AliSimulator");
1039   AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
1040   manager->SetInputStream(0, "galice.root");
1041   manager->SetOutputFile("H2Dfile");
1042   
1043   /* AliDigitizer* dig =*/ CreateDigitizer(manager);
1044   manager->Exec("");
1045   delete manager;
1046 }
1047
1048 //____________________________________________________________________
1049 void 
1050 AliFMD::Hits2SDigits() 
1051 {
1052   // Create AliFMDSDigit's from AliFMDHit's.  This is done by creating
1053   // an AliFMDSDigitizer object, and executing it. 
1054   // 
1055   AliFMDSDigitizer* digitizer = new AliFMDSDigitizer("galice.root");
1056   digitizer->Exec("");
1057   delete digitizer;
1058 }
1059
1060   
1061 //____________________________________________________________________
1062 AliDigitizer* 
1063 AliFMD::CreateDigitizer(AliRunDigitizer* manager) const
1064 {
1065   // Create a digitizer object 
1066   AliFMDDigitizer* digitizer = new AliFMDDigitizer(manager);
1067   return digitizer;
1068 }
1069
1070 //====================================================================
1071 //
1072 // Raw data simulation 
1073 //
1074 //__________________________________________________________________
1075 void 
1076 AliFMD::Digits2Raw() 
1077 {
1078   // Turn digits into raw data. 
1079   // 
1080   // This uses the class AliFMDRawWriter to do the job.   Please refer
1081   // to that class for more information. 
1082   AliFMDRawWriter writer(this);
1083   writer.Exec();
1084 }
1085
1086
1087 //====================================================================
1088 //
1089 // Utility 
1090 //
1091 //__________________________________________________________________
1092 void 
1093 AliFMD::Browse(TBrowser* b) 
1094 {
1095   // Browse this object. 
1096   //
1097   AliDebug(30, "\tBrowsing the FMD");
1098   AliDetector::Browse(b);
1099   b->Add(AliFMDGeometry::Instance());
1100 }
1101
1102 //____________________________________________________________________  
1103 void
1104 AliFMD::AddAlignableVolumes() const
1105 {
1106   //
1107   // Create entries for alignable volumes associating the symbolic volume
1108   // name with the corresponding volume path. Needs to be syncronized with
1109   // eventual changes in the geometry.
1110   // 
1111   // This code was made by Raffaele Grosso <rgrosso@mail.cern.ch>.  I
1112   // (cholm) will probably want to change it.   For one, I think it
1113   // should be the job of the geometry manager to deal with this. 
1114   for(size_t f = 1; f <= 3; f++){ // Detector 1,2,3
1115     for(size_t tb =  0; tb <2 ; tb++){ // Top/Bottom 
1116       char     stb = tb == 0 ? 'T' : 'B';
1117       unsigned min = tb == 0 ? 0   : 5;
1118
1119       TString halfVol(Form("/ALIC_1/F%dM%c_%d", f, stb, f));
1120       TString halfSym(halfVol);
1121       if(!gGeoManager->SetAlignableEntry(halfSym.Data(),halfVol.Data()))
1122         AliFatal(Form("Alignable entry %s not created. "
1123                       "Volume path %s not valid", 
1124                       halfSym.Data(),halfVol.Data()));
1125       for(size_t io = 0; io < 2; io++){ // inner, outer 
1126         if (f==1 && io==1) continue; // Only one ring in FMD1 
1127         if(tb == 1 && io==1) min=10;
1128         char     sio = (io == 0 ? 'I' : 'O');
1129         unsigned nio = (io == 0 ? 3   : 9);
1130         unsigned max = (io == 0 ? 5   : 10) + min;
1131         
1132         for(size_t i = min; i < max; i++) { // Modules
1133           TString modVol(Form("%s/F%c%cV_7%d/F%cSE_%d", halfVol.Data(), 
1134                               sio, stb, nio, sio, i));
1135           TString modSym(modVol);
1136           if(!gGeoManager->SetAlignableEntry(modSym.Data(),modVol.Data()))
1137             AliFatal(Form("Alignable entry %s not created. "
1138                           "Volume path %s not valid", 
1139                           modSym.Data(), modVol.Data()));
1140         }
1141       }
1142     }
1143   }
1144
1145 }
1146 //___________________________________________________________________
1147 //
1148 // EOF
1149 //