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