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