1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //____________________________________________________________________
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
25 // This is the base class for all FMD manager classes.
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
32 // +----------+ +----------+
33 // | AliFMDv1 | | AliFMDv1 |
34 // +----------+ +----------+
36 // +----+--------------+
38 // | +------------+ 1 +---------------+
39 // | +- | AliFMDRing |<>--| AliFMDPolygon |
40 // V 2 | +------------+ +---------------+
43 // +--------+<>--+ V 1..2
44 // 3 | +-------------------+
45 // +-| AliFMDSubDetector |
46 // +-------------------+
49 // +-------------+-------------+
51 // +---------+ +---------+ +---------+
52 // | AliFMD1 | | AliFMD2 | | AliFMD3 |
53 // +---------+ +---------+ +---------+
57 // This defines the interface for the various parts of AliROOT that
58 // uses the FMD, like AliFMDDigitizer, AliFMDReconstructor, and so
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.
67 // It uses the objects of class AliFMDSubDetector to do the various
68 // stuff for FMD1, 2, and 3
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.
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.
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.
93 // The responsible person for this module is Alla Maevskaia
94 // <Alla.Maevskaia@cern.ch>.
96 // Many modifications by Christian Holm Christensen <cholm@nbi.dk>
99 // These files are not in the same directory, so there's no reason to
100 // ask the preprocessor to search in the current directory for these
101 // files by including them with `#include "..."'
102 #include <TClonesArray.h> // ROOT_TClonesArray
103 #include <TGeometry.h> // ROOT_TGeomtry
104 #include <TNode.h> // ROOT_TNode
105 #include <TTUBE.h> // ROOT_TTUBE
106 #include <TTree.h> // ROOT_TTree
107 #include <TVirtualMC.h> // ROOT_TVirtualMC
108 #include <TBrowser.h> // ROOT_TBrowser
109 #include <TMath.h> // ROOT_TMath
111 #include <AliRunDigitizer.h> // ALIRUNDIGITIZER_H
112 #include <AliLoader.h> // ALILOADER_H
113 #include <AliRun.h> // ALIRUN_H
114 #include <AliMC.h> // ALIMC_H
115 #include <AliLog.h> // ALILOG_H
116 #include <AliMagF.h> // ALIMAGF_H
117 #include "AliFMD.h" // ALIFMD_H
118 #include "AliFMDDigit.h" // ALIFMDDIGIG_H
119 #include "AliFMDHit.h" // ALIFMDHIT_H
120 #include "AliFMDDigitizer.h" // ALIFMDDIGITIZER_H
121 #include "AliFMD1.h" // ALIFMD1_H
122 #include "AliFMD2.h" // ALIFMD2_H
123 #include "AliFMD3.h" // ALIFMD3_H
124 #include "AliFMDRawWriter.h" // ALIFMDRAWWRITER_H
126 //____________________________________________________________________
129 //____________________________________________________________________
130 const Char_t* AliFMD::fgkShortLegName = "FSSL";
131 const Char_t* AliFMD::fgkLongLegName = "FSLL";
133 //____________________________________________________________________
142 fPrintboardRotationId(0),
143 fIdentityRotationId(0),
153 fAltroChannelSize(0),
157 // Default constructor for class AliFMD
159 AliDebug(0, "\tDefault CTOR");
165 //____________________________________________________________________
166 AliFMD::AliFMD(const AliFMD& other)
167 : AliDetector(other),
168 fInner(other.fInner),
169 fOuter(other.fOuter),
173 fSDigits(other.fSDigits),
174 fNsdigits(other.fNsdigits),
175 fPrintboardRotationId(other.fPrintboardRotationId),
176 fIdentityRotationId(other.fIdentityRotationId),
177 fShortLegId(other.fShortLegId),
178 fLongLegId(other.fLongLegId),
179 fLegLength(other.fLegLength),
180 fLegRadius(other.fLegRadius),
181 fModuleSpacing(other.fModuleSpacing),
182 fSiDensity(other.fSiDensity),
183 fSiThickness(other.fSiThickness),
184 fSiDeDxMip(other.fSiDeDxMip),
185 fVA1MipRange(other.fVA1MipRange),
186 fAltroChannelSize(other.fAltroChannelSize),
187 fSampleRate(other.fSampleRate)
192 //____________________________________________________________________
193 AliFMD::AliFMD(const char *name, const char *title, bool detailed)
194 : AliDetector (name, title),
202 fPrintboardRotationId(0),
203 fIdentityRotationId(0),
213 fAltroChannelSize(0),
217 // Standard constructor for Forward Multiplicity Detector
219 AliDebug(0, "\tStandard CTOR");
221 // Initialise Hit array
223 gAlice->GetMCApp()->AddHitList(fHits);
225 // (S)Digits for the detectors disk
229 // CHC: What is this?
231 SetMarkerColor(kRed);
232 SetLineColor(kYellow);
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();
242 // Specify parameters of sub-volume managers
243 fFMD1->SetInner(fInner);
246 fFMD2->SetInner(fInner);
247 fFMD2->SetOuter(fOuter);
249 fFMD3->SetInner(fInner);
250 fFMD3->SetOuter(fOuter);
259 SetAltroChannelSize();
262 fInner->SetLowR(4.3);
263 fInner->SetHighR(17.2);
264 fInner->SetWaferRadius(13.4/2);
265 fInner->SetTheta(36/2);
266 fInner->SetNStrips(512);
267 fInner->SetSiThickness(fSiThickness);
268 fInner->SetPrintboardThickness(.11);
269 fInner->SetBondingWidth(.5);
271 fOuter->SetLowR(15.6);
272 fOuter->SetHighR(28.0);
273 fOuter->SetWaferRadius(13.4/2);
274 fOuter->SetTheta(18/2);
275 fOuter->SetNStrips(256);
276 fOuter->SetSiThickness(fSiThickness);
277 fOuter->SetPrintboardThickness(.1);
278 fOuter->SetBondingWidth(.5);
281 fFMD1->SetHoneycombThickness(1);
282 fFMD1->SetInnerZ(340.0);
284 fFMD2->SetHoneycombThickness(1);
285 fFMD2->SetInnerZ(83.4);
286 fFMD2->SetOuterZ(75.2);
288 fFMD3->SetHoneycombThickness(1);
289 fFMD3->SetInnerZ(-62.8);
290 fFMD3->SetOuterZ(-75.2);
293 //____________________________________________________________________
296 // Destructor for base class AliFMD
314 //____________________________________________________________________
316 AliFMD::operator=(const AliFMD& other)
318 AliDetector::operator=(other);
319 fInner = other.fInner;
320 fOuter = other.fOuter;
324 fSDigits = other.fSDigits;
325 fNsdigits = other.fNsdigits;
326 fSiDensity = other.fSiDensity;
327 fPrintboardRotationId = other.fPrintboardRotationId;
328 fIdentityRotationId = other.fIdentityRotationId;
329 fShortLegId = other.fShortLegId;
330 fLongLegId = other.fLongLegId;
331 fLegLength = other.fLegLength;
332 fLegRadius = other.fLegRadius;
333 fModuleSpacing = other.fModuleSpacing;
334 fSiDensity = other.fSiDensity;
335 fSiThickness = other.fSiThickness;
336 fVA1MipRange = other.fVA1MipRange;
337 fAltroChannelSize = other.fAltroChannelSize;
338 fSampleRate = other.fSampleRate;
343 //====================================================================
345 // GEometry ANd Traking
347 //____________________________________________________________________
349 AliFMD::CreateGeometry()
352 // Create the geometry of Forward Multiplicity Detector. The actual
353 // construction of the geometry is delegated to the class AliFMDRing
354 // and AliFMDSubDetector and the relevant derived classes.
356 // The flow of this member function is:
358 // FOR rings fInner and fOuter DO
359 // AliFMDRing::Init();
362 // Set up hybrud card support (leg) volume shapes
364 // FOR rings fInner and fOuter DO
365 // AliFMDRing::SetupGeometry();
368 // FOR subdetectors fFMD1, fFMD2, and fFMD3 DO
369 // AliFMDSubDetector::SetupGeomtry();
372 // FOR subdetectors fFMD1, fFMD2, and fFMD3 DO
373 // AliFMDSubDetector::Geomtry();
377 // DebugGuard guard("AliFMD::CreateGeometry");
378 AliDebug(10, "\tCreating geometry");
385 Int_t airId = (*fIdtmed)[kAirId];
386 Int_t alId = (*fIdtmed)[kAlId];
387 Int_t cId = (*fIdtmed)[kCarbonId];
388 Int_t siId = (*fIdtmed)[kSiId];
389 Int_t pcbId = (*fIdtmed)[kPcbId];
390 Int_t plaId = (*fIdtmed)[kPlasticId];
391 Int_t pbId = fPrintboardRotationId;
392 Int_t idId = fIdentityRotationId;
394 par[0] = fLegRadius - .1;
396 par[2] = fLegLength / 2;
397 fShortLegId = gMC->Gsvolu(fgkShortLegName,"TUBE", plaId, par, 3);
399 par[2] += fModuleSpacing / 2;
400 fLongLegId = gMC->Gsvolu(fgkLongLegName,"TUBE", plaId, par, 3);
402 fInner->SetupGeometry(airId, siId, pcbId, pbId, idId);
403 fOuter->SetupGeometry(airId, siId, pcbId, pbId, idId);
405 fFMD1->SetupGeometry(airId, alId, cId);
406 fFMD2->SetupGeometry(airId, alId, cId);
407 fFMD3->SetupGeometry(airId, alId, cId);
409 fFMD1->Geometry("ALIC", pbId, idId);
410 fFMD2->Geometry("ALIC", pbId, idId);
411 fFMD3->Geometry("ALIC", pbId, idId);
414 //____________________________________________________________________
415 void AliFMD::CreateMaterials()
417 // Register various materials and tracking mediums with the
420 // Currently defined materials and mediums are
422 // FMD Air Normal air
423 // FMD Si Active silicon of sensors
424 // FMD Carbon Normal carbon used in support, etc.
425 // FMD Kapton Carbon used in Honeycomb
426 // FMD PCB Printed circuit board material
427 // FMD Plastic Material for support legs
429 // Also defined are two rotation matricies.
431 // DebugGuard guard("AliFMD::CreateMaterials");
432 AliDebug(10, "\tCreating materials");
436 Double_t density = 0;
437 Double_t radiationLength = 0;
438 Double_t absorbtionLength = 999;
439 Int_t fieldType = gAlice->Field()->Integ(); // Field type
440 Double_t maxField = gAlice->Field()->Max(); // Field max.
441 Double_t maxBending = 0; // Max Angle
442 Double_t maxStepSize = 0.001; // Max step size
443 Double_t maxEnergyLoss = 1; // Max Delta E
444 Double_t precision = 0.001; // Precision
445 Double_t minStepSize = 0.001; // Minimum step size
450 density = fSiDensity;
451 radiationLength = 9.36;
457 AliMaterial(id, "FMD Si$", a, z, density, radiationLength, absorbtionLength);
458 AliMedium(kSiId, "FMD Si$",id,1,fieldType,maxField,maxBending,
459 maxStepSize,maxEnergyLoss,precision,minStepSize);
466 radiationLength = 18.8;
472 AliMaterial(id, "FMD Carbon$", a, z, density, radiationLength,
474 AliMedium(kCarbonId, "FMD Carbon$",id,0,fieldType,maxField,maxBending,
475 maxStepSize,maxEnergyLoss,precision,minStepSize);
481 radiationLength = 8.9;
483 AliMaterial(id, "FMD Aluminum$", a, z, density, radiationLength,
485 AliMedium(kAlId, "FMD Aluminum$", id, 0, fieldType, maxField, maxBending,
486 maxStepSize, maxEnergyLoss, precision, minStepSize);
491 Float_t as[] = { 12.0107, 14.0067, 15.9994,
492 1.00794, 28.0855, 107.8682 };
493 Float_t zs[] = { 6., 7., 8.,
495 Float_t ws[] = { 0.039730642, 0.001396798, 0.01169634,
496 0.004367771, 0.844665, 0.09814344903 };
503 AliMixture(id, "FMD Si Chip$", as, zs, density, 6, ws);
504 AliMedium(kSiChipId, "FMD Si Chip$", id, 0, fieldType, maxField,
505 maxBending, maxStepSize, maxEnergyLoss, precision, minStepSize);
511 Float_t as[] = { 1.00794, 12.0107, 14.010, 15.9994};
512 Float_t zs[] = { 1., 6., 7., 8.};
513 Float_t ws[] = { 0.026362, 0.69113, 0.07327, 0.209235};
520 AliMixture(id, "FMD Kaption$", as, zs, density, 4, ws);
521 AliMedium(kAlId, "FMD Kaption$",id,0,fieldType,maxField,maxBending,
522 maxStepSize,maxEnergyLoss,precision,minStepSize);
528 Float_t as[] = { 12.0107, 14.0067, 15.9994, 39.948 };
529 Float_t zs[] = { 6., 7., 8., 18. };
530 Float_t ws[] = { 0.000124, 0.755267, 0.231781, 0.012827 };
537 AliMixture(id, "FMD Air$", as, zs, density, 4, ws);
538 AliMedium(kAirId, "FMD Air$", id,0,fieldType,maxField,maxBending,
539 maxStepSize,maxEnergyLoss,precision,minStepSize);
544 Float_t zs[] = { 14., 20., 13., 12.,
548 Float_t as[] = { 28.0855, 40.078, 26.981538, 24.305,
549 10.811, 47.867, 22.98977, 39.0983,
550 55.845, 18.9984, 15.9994, 12.0107,
552 Float_t ws[] = { 0.15144894, 0.08147477, 0.04128158, 0.00904554,
553 0.01397570, 0.00287685, 0.00445114, 0.00498089,
554 0.00209828, 0.00420000, 0.36043788, 0.27529426,
555 0.01415852, 0.03427566};
562 AliMixture(id, "FMD PCB$", as, zs, density, 14, ws);
563 AliMedium(kPcbId, "FMD PCB$", id,0,fieldType,maxField,maxBending,
564 maxStepSize,maxEnergyLoss,precision,minStepSize);
569 Float_t as[] = { 1.01, 12.01 };
570 Float_t zs[] = { 1., 6. };
571 Float_t ws[] = { 1., 1. };
578 AliMixture(id, "FMD Plastic$", as, zs, density, -2, ws);
579 AliMedium(kPlasticId, "FMD Plastic$", id,0,fieldType,maxField,maxBending,
580 maxStepSize,maxEnergyLoss,precision,minStepSize);
582 AliMatrix(fPrintboardRotationId, 90, 90, 0, 90, 90, 0);
583 AliMatrix(fIdentityRotationId, 90, 0, 90, 90, 0, 0);
586 //____________________________________________________________________
591 // Initialis the FMD after it has been built
595 cout << "\n" << ClassName() << ": " << flush;
596 for (i = 0; i < 35; i++) cout << "*";
597 cout << " FMD_INIT ";
598 for (i = 0; i < 35; i++) cout << "*";
599 cout << "\n" << ClassName() << ": " << flush;
601 // Here the FMD initialisation code (if any!)
602 for (i = 0; i < 80; i++) cout << "*";
609 //====================================================================
611 // Graphics and event display
613 //____________________________________________________________________
615 AliFMD::BuildGeometry()
618 // Build simple ROOT TNode geometry for event display
620 // Build a simplified geometry of the FMD used for event display
622 // The actual building of the TNodes is done by
623 // AliFMDSubDetector::SimpleGeometry.
624 AliDebug(10, "\tCreating a simplified geometry");
626 TNode* top = gAlice->GetGeometry()->GetNode("alice");
628 fFMD1->SimpleGeometry(fNodes, top, GetLineColor(), 0);
629 fFMD2->SimpleGeometry(fNodes, top, GetLineColor(), 0);
630 fFMD3->SimpleGeometry(fNodes, top, GetLineColor(), 0);
633 //____________________________________________________________________
635 AliFMD::DrawDetector()
638 // Draw a shaded view of the Forward multiplicity detector
640 // DebugGuard guard("AliFMD::DrawDetector");
641 AliDebug(10, "\tDraw detector");
643 //Set ALIC mother transparent
644 gMC->Gsatt("ALIC","SEEN",0);
645 gMC->Gsatt(fgkShortLegName,"SEEN",1);
646 gMC->Gsatt(fgkLongLegName,"SEEN",1);
648 //Set volumes visible
656 gMC->Gdopt("hide", "on");
657 gMC->Gdopt("shad", "on");
658 gMC->Gsatt("*", "fill", 7);
659 gMC->SetClipBox(".");
660 gMC->SetClipBox("*", 0, 1000, -1000, 1000, -1000, 1000);
662 gMC->Gdraw("alic", 40, 30, 0, 12, 12, .055, .055);
663 gMC->Gdhead(1111, "Forward Multiplicity Detector");
664 gMC->Gdman(16, 10, "MAN");
665 gMC->Gdopt("hide", "off");
668 //____________________________________________________________________
670 AliFMD::DistanceToPrimitive(Int_t, Int_t)
673 // Calculate the distance from the mouse to the FMD on the screen
679 //====================================================================
681 // Hit and Digit managment
683 //____________________________________________________________________
685 AliFMD::MakeBranch(Option_t * option)
687 // Create Tree branches for the FMD.
691 // H Make a branch of TClonesArray of AliFMDHit's
692 // D Make a branch of TClonesArray of AliFMDDigit's
693 // S Make a branch of TClonesArray of AliFMDSDigit's
695 const Int_t kBufferSize = 16000;
696 TString branchname(GetName());
699 if (opt.Contains("H", TString::kIgnoreCase)) {
701 AliDetector::MakeBranch(option);
703 if (opt.Contains("D", TString::kIgnoreCase)) {
705 MakeBranchInTree(fLoader->TreeD(), branchname.Data(),
706 &fDigits, kBufferSize, 0);
708 if (opt.Contains("S", TString::kIgnoreCase)) {
710 MakeBranchInTree(fLoader->TreeS(), branchname.Data(),
711 &fSDigits, kBufferSize, 0);
715 //____________________________________________________________________
717 AliFMD::SetTreeAddress()
719 // Set branch address for the Hits, Digits, and SDigits Tree.
720 if (fLoader->TreeH()) HitsArray();
721 AliDetector::SetTreeAddress();
723 TTree *treeD = fLoader->TreeD();
726 TBranch* branch = treeD->GetBranch ("FMD");
727 if (branch) branch->SetAddress(&fDigits);
730 TTree *treeS = fLoader->TreeS();
733 TBranch* branch = treeS->GetBranch ("FMD");
734 if (branch) branch->SetAddress(&fSDigits);
740 //____________________________________________________________________
742 AliFMD::SetHitsAddressBranch(TBranch *b)
744 // Set the TClonesArray to read hits into.
745 b->SetAddress(&fHits);
748 //____________________________________________________________________
750 AliFMD::AddHit(Int_t track, Int_t *vol, Float_t *hits)
752 // Add a hit to the hits tree
754 // The information of the two arrays are decoded as
758 // ivol[0] [UShort_t ] Detector #
759 // ivol[1] [Char_t ] Ring ID
760 // ivol[2] [UShort_t ] Sector #
761 // ivol[3] [UShort_t ] Strip #
762 // hits[0] [Float_t ] Track's X-coordinate at hit
763 // hits[1] [Float_t ] Track's Y-coordinate at hit
764 // hits[3] [Float_t ] Track's Z-coordinate at hit
765 // hits[4] [Float_t ] X-component of track's momentum
766 // hits[5] [Float_t ] Y-component of track's momentum
767 // hits[6] [Float_t ] Z-component of track's momentum
768 // hits[7] [Float_t ] Energy deposited by track
769 // hits[8] [Int_t ] Track's particle Id #
770 // hits[9] [Float_t ] Time when the track hit
774 UShort_t(vol[0]), // Detector #
775 Char_t(vol[1]), // Ring ID
776 UShort_t(vol[2]), // Sector #
777 UShort_t(vol[3]), // Strip #
784 hits[6], // Energy loss
785 Int_t(hits[7]), // PDG
789 //____________________________________________________________________
791 AliFMD::AddHit(Int_t track,
807 // Add a hit to the list
812 // detector Detector # (1, 2, or 3)
813 // ring Ring ID ('I' or 'O')
814 // sector Sector # (For inner/outer rings: 0-19/0-39)
815 // strip Strip # (For inner/outer rings: 0-511/0-255)
816 // x Track's X-coordinate at hit
817 // y Track's Y-coordinate at hit
818 // z Track's Z-coordinate at hit
819 // px X-component of track's momentum
820 // py Y-component of track's momentum
821 // pz Z-component of track's momentum
822 // edep Energy deposited by track
823 // pdg Track's particle Id #
824 // t Time when the track hit
826 TClonesArray& a = *(HitsArray());
827 // Search through the list of already registered hits, and see if we
828 // find a hit with the same parameters. If we do, then don't create
829 // a new hit, but rather update the energy deposited in the hit.
830 // This is done, so that a FLUKA based simulation will get the
831 // number of hits right, not just the enerrgy deposition.
832 for (Int_t i = 0; i < fNhits; i++) {
833 if (!a.At(i)) continue;
834 AliFMDHit* hit = static_cast<AliFMDHit*>(a.At(i));
835 if (hit->Detector() == detector
836 && hit->Ring() == ring
837 && hit->Sector() == sector
838 && hit->Strip() == strip
839 && hit->Track() == track) {
840 Warning("AddHit", "already had a hit in FMD%d%c[%2d,%3d] for track # %d,"
841 " adding energy (%f) to that hit (%f) -> %f",
842 detector, ring, sector, strip, track, edep, hit->Edep(),
844 hit->SetEdep(hit->Edep() + edep);
848 // If hit wasn't already registered, do so know.
849 new (a[fNhits]) AliFMDHit(fIshunt, track, detector, ring, sector, strip,
850 x, y, z, px, py, pz, edep, pdg, t);
854 //____________________________________________________________________
856 AliFMD::AddDigit(Int_t* digits)
858 // Add a digit to the Digit tree
862 // digits[0] [UShort_t] Detector #
863 // digits[1] [Char_t] Ring ID
864 // digits[2] [UShort_t] Sector #
865 // digits[3] [UShort_t] Strip #
866 // digits[4] [UShort_t] ADC Count
867 // digits[5] [Short_t] ADC Count, -1 if not used
868 // digits[6] [Short_t] ADC Count, -1 if not used
870 AddDigit(UShort_t(digits[0]), // Detector #
871 Char_t(digits[1]), // Ring ID
872 UShort_t(digits[2]), // Sector #
873 UShort_t(digits[3]), // Strip #
874 UShort_t(digits[4]), // ADC Count1
875 Short_t(digits[5]), // ADC Count2
876 Short_t(digits[6])); // ADC Count3
879 //____________________________________________________________________
881 AliFMD::AddDigit(UShort_t detector,
889 // add a real digit - as coming from data
893 // detector Detector # (1, 2, or 3)
894 // ring Ring ID ('I' or 'O')
895 // sector Sector # (For inner/outer rings: 0-19/0-39)
896 // strip Strip # (For inner/outer rings: 0-511/0-255)
897 // count1 ADC count (a 10-bit word)
898 // count2 ADC count (a 10-bit word), or -1 if not used
899 // count3 ADC count (a 10-bit word), or -1 if not used
900 TClonesArray& a = *(DigitsArray());
903 AliFMDDigit(detector, ring, sector, strip, count1, count2, count3);
906 //____________________________________________________________________
908 AliFMD::AddSDigit(Int_t* digits)
910 // Add a digit to the SDigit tree
914 // digits[0] [UShort_t] Detector #
915 // digits[1] [Char_t] Ring ID
916 // digits[2] [UShort_t] Sector #
917 // digits[3] [UShort_t] Strip #
918 // digits[4] [Float_t] Total energy deposited
919 // digits[5] [UShort_t] ADC Count
920 // digits[6] [Short_t] ADC Count, -1 if not used
921 // digits[7] [Short_t] ADC Count, -1 if not used
923 AddSDigit(UShort_t(digits[0]), // Detector #
924 Char_t(digits[1]), // Ring ID
925 UShort_t(digits[2]), // Sector #
926 UShort_t(digits[3]), // Strip #
927 Float_t(digits[4]), // Edep
928 UShort_t(digits[5]), // ADC Count1
929 Short_t(digits[6]), // ADC Count2
930 Short_t(digits[7])); // ADC Count3
933 //____________________________________________________________________
935 AliFMD::AddSDigit(UShort_t detector,
944 // add a summable digit
948 // detector Detector # (1, 2, or 3)
949 // ring Ring ID ('I' or 'O')
950 // sector Sector # (For inner/outer rings: 0-19/0-39)
951 // strip Strip # (For inner/outer rings: 0-511/0-255)
952 // edep Total energy deposited
953 // count1 ADC count (a 10-bit word)
954 // count2 ADC count (a 10-bit word), or -1 if not used
955 // count3 ADC count (a 10-bit word), or -1 if not used
957 TClonesArray& a = *(SDigitsArray());
960 AliFMDSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
963 //____________________________________________________________________
965 AliFMD::ResetSDigits()
968 // Reset number of digits and the digits array for this detector
971 if (fSDigits) fSDigits->Clear();
975 //____________________________________________________________________
979 // Initialize hit array if not already, and return pointer to it.
981 fHits = new TClonesArray("AliFMDHit", 1000);
987 //____________________________________________________________________
989 AliFMD::DigitsArray()
991 // Initialize digit array if not already, and return pointer to it.
993 fDigits = new TClonesArray("AliFMDDigit", 1000);
999 //____________________________________________________________________
1001 AliFMD::SDigitsArray()
1003 // Initialize digit array if not already, and return pointer to it.
1005 fSDigits = new TClonesArray("AliFMDSDigit", 1000);
1011 //====================================================================
1015 //____________________________________________________________________
1017 AliFMD::Hits2Digits()
1019 // Create AliFMDDigit's from AliFMDHit's. This is done by making a
1020 // AliFMDDigitizer, and executing that code.
1022 AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
1023 manager->SetInputStream(0, "galice.root");
1024 manager->SetOutputFile("H2Dfile");
1026 /* AliDigitizer* dig =*/ CreateDigitizer(manager);
1030 //____________________________________________________________________
1032 AliFMD::Hits2SDigits()
1034 // Create AliFMDSDigit's from AliFMDHit's. This is done by creating
1035 // an AliFMDSDigitizer object, and executing it.
1037 AliFMDSDigitizer* digitizer = new AliFMDSDigitizer("galice.root");
1038 digitizer->SetSampleRate(fSampleRate);
1039 digitizer->SetVA1MipRange(fVA1MipRange);
1040 digitizer->SetAltroChannelSize(fAltroChannelSize);
1041 digitizer->Exec("");
1045 //____________________________________________________________________
1047 AliFMD::CreateDigitizer(AliRunDigitizer* manager) const
1049 // Create a digitizer object
1050 AliFMDDigitizer* digitizer = new AliFMDDigitizer(manager);
1051 digitizer->SetSampleRate(fSampleRate);
1052 digitizer->SetVA1MipRange(fVA1MipRange);
1053 digitizer->SetAltroChannelSize(fAltroChannelSize);
1057 //====================================================================
1059 // Raw data simulation
1061 //__________________________________________________________________
1063 AliFMD::Digits2Raw()
1065 // Turn digits into raw data.
1067 // This uses the class AliFMDRawWriter to do the job. Please refer
1068 // to that class for more information.
1069 AliFMDRawWriter writer(this);
1070 writer.SetSampleRate(fSampleRate);
1074 //==================================================================
1076 // Various setter functions for the common paramters
1079 //__________________________________________________________________
1081 AliFMD::SetLegLength(Double_t length)
1083 // Set lenght of plastic legs that hold the hybrid (print board and
1084 // silicon sensor) onto the honeycomp support
1086 AliDebug(10, Form("\tLeg length set to %lf cm", length));
1087 fLegLength = length;
1088 fInner->SetLegLength(fLegLength);
1089 fOuter->SetLegLength(fLegLength);
1092 //__________________________________________________________________
1094 AliFMD::SetLegOffset(Double_t offset)
1096 // Set offset from edge of hybrid to plastic legs that hold the
1097 // hybrid (print board and silicon sensor) onto the honeycomp
1100 AliDebug(10, Form("\tLeg offset set to %lf cm", offset));
1101 fInner->SetLegOffset(offset);
1102 fOuter->SetLegOffset(offset);
1105 //__________________________________________________________________
1107 AliFMD::SetLegRadius(Double_t radius)
1109 // Set the diameter of the plastic legs that hold the hybrid (print
1110 // board and silicon sensor) onto the honeycomp support
1112 AliDebug(10, Form("\tLeg radius set to %lf cm", radius));
1113 fLegRadius = radius;
1114 fInner->SetLegRadius(fLegRadius);
1115 fOuter->SetLegRadius(fLegRadius);
1118 //__________________________________________________________________
1120 AliFMD::SetModuleSpacing(Double_t spacing)
1122 // Set the distance between the front and back sensor modules
1123 // (module staggering).
1125 AliDebug(10, Form("\tModule spacing set to %lf cm", spacing));
1126 fModuleSpacing = spacing;
1127 fInner->SetModuleSpacing(fModuleSpacing);
1128 fOuter->SetModuleSpacing(fModuleSpacing);
1131 //====================================================================
1135 //__________________________________________________________________
1137 AliFMD::Browse(TBrowser* b)
1139 // Browse this object.
1141 AliDebug(30, "\tBrowsing the FMD");
1142 AliDetector::Browse(b);
1143 if (fInner) b->Add(fInner, "Inner Ring");
1144 if (fOuter) b->Add(fOuter, "Outer Ring");
1145 if (fFMD1) b->Add(fFMD1, "FMD1 SubDetector");
1146 if (fFMD2) b->Add(fFMD2, "FMD2 SubDetector");
1147 if (fFMD3) b->Add(fFMD3, "FMD3 SubDetector");
1151 //___________________________________________________________________