]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMD.cxx
Bugfix.
[u/mrichter/AliRoot.git] / FMD / AliFMD.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //____________________________________________________________________
19 //                                                                          
20 // Forward Multiplicity Detector based on Silicon wafers. This class
21 // contains the base procedures for the Forward Multiplicity detector
22 // Detector consists of 5 Si volumes covered pseudorapidity interval
23 // from 1.7 to 5.1.
24 //                                                       
25 // This is the base class for all FMD manager classes. 
26 //                    
27 // The actual code is done by various separate classes.   Below is
28 // diagram showing the relationship between the various FMD classes
29 // that handles the geometry 
30 //
31 //
32 //       +----------+   +----------+   
33 //       | AliFMDv1 |   | AliFMDv1 |   
34 //       +----------+   +----------+   
35 //            |              |
36 //       +----+--------------+
37 //       |
38 //       |           +------------+ 1  +---------------+
39 //       |        +- | AliFMDRing |<>--| AliFMDPolygon | 
40 //       V     2  |  +------------+    +---------------+   
41 //  +--------+<>--+        |
42 //  | AliFMD |             ^                       
43 //  +--------+<>--+        V 1..2                     
44 //             3  | +-------------------+ 
45 //                +-| AliFMDSubDetector | 
46 //                  +-------------------+
47 //                           ^              
48 //                           |
49 //             +-------------+-------------+
50 //             |             |             |          
51 //        +---------+   +---------+   +---------+
52 //        | AliFMD1 |   | AliFMD2 |   | AliFMD3 |
53 //        +---------+   +---------+   +---------+
54 //      
55 //
56 // *  AliFMD 
57 //    This defines the interface for the various parts of AliROOT that
58 //    uses the FMD, like AliFMDDigitizer, AliFMDReconstructor, and so
59 //    on. 
60 //
61 // *  AliFMDv1 
62 //    This is a concrete implementation of the AliFMD interface. 
63 //    It is the responsibility of this class to create the FMD
64 //    geometry, process hits in the FMD, and serve hits and digits to
65 //    the various clients. 
66 //  
67 //    It uses the objects of class AliFMDSubDetector to do the various
68 //    stuff for FMD1, 2, and 3 
69 //
70 // *  AliFMDRing 
71 //    This class contains all stuff needed to do with a ring.  It's
72 //    used by the AliFMDSubDetector objects to instantise inner and
73 //    outer rings.  The AliFMDRing objects are shared by the
74 //    AliFMDSubDetector objects, and owned by the AliFMDv1 object. 
75 //
76 // *  AliFMDPolygon 
77 //    The code I lifted from TGeoPolygon to help with the geometry of
78 //    the modules, as well as to decide wether a hit is actually with
79 //    in the real module shape.  The point is, that the shape of the
80 //    various ring modules are really polygons (much like the lid of a
81 //    coffin), but it's segmented at constant radius.  That is very
82 //    hard to implement using GEANT 3.21 shapes, so instead the
83 //    modules are implemented as TUBS (tube sections), and in the step
84 //    procedure we do the test whether the track was inside the real
85 //    shape of the module.  
86 //
87 // *  AliFMD1, AliFMD2, and AliFMD3 
88 //    These are specialisation of AliFMDSubDetector, that contains the
89 //    particularities of each of the sub-detector system.  It is
90 //    envisioned that the classes should also define the support
91 //    volumes and material for each of the detectors.                          
92 //                                                                          
93 // The responsible person for this module is Alla Maevskaia
94 // <Alla.Maevskaia@cern.ch>.
95 //
96 // Many modifications by Christian Holm Christensen <cholm@nbi.dk>
97 //
98
99 // These files are not in the same directory, so there's no reason to
100 // ask the preprocessor to search in the current directory for these
101 // files by including them with `#include "..."' 
102 #include <TClonesArray.h>       // ROOT_TClonesArray
103 #include <TGeometry.h>          // ROOT_TGeomtry
104 #include <TNode.h>              // ROOT_TNode
105 #include <TTUBE.h>              // ROOT_TTUBE
106 #include <TTree.h>              // ROOT_TTree
107 #include <TVirtualMC.h>         // ROOT_TVirtualMC
108 #include <TBrowser.h>           // ROOT_TBrowser
109 #include <TMath.h>              // ROOT_TMath
110
111 #include <AliRunDigitizer.h>    // ALIRUNDIGITIZER_H
112 #include <AliLoader.h>          // ALILOADER_H
113 #include <AliRun.h>             // ALIRUN_H
114 #include <AliMC.h>              // ALIMC_H
115 #include <AliLog.h>             // ALILOG_H
116 #include <AliMagF.h>            // ALIMAGF_H
117 #include "AliFMD.h"             // ALIFMD_H
118 #include "AliFMDDigit.h"        // ALIFMDDIGIG_H
119 #include "AliFMDHit.h"          // ALIFMDHIT_H
120 #include "AliFMDDigitizer.h"    // ALIFMDDIGITIZER_H
121 #include "AliFMD1.h"            // ALIFMD1_H
122 #include "AliFMD2.h"            // ALIFMD2_H
123 #include "AliFMD3.h"            // ALIFMD3_H
124 #include "AliFMDRawWriter.h"    // ALIFMDRAWWRITER_H
125
126 //____________________________________________________________________
127 ClassImp(AliFMD);
128
129 //____________________________________________________________________
130 const Char_t* AliFMD::fgkShortLegName = "FSSL";
131 const Char_t* AliFMD::fgkLongLegName  = "FSLL";
132
133 //____________________________________________________________________
134 AliFMD::AliFMD()
135   : fInner(0), 
136     fOuter(0),
137     fFMD1(0),
138     fFMD2(0), 
139     fFMD3(0), 
140     fSDigits(0), 
141     fNsdigits(0),
142     fPrintboardRotationId(0),
143     fIdentityRotationId(0),
144     fShortLegId(0),
145     fLongLegId(0),
146     fLegLength(0),
147     fLegRadius(0),
148     fModuleSpacing(0), 
149     fSiDensity(0),
150     fSiThickness(0),
151     fSiDeDxMip(1.664),
152     fVA1MipRange(0),
153     fAltroChannelSize(0),
154     fSampleRate(0)
155 {
156   //
157   // Default constructor for class AliFMD
158   //
159   AliDebug(0, "\tDefault CTOR");
160   fHits     = 0;
161   fDigits   = 0;
162   fIshunt   = 0;
163 }
164
165 //____________________________________________________________________
166 AliFMD::AliFMD(const AliFMD& other)
167   : AliDetector(other),
168     fInner(other.fInner), 
169     fOuter(other.fOuter),
170     fFMD1(other.fFMD1),
171     fFMD2(other.fFMD2), 
172     fFMD3(other.fFMD3), 
173     fSDigits(other.fSDigits), 
174     fNsdigits(other.fNsdigits),
175     fPrintboardRotationId(other.fPrintboardRotationId),
176     fIdentityRotationId(other.fIdentityRotationId),
177     fShortLegId(other.fShortLegId),
178     fLongLegId(other.fLongLegId),
179     fLegLength(other.fLegLength),
180     fLegRadius(other.fLegRadius),
181     fModuleSpacing(other.fModuleSpacing), 
182     fSiDensity(other.fSiDensity),
183     fSiThickness(other.fSiThickness),
184     fSiDeDxMip(other.fSiDeDxMip),
185     fVA1MipRange(other.fVA1MipRange),
186     fAltroChannelSize(other.fAltroChannelSize),
187     fSampleRate(other.fSampleRate)
188 {
189   // Copy constructor 
190 }
191
192 //____________________________________________________________________
193 AliFMD::AliFMD(const char *name, const char *title, bool detailed)
194   : AliDetector (name, title),
195     fInner(0), 
196     fOuter(0),
197     fFMD1(0),
198     fFMD2(0), 
199     fFMD3(0),
200     fSDigits(0),
201     fNsdigits(0),
202     fPrintboardRotationId(0),
203     fIdentityRotationId(0),
204     fShortLegId(0),
205     fLongLegId(0),
206     fLegLength(0),
207     fLegRadius(0),
208     fModuleSpacing(0), 
209     fSiDensity(0),
210     fSiThickness(0),
211     fSiDeDxMip(1.664),
212     fVA1MipRange(0),
213     fAltroChannelSize(0),
214     fSampleRate(0)
215 {
216   //
217   // Standard constructor for Forward Multiplicity Detector
218   //
219   AliDebug(0, "\tStandard CTOR");
220
221   // Initialise Hit array
222   HitsArray();
223   gAlice->GetMCApp()->AddHitList(fHits);
224
225   // (S)Digits for the detectors disk
226   DigitsArray();
227   SDigitsArray();
228   
229   // CHC: What is this?
230   fIshunt = 0;
231   SetMarkerColor(kRed);
232   SetLineColor(kYellow);
233   SetSiDensity();
234
235   // Create sub-volume managers 
236   fInner = new AliFMDRing('I', detailed);
237   fOuter = new AliFMDRing('O', detailed);
238   fFMD1  = new AliFMD1();
239   fFMD2  = new AliFMD2();
240   fFMD3  = new AliFMD3();
241
242   // Specify parameters of sub-volume managers 
243   fFMD1->SetInner(fInner);
244   fFMD1->SetOuter(0);
245
246   fFMD2->SetInner(fInner);
247   fFMD2->SetOuter(fOuter);
248   
249   fFMD3->SetInner(fInner);
250   fFMD3->SetOuter(fOuter);
251
252   SetLegLength();
253   SetLegRadius();
254   SetLegOffset();
255   SetModuleSpacing();
256   SetSiThickness();
257   SetSiDensity();
258   SetVA1MipRange();
259   SetAltroChannelSize();
260   SetSampleRate();
261   
262   fInner->SetLowR(4.3);
263   fInner->SetHighR(17.2);
264   fInner->SetWaferRadius(13.4/2);
265   fInner->SetTheta(36/2);
266   fInner->SetNStrips(512);
267   fInner->SetSiThickness(fSiThickness);
268   fInner->SetPrintboardThickness(.11);
269   fInner->SetBondingWidth(.5);
270
271   fOuter->SetLowR(15.6);
272   fOuter->SetHighR(28.0);
273   fOuter->SetWaferRadius(13.4/2);
274   fOuter->SetTheta(18/2);
275   fOuter->SetNStrips(256);
276   fOuter->SetSiThickness(fSiThickness);
277   fOuter->SetPrintboardThickness(.1);
278   fOuter->SetBondingWidth(.5);
279   
280   
281   fFMD1->SetHoneycombThickness(1);
282   fFMD1->SetInnerZ(340.0);
283   
284   fFMD2->SetHoneycombThickness(1);
285   fFMD2->SetInnerZ(83.4);
286   fFMD2->SetOuterZ(75.2);
287
288   fFMD3->SetHoneycombThickness(1);
289   fFMD3->SetInnerZ(-62.8);
290   fFMD3->SetOuterZ(-75.2);
291 }
292
293 //____________________________________________________________________
294 AliFMD::~AliFMD ()
295 {
296   // Destructor for base class AliFMD
297   if (fHits) {
298     fHits->Delete();
299     delete fHits;
300     fHits = 0;
301   }
302   if (fDigits) {
303     fDigits->Delete();
304     delete fDigits;
305     fDigits = 0;
306   }
307   if (fSDigits) {
308     fSDigits->Delete();
309     delete fSDigits;
310     fSDigits = 0;
311   }
312 }
313
314 //____________________________________________________________________
315 AliFMD&
316 AliFMD::operator=(const AliFMD& other)
317 {
318   AliDetector::operator=(other);
319   fInner                = other.fInner; 
320   fOuter                = other.fOuter;
321   fFMD1                 = other.fFMD1;
322   fFMD2                 = other.fFMD2; 
323   fFMD3                 = other.fFMD3; 
324   fSDigits              = other.fSDigits; 
325   fNsdigits             = other.fNsdigits;
326   fSiDensity            = other.fSiDensity;
327   fPrintboardRotationId = other.fPrintboardRotationId;
328   fIdentityRotationId   = other.fIdentityRotationId;
329   fShortLegId           = other.fShortLegId;
330   fLongLegId            = other.fLongLegId;
331   fLegLength            = other.fLegLength;
332   fLegRadius            = other.fLegRadius;
333   fModuleSpacing        = other.fModuleSpacing; 
334   fSiDensity            = other.fSiDensity;
335   fSiThickness          = other.fSiThickness;
336   fVA1MipRange          = other.fVA1MipRange;
337   fAltroChannelSize     = other.fAltroChannelSize;
338   fSampleRate           = other.fSampleRate;
339
340   return *this;
341 }
342
343 //====================================================================
344 //
345 // GEometry ANd Traking
346 //
347 //____________________________________________________________________
348 void 
349 AliFMD::CreateGeometry()
350 {
351   //
352   // Create the geometry of Forward Multiplicity Detector.  The actual
353   // construction of the geometry is delegated to the class AliFMDRing
354   // and AliFMDSubDetector and the relevant derived classes. 
355   //
356   // The flow of this member function is:
357   //
358   //   FOR rings fInner and fOuter DO  
359   //     AliFMDRing::Init();
360   //   END FOR
361   // 
362   //   Set up hybrud card support (leg) volume shapes  
363   // 
364   //   FOR rings fInner and fOuter DO  
365   //     AliFMDRing::SetupGeometry();
366   //   END FOR
367   // 
368   //   FOR subdetectors fFMD1, fFMD2, and fFMD3 DO 
369   //     AliFMDSubDetector::SetupGeomtry();
370   //   END FOR
371   // 
372   //   FOR subdetectors fFMD1, fFMD2, and fFMD3 DO 
373   //     AliFMDSubDetector::Geomtry();
374   //   END FOR
375   //   
376
377   // DebugGuard guard("AliFMD::CreateGeometry");
378   AliDebug(10, "\tCreating geometry");
379
380   fInner->Init();
381   fOuter->Init();
382
383   Double_t par[3];
384
385   Int_t airId = (*fIdtmed)[kAirId];
386   Int_t alId  = (*fIdtmed)[kAlId];
387   Int_t cId   = (*fIdtmed)[kCarbonId];
388   Int_t siId  = (*fIdtmed)[kSiId];
389   Int_t pcbId = (*fIdtmed)[kPcbId];
390   Int_t plaId = (*fIdtmed)[kPlasticId];
391   Int_t pbId  = fPrintboardRotationId;
392   Int_t idId  = fIdentityRotationId;
393   
394   par[0]      =  fLegRadius - .1;
395   par[1]      =  fLegRadius;
396   par[2]      =  fLegLength / 2;
397   fShortLegId =  gMC->Gsvolu(fgkShortLegName,"TUBE", plaId, par, 3);
398   
399   par[2]      += fModuleSpacing / 2;
400   fLongLegId  =  gMC->Gsvolu(fgkLongLegName,"TUBE", plaId, par, 3);
401
402   fInner->SetupGeometry(airId, siId, pcbId, pbId, idId);
403   fOuter->SetupGeometry(airId, siId, pcbId, pbId, idId);
404
405   fFMD1->SetupGeometry(airId, alId, cId);
406   fFMD2->SetupGeometry(airId, alId, cId);
407   fFMD3->SetupGeometry(airId, alId, cId);
408   
409   fFMD1->Geometry("ALIC", pbId, idId);
410   fFMD2->Geometry("ALIC", pbId, idId);
411   fFMD3->Geometry("ALIC", pbId, idId);    
412 }    
413
414 //____________________________________________________________________
415 void AliFMD::CreateMaterials() 
416 {
417   // Register various materials and tracking mediums with the
418   // backend.   
419   // 
420   // Currently defined materials and mediums are 
421   // 
422   //    FMD Air         Normal air 
423   //    FMD Si          Active silicon of sensors 
424   //    FMD Carbon      Normal carbon used in support, etc. 
425   //    FMD Kapton      Carbon used in Honeycomb
426   //    FMD PCB         Printed circuit board material 
427   //    FMD Plastic     Material for support legs 
428   // 
429   // Also defined are two rotation matricies. 
430   //
431   // DebugGuard guard("AliFMD::CreateMaterials");
432   AliDebug(10, "\tCreating materials");
433   Int_t    id;
434   Double_t a                = 0;
435   Double_t z                = 0;
436   Double_t density          = 0;
437   Double_t radiationLength  = 0;
438   Double_t absorbtionLength = 999;
439   Int_t    fieldType        = gAlice->Field()->Integ();     // Field type 
440   Double_t maxField         = gAlice->Field()->Max();     // Field max.
441   Double_t maxBending       = 0;     // Max Angle
442   Double_t maxStepSize      = 0.001; // Max step size 
443   Double_t maxEnergyLoss    = 1;     // Max Delta E
444   Double_t precision        = 0.001; // Precision
445   Double_t minStepSize      = 0.001; // Minimum step size 
446  
447   // Silicon 
448   a                = 28.0855;
449   z                = 14.;
450   density          = fSiDensity;
451   radiationLength  = 9.36;
452   maxBending       = 1;
453   maxStepSize      = .001;
454   precision        = .001;
455   minStepSize      = .001;
456   id               = kSiId;
457   AliMaterial(id, "FMD Si$", a, z, density, radiationLength, absorbtionLength);
458   AliMedium(kSiId, "FMD Si$",id,1,fieldType,maxField,maxBending,
459             maxStepSize,maxEnergyLoss,precision,minStepSize);
460
461
462   // Carbon 
463   a                = 12.011;
464   z                = 6.;
465   density          = 2.265;
466   radiationLength  = 18.8;
467   maxBending       = 10;
468   maxStepSize      = .01;
469   precision        = .003;
470   minStepSize      = .003;
471   id               = kCarbonId;
472   AliMaterial(id, "FMD Carbon$", a, z, density, radiationLength, 
473               absorbtionLength);
474   AliMedium(kCarbonId, "FMD Carbon$",id,0,fieldType,maxField,maxBending,
475             maxStepSize,maxEnergyLoss,precision,minStepSize);
476
477   // Aluminum
478   a                = 26.981539;
479   z                = 13.;
480   density          = 2.7;
481   radiationLength  = 8.9;
482   id               = kAlId;
483   AliMaterial(id, "FMD Aluminum$", a, z, density, radiationLength, 
484               absorbtionLength);
485   AliMedium(kAlId, "FMD Aluminum$", id, 0, fieldType, maxField, maxBending,
486             maxStepSize, maxEnergyLoss, precision, minStepSize);
487   
488   
489   // Silicon chip 
490   {
491     Float_t as[] = { 12.0107,      14.0067,      15.9994,
492                      1.00794,      28.0855,     107.8682 };
493     Float_t zs[] = {  6.,           7.,           8.,
494                       1.,          14.,          47. };
495     Float_t ws[] = {  0.039730642,  0.001396798,  0.01169634,
496                       0.004367771,  0.844665,     0.09814344903 };
497     density = 2.36436;
498     maxBending       = 10;
499     maxStepSize      = .01;
500     precision        = .003;
501     minStepSize      = .003;
502     id = kSiChipId;
503     AliMixture(id, "FMD Si Chip$", as, zs, density, 6, ws);
504     AliMedium(kSiChipId, "FMD Si Chip$", id, 0, fieldType, maxField, 
505               maxBending, maxStepSize, maxEnergyLoss, precision, minStepSize);
506   }
507   
508 #if 0
509   // Kaption
510   {
511     Float_t as[] = { 1.00794,  12.0107,  14.010,   15.9994};
512     Float_t zs[] = { 1.,        6.,       7.,       8.};
513     Float_t ws[] = { 0.026362,  0.69113,  0.07327,  0.209235};
514     density          = 1.42;
515     maxBending       = 1;
516     maxStepSize      = .001;
517     precision        = .001;
518     minStepSize      = .001;
519     id               = KaptionId;
520     AliMixture(id, "FMD Kaption$", as, zs, density, 4, ws);
521     AliMedium(kAlId, "FMD Kaption$",id,0,fieldType,maxField,maxBending,
522               maxStepSize,maxEnergyLoss,precision,minStepSize);
523   }
524 #endif
525
526   // Air
527   {
528     Float_t as[] = { 12.0107, 14.0067,   15.9994,  39.948 };
529     Float_t zs[] = {  6.,      7.,       8.,       18. };
530     Float_t ws[] = { 0.000124, 0.755267, 0.231781, 0.012827 }; 
531     density      = .00120479;
532     maxBending   = 1;
533     maxStepSize  = .001;
534     precision    = .001;
535     minStepSize  = .001;
536     id           = kAirId;
537     AliMixture(id, "FMD Air$", as, zs, density, 4, ws);
538     AliMedium(kAirId, "FMD Air$", id,0,fieldType,maxField,maxBending,
539               maxStepSize,maxEnergyLoss,precision,minStepSize);
540   }
541   
542   // PCB
543   {
544     Float_t zs[] = { 14.,         20.,         13.,         12.,
545                       5.,         22.,         11.,         19.,
546                      26.,          9.,          8.,          6.,
547                       7.,          1.};
548     Float_t as[] = { 28.0855,     40.078,      26.981538,   24.305, 
549                      10.811,      47.867,      22.98977,    39.0983,
550                      55.845,      18.9984,     15.9994,     12.0107,
551                      14.0067,      1.00794};
552     Float_t ws[] = {  0.15144894,  0.08147477,  0.04128158,  0.00904554, 
553                       0.01397570,  0.00287685,  0.00445114,  0.00498089,
554                       0.00209828,  0.00420000,  0.36043788,  0.27529426,
555                       0.01415852,  0.03427566};
556     density      = 1.8;
557     maxBending   = 1;
558     maxStepSize  = .001;
559     precision    = .001;
560     minStepSize  = .001;
561     id           = kPcbId;
562     AliMixture(id, "FMD PCB$", as, zs, density, 14, ws);
563     AliMedium(kPcbId, "FMD PCB$", id,0,fieldType,maxField,maxBending,
564               maxStepSize,maxEnergyLoss,precision,minStepSize);
565   }
566   
567   // Plastic 
568   {
569     Float_t as[] = { 1.01, 12.01 };
570     Float_t zs[] = { 1.,   6.    };
571     Float_t ws[] = { 1.,   1.    };
572     density      = 1.03;
573     maxBending   = 10;
574     maxStepSize  = .01;
575     precision    = .003;
576     minStepSize  = .003;
577     id           = kPlasticId;
578     AliMixture(id, "FMD Plastic$", as, zs, density, -2, ws);
579     AliMedium(kPlasticId, "FMD Plastic$", id,0,fieldType,maxField,maxBending,
580                 maxStepSize,maxEnergyLoss,precision,minStepSize);
581   }
582   AliMatrix(fPrintboardRotationId, 90, 90, 0, 90, 90, 0);
583   AliMatrix(fIdentityRotationId, 90, 0, 90, 90, 0, 0);
584 }
585
586 //____________________________________________________________________
587 void  
588 AliFMD::Init()
589 {
590   //
591   // Initialis the FMD after it has been built
592   Int_t i;
593   //
594   if (fDebug) {
595     cout << "\n" << ClassName() << ": " << flush;
596     for (i = 0; i < 35; i++) cout << "*";
597     cout << " FMD_INIT ";
598     for (i = 0; i < 35; i++) cout << "*";
599     cout << "\n" << ClassName() << ": " << flush;
600     //
601     // Here the FMD initialisation code (if any!)
602     for (i = 0; i < 80; i++) cout << "*";
603     cout << endl;
604   }
605   //
606   //
607 }
608
609 //====================================================================
610 //
611 // Graphics and event display
612 //
613 //____________________________________________________________________
614 void 
615 AliFMD::BuildGeometry()
616 {
617   //
618   // Build simple ROOT TNode geometry for event display
619   //
620   // Build a simplified geometry of the FMD used for event display  
621   // 
622   // The actual building of the TNodes is done by
623   // AliFMDSubDetector::SimpleGeometry. 
624   AliDebug(10, "\tCreating a simplified geometry");
625
626   TNode* top = gAlice->GetGeometry()->GetNode("alice");
627   
628   fFMD1->SimpleGeometry(fNodes, top, GetLineColor(), 0);
629   fFMD2->SimpleGeometry(fNodes, top, GetLineColor(), 0);
630   fFMD3->SimpleGeometry(fNodes, top, GetLineColor(), 0);
631 }
632
633 //____________________________________________________________________
634 void 
635 AliFMD::DrawDetector()
636 {
637   //
638   // Draw a shaded view of the Forward multiplicity detector
639   //
640   // DebugGuard guard("AliFMD::DrawDetector");
641   AliDebug(10, "\tDraw detector");
642   
643   //Set ALIC mother transparent
644   gMC->Gsatt("ALIC","SEEN",0);
645   gMC->Gsatt(fgkShortLegName,"SEEN",1);
646   gMC->Gsatt(fgkLongLegName,"SEEN",1);
647
648   //Set volumes visible
649   fFMD1->Gsatt();
650   fFMD2->Gsatt();
651   fFMD3->Gsatt();
652   fInner->Gsatt();
653   fOuter->Gsatt();
654
655   //
656   gMC->Gdopt("hide", "on");
657   gMC->Gdopt("shad", "on");
658   gMC->Gsatt("*", "fill", 7);
659   gMC->SetClipBox(".");
660   gMC->SetClipBox("*", 0, 1000, -1000, 1000, -1000, 1000);
661   gMC->DefaultRange();
662   gMC->Gdraw("alic", 40, 30, 0, 12, 12, .055, .055);
663   gMC->Gdhead(1111, "Forward Multiplicity Detector");
664   gMC->Gdman(16, 10, "MAN");
665   gMC->Gdopt("hide", "off");
666 }
667
668 //____________________________________________________________________
669 Int_t 
670 AliFMD::DistanceToPrimitive(Int_t, Int_t)
671 {
672   //
673   // Calculate the distance from the mouse to the FMD on the screen
674   // Dummy routine
675   //
676   return 9999;
677 }
678
679 //====================================================================
680 //
681 // Hit and Digit managment 
682 //
683 //____________________________________________________________________
684 void 
685 AliFMD::MakeBranch(Option_t * option)
686 {
687   // Create Tree branches for the FMD.
688   //
689   // Options:
690   //
691   //    H          Make a branch of TClonesArray of AliFMDHit's
692   //    D          Make a branch of TClonesArray of AliFMDDigit's
693   //    S          Make a branch of TClonesArray of AliFMDSDigit's
694   // 
695   const Int_t kBufferSize = 16000;
696   TString branchname(GetName());
697   TString opt(option);
698   
699   if (opt.Contains("H", TString::kIgnoreCase)) {
700     HitsArray();
701     AliDetector::MakeBranch(option); 
702   }
703   if (opt.Contains("D", TString::kIgnoreCase)) { 
704     DigitsArray();
705     MakeBranchInTree(fLoader->TreeD(), branchname.Data(),
706                      &fDigits, kBufferSize, 0);
707   }
708   if (opt.Contains("S", TString::kIgnoreCase)) { 
709     SDigitsArray();
710     MakeBranchInTree(fLoader->TreeS(), branchname.Data(),
711                      &fSDigits, kBufferSize, 0);
712   }
713 }
714
715 //____________________________________________________________________
716 void 
717 AliFMD::SetTreeAddress()
718 {
719   // Set branch address for the Hits, Digits, and SDigits Tree.
720   if (fLoader->TreeH()) HitsArray();
721   AliDetector::SetTreeAddress();
722
723   TTree *treeD = fLoader->TreeD();
724   if (treeD) {
725     DigitsArray();
726     TBranch* branch = treeD->GetBranch ("FMD");
727     if (branch) branch->SetAddress(&fDigits);
728   }
729
730   TTree *treeS = fLoader->TreeS();
731   if (treeS) {
732     SDigitsArray();
733     TBranch* branch = treeS->GetBranch ("FMD");
734     if (branch) branch->SetAddress(&fSDigits);
735   }
736 }
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   AddHit(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 void 
791 AliFMD::AddHit(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 {
806   //
807   // Add a hit to the list
808   //
809   // Parameters:
810   // 
811   //    track     Track #
812   //    detector  Detector # (1, 2, or 3)                      
813   //    ring      Ring ID ('I' or 'O')
814   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
815   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
816   //    x         Track's X-coordinate at hit
817   //    y         Track's Y-coordinate at hit
818   //    z         Track's Z-coordinate at hit
819   //    px        X-component of track's momentum 
820   //    py        Y-component of track's momentum
821   //    pz        Z-component of track's momentum
822   //    edep      Energy deposited by track
823   //    pdg       Track's particle Id #
824   //    t         Time when the track hit 
825   // 
826   TClonesArray& a = *(HitsArray());
827   // Search through the list of already registered hits, and see if we
828   // find a hit with the same parameters.  If we do, then don't create
829   // a new hit, but rather update the energy deposited in the hit.
830   // This is done, so that a FLUKA based simulation will get the
831   // number of hits right, not just the enerrgy deposition. 
832   for (Int_t i = 0; i < fNhits; i++) {
833     if (!a.At(i)) continue;
834     AliFMDHit* hit = static_cast<AliFMDHit*>(a.At(i));
835     if (hit->Detector() == detector 
836         && hit->Ring() == ring
837         && hit->Sector() == sector 
838         && hit->Strip() == strip
839         && hit->Track() == track) {
840       Warning("AddHit", "already had a hit in FMD%d%c[%2d,%3d] for track # %d,"
841               " adding energy (%f) to that hit (%f) -> %f", 
842               detector, ring, sector, strip, track, edep, hit->Edep(),
843               hit->Edep() + edep);
844       hit->SetEdep(hit->Edep() + edep);
845       return;
846     }
847   }
848   // If hit wasn't already registered, do so know. 
849   new (a[fNhits]) AliFMDHit(fIshunt, track, detector, ring, sector, strip, 
850                             x, y, z, px, py, pz, edep, pdg, t);
851   fNhits++;
852 }
853
854 //____________________________________________________________________
855 void 
856 AliFMD::AddDigit(Int_t* digits)
857 {
858   // Add a digit to the Digit tree 
859   // 
860   // Paramters 
861   //
862   //    digits[0]  [UShort_t] Detector #
863   //    digits[1]  [Char_t]   Ring ID
864   //    digits[2]  [UShort_t] Sector #
865   //    digits[3]  [UShort_t] Strip #
866   //    digits[4]  [UShort_t] ADC Count 
867   //    digits[5]  [Short_t]  ADC Count, -1 if not used
868   //    digits[6]  [Short_t]  ADC Count, -1 if not used 
869   // 
870   AddDigit(UShort_t(digits[0]),  // Detector #
871            Char_t(digits[1]),    // Ring ID
872            UShort_t(digits[2]),  // Sector #
873            UShort_t(digits[3]),  // Strip #
874            UShort_t(digits[4]),  // ADC Count1 
875            Short_t(digits[5]),   // ADC Count2 
876            Short_t(digits[6]));  // ADC Count3 
877 }
878
879 //____________________________________________________________________
880 void 
881 AliFMD::AddDigit(UShort_t detector, 
882                  Char_t   ring, 
883                  UShort_t sector, 
884                  UShort_t strip, 
885                  UShort_t count1, 
886                  Short_t  count2,
887                  Short_t  count3)
888 {
889   // add a real digit - as coming from data
890   // 
891   // Parameters 
892   //
893   //    detector  Detector # (1, 2, or 3)                      
894   //    ring      Ring ID ('I' or 'O')
895   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
896   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
897   //    count1    ADC count (a 10-bit word)
898   //    count2    ADC count (a 10-bit word), or -1 if not used
899   //    count3    ADC count (a 10-bit word), or -1 if not used
900   TClonesArray& a = *(DigitsArray());
901   
902   new (a[fNdigits++]) 
903     AliFMDDigit(detector, ring, sector, strip, count1, count2, count3);
904 }
905
906 //____________________________________________________________________
907 void 
908 AliFMD::AddSDigit(Int_t* digits)
909 {
910   // Add a digit to the SDigit tree 
911   // 
912   // Paramters 
913   //
914   //    digits[0]  [UShort_t] Detector #
915   //    digits[1]  [Char_t]   Ring ID
916   //    digits[2]  [UShort_t] Sector #
917   //    digits[3]  [UShort_t] Strip #
918   //    digits[4]  [Float_t]  Total energy deposited 
919   //    digits[5]  [UShort_t] ADC Count 
920   //    digits[6]  [Short_t]  ADC Count, -1 if not used
921   //    digits[7]  [Short_t]  ADC Count, -1 if not used 
922   // 
923   AddSDigit(UShort_t(digits[0]),  // Detector #
924             Char_t(digits[1]),    // Ring ID
925             UShort_t(digits[2]),  // Sector #
926             UShort_t(digits[3]),  // Strip #
927             Float_t(digits[4]),   // Edep
928             UShort_t(digits[5]),  // ADC Count1 
929             Short_t(digits[6]),   // ADC Count2 
930             Short_t(digits[7]));  // ADC Count3 
931 }
932
933 //____________________________________________________________________
934 void 
935 AliFMD::AddSDigit(UShort_t detector, 
936                   Char_t   ring, 
937                   UShort_t sector, 
938                   UShort_t strip, 
939                   Float_t  edep,
940                   UShort_t count1, 
941                   Short_t  count2,
942                   Short_t  count3)
943 {
944   // add a summable digit
945   // 
946   // Parameters 
947   //
948   //    detector  Detector # (1, 2, or 3)                      
949   //    ring      Ring ID ('I' or 'O')
950   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
951   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
952   //    edep      Total energy deposited
953   //    count1    ADC count (a 10-bit word)
954   //    count2    ADC count (a 10-bit word), or -1 if not used
955   //    count3    ADC count (a 10-bit word), or -1 if not used
956   //
957   TClonesArray& a = *(SDigitsArray());
958   
959   new (a[fNsdigits++]) 
960     AliFMDSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
961 }
962
963 //____________________________________________________________________
964 void 
965 AliFMD::ResetSDigits()
966 {
967   //
968   // Reset number of digits and the digits array for this detector
969   //
970   fNsdigits   = 0;
971   if (fSDigits) fSDigits->Clear();
972 }
973
974
975 //____________________________________________________________________
976 TClonesArray*
977 AliFMD::HitsArray() 
978 {
979   // Initialize hit array if not already, and return pointer to it. 
980   if (!fHits) { 
981     fHits = new TClonesArray("AliFMDHit", 1000);
982     fNhits = 0;
983   }
984   return fHits;
985 }
986
987 //____________________________________________________________________
988 TClonesArray*
989 AliFMD::DigitsArray() 
990 {
991   // Initialize digit array if not already, and return pointer to it. 
992   if (!fDigits) { 
993     fDigits = new TClonesArray("AliFMDDigit", 1000);
994     fNdigits = 0;
995   }
996   return fDigits;
997 }
998
999 //____________________________________________________________________
1000 TClonesArray*
1001 AliFMD::SDigitsArray() 
1002 {
1003   // Initialize digit array if not already, and return pointer to it. 
1004   if (!fSDigits) { 
1005     fSDigits = new TClonesArray("AliFMDSDigit", 1000);
1006     fNsdigits = 0;
1007   }
1008   return fSDigits;
1009 }
1010
1011 //====================================================================
1012 //
1013 // Digitization 
1014 //
1015 //____________________________________________________________________
1016 void 
1017 AliFMD::Hits2Digits() 
1018 {
1019   // Create AliFMDDigit's from AliFMDHit's.  This is done by making a
1020   // AliFMDDigitizer, and executing that code.
1021   // 
1022   AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
1023   manager->SetInputStream(0, "galice.root");
1024   manager->SetOutputFile("H2Dfile");
1025   
1026   /* AliDigitizer* dig =*/ CreateDigitizer(manager);
1027   manager->Exec("");
1028 }
1029
1030 //____________________________________________________________________
1031 void 
1032 AliFMD::Hits2SDigits() 
1033 {
1034   // Create AliFMDSDigit's from AliFMDHit's.  This is done by creating
1035   // an AliFMDSDigitizer object, and executing it. 
1036   // 
1037   AliFMDSDigitizer* digitizer = new AliFMDSDigitizer("galice.root");
1038   digitizer->SetSampleRate(fSampleRate);
1039   digitizer->SetVA1MipRange(fVA1MipRange);
1040   digitizer->SetAltroChannelSize(fAltroChannelSize);
1041   digitizer->Exec("");
1042 }
1043
1044   
1045 //____________________________________________________________________
1046 AliDigitizer* 
1047 AliFMD::CreateDigitizer(AliRunDigitizer* manager) const
1048 {
1049   // Create a digitizer object 
1050   AliFMDDigitizer* digitizer = new AliFMDDigitizer(manager);
1051   digitizer->SetSampleRate(fSampleRate);
1052   digitizer->SetVA1MipRange(fVA1MipRange);
1053   digitizer->SetAltroChannelSize(fAltroChannelSize);
1054   return digitizer;
1055 }
1056
1057 //====================================================================
1058 //
1059 // Raw data simulation 
1060 //
1061 //__________________________________________________________________
1062 void 
1063 AliFMD::Digits2Raw() 
1064 {
1065   // Turn digits into raw data. 
1066   // 
1067   // This uses the class AliFMDRawWriter to do the job.   Please refer
1068   // to that class for more information. 
1069   AliFMDRawWriter writer(this);
1070   writer.SetSampleRate(fSampleRate);
1071   writer.Exec();
1072 }
1073
1074 //==================================================================
1075 //
1076 // Various setter functions for the common paramters 
1077 //
1078
1079 //__________________________________________________________________
1080 void 
1081 AliFMD::SetLegLength(Double_t length) 
1082 {
1083   // Set lenght of plastic legs that hold the hybrid (print board and
1084   // silicon sensor) onto the honeycomp support
1085   //
1086   AliDebug(10, Form("\tLeg length set to %lf cm", length));
1087   fLegLength = length;
1088   fInner->SetLegLength(fLegLength);
1089   fOuter->SetLegLength(fLegLength);
1090 }
1091
1092 //__________________________________________________________________
1093 void 
1094 AliFMD::SetLegOffset(Double_t offset) 
1095 {
1096   // Set offset from edge of hybrid to plastic legs that hold the
1097   // hybrid (print board and silicon sensor) onto the honeycomp
1098   // support 
1099   //
1100   AliDebug(10, Form("\tLeg offset set to %lf cm", offset));
1101   fInner->SetLegOffset(offset);
1102   fOuter->SetLegOffset(offset);
1103 }
1104
1105 //__________________________________________________________________
1106 void 
1107 AliFMD::SetLegRadius(Double_t radius) 
1108 {
1109   // Set the diameter of the plastic legs that hold the hybrid (print
1110   // board and silicon sensor) onto the honeycomp support
1111   //
1112   AliDebug(10, Form("\tLeg radius set to %lf cm", radius));
1113   fLegRadius = radius;
1114   fInner->SetLegRadius(fLegRadius);
1115   fOuter->SetLegRadius(fLegRadius);
1116 }
1117
1118 //__________________________________________________________________
1119 void 
1120 AliFMD::SetModuleSpacing(Double_t spacing) 
1121 {
1122   // Set the distance between the front and back sensor modules
1123   // (module staggering). 
1124   //
1125   AliDebug(10, Form("\tModule spacing set to %lf cm", spacing));  
1126   fModuleSpacing = spacing;
1127   fInner->SetModuleSpacing(fModuleSpacing);
1128   fOuter->SetModuleSpacing(fModuleSpacing);
1129 }
1130
1131 //====================================================================
1132 //
1133 // Utility 
1134 //
1135 //__________________________________________________________________
1136 void 
1137 AliFMD::Browse(TBrowser* b) 
1138 {
1139   // Browse this object. 
1140   //
1141   AliDebug(30, "\tBrowsing the FMD");
1142   AliDetector::Browse(b);
1143   if (fInner) b->Add(fInner, "Inner Ring");
1144   if (fOuter) b->Add(fOuter, "Outer Ring");
1145   if (fFMD1)  b->Add(fFMD1,  "FMD1 SubDetector");
1146   if (fFMD2)  b->Add(fFMD2,  "FMD2 SubDetector");
1147   if (fFMD3)  b->Add(fFMD3,  "FMD3 SubDetector");
1148 }
1149
1150
1151 //___________________________________________________________________
1152 //
1153 // EOF
1154 //