]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMD.cxx
Fixed some coding convention violations
[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 #ifndef ROOT_TClonesArray
100 #include <TClonesArray.h>
101 #endif
102 #ifndef ROOT_TGeomtry
103 # include <TGeometry.h>
104 #endif
105 #ifndef ROOT_TNode
106 # include <TNode.h>
107 #endif
108 #ifndef ROOT_TTUBE
109 # include <TTUBE.h>
110 #endif
111 #ifndef ROOT_TTree
112 # include <TTree.h>
113 #endif
114 #ifndef ROOT_TVirtualMC
115 # include <TVirtualMC.h>
116 #endif
117 #ifndef ROOT_TBrowser
118 # include <TBrowser.h>
119 #endif
120 #ifndef ROOT_TMath
121 # include <TMath.h>
122 #endif
123
124 #ifndef ALIRUNDIGITIZER_H
125 # include "AliRunDigitizer.h"
126 #endif
127 #ifndef ALILOADER_H
128 # include "AliLoader.h"
129 #endif
130 #ifndef ALIRUN_H
131 # include "AliRun.h"
132 #endif
133 #ifndef ALIMC_H
134 # include "AliMC.h"
135 #endif
136 #ifndef ALILOG_H
137 # include "AliLog.h"
138 #endif
139 #ifndef ALIMAGF_H
140 # include "AliMagF.h"
141 #endif
142 #ifndef ALIFMD_H
143 # include "AliFMD.h"
144 #endif
145 #ifndef ALIFMDDIGIG_H
146 # include "AliFMDDigit.h"
147 #endif
148 #ifndef ALIFMDHIT_H
149 # include "AliFMDHit.h"
150 #endif
151 #ifndef ALIFMDDIGITIZER_H
152 # include "AliFMDDigitizer.h"
153 #endif
154 #ifndef ALIFMD1_H
155 # include "AliFMD1.h"
156 #endif
157 #ifndef ALIFMD2_H
158 # include "AliFMD2.h"
159 #endif
160 #ifndef ALIFMD3_H
161 # include "AliFMD3.h"
162 #endif
163 #ifndef ALIALTROBUFFER_H
164 # include "AliAltroBuffer.h"
165 #endif
166
167 //____________________________________________________________________
168 ClassImp(AliFMD);
169
170 //____________________________________________________________________
171 AliFMD::AliFMD()
172   : fInner(0), 
173     fOuter(0),
174     fFMD1(0),
175     fFMD2(0), 
176     fFMD3(0), 
177     fSDigits(0), 
178     fNsdigits(0),
179     fSiDensity(0),
180     fPrintboardRotationId(0),
181     fIdentityRotationId(0),
182     fShortLegId(0),
183     fLongLegId(0),
184     fLegLength(0),
185     fLegRadius(0),
186     fModuleSpacing(0)
187 {
188   //
189   // Default constructor for class AliFMD
190   //
191   AliDebug(0, "Default CTOR");
192   fHits     = 0;
193   fDigits   = 0;
194   fIshunt   = 0;
195 }
196
197 //____________________________________________________________________
198 AliFMD::AliFMD(const char *name, const char *title, bool detailed)
199   : AliDetector (name, title),
200     fInner(0), 
201     fOuter(0),
202     fFMD1(0),
203     fFMD2(0), 
204     fFMD3(0),
205     fSDigits(0),
206     fNsdigits(0),
207     fSiDensity(0),
208     fPrintboardRotationId(0),
209     fIdentityRotationId(0),
210     fShortLegId(0),
211     fLongLegId(0),
212     fLegLength(0),
213     fLegRadius(0),
214     fModuleSpacing(0)
215 {
216   //
217   // Standard constructor for Forward Multiplicity Detector
218   //
219   AliDebug(0, "Standard 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   
257   fInner->SetLowR(4.3);
258   fInner->SetHighR(17.2);
259   fInner->SetWaferRadius(13.4/2);
260   fInner->SetTheta(36/2);
261   fInner->SetNStrips(512);
262   fInner->SetSiThickness(.03);
263   fInner->SetPrintboardThickness(.11);
264   fInner->SetBondingWidth(.5);
265
266   fOuter->SetLowR(15.6);
267   fOuter->SetHighR(28.0);
268   fOuter->SetWaferRadius(13.4/2);
269   fOuter->SetTheta(18/2);
270   fOuter->SetNStrips( 256);
271   fOuter->SetSiThickness(.03);
272   fOuter->SetPrintboardThickness(.1);
273   fOuter->SetBondingWidth(.5);
274   
275   
276   fFMD1->SetHoneycombThickness(1);
277   fFMD1->SetInnerZ(340.0);
278   
279   fFMD2->SetHoneycombThickness(1);
280   fFMD2->SetInnerZ(83.4);
281   fFMD2->SetOuterZ(75.2);
282
283   fFMD3->SetHoneycombThickness(1);
284   fFMD3->SetInnerZ(-62.8);
285   fFMD3->SetOuterZ(-75.2);
286 }
287
288 //____________________________________________________________________
289 AliFMD::~AliFMD ()
290 {
291   // Destructor for base class AliFMD
292   if (fHits) {
293     fHits->Delete();
294     delete fHits;
295     fHits = 0;
296   }
297   if (fDigits) {
298     fDigits->Delete();
299     delete fDigits;
300     fDigits = 0;
301   }
302   if (fSDigits) {
303     fSDigits->Delete();
304     delete fSDigits;
305     fSDigits = 0;
306   }
307 }
308
309 //====================================================================
310 //
311 // GEometry ANd Traking
312 //
313 //____________________________________________________________________
314 void 
315 AliFMD::CreateGeometry()
316 {
317   //
318   // Create the geometry of Forward Multiplicity Detector.  The actual
319   // construction of the geometry is delegated to the class AliFMDRing
320   // and AliFMDSubDetector and the relevant derived classes. 
321   //
322   // The flow of this member function is:
323   //
324   //   FOR rings fInner and fOuter DO  
325   //     AliFMDRing::Init();
326   //   END FOR
327   // 
328   //   Set up hybrud card support (leg) volume shapes  
329   // 
330   //   FOR rings fInner and fOuter DO  
331   //     AliFMDRing::SetupGeometry();
332   //   END FOR
333   // 
334   //   FOR subdetectors fFMD1, fFMD2, and fFMD3 DO 
335   //     AliFMDSubDetector::SetupGeomtry();
336   //   END FOR
337   // 
338   //   FOR subdetectors fFMD1, fFMD2, and fFMD3 DO 
339   //     AliFMDSubDetector::Geomtry();
340   //   END FOR
341   //   
342
343   // DebugGuard guard("AliFMD::CreateGeometry");
344   AliDebug(10, "Creating geometry");
345
346   fInner->Init();
347   fOuter->Init();
348
349   TString name;
350   Double_t par[3];
351
352   par[0]      =  fLegRadius - .1;
353   par[1]      =  fLegRadius;
354   par[2]      =  fLegLength / 2;
355   name        =  "FSL";
356   fShortLegId =  gMC->Gsvolu(name.Data(),"TUBE",(*fIdtmed)[kPlasticId],par,3);
357   
358   par[2]      += fModuleSpacing / 2;
359   name        = "FLL";
360   fLongLegId  =  gMC->Gsvolu(name.Data(),"TUBE",(*fIdtmed)[kPlasticId],par,3);
361
362   fInner->SetupGeometry((*fIdtmed)[kAirId], 
363                         (*fIdtmed)[kSiId], 
364                         (*fIdtmed)[kPcbId], 
365                         fPrintboardRotationId, 
366                         fIdentityRotationId);
367   fOuter->SetupGeometry((*fIdtmed)[kAirId], 
368                         (*fIdtmed)[kSiId], 
369                         (*fIdtmed)[kPcbId], 
370                         fPrintboardRotationId, 
371                         fIdentityRotationId);
372
373   fFMD1->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
374   fFMD2->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
375   fFMD3->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
376   
377   fFMD1->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);
378   fFMD2->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);
379   fFMD3->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);    
380 }    
381
382 //____________________________________________________________________
383 void AliFMD::CreateMaterials() 
384 {
385   // Register various materials and tracking mediums with the
386   // backend.   
387   // 
388   // Currently defined materials and mediums are 
389   // 
390   //    FMD Air         Normal air 
391   //    FMD Si          Active silicon of sensors 
392   //    FMD Carbon      Normal carbon used in support, etc. 
393   //    FMD Kapton      Carbon used in Honeycomb
394   //    FMD PCB         Printed circuit board material 
395   //    FMD Plastic     Material for support legs 
396   // 
397   // Also defined are two rotation matricies. 
398   //
399   // DebugGuard guard("AliFMD::CreateMaterials");
400   AliDebug(10, "Creating materials");
401   Int_t    id;
402   Double_t a                = 0;
403   Double_t z                = 0;
404   Double_t density          = 0;
405   Double_t radiationLength  = 0;
406   Double_t absorbtionLength = 999;
407   Int_t    fieldType        = gAlice->Field()->Integ();     // Field type 
408   Double_t maxField         = gAlice->Field()->Max();     // Field max.
409   Double_t maxBending       = 0;     // Max Angle
410   Double_t maxStepSize      = 0.001; // Max step size 
411   Double_t maxEnergyLoss    = 1;     // Max Delta E
412   Double_t precision        = 0.001; // Precision
413   Double_t minStepSize      = 0.001; // Minimum step size 
414  
415   // Silicon 
416   a                = 28.0855;
417   z                = 14.;
418   density          = fSiDensity;
419   radiationLength  = 9.36;
420   maxBending       = 1;
421   maxStepSize      = .001;
422   precision        = .001;
423   minStepSize      = .001;
424   id               = kSiId;
425   AliMaterial(id, "FMD Si$", a, z, density, radiationLength, absorbtionLength);
426   AliMedium(kSiId, "FMD Si$",id,1,fieldType,maxField,maxBending,
427             maxStepSize,maxEnergyLoss,precision,minStepSize);
428
429
430   // Carbon 
431   a                = 12.011;
432   z                = 6.;
433   density          = 2.265;
434   radiationLength  = 18.8;
435   maxBending       = 10;
436   maxStepSize      = .01;
437   precision        = .003;
438   minStepSize      = .003;
439   id               = kCarbonId;
440   AliMaterial(id, "FMD Carbon$", a, z, density, radiationLength, 
441               absorbtionLength);
442   AliMedium(kCarbonId, "FMD Carbon$",id,0,fieldType,maxField,maxBending,
443             maxStepSize,maxEnergyLoss,precision,minStepSize);
444
445   // Silicon chip 
446   {
447     Float_t as[] = { 12.0107,      14.0067,      15.9994,
448                      1.00794,      28.0855,     107.8682 };
449     Float_t zs[] = {  6.,           7.,           8.,
450                       1.,          14.,          47. };
451     Float_t ws[] = {  0.039730642,  0.001396798,  0.01169634,
452                       0.004367771,  0.844665,     0.09814344903 };
453     density = 2.36436;
454     maxBending       = 10;
455     maxStepSize      = .01;
456     precision        = .003;
457     minStepSize      = .003;
458     id = kSiChipId;
459     AliMixture(id, "FMD Si Chip$", as, zs, density, 6, ws);
460     AliMedium(kSiChipId, "FMD Si Chip$", id, 0, fieldType, maxField, 
461               maxBending, maxStepSize, maxEnergyLoss, precision, minStepSize);
462   }
463   
464
465   // Kaption 
466   {
467     Float_t as[] = { 1.00794,  12.0107,  14.010,   15.9994};
468     Float_t zs[] = { 1.,        6.,       7.,       8.};
469     Float_t ws[] = { 0.026362,  0.69113,  0.07327,  0.209235};
470     density          = 1.42;
471     maxBending       = 1;
472     maxStepSize      = .001;
473     precision        = .001;
474     minStepSize      = .001;
475     id               = kKaptionId;
476     AliMixture(id, "FMD Kaption$", as, zs, density, 4, ws);
477     AliMedium(kKaptionId, "FMD Kaption$",id,0,fieldType,maxField,maxBending,
478               maxStepSize,maxEnergyLoss,precision,minStepSize);
479   }
480   
481   // Air
482   {
483     Float_t as[] = { 12.0107, 14.0067,   15.9994,  39.948 };
484     Float_t zs[] = {  6.,      7.,       8.,       18. };
485     Float_t ws[] = { 0.000124, 0.755267, 0.231781, 0.012827 }; 
486     density      = .00120479;
487     maxBending   = 1;
488     maxStepSize  = .001;
489     precision    = .001;
490     minStepSize  = .001;
491     id           = kAirId;
492     AliMixture(id, "FMD Air$", as, zs, density, 4, ws);
493     AliMedium(kAirId, "FMD Air$", id,0,fieldType,maxField,maxBending,
494               maxStepSize,maxEnergyLoss,precision,minStepSize);
495   }
496   
497   // PCB
498   {
499     Float_t zs[] = { 14.,         20.,         13.,         12.,
500                       5.,         22.,         11.,         19.,
501                      26.,          9.,          8.,          6.,
502                       7.,          1.};
503     Float_t as[] = { 28.0855,     40.078,      26.981538,   24.305, 
504                      10.811,      47.867,      22.98977,    39.0983,
505                      55.845,      18.9984,     15.9994,     12.0107,
506                      14.0067,      1.00794};
507     Float_t ws[] = {  0.15144894,  0.08147477,  0.04128158,  0.00904554, 
508                       0.01397570,  0.00287685,  0.00445114,  0.00498089,
509                       0.00209828,  0.00420000,  0.36043788,  0.27529426,
510                       0.01415852,  0.03427566};
511     density      = 1.8;
512     maxBending   = 1;
513     maxStepSize  = .001;
514     precision    = .001;
515     minStepSize  = .001;
516     id           = kPcbId;
517     AliMixture(id, "FMD PCB$", as, zs, density, 14, ws);
518     AliMedium(kPcbId, "FMD PCB$", id,1,fieldType,maxField,maxBending,
519               maxStepSize,maxEnergyLoss,precision,minStepSize);
520   }
521   
522   // Plastic 
523   {
524     Float_t as[] = { 1.01, 12.01 };
525     Float_t zs[] = { 1.,   6.    };
526     Float_t ws[] = { 1.,   1.    };
527     density      = 1.03;
528     maxBending   = 10;
529     maxStepSize  = .01;
530     precision    = .003;
531     minStepSize  = .003;
532     id           = kPlasticId;
533     AliMixture(id, "FMD Plastic$", as, zs, density, -2, ws);
534     AliMedium(kPlasticId, "FMD Plastic$", id,0,fieldType,maxField,maxBending,
535                 maxStepSize,maxEnergyLoss,precision,minStepSize);
536   }
537   AliMatrix(fPrintboardRotationId, 90, 90, 0, 90, 90, 0);
538   AliMatrix(fIdentityRotationId, 90, 0, 90, 90, 0, 0);
539 }
540
541 //____________________________________________________________________
542 void  
543 AliFMD::Init()
544 {
545   //
546   // Initialis the FMD after it has been built
547   Int_t i;
548   //
549   if (fDebug) {
550     std::cout << "\n" << ClassName() << ": " << std::flush;
551     for (i = 0; i < 35; i++) std::cout << "*";
552     std::cout << " FMD_INIT ";
553     for (i = 0; i < 35; i++) std::cout << "*";
554     std::cout << "\n" << ClassName() << ": " << std::flush;
555     //
556     // Here the FMD initialisation code (if any!)
557     for (i = 0; i < 80; i++) std::cout << "*";
558     std::cout << std::endl;
559   }
560   //
561   //
562 }
563
564 //====================================================================
565 //
566 // Graphics and event display
567 //
568 //____________________________________________________________________
569 void 
570 AliFMD::BuildGeometry()
571 {
572   //
573   // Build simple ROOT TNode geometry for event display
574   //
575   // Build a simplified geometry of the FMD used for event display  
576   // 
577   // The actual building of the TNodes is done by
578   // AliFMDSubDetector::SimpleGeometry. 
579   AliDebug(10, "Creating a simplified geometry");
580
581   TNode* top = gAlice->GetGeometry()->GetNode("alice");
582   
583   fFMD1->SimpleGeometry(fNodes, top, GetLineColor(), 0);
584   fFMD2->SimpleGeometry(fNodes, top, GetLineColor(), 0);
585   fFMD3->SimpleGeometry(fNodes, top, GetLineColor(), 0);
586 }
587
588 //____________________________________________________________________
589 void 
590 AliFMD::DrawDetector()
591 {
592   //
593   // Draw a shaded view of the Forward multiplicity detector
594   //
595   // DebugGuard guard("AliFMD::DrawDetector");
596   AliDebug(10, "Draw detector");
597   
598   //Set ALIC mother transparent
599   gMC->Gsatt("ALIC","SEEN",0);
600
601   //Set volumes visible
602   fFMD1->Gsatt();
603   fFMD2->Gsatt();
604   fFMD3->Gsatt();
605   fInner->Gsatt();
606   fOuter->Gsatt();
607
608   //
609   gMC->Gdopt("hide", "on");
610   gMC->Gdopt("shad", "on");
611   gMC->Gsatt("*", "fill", 7);
612   gMC->SetClipBox(".");
613   gMC->SetClipBox("*", 0, 1000, -1000, 1000, -1000, 1000);
614   gMC->DefaultRange();
615   gMC->Gdraw("alic", 40, 30, 0, 12, 12, .055, .055);
616   gMC->Gdhead(1111, "Forward Multiplicity Detector");
617   gMC->Gdman(16, 10, "MAN");
618   gMC->Gdopt("hide", "off");
619 }
620
621 //____________________________________________________________________
622 const Int_t 
623 AliFMD::DistanceToPrimitive(Int_t, Int_t)
624 {
625   //
626   // Calculate the distance from the mouse to the FMD on the screen
627   // Dummy routine
628   //
629   return 9999;
630 }
631
632 //====================================================================
633 //
634 // Hit and Digit managment 
635 //
636 //____________________________________________________________________
637 void 
638 AliFMD::MakeBranch(Option_t * option)
639 {
640   // Create Tree branches for the FMD.
641   //
642   // Options:
643   //
644   //    H          Make a branch of TClonesArray of AliFMDHit's
645   //    D          Make a branch of TClonesArray of AliFMDDigit's
646   //    S          Make a branch of TClonesArray of AliFMDSDigit's
647   // 
648   const Int_t kBufferSize = 16000;
649   TString branchname(GetName());
650   TString opt(option);
651   
652   if (opt.Contains("H", TString::kIgnoreCase)) {
653     HitsArray();
654     AliDetector::MakeBranch(option); 
655   }
656   if (opt.Contains("D", TString::kIgnoreCase)) { 
657     DigitsArray();
658     MakeBranchInTree(fLoader->TreeD(), branchname.Data(),
659                      &fDigits, kBufferSize, 0);
660   }
661   if (opt.Contains("S", TString::kIgnoreCase)) { 
662     SDigitsArray();
663     MakeBranchInTree(fLoader->TreeS(), branchname.Data(),
664                      &fSDigits, kBufferSize, 0);
665   }
666 }
667
668 //____________________________________________________________________
669 void 
670 AliFMD::SetTreeAddress()
671 {
672   // Set branch address for the Hits, Digits, and SDigits Tree.
673   if (fLoader->TreeH()) HitsArray();
674   AliDetector::SetTreeAddress();
675
676   TTree *treeD = fLoader->TreeD();
677   if (treeD) {
678     DigitsArray();
679     TBranch* branch = treeD->GetBranch ("FMD");
680     if (branch) branch->SetAddress(&fDigits);
681   }
682
683   TTree *treeS = fLoader->TreeS();
684   if (treeS) {
685     SDigitsArray();
686     TBranch* branch = treeS->GetBranch ("FMD");
687     if (branch) branch->SetAddress(&fSDigits);
688   }
689 }
690
691
692
693 //____________________________________________________________________
694 void 
695 AliFMD::SetHitsAddressBranch(TBranch *b)
696 {
697   // Set the TClonesArray to read hits into. 
698   b->SetAddress(&fHits);
699 }
700
701 //____________________________________________________________________
702 void 
703 AliFMD::AddHit(Int_t track, Int_t *vol, Float_t *hits) 
704 {
705   // Add a hit to the hits tree 
706   // 
707   // The information of the two arrays are decoded as 
708   // 
709   // Parameters
710   //    track                Track #
711   //    ivol[0]  [UShort_t ] Detector # 
712   //    ivol[1]  [Char_t   ] Ring ID 
713   //    ivol[2]  [UShort_t ] Sector #
714   //    ivol[3]  [UShort_t ] Strip # 
715   //    hits[0]  [Float_t  ] Track's X-coordinate at hit 
716   //    hits[1]  [Float_t  ] Track's Y-coordinate at hit
717   //    hits[3]  [Float_t  ] Track's Z-coordinate at hit
718   //    hits[4]  [Float_t  ] X-component of track's momentum             
719   //    hits[5]  [Float_t  ] Y-component of track's momentum             
720   //    hits[6]  [Float_t  ] Z-component of track's momentum            
721   //    hits[7]  [Float_t  ] Energy deposited by track                  
722   //    hits[8]  [Int_t    ] Track's particle Id # 
723   //    hits[9]  [Float_t  ] Time when the track hit
724   // 
725   // 
726   AddHit(track, 
727          UShort_t(vol[0]),  // Detector # 
728          Char_t(vol[1]),    // Ring ID
729          UShort_t(vol[2]),  // Sector # 
730          UShort_t(vol[3]),  // Strip # 
731          hits[0],           // X
732          hits[1],           // Y
733          hits[2],           // Z
734          hits[3],           // Px
735          hits[4],           // Py
736          hits[5],           // Pz
737          hits[6],           // Energy loss 
738          Int_t(hits[7]),    // PDG 
739          hits[8]);          // Time
740 }
741
742 //____________________________________________________________________
743 void 
744 AliFMD::AddHit(Int_t    track, 
745                UShort_t detector, 
746                Char_t   ring, 
747                UShort_t sector, 
748                UShort_t strip, 
749                Float_t  x, 
750                Float_t  y, 
751                Float_t  z,
752                Float_t  px, 
753                Float_t  py, 
754                Float_t  pz,
755                Float_t  edep,
756                Int_t    pdg,
757                Float_t  t)
758 {
759   //
760   // Add a hit to the list
761   //
762   // Parameters:
763   // 
764   //    track     Track #
765   //    detector  Detector # (1, 2, or 3)                      
766   //    ring      Ring ID ('I' or 'O')
767   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
768   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
769   //    x         Track's X-coordinate at hit
770   //    y         Track's Y-coordinate at hit
771   //    z         Track's Z-coordinate at hit
772   //    px        X-component of track's momentum 
773   //    py        Y-component of track's momentum
774   //    pz        Z-component of track's momentum
775   //    edep      Energy deposited by track
776   //    pdg       Track's particle Id #
777   //    t         Time when the track hit 
778   // 
779   TClonesArray& a = *(HitsArray());
780   // Search through the list of already registered hits, and see if we
781   // find a hit with the same parameters.  If we do, then don't create
782   // a new hit, but rather update the energy deposited in the hit.
783   // This is done, so that a FLUKA based simulation will get the
784   // number of hits right, not just the enerrgy deposition. 
785   for (Int_t i = 0; i < fNhits; i++) {
786     if (!a.At(i)) continue;
787     AliFMDHit* hit = static_cast<AliFMDHit*>(a.At(i));
788     if (hit->Detector() == detector 
789         && hit->Ring() == ring
790         && hit->Sector() == sector 
791         && hit->Strip() == strip
792         && hit->Track() == track) {
793       Warning("AddHit", "already had a hit in FMD%d%c[%2d,%3d] for track # %d,"
794               " adding energy (%f) to that hit (%f) -> %f", 
795               detector, ring, sector, strip, track, edep, hit->Edep(),
796               hit->Edep() + edep);
797       hit->SetEdep(hit->Edep() + edep);
798       return;
799     }
800   }
801   // If hit wasn't already registered, do so know. 
802   new (a[fNhits]) AliFMDHit(fIshunt, track, detector, ring, sector, strip, 
803                             x, y, z, px, py, pz, edep, pdg, t);
804   fNhits++;
805 }
806
807 //____________________________________________________________________
808 void 
809 AliFMD::AddDigit(Int_t* digits)
810 {
811   // Add a digit to the Digit tree 
812   // 
813   // Paramters 
814   //
815   //    digits[0]  [UShort_t] Detector #
816   //    digits[1]  [Char_t]   Ring ID
817   //    digits[2]  [UShort_t] Sector #
818   //    digits[3]  [UShort_t] Strip #
819   //    digits[4]  [UShort_t] ADC Count 
820   //    digits[5]  [Short_t]  ADC Count, -1 if not used
821   //    digits[6]  [Short_t]  ADC Count, -1 if not used 
822   // 
823   AddDigit(UShort_t(digits[0]),  // Detector #
824            Char_t(digits[1]),    // Ring ID
825            UShort_t(digits[2]),  // Sector #
826            UShort_t(digits[3]),  // Strip #
827            UShort_t(digits[4]),  // ADC Count1 
828            Short_t(digits[5]),   // ADC Count2 
829            Short_t(digits[6]));  // ADC Count3 
830 }
831
832 //____________________________________________________________________
833 void 
834 AliFMD::AddDigit(UShort_t detector, 
835                  Char_t   ring, 
836                  UShort_t sector, 
837                  UShort_t strip, 
838                  UShort_t count1, 
839                  Short_t  count2,
840                  Short_t  count3)
841 {
842   // add a real digit - as coming from data
843   // 
844   // Parameters 
845   //
846   //    detector  Detector # (1, 2, or 3)                      
847   //    ring      Ring ID ('I' or 'O')
848   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
849   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
850   //    count1    ADC count (a 10-bit word)
851   //    count2    ADC count (a 10-bit word), or -1 if not used
852   //    count3    ADC count (a 10-bit word), or -1 if not used
853   TClonesArray& a = *(DigitsArray());
854   
855   new (a[fNdigits++]) 
856     AliFMDDigit(detector, ring, sector, strip, count1, count2, count3);
857 }
858
859 //____________________________________________________________________
860 void 
861 AliFMD::AddSDigit(Int_t* digits)
862 {
863   // Add a digit to the SDigit tree 
864   // 
865   // Paramters 
866   //
867   //    digits[0]  [UShort_t] Detector #
868   //    digits[1]  [Char_t]   Ring ID
869   //    digits[2]  [UShort_t] Sector #
870   //    digits[3]  [UShort_t] Strip #
871   //    digits[4]  [Float_t]  Total energy deposited 
872   //    digits[5]  [UShort_t] ADC Count 
873   //    digits[6]  [Short_t]  ADC Count, -1 if not used
874   //    digits[7]  [Short_t]  ADC Count, -1 if not used 
875   // 
876   AddSDigit(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             Float_t(digits[4]),   // Edep
881             UShort_t(digits[5]),  // ADC Count1 
882             Short_t(digits[6]),   // ADC Count2 
883             Short_t(digits[7]));  // ADC Count3 
884 }
885
886 //____________________________________________________________________
887 void 
888 AliFMD::AddSDigit(UShort_t detector, 
889                   Char_t   ring, 
890                   UShort_t sector, 
891                   UShort_t strip, 
892                   Float_t  edep,
893                   UShort_t count1, 
894                   Short_t  count2,
895                   Short_t  count3)
896 {
897   // add a summable digit
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   //    edep      Total energy deposited
906   //    count1    ADC count (a 10-bit word)
907   //    count2    ADC count (a 10-bit word), or -1 if not used
908   //    count3    ADC count (a 10-bit word), or -1 if not used
909   //
910   TClonesArray& a = *(SDigitsArray());
911   
912   new (a[fNsdigits++]) 
913     AliFMDSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
914 }
915
916 //____________________________________________________________________
917 void 
918 AliFMD::ResetSDigits()
919 {
920   //
921   // Reset number of digits and the digits array for this detector
922   //
923   fNsdigits   = 0;
924   if (fSDigits) fSDigits->Clear();
925 }
926
927
928 //____________________________________________________________________
929 TClonesArray*
930 AliFMD::HitsArray() 
931 {
932   // Initialize hit array if not already, and return pointer to it. 
933   if (!fHits) { 
934     fHits = new TClonesArray("AliFMDHit", 1000);
935     fNhits = 0;
936   }
937   return fHits;
938 }
939
940 //____________________________________________________________________
941 TClonesArray*
942 AliFMD::DigitsArray() 
943 {
944   // Initialize digit array if not already, and return pointer to it. 
945   if (!fDigits) { 
946     fDigits = new TClonesArray("AliFMDDigit", 1000);
947     fNdigits = 0;
948   }
949   return fDigits;
950 }
951
952 //____________________________________________________________________
953 TClonesArray*
954 AliFMD::SDigitsArray() 
955 {
956   // Initialize digit array if not already, and return pointer to it. 
957   if (!fSDigits) { 
958     fSDigits = new TClonesArray("AliFMDSDigit", 1000);
959     fNsdigits = 0;
960   }
961   return fSDigits;
962 }
963
964 //====================================================================
965 //
966 // Digitization 
967 //
968 //____________________________________________________________________
969 void 
970 AliFMD::Hits2Digits() 
971 {
972   // Create AliFMDDigit's from AliFMDHit's.  This is done by making a
973   // AliFMDDigitizer, and executing that code.
974   // 
975   AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
976   manager->SetInputStream(0, "galice.root");
977   manager->SetOutputFile("H2Dfile");
978   
979   /* AliDigitizer* dig =*/ CreateDigitizer(manager);
980   manager->Exec("");
981 }
982
983 //____________________________________________________________________
984 void 
985 AliFMD::Hits2SDigits() 
986 {
987   // Create AliFMDSDigit's from AliFMDHit's.  This is done by creating
988   // an AliFMDSDigitizer object, and executing it. 
989   // 
990   AliDigitizer* sdig = new AliFMDSDigitizer("galice.root");
991   sdig->Exec("");
992 }
993
994   
995 //____________________________________________________________________
996 AliDigitizer* 
997 AliFMD::CreateDigitizer(AliRunDigitizer* manager) const
998 {
999   // Create a digitizer object 
1000   return new AliFMDDigitizer(manager);
1001 }
1002
1003 //====================================================================
1004 //
1005 // Raw data simulation 
1006 //
1007 //__________________________________________________________________
1008 void 
1009 AliFMD::Digits2Raw() 
1010 {
1011   // Turn digits into raw data. 
1012   // 
1013   // Digits are read from the Digit branch, and processed to make
1014   // three DDL files, one for each of the sub-detectors FMD1, FMD2,
1015   // and FMD3. 
1016   //
1017   // The raw data files consists of a header, followed by ALTRO
1018   // formatted blocks.  
1019   // 
1020   //          +-------------+
1021   //          | Header      |
1022   //          +-------------+
1023   //          | ALTRO Block |
1024   //          | ...         |
1025   //          +-------------+
1026   //          DDL file 
1027   // 
1028   // An ALTRO formatted block, in the FMD context, consists of a
1029   // number of counts followed by a trailer. 
1030   // 
1031   //          +------------------+
1032   //          | Count            |
1033   //          | ...              |
1034   //          | possible fillers |
1035   //          +------------------+
1036   //          | Trailer          |
1037   //          +------------------+
1038   //          ALTRO block 
1039   // 
1040   // The counts are listed backwards, that is, starting with the
1041   // latest count, and ending in the first. 
1042   // 
1043   // Each count consist of 1 or more ADC samples of the VA1_ALICE
1044   // pre-amp. signal.  Just how many samples are used depends on
1045   // whether the ALTRO over samples the pre-amp.  Each sample is a
1046   // 10-bit word, and the samples are grouped into 40-bit blocks 
1047   //
1048   //          +------------------------------------+
1049   //          |  S(n)   | S(n-1) | S(n-2) | S(n-3) |
1050   //          |  ...    | ...    | ...    | ...    |
1051   //          |  S(2)   | S(1)   | AA     | AA     |
1052   //          +------------------------------------+
1053   //          Counts + possible filler 
1054   //
1055   // The trailer of the number of words of signales, the starting
1056   // strip number, the sector number, and the ring ID; each 10-bit
1057   // words,  packed into 40-bits. 
1058   // 
1059   //          +------------------------------------+
1060   //          | # words | start  | sector | ring   |
1061   //          +------------------------------------+
1062   //          Trailer
1063   // 
1064   // Note, that this method assumes that the digits are ordered. 
1065   //
1066   AliFMD* fmd = static_cast<AliFMD*>(gAlice->GetDetector(GetName()));
1067   fLoader->LoadDigits();
1068   TTree* digitTree = fLoader->TreeD();
1069   if (!digitTree) {
1070     Error("Digits2Raw", "no digit tree");
1071     return;
1072   }
1073   
1074   TClonesArray* digits = new TClonesArray("AliFMDDigit", 1000);
1075   fmd->SetTreeAddress();
1076   TBranch* digitBranch = digitTree->GetBranch(GetName());
1077   if (!digitBranch) {
1078     Error("Digits2Raw", "no branch for %s", GetName());
1079     return;
1080   }
1081   digitBranch->SetAddress(&digits);
1082   
1083   Int_t nEvents = Int_t(digitTree->GetEntries());
1084   for (Int_t event = 0; event < nEvents; event++) {
1085     fmd->ResetDigits();
1086     digitTree->GetEvent(event);
1087     
1088     Int_t nDigits = digits->GetEntries();
1089     if (nDigits < 1) continue;
1090
1091
1092     UShort_t prevDetector = 0;
1093     Char_t   prevRing     = '\0';
1094     UShort_t prevSector   = 0;
1095     // UShort_t prevStrip    = 0;
1096
1097     // The first seen strip number for a channel 
1098     UShort_t startStrip   = 0;
1099     
1100     // Which channel number in the ALTRO channel we're at 
1101     UShort_t offset       = 0;
1102
1103     // How many times the ALTRO Samples one VA1_ALICE channel 
1104     Int_t sampleRate = 1;
1105
1106     // A buffer to hold 1 ALTRO channel - Normally, one ALTRO channel
1107     // holds 128 VA1_ALICE channels, sampled at a rate of `sampleRate' 
1108     TArrayI channel(128 * sampleRate);
1109     
1110     // The Altro buffer 
1111     AliAltroBuffer* altro = 0;
1112     
1113     // Loop over the digits in the event.  Note, that we assume the
1114     // the digits are in order in the branch.   If they were not, we'd
1115     // have to cache all channels before we could write the data to
1116     // the ALTRO buffer, or we'd have to set up a map of the digits. 
1117     for (Int_t i = 0; i < nDigits; i++) {
1118       // Get the digit
1119       AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
1120
1121       UShort_t det    = digit->Detector();
1122       Char_t   ring   = digit->Ring();
1123       UShort_t sector = digit->Sector();
1124       UShort_t strip  = digit->Strip();
1125       if (det != prevDetector) {
1126         AliDebug(10, Form("FMD: New DDL, was %d, now %d",
1127                           kBaseDDL + prevDetector - 1,
1128                           kBaseDDL + det - 1));
1129         // If an altro exists, delete the object, flushing the data to
1130         // disk, and closing the file. 
1131         if (altro) { 
1132           // When the first argument is false, we write the real
1133           // header. 
1134           AliDebug(10, Form("New altro: Write channel at %d Strip: %d "
1135                             "Sector: %d  Ring: %d", 
1136                             i, startStrip, prevSector, prevRing));
1137           // TPC to FMD translations 
1138           // 
1139           //    TPC                FMD
1140           //    ----------+-----------
1141           //    pad       |      strip
1142           //    row       |     sector
1143           //    sector    |       ring
1144           // 
1145           altro->WriteChannel(Int_t(startStrip), 
1146                               Int_t(prevSector), 
1147                               Int_t((prevRing == 'I' ? 0 : 1)), 
1148                               channel.fN, channel.fArray, 0);
1149           altro->Flush();
1150           altro->WriteDataHeader(kFALSE, kFALSE);
1151           delete altro;
1152           altro = 0;
1153         }
1154
1155         prevDetector = det;
1156         // Need to open a new DDL! 
1157         Int_t ddlId = kBaseDDL + det - 1;
1158         TString filename(Form("%s_%d.ddl", GetName(),  ddlId));
1159
1160         AliDebug(10, Form("New altro buffer with DDL file %s", 
1161                           filename.Data()));
1162         AliDebug(10, Form("New altro at %d", i));
1163         // Create a new altro buffer - a `1' as the second argument
1164         // means `write mode' 
1165         altro = new AliAltroBuffer(filename.Data(), 1);
1166         
1167         // Write a dummy (first argument is true) header to the DDL
1168         // file - later on, when we close the file, we write the real
1169         // header
1170         altro->WriteDataHeader(kTRUE, kFALSE);
1171
1172         // Figure out the sample rate 
1173         if (digit->Count2() > 0) sampleRate = 2;
1174         if (digit->Count3() > 0) sampleRate = 3;
1175
1176         channel.Set(128 * sampleRate);
1177         offset     = 0;
1178         prevRing   = ring;
1179         prevSector = sector;
1180         startStrip = strip;
1181       }
1182       else if (offset == 128                        
1183                || digit->Ring() != prevRing 
1184                || digit->Sector() != prevSector) {
1185         // Force a new Altro channel
1186         AliDebug(10, Form("Flushing channel to disk because %s",
1187                           (offset == 128 ? "channel is full" :
1188                            (ring != prevRing ? "new ring up" :
1189                             "new sector up"))));
1190         AliDebug(10, Form("New Channel: Write channel at %d Strip: %d "
1191                           "Sector: %d  Ring: %d", 
1192                           i, startStrip, prevSector, prevRing));
1193         altro->WriteChannel(Int_t(startStrip), 
1194                             Int_t(prevSector), 
1195                             Int_t((prevRing == 'I' ? 0 : 1)), 
1196                             channel.fN, channel.fArray, 0);
1197         // Reset and update channel variables 
1198         channel.Reset(0);
1199         offset     = 0; 
1200         startStrip = strip;
1201         prevRing   = ring;
1202         prevSector = sector;
1203       }
1204
1205       // Store the counts of the ADC in the channel buffer 
1206       channel[offset * sampleRate] = digit->Count1();
1207       if (sampleRate > 1) 
1208         channel[offset * sampleRate + 1] = digit->Count2();
1209       if (sampleRate > 2) 
1210         channel[offset * sampleRate + 2] = digit->Count3();
1211       offset++;
1212     }
1213     // Finally, we need to close the final ALTRO buffer if it wasn't
1214     // already 
1215     if (altro) {
1216       altro->Flush();
1217       altro->WriteDataHeader(kFALSE, kFALSE);
1218       delete altro;
1219     }
1220   }
1221   fLoader->UnloadDigits();
1222 }
1223
1224 //==================================================================
1225 //
1226 // Various setter functions for the common paramters 
1227 //
1228
1229 //__________________________________________________________________
1230 void 
1231 AliFMD::SetLegLength(Double_t length) 
1232 {
1233   // Set lenght of plastic legs that hold the hybrid (print board and
1234   // silicon sensor) onto the honeycomp support
1235   //
1236   // DebugGuard guard("AliFMD::SetLegLength");
1237   AliDebug(10, "AliFMD::SetLegLength");
1238   fLegLength = length;
1239   fInner->SetLegLength(fLegLength);
1240   fOuter->SetLegLength(fLegLength);
1241 }
1242
1243 //__________________________________________________________________
1244 void 
1245 AliFMD::SetLegOffset(Double_t offset) 
1246 {
1247   // Set offset from edge of hybrid to plastic legs that hold the
1248   // hybrid (print board and silicon sensor) onto the honeycomp
1249   // support 
1250   //
1251   // DebugGuard guard("AliFMD::SetLegOffset");
1252   AliDebug(10, "AliFMD::SetLegOffset");
1253   fInner->SetLegOffset(offset);
1254   fOuter->SetLegOffset(offset);
1255 }
1256
1257 //__________________________________________________________________
1258 void 
1259 AliFMD::SetLegRadius(Double_t radius) 
1260 {
1261   // Set the diameter of the plastic legs that hold the hybrid (print
1262   // board and silicon sensor) onto the honeycomp support
1263   //
1264   // DebugGuard guard("AliFMD::SetLegRadius");
1265   AliDebug(10, "AliFMD::SetLegRadius");
1266   fLegRadius = radius;
1267   fInner->SetLegRadius(fLegRadius);
1268   fOuter->SetLegRadius(fLegRadius);
1269 }
1270
1271 //__________________________________________________________________
1272 void 
1273 AliFMD::SetModuleSpacing(Double_t spacing) 
1274 {
1275   // Set the distance between the front and back sensor modules
1276   // (module staggering). 
1277   //
1278   // DebugGuard guard("AliFMD::SetModuleSpacing");
1279   AliDebug(10, "AliFMD::SetModuleSpacing");  
1280   fModuleSpacing = spacing;
1281   fInner->SetModuleSpacing(fModuleSpacing);
1282   fOuter->SetModuleSpacing(fModuleSpacing);
1283 }
1284
1285 //====================================================================
1286 //
1287 // Utility 
1288 //
1289 //__________________________________________________________________
1290 void 
1291 AliFMD::Browse(TBrowser* b) 
1292 {
1293   // Browse this object. 
1294   //
1295   AliDebug(10, "AliFMD::Browse");
1296   AliDetector::Browse(b);
1297   if (fInner) b->Add(fInner, "Inner Ring");
1298   if (fOuter) b->Add(fOuter, "Outer Ring");
1299   if (fFMD1)  b->Add(fFMD1,  "FMD1 SubDetector");
1300   if (fFMD2)  b->Add(fFMD2,  "FMD2 SubDetector");
1301   if (fFMD3)  b->Add(fFMD3,  "FMD3 SubDetector");
1302 }
1303
1304
1305 //___________________________________________________________________
1306 //
1307 // EOF
1308 //