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