]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMD.cxx
Additional protection and warning message.
[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 3 sub-detectors FMD1, FMD2, and FMD3, each of
23 // which has 1 or 2 rings of silicon sensors. 
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 simulation
30 //
31 //
32 //       +----------+   +----------+   
33 //       | AliFMDv1 |   | AliFMDv0 |   
34 //       +----------+   +----------+   
35 //            |              |                    +-----------------+
36 //       +----+--------------+                 +--| AliFMDDigitizer |
37 //       |                                     |  +-----------------+
38 //       |           +---------------------+   |
39 //       |        +--| AliFMDBaseDigitizer |<--+
40 //       V     1  |  +---------------------+   |
41 //  +--------+<>--+                            |  +------------------+
42 //  | AliFMD |                                 +--| AliFMDSDigitizer |    
43 //  +--------+<>--+                               +------------------+       
44 //             1  |  +-----------------+ 
45 //                +--| AliFMDSimulator |
46 //                   +-----------------+
47 //                            ^              
48 //                            |
49 //              +-------------+-------------+
50 //              |                           |         
51 //     +--------------------+   +-------------------+
52 //     | AliFMDGeoSimulator |   | AliFMDG3Simulator | 
53 //     +--------------------+   +-------------------+
54 //               ^                         ^
55 //               |                         |
56 //  +-----------------------+   +----------------------+
57 //  | AliFMDGeoOldSimulator |   | AliFMDG3OldSimulator |
58 //  +-----------------------+   +----------------------+
59 //      
60 //
61 // *  AliFMD 
62 //    This defines the interface for the various parts of AliROOT that
63 //    uses the FMD, like AliFMDSimulator, AliFMDDigitizer, 
64 //    AliFMDReconstructor, and so on. 
65 //
66 // *  AliFMDv0
67 //    This is a concrete implementation of the AliFMD interface. 
68 //    It is the responsibility of this class to create the FMD
69 //    geometry.
70 //
71 // *  AliFMDv1 
72 //    This is a concrete implementation of the AliFMD interface. 
73 //    It is the responsibility of this class to create the FMD
74 //    geometry, process hits in the FMD, and serve hits and digits to
75 //    the various clients. 
76 //  
77 // *  AliFMDSimulator
78 //    This is the base class for the FMD simulation tasks.   The
79 //    simulator tasks are responsible to implment the geoemtry, and
80 //    process hits. 
81 //                                                                          
82 // *  AliFMDGeoSimulator
83 //    This is a concrete implementation of the AliFMDSimulator that
84 //    uses the TGeo classes directly only. 
85 //
86 // *  AliFMDG3Simulator
87 //    This is a concrete implementation of the AliFMDSimulator that
88 //    uses the TVirtualMC interface with GEANT 3.21-like messages.
89 //
90
91 // These files are not in the same directory, so there's no reason to
92 // ask the preprocessor to search in the current directory for these
93 // files by including them with `#include "..."' 
94 #include <math.h>               // __CMATH__
95 #include <TClonesArray.h>       // ROOT_TClonesArray
96 #include <TGeometry.h>          // ROOT_TGeomtry
97 #include <TNode.h>              // ROOT_TNode
98 #include <TXTRU.h>              // ROOT_TXTRU
99 #include <TRotMatrix.h>         // ROOT_TRotMatrix
100 #include <TTUBE.h>              // ROOT_TTUBE
101 #include <TTree.h>              // ROOT_TTree
102 #include <TBrowser.h>           // ROOT_TBrowser
103 #include <TMath.h>              // ROOT_TMath
104 #include <TVirtualMC.h>         // ROOT_TVirtualMC
105
106 #include <AliRunDigitizer.h>    // ALIRUNDIGITIZER_H
107 #include <AliLoader.h>          // ALILOADER_H
108 #include <AliRun.h>             // ALIRUN_H
109 #include <AliMC.h>              // ALIMC_H
110 #include "AliMagF.h"            // ALIMAGF_H
111 #include <AliLog.h>             // ALILOG_H
112 #include "AliFMD.h"             // ALIFMD_H
113 #include "AliFMDDigit.h"        // ALIFMDDIGIG_H
114 #include "AliFMDHit.h"          // ALIFMDHIT_H
115 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
116 #include "AliFMDDetector.h"     // ALIFMDDETECTOR_H
117 #include "AliFMDRing.h"         // ALIFMDRING_H
118 #include "AliFMDDigitizer.h"    // ALIFMDDIGITIZER_H
119 #ifdef USE_PRE_MOVE
120 #include "AliFMDSimulator.h"    // ALIFMDSIMULATOR_H
121 #include "AliFMDG3Simulator.h"  // ALIFMDG3SIMULATOR_H
122 #include "AliFMDGeoSimulator.h" // ALIFMDGEOSIMULATOR_H
123 #include "AliFMDG3OldSimulator.h"       // ALIFMDG3OLDSIMULATOR_H
124 #include "AliFMDGeoOldSimulator.h"      // ALIFMDGEOOLDSIMULATOR_H
125 #else
126 #include "AliFMDGeometryBuilderSimple.h"
127 #endif
128 #include "AliFMDRawWriter.h"    // ALIFMDRAWWRITER_H
129 #include <TVector2.h>
130
131 //____________________________________________________________________
132 ClassImp(AliFMD)
133 #if 0
134   ; // This is to keep Emacs from indenting the next line 
135 #endif 
136
137 //____________________________________________________________________
138 AliFMD::AliFMD()
139   : AliDetector(),
140     fSDigits(0), 
141     fNsdigits(0),
142     fDetailed(kTRUE),
143 #ifdef USE_PRE_MOVE
144     fSimulator(0), 
145 #endif
146     fBad(0)
147 {
148   //
149   // Default constructor for class AliFMD
150   //
151   AliDebug(10, "\tDefault CTOR");
152   fHits        = 0;
153   fDigits      = 0;
154   fIshunt      = 0;
155   fUseOld      = kFALSE;
156   fUseAssembly = kTRUE;
157   fBad         = new TClonesArray("AliFMDHit");
158 }
159
160 //____________________________________________________________________
161 AliFMD::AliFMD(const AliFMD& other)
162   : AliDetector(other),
163     fSDigits(other.fSDigits), 
164     fNsdigits(other.fNsdigits),
165     fDetailed(other.fDetailed),
166 #ifdef USE_PRE_MOVE
167     fSimulator(other.fSimulator),
168 #endif
169     fBad(other.fBad)
170 {
171   // Copy constructor 
172   fUseOld      = other.fUseOld;
173   fUseAssembly = other.fUseAssembly;
174 }
175
176 //____________________________________________________________________
177 AliFMD::AliFMD(const char *name, const char *title)
178   : AliDetector (name, title),
179     fSDigits(0),
180     fNsdigits(0),
181     fDetailed(kTRUE),
182 #ifdef USE_PRE_MOVE
183     fSimulator(0), 
184 #endif
185     fBad(0)
186 {
187   //
188   // Standard constructor for Forward Multiplicity Detector
189   //
190   AliDebug(10, "\tStandard CTOR");
191   fUseOld      = kFALSE;
192   fUseAssembly = kFALSE;
193   fBad         = new TClonesArray("AliFMDHit");
194   
195   // Initialise Hit array
196   HitsArray();
197   gAlice->GetMCApp()->AddHitList(fHits);
198
199   // (S)Digits for the detectors disk
200   DigitsArray();
201   SDigitsArray();
202   
203   // CHC: What is this?
204   fIshunt = 0;
205   SetMarkerColor(kRed);
206   SetLineColor(kYellow);
207 }
208
209 //____________________________________________________________________
210 AliFMD::~AliFMD ()
211 {
212   // Destructor for base class AliFMD
213   if (fHits) {
214     fHits->Delete();
215     delete fHits;
216     fHits = 0;
217   }
218   if (fDigits) {
219     fDigits->Delete();
220     delete fDigits;
221     fDigits = 0;
222   }
223   if (fSDigits) {
224     fSDigits->Delete();
225     delete fSDigits;
226     fSDigits = 0;
227   }
228   if (fBad) {
229     fBad->Delete();
230     delete fBad;
231     fBad = 0;
232   }
233 }
234
235 //____________________________________________________________________
236 AliFMD&
237 AliFMD::operator=(const AliFMD& other)
238 {
239   // Assignment operator
240   AliDetector::operator=(other);
241   fSDigits              = other.fSDigits; 
242   fNsdigits             = other.fNsdigits;
243   fDetailed             = other.fDetailed;
244 #ifdef USE_PRE_MOVE
245   fSimulator            = other.fSimulator;
246 #endif
247   fBad                  = other.fBad;
248   return *this;
249 }
250
251 //====================================================================
252 //
253 // GEometry ANd Traking
254 //
255 //____________________________________________________________________
256 void 
257 AliFMD::CreateGeometry()
258 {
259   //
260   // Create the geometry of Forward Multiplicity Detector.  The actual
261   // construction of the geometry is delegated to the class AliFMDRing
262   // and AliFMDSubDetector and the relevant derived classes. 
263   //
264   // The flow of this member function is:
265   //
266   //   FOR rings fInner and fOuter DO  
267   //     AliFMDRing::Init();
268   //   END FOR
269   // 
270   //   Set up hybrud card support (leg) volume shapes  
271   // 
272   //   FOR rings fInner and fOuter DO  
273   //     AliFMDRing::SetupGeometry();
274   //   END FOR
275   // 
276   //   FOR subdetectors fFMD1, fFMD2, and fFMD3 DO 
277   //     AliFMDSubDetector::SetupGeomtry();
278   //   END FOR
279   // 
280   //   FOR subdetectors fFMD1, fFMD2, and fFMD3 DO 
281   //     AliFMDSubDetector::Geomtry();
282   //   END FOR
283   //
284 #ifndef USE_PRE_MOVE
285   AliFMDGeometry*  fmd = AliFMDGeometry::Instance();
286   if (fUseOld) fmd->SetBuilder(new AliFMDGeometryBuilderSimple(fDetailed));
287   fmd->SetDetailed(fDetailed);
288   fmd->UseAssembly(fUseAssembly);
289   fmd->Build();
290 #else
291   if (!fSimulator) {
292     AliFatal("Simulator object not made yet!");
293     return;
294   }
295   fSimulator->DefineGeometry();
296 #endif
297 }    
298
299 //____________________________________________________________________
300 void AliFMD::CreateMaterials() 
301 {
302   // Define the materials and tracking mediums needed by the FMD
303   // simulation.   These mediums are made by sending the messages
304   // AliMaterial, AliMixture, and AliMedium to the passed AliModule
305   // object module.   The defined mediums are 
306   // 
307   //    FMD Si$         Silicon (active medium in sensors)
308   //    FMD C$          Carbon fibre (support cone for FMD3 and vacuum pipe)
309   //    FMD Al$         Aluminium (honeycomb support plates)
310   //    FMD PCB$        Printed Circuit Board (FEE board with VA1_3)
311   //    FMD Chip$       Electronics chips (currently not used)
312   //    FMD Air$        Air (Air in the FMD)
313   //    FMD Plastic$    Plastic (Support legs for the hybrid cards)
314   //
315   // Pointers to TGeoMedium objects are retrived from the TGeoManager
316   // singleton.  These pointers are later used when setting up the
317   // geometry 
318   AliDebug(10, "\tCreating materials");
319   // Get pointer to geometry singleton object. 
320   AliFMDGeometry* geometry = AliFMDGeometry::Instance();
321   geometry->Init();
322 #if 0
323   if (gGeoManager && gGeoManager->GetMedium("FMD Si$")) {
324     // We need to figure out the some stuff about the geometry
325     fmd->ExtractGeomInfo();
326     return;
327   }
328 #endif  
329 #ifndef USE_PRE_MOVE     
330   Int_t    id;
331   Double_t a                = 0;
332   Double_t z                = 0;
333   Double_t density          = 0;
334   Double_t radiationLength  = 0;
335   Double_t absorbtionLength = 999;
336   Int_t    fieldType        = gAlice->Field()->Integ();     // Field type 
337   Double_t maxField         = gAlice->Field()->Max();     // Field max.
338   Double_t maxBending       = 0;     // Max Angle
339   Double_t maxStepSize      = 0.001; // Max step size 
340   Double_t maxEnergyLoss    = 1;     // Max Delta E
341   Double_t precision        = 0.001; // Precision
342   Double_t minStepSize      = 0.001; // Minimum step size 
343  
344   // Silicon 
345   a                = 28.0855;
346   z                = 14.;
347   density          = geometry->GetSiDensity();
348   radiationLength  = 9.36;
349   maxBending       = 1;
350   maxStepSize      = .001;
351   precision        = .001;
352   minStepSize      = .001;
353   id               = kSiId;
354   AliMaterial(id, "Si$", a, z, density, radiationLength, absorbtionLength);
355   AliMedium(kSiId, "Si$", id,1,fieldType,maxField,maxBending,
356             maxStepSize,maxEnergyLoss,precision,minStepSize);
357   
358
359   // Carbon 
360   a                = 12.011;
361   z                = 6.;
362   density          = 2.265;
363   radiationLength  = 18.8;
364   maxBending       = 10;
365   maxStepSize      = .01;
366   precision        = .003;
367   minStepSize      = .003;
368   id               = kCarbonId;
369   AliMaterial(id, "Carbon$", a, z, density, radiationLength, absorbtionLength);
370   AliMedium(kCarbonId, "Carbon$", id,0,fieldType,maxField,maxBending,
371                     maxStepSize,maxEnergyLoss,precision,minStepSize);
372
373   // Aluminum
374   a                = 26.981539;
375   z                = 13.;
376   density          = 2.7;
377   radiationLength  = 8.9;
378   id               = kAlId;
379   AliMaterial(id, "Aluminum$",a,z, density, radiationLength, absorbtionLength);
380   AliMedium(kAlId, "Aluminum$", id, 0, fieldType, maxField, maxBending,
381             maxStepSize, maxEnergyLoss, precision, minStepSize);
382   
383   
384   // Copper 
385   a                = 63.546;
386   z                = 29;
387   density          =  8.96;
388   radiationLength  =  1.43;
389   id               = kCopperId;
390   AliMaterial(id, "Copper$", 
391                       a, z, density, radiationLength, absorbtionLength);
392   AliMedium(kCopperId, "Copper$", id, 0, fieldType, maxField, maxBending,
393             maxStepSize, maxEnergyLoss, precision, minStepSize);
394   
395
396   // Silicon chip 
397   {
398     Float_t as[] = { 12.0107,      14.0067,      15.9994,
399                       1.00794,     28.0855,     107.8682 };
400     Float_t zs[] = {  6.,           7.,           8.,
401                       1.,          14.,          47. };
402     Float_t ws[] = {  0.039730642,  0.001396798,  0.01169634,
403                       0.004367771,  0.844665,     0.09814344903 };
404     density          = 2.36436;
405     maxBending       = 10;
406     maxStepSize      = .01;
407     precision        = .003;
408     minStepSize      = .003;
409     id               = kSiChipId;
410     AliMixture(id, "Si Chip$", as, zs, density, 6, ws);
411     AliMedium(kSiChipId, "Si Chip$",  id, 0, fieldType, maxField, maxBending, 
412               maxStepSize, maxEnergyLoss, precision, minStepSize);
413   }
414   
415   // Kaption
416   {
417     Float_t as[] = { 1.00794,  12.0107,  14.010,   15.9994};
418     Float_t zs[] = { 1.,        6.,       7.,       8.};
419     Float_t ws[] = { 0.026362,  0.69113,  0.07327,  0.209235};
420     density          = 1.42;
421     maxBending       = 1;
422     maxStepSize      = .001;
423     precision        = .001;
424     minStepSize      = .001;
425     id               = kKaptonId;
426     AliMixture(id, "Kaption$", as, zs, density, 4, ws);
427     AliMedium(kKaptonId, "Kaption$", id,0,fieldType,maxField,maxBending,
428               maxStepSize,maxEnergyLoss,precision,minStepSize);
429   }
430
431   // Air
432   {
433     Float_t as[] = { 12.0107, 14.0067,   15.9994,  39.948 };
434     Float_t zs[] = {  6.,      7.,       8.,       18. };
435     Float_t ws[] = { 0.000124, 0.755267, 0.231781, 0.012827 }; 
436     density      = .00120479;
437     maxBending   = 1;
438     maxStepSize  = .001;
439     precision    = .001;
440     minStepSize  = .001;
441     id           = kAirId;
442     AliMixture(id, "Air$", as, zs, density, 4, ws);
443     AliMedium(kAirId, "Air$", id,0,fieldType,maxField,maxBending,
444               maxStepSize,maxEnergyLoss,precision,minStepSize);
445   }
446   
447   // PCB
448   {
449     Float_t zs[] = { 14.,         20.,         13.,         12.,
450                       5.,         22.,         11.,         19.,
451                      26.,          9.,          8.,          6.,
452                       7.,          1.};
453     Float_t as[] = { 28.0855,     40.078,      26.981538,   24.305, 
454                      10.811,      47.867,      22.98977,    39.0983,
455                      55.845,      18.9984,     15.9994,     12.0107,
456                      14.0067,      1.00794};
457     Float_t ws[] = {  0.15144894,  0.08147477,  0.04128158,  0.00904554, 
458                       0.01397570,  0.00287685,  0.00445114,  0.00498089,
459                       0.00209828,  0.00420000,  0.36043788,  0.27529426,
460                       0.01415852,  0.03427566};
461     density      = 1.8;
462     maxBending   = 1;
463     maxStepSize  = .001;
464     precision    = .001;
465     minStepSize  = .001;
466     id           = kPcbId;
467     AliMixture(id, "PCB$", as, zs, density, 14, ws);
468     AliMedium(kPcbId, "PCB$", id,0,fieldType,maxField,maxBending,
469               maxStepSize,maxEnergyLoss,precision,minStepSize);
470   }
471   
472   // Plastic 
473   {
474     Float_t as[] = { 1.01, 12.01 };
475     Float_t zs[] = { 1.,   6.    };
476     Float_t ws[] = { 1.,   1.    };
477     density      = 1.03;
478     maxBending   = 10;
479     maxStepSize  = .01;
480     precision    = .003;
481     minStepSize  = .003;
482     id           = kPlasticId;
483     AliMixture(id, "Plastic$", as, zs, density, -2, ws);
484     AliMedium(kPlasticId, "Plastic$", id,0,fieldType,maxField,maxBending,
485               maxStepSize,maxEnergyLoss,precision,minStepSize);
486   }
487 #else
488   AliDebug(10, "\tCreating materials");
489
490   if (fSimulator) {
491     AliFatal("Simulator object already instantised!");
492     return;
493   }
494   TVirtualMC* mc = TVirtualMC::GetMC();
495
496   Bool_t geo = mc->IsRootGeometrySupported();
497   if (geo) {
498     if (fUseOld) 
499       fSimulator = new AliFMDGeoOldSimulator(this, fDetailed);
500     else 
501       fSimulator = new AliFMDGeoSimulator(this, fDetailed);
502   }
503   else {
504     if (fUseOld) 
505       fSimulator = new AliFMDG3OldSimulator(this, fDetailed);
506     else    
507       fSimulator = new AliFMDG3Simulator(this, fDetailed);
508   }
509   AliDebug(1, Form("using a %s as simulation backend", 
510                    fSimulator->IsA()->GetName()));
511   fSimulator->SetDetailed(fDetailed);
512   fSimulator->UseAssembly(fUseAssembly);
513   fSimulator->DefineMaterials();
514 #endif
515 }
516
517 //____________________________________________________________________
518 void  
519 AliFMD::Init()
520 {}
521
522 //____________________________________________________________________
523 void
524 AliFMD::FinishEvent()
525 {
526 #ifndef USE_PRE_MOVE
527   if (fBad && fBad->GetEntries() > 0) {
528     AliWarning((Form("EndEvent", "got %d 'bad' hits", fBad->GetEntries())));
529     TIter next(fBad);
530     AliFMDHit* hit;
531     while ((hit = static_cast<AliFMDHit*>(next()))) 
532       hit->Print("D");
533     fBad->Clear();
534   }
535 #else
536   if (fSimulator) fSimulator->EndEvent();
537 #endif
538 }
539
540
541 //====================================================================
542 //
543 // Graphics and event display
544 //
545 //____________________________________________________________________
546 void 
547 AliFMD::BuildGeometry()
548 {
549   //
550   // Build simple ROOT TNode geometry for event display
551   //
552   // Build a simplified geometry of the FMD used for event display  
553   // 
554   // The actual building of the TNodes is done by
555   // AliFMDSubDetector::SimpleGeometry. 
556   AliDebug(10, "\tCreating a simplified geometry");
557
558   AliFMDGeometry* fmd = AliFMDGeometry::Instance();
559   
560   static TXTRU*     innerShape = 0;
561   static TXTRU*     outerShape = 0;
562   static TObjArray* innerRot   = 0;
563   static TObjArray* outerRot   = 0;
564
565   if (!innerShape || !outerShape) {
566     // Make the shapes for the modules 
567     for (Int_t i = 0; i < 2; i++) {
568       AliFMDRing* r = 0;
569       switch (i) {
570       case 0: r = fmd->GetRing('I'); break;
571       case 1: r = fmd->GetRing('O'); break;
572       }
573       if (!r) {
574         AliError(Form("no ring found for i=%d", i));
575         return;
576       }
577       Double_t    siThick  = r->GetSiThickness();
578       const Int_t nv       = r->GetNVerticies();
579       Double_t    theta    = r->GetTheta();
580       Int_t       nmod     = r->GetNModules();
581       
582       TXTRU* shape = new TXTRU(r->GetName(), r->GetTitle(), "void", nv, 2);
583       for (Int_t j = 0; j < nv; j++) {
584         TVector2* vv = r->GetVertex(nv - 1 - j);
585         shape->DefineVertex(j, vv->X(), vv->Y());
586       }
587       shape->DefineSection(0, -siThick / 2, 1, 0, 0);
588       shape->DefineSection(1, +siThick / 2, 1, 0, 0);
589       shape->SetLineColor(GetLineColor());
590       
591       TObjArray* rots = new TObjArray(nmod);
592       for (Int_t j = 0; j < nmod; j++) {
593         Double_t th = (j + .5) * theta * 2;
594         TString name(Form("FMD_ring_%c_rot_%02d", r->GetId(), j));
595         TString title(Form("FMD Ring %c Rotation # %d", r->GetId(), j));
596         TRotMatrix* rot = new TRotMatrix(name.Data(), title.Data(),
597                                          90, th, 90, fmod(90+th,360), 0, 0);
598         rots->AddAt(rot, j);
599       }
600       
601       switch (r->GetId()) {
602       case 'i':
603       case 'I': innerShape = shape; innerRot = rots; break;
604       case 'o':
605       case 'O': outerShape = shape; outerRot = rots; break;
606       }
607     }
608   }
609   
610   TNode* top = gAlice->GetGeometry()->GetNode("alice");
611   
612   for (Int_t i = 1; i <= 3; i++) {
613     AliFMDDetector* det = fmd->GetDetector(i);
614     if (!det) {
615       Warning("BuildGeometry", "FMD%d seems to be disabled", i);
616       continue;
617     }
618     Double_t w  = 0;
619     Double_t rh = det->GetRing('I')->GetHighR();
620     Char_t   id = 'I';
621     if (det->GetRing('O')) {
622       w  = TMath::Abs(det->GetRingZ('O') - det->GetRingZ('I'));
623       id = (TMath::Abs(det->GetRingZ('O')) 
624             > TMath::Abs(det->GetRingZ('I')) ? 'O' : 'I');
625       rh = det->GetRing('O')->GetHighR();
626     }
627     w += (det->GetRing(id)->GetModuleSpacing() +
628           det->GetRing(id)->GetSiThickness());
629     TShape* shape = new TTUBE(det->GetName(), det->GetTitle(), "void",
630                               det->GetRing('I')->GetLowR(), rh, w / 2);
631     Double_t z = (det->GetRingZ('I') - w / 2);
632     if (z > 0) z += det->GetRing(id)->GetModuleSpacing();
633     top->cd();
634     TNode* node = new TNode(det->GetName(), det->GetTitle(), shape, 
635                             0, 0, z, 0);
636     fNodes->Add(node);
637     
638     for (Int_t j = 0; j < 2; j++) {
639       AliFMDRing* r      = 0;
640       TShape*     rshape = 0;
641       TObjArray*  rots   = 0;
642       switch (j) {
643       case 0: 
644         r = det->GetRing('I'); rshape = innerShape; rots = innerRot; break;
645       case 1: 
646         r = det->GetRing('O'); rshape = outerShape; rots = outerRot; break;
647       }
648       if (!r) continue;
649       
650       Double_t    siThick  = r->GetSiThickness();
651       Int_t       nmod     = r->GetNModules();
652       Double_t    modspace = r->GetModuleSpacing();
653       Double_t    rz       = - (z - det->GetRingZ(r->GetId()));
654       
655       for (Int_t k = 0; k < nmod; k++) {
656         node->cd();
657         Double_t    offz    = (k % 2 == 1 ? modspace : 0);
658         TRotMatrix* rot     = static_cast<TRotMatrix*>(rots->At(k));
659         TString name(Form("%s%c_module_%02d", det->GetName(), r->GetId(),k));
660         TString title(Form("%s%c Module %d", det->GetName(), r->GetId(),k));
661         TNode* mnod = new TNode(name.Data(), title.Data(), rshape, 
662                                 0, 0, rz - siThick / 2 
663                                 + TMath::Sign(offz,z), rot);
664         mnod->SetLineColor(GetLineColor());
665         fNodes->Add(mnod);
666       } // for (Int_t k = 0 ; ...)
667     } // for (Int_t j = 0 ; ...)
668   } // for (Int_t i = 1 ; ...)
669 }
670
671 //____________________________________________________________________
672 void 
673 AliFMD::DrawDetector()
674 {
675   //
676   // Draw a shaded view of the Forward multiplicity detector
677   //
678   // DebugGuard guard("AliFMD::DrawDetector");
679   AliDebug(10, "\tDraw detector");
680   
681 #if 0
682   //Set ALIC mother transparent
683   gMC->Gsatt("ALIC","SEEN",0);
684   //
685   gMC->Gdopt("hide", "on");
686   gMC->Gdopt("shad", "on");
687   gMC->Gsatt("*", "fill", 7);
688   gMC->SetClipBox(".");
689   gMC->SetClipBox("*", 0, 1000, -1000, 1000, -1000, 1000);
690   gMC->DefaultRange();
691   gMC->Gdraw("alic", 40, 30, 0, 12, 12, .055, .055);
692   gMC->Gdhead(1111, "Forward Multiplicity Detector");
693   gMC->Gdman(16, 10, "MAN");
694   gMC->Gdopt("hide", "off");
695 #endif
696 }
697
698 //____________________________________________________________________
699 Int_t 
700 AliFMD::DistanceToPrimitive(Int_t, Int_t)
701 {
702   //
703   // Calculate the distance from the mouse to the FMD on the screen
704   // Dummy routine
705   //
706   return 9999;
707 }
708
709 //====================================================================
710 //
711 // Hit and Digit managment 
712 //
713 //____________________________________________________________________
714 void 
715 AliFMD::MakeBranch(Option_t * option)
716 {
717   // Create Tree branches for the FMD.
718   //
719   // Options:
720   //
721   //    H          Make a branch of TClonesArray of AliFMDHit's
722   //    D          Make a branch of TClonesArray of AliFMDDigit's
723   //    S          Make a branch of TClonesArray of AliFMDSDigit's
724   // 
725   const Int_t kBufferSize = 16000;
726   TString branchname(GetName());
727   TString opt(option);
728   
729   if (opt.Contains("H", TString::kIgnoreCase)) {
730     HitsArray();
731     AliDetector::MakeBranch(option); 
732   }
733   if (opt.Contains("D", TString::kIgnoreCase)) { 
734     DigitsArray();
735     MakeBranchInTree(fLoader->TreeD(), branchname.Data(),
736                      &fDigits, kBufferSize, 0);
737   }
738   if (opt.Contains("S", TString::kIgnoreCase)) { 
739     SDigitsArray();
740     MakeBranchInTree(fLoader->TreeS(), branchname.Data(),
741                      &fSDigits, kBufferSize, 0);
742   }
743 }
744
745 //____________________________________________________________________
746 void 
747 AliFMD::SetTreeAddress()
748 {
749   // Set branch address for the Hits, Digits, and SDigits Tree.
750   if (fLoader->TreeH()) HitsArray();
751   AliDetector::SetTreeAddress();
752
753   TTree *treeD = fLoader->TreeD();
754   if (treeD) {
755     DigitsArray();
756     TBranch* branch = treeD->GetBranch ("FMD");
757     if (branch) branch->SetAddress(&fDigits);
758   }
759
760   TTree *treeS = fLoader->TreeS();
761   if (treeS) {
762     SDigitsArray();
763     TBranch* branch = treeS->GetBranch ("FMD");
764     if (branch) branch->SetAddress(&fSDigits);
765   }
766 }
767
768
769
770 //____________________________________________________________________
771 void 
772 AliFMD::SetHitsAddressBranch(TBranch *b)
773 {
774   // Set the TClonesArray to read hits into. 
775   b->SetAddress(&fHits);
776 }
777
778 //____________________________________________________________________
779 void 
780 AliFMD::AddHit(Int_t track, Int_t *vol, Float_t *hits) 
781 {
782   // Add a hit to the hits tree 
783   // 
784   // The information of the two arrays are decoded as 
785   // 
786   // Parameters
787   //    track                Track #
788   //    ivol[0]  [UShort_t ] Detector # 
789   //    ivol[1]  [Char_t   ] Ring ID 
790   //    ivol[2]  [UShort_t ] Sector #
791   //    ivol[3]  [UShort_t ] Strip # 
792   //    hits[0]  [Float_t  ] Track's X-coordinate at hit 
793   //    hits[1]  [Float_t  ] Track's Y-coordinate at hit
794   //    hits[3]  [Float_t  ] Track's Z-coordinate at hit
795   //    hits[4]  [Float_t  ] X-component of track's momentum             
796   //    hits[5]  [Float_t  ] Y-component of track's momentum             
797   //    hits[6]  [Float_t  ] Z-component of track's momentum            
798   //    hits[7]  [Float_t  ] Energy deposited by track                  
799   //    hits[8]  [Int_t    ] Track's particle Id # 
800   //    hits[9]  [Float_t  ] Time when the track hit
801   // 
802   // 
803   AddHitByFields(track, 
804                  UShort_t(vol[0]),  // Detector # 
805                  Char_t(vol[1]),    // Ring ID
806                  UShort_t(vol[2]),  // Sector # 
807                  UShort_t(vol[3]),  // Strip # 
808                  hits[0],           // X
809                  hits[1],           // Y
810                  hits[2],           // Z
811                  hits[3],           // Px
812                  hits[4],           // Py
813                  hits[5],           // Pz
814                  hits[6],           // Energy loss 
815                  Int_t(hits[7]),    // PDG 
816                  hits[8]);          // Time
817 }
818
819 //____________________________________________________________________
820 AliFMDHit*
821 AliFMD::AddHitByFields(Int_t    track, 
822                        UShort_t detector, 
823                        Char_t   ring, 
824                        UShort_t sector, 
825                        UShort_t strip, 
826                        Float_t  x, 
827                        Float_t  y, 
828                        Float_t  z,
829                        Float_t  px, 
830                        Float_t  py, 
831                        Float_t  pz,
832                        Float_t  edep,
833                        Int_t    pdg,
834                        Float_t  t, 
835                        Float_t  l, 
836                        Bool_t   stop)
837 {
838   //
839   // Add a hit to the list
840   //
841   // Parameters:
842   // 
843   //    track     Track #
844   //    detector  Detector # (1, 2, or 3)                      
845   //    ring      Ring ID ('I' or 'O')
846   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
847   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
848   //    x         Track's X-coordinate at hit
849   //    y         Track's Y-coordinate at hit
850   //    z         Track's Z-coordinate at hit
851   //    px        X-component of track's momentum 
852   //    py        Y-component of track's momentum
853   //    pz        Z-component of track's momentum
854   //    edep      Energy deposited by track
855   //    pdg       Track's particle Id #
856   //    t         Time when the track hit 
857   //    l         Track length through the material. 
858   //    stop      Whether track was stopped or disappeared
859   // 
860   TClonesArray& a = *(HitsArray());
861   // Search through the list of already registered hits, and see if we
862   // find a hit with the same parameters.  If we do, then don't create
863   // a new hit, but rather update the energy deposited in the hit.
864   // This is done, so that a FLUKA based simulation will get the
865   // number of hits right, not just the enerrgy deposition. 
866   AliFMDHit* hit = 0;
867   for (Int_t i = 0; i < fNhits; i++) {
868     if (!a.At(i)) continue;
869     hit = static_cast<AliFMDHit*>(a.At(i));
870     if (hit->Detector() == detector 
871         && hit->Ring() == ring
872         && hit->Sector() == sector 
873         && hit->Strip() == strip
874         && hit->Track() == track) {
875       AliDebug(1, Form("already had a hit in FMD%d%c[%2d,%3d] for track # %d,"
876                        " adding energy (%f) to that hit (%f) -> %f", 
877                        detector, ring, sector, strip, track, edep, hit->Edep(),
878                        hit->Edep() + edep));
879       hit->SetEdep(hit->Edep() + edep);
880       return hit;
881     }
882   }
883   // If hit wasn't already registered, do so know. 
884   hit = new (a[fNhits]) AliFMDHit(fIshunt, track, detector, ring, sector, 
885                                   strip, x, y, z, px, py, pz, edep, pdg, t, 
886                                   l, stop);
887   fNhits++;
888   return hit;
889 }
890
891 //____________________________________________________________________
892 void 
893 AliFMD::AddDigit(Int_t* digits, Int_t*)
894 {
895   // Add a digit to the Digit tree 
896   // 
897   // Paramters 
898   //
899   //    digits[0]  [UShort_t] Detector #
900   //    digits[1]  [Char_t]   Ring ID
901   //    digits[2]  [UShort_t] Sector #
902   //    digits[3]  [UShort_t] Strip #
903   //    digits[4]  [UShort_t] ADC Count 
904   //    digits[5]  [Short_t]  ADC Count, -1 if not used
905   //    digits[6]  [Short_t]  ADC Count, -1 if not used 
906   // 
907   AddDigitByFields(UShort_t(digits[0]),  // Detector #
908                    Char_t(digits[1]),    // Ring ID
909                    UShort_t(digits[2]),  // Sector #
910                    UShort_t(digits[3]),  // Strip #
911                    UShort_t(digits[4]),  // ADC Count1 
912                    Short_t(digits[5]),   // ADC Count2 
913                    Short_t(digits[6]));  // ADC Count3 
914 }
915
916 //____________________________________________________________________
917 void 
918 AliFMD::AddDigitByFields(UShort_t detector, 
919                          Char_t   ring, 
920                          UShort_t sector, 
921                          UShort_t strip, 
922                          UShort_t count1, 
923                          Short_t  count2,
924                          Short_t  count3)
925 {
926   // add a real digit - as coming from data
927   // 
928   // Parameters 
929   //
930   //    detector  Detector # (1, 2, or 3)                      
931   //    ring      Ring ID ('I' or 'O')
932   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
933   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
934   //    count1    ADC count (a 10-bit word)
935   //    count2    ADC count (a 10-bit word), or -1 if not used
936   //    count3    ADC count (a 10-bit word), or -1 if not used
937   TClonesArray& a = *(DigitsArray());
938   
939   new (a[fNdigits++]) 
940     AliFMDDigit(detector, ring, sector, strip, count1, count2, count3);
941 }
942
943 //____________________________________________________________________
944 void 
945 AliFMD::AddSDigit(Int_t* digits)
946 {
947   // Add a digit to the SDigit tree 
948   // 
949   // Paramters 
950   //
951   //    digits[0]  [UShort_t] Detector #
952   //    digits[1]  [Char_t]   Ring ID
953   //    digits[2]  [UShort_t] Sector #
954   //    digits[3]  [UShort_t] Strip #
955   //    digits[4]  [Float_t]  Total energy deposited 
956   //    digits[5]  [UShort_t] ADC Count 
957   //    digits[6]  [Short_t]  ADC Count, -1 if not used
958   //    digits[7]  [Short_t]  ADC Count, -1 if not used 
959   // 
960   AddSDigitByFields(UShort_t(digits[0]),  // Detector #
961                     Char_t(digits[1]),    // Ring ID
962                     UShort_t(digits[2]),  // Sector #
963                     UShort_t(digits[3]),  // Strip #
964                     Float_t(digits[4]),   // Edep
965                     UShort_t(digits[5]),  // ADC Count1 
966                     Short_t(digits[6]),   // ADC Count2 
967                     Short_t(digits[7]));  // ADC Count3 
968 }
969
970 //____________________________________________________________________
971 void 
972 AliFMD::AddSDigitByFields(UShort_t detector, 
973                           Char_t   ring, 
974                           UShort_t sector, 
975                           UShort_t strip, 
976                           Float_t  edep,
977                           UShort_t count1, 
978                           Short_t  count2,
979                           Short_t  count3)
980 {
981   // add a summable digit
982   // 
983   // Parameters 
984   //
985   //    detector  Detector # (1, 2, or 3)                      
986   //    ring      Ring ID ('I' or 'O')
987   //    sector    Sector # (For inner/outer rings: 0-19/0-39)
988   //    strip     Strip # (For inner/outer rings: 0-511/0-255)
989   //    edep      Total energy deposited
990   //    count1    ADC count (a 10-bit word)
991   //    count2    ADC count (a 10-bit word), or -1 if not used
992   //    count3    ADC count (a 10-bit word), or -1 if not used
993   //
994   TClonesArray& a = *(SDigitsArray());
995   
996   new (a[fNsdigits++]) 
997     AliFMDSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
998 }
999
1000 //____________________________________________________________________
1001 void 
1002 AliFMD::ResetSDigits()
1003 {
1004   //
1005   // Reset number of digits and the digits array for this detector
1006   //
1007   fNsdigits   = 0;
1008   if (fSDigits) fSDigits->Clear();
1009 }
1010
1011
1012 //____________________________________________________________________
1013 TClonesArray*
1014 AliFMD::HitsArray() 
1015 {
1016   // Initialize hit array if not already, and return pointer to it. 
1017   if (!fHits) { 
1018     fHits = new TClonesArray("AliFMDHit", 1000);
1019     fNhits = 0;
1020   }
1021   return fHits;
1022 }
1023
1024 //____________________________________________________________________
1025 TClonesArray*
1026 AliFMD::DigitsArray() 
1027 {
1028   // Initialize digit array if not already, and return pointer to it. 
1029   if (!fDigits) { 
1030     fDigits = new TClonesArray("AliFMDDigit", 1000);
1031     fNdigits = 0;
1032   }
1033   return fDigits;
1034 }
1035
1036 //____________________________________________________________________
1037 TClonesArray*
1038 AliFMD::SDigitsArray() 
1039 {
1040   // Initialize digit array if not already, and return pointer to it. 
1041   if (!fSDigits) { 
1042     fSDigits = new TClonesArray("AliFMDSDigit", 1000);
1043     fNsdigits = 0;
1044   }
1045   return fSDigits;
1046 }
1047
1048 //====================================================================
1049 //
1050 // Digitization 
1051 //
1052 //____________________________________________________________________
1053 void 
1054 AliFMD::Hits2Digits() 
1055 {
1056   // Create AliFMDDigit's from AliFMDHit's.  This is done by making a
1057   // AliFMDDigitizer, and executing that code.
1058   // 
1059   Warning("Hits2Digits", "Try not to use this method.\n"
1060           "Instead, use AliSimulator");
1061   AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
1062   manager->SetInputStream(0, "galice.root");
1063   manager->SetOutputFile("H2Dfile");
1064   
1065   /* AliDigitizer* dig =*/ CreateDigitizer(manager);
1066   manager->Exec("");
1067   delete manager;
1068 }
1069
1070 //____________________________________________________________________
1071 void 
1072 AliFMD::Hits2SDigits() 
1073 {
1074   // Create AliFMDSDigit's from AliFMDHit's.  This is done by creating
1075   // an AliFMDSDigitizer object, and executing it. 
1076   // 
1077   AliFMDSDigitizer* digitizer = new AliFMDSDigitizer("galice.root");
1078   digitizer->Exec("");
1079   delete digitizer;
1080 }
1081
1082   
1083 //____________________________________________________________________
1084 AliDigitizer* 
1085 AliFMD::CreateDigitizer(AliRunDigitizer* manager) const
1086 {
1087   // Create a digitizer object 
1088   AliFMDDigitizer* digitizer = new AliFMDDigitizer(manager);
1089   return digitizer;
1090 }
1091
1092 //====================================================================
1093 //
1094 // Raw data simulation 
1095 //
1096 //__________________________________________________________________
1097 void 
1098 AliFMD::Digits2Raw() 
1099 {
1100   // Turn digits into raw data. 
1101   // 
1102   // This uses the class AliFMDRawWriter to do the job.   Please refer
1103   // to that class for more information. 
1104   AliFMDRawWriter writer(this);
1105   writer.Exec();
1106 }
1107
1108
1109 //====================================================================
1110 //
1111 // Utility 
1112 //
1113 //__________________________________________________________________
1114 void 
1115 AliFMD::Browse(TBrowser* b) 
1116 {
1117   // Browse this object. 
1118   //
1119   AliDebug(30, "\tBrowsing the FMD");
1120   AliDetector::Browse(b);
1121 #ifdef USE_PRE_MOVE
1122   if (fSimulator) b->Add(fSimulator);
1123 #endif
1124   b->Add(AliFMDGeometry::Instance());
1125 }
1126
1127 //___________________________________________________________________
1128 //
1129 // EOF
1130 //