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 // The actual code is done by various separate classes. Below is
26 // diagram showing the relationship between the various FMD classes
27 // that handles the geometry
30 // +----------+ +----------+
31 // | AliFMDv1 | | AliFMDv1 |
32 // +----------+ +----------+
34 // +----+--------------+
36 // | +------------+ 1 +---------------+
37 // | +- | AliFMDRing |<>--| AliFMDPolygon |
38 // V 2 | +------------+ +---------------+
41 // +--------+<>--+ V 1..2
42 // 3 | +-------------------+
43 // +-| AliFMDSubDetector |
44 // +-------------------+
47 // +-------------+-------------+
49 // +---------+ +---------+ +---------+
50 // | AliFMD1 | | AliFMD2 | | AliFMD3 |
51 // +---------+ +---------+ +---------+
55 // This defines the interface for the various parts of AliROOT that
56 // uses the FMD, like AliFMDDigitizer, AliFMDReconstructor, and so
60 // This is a concrete implementation of the AliFMD interface.
61 // It is the responsibility of this class to create the FMD
62 // geometry, process hits in the FMD, and serve hits and digits to
63 // the various clients.
65 // It uses the objects of class AliFMDSubDetector to do the various
66 // stuff for FMD1, 2, and 3
69 // This class contains all stuff needed to do with a ring. It's
70 // used by the AliFMDSubDetector objects to instantise inner and
71 // outer rings. The AliFMDRing objects are shared by the
72 // AliFMDSubDetector objects, and owned by the AliFMDv1 object.
75 // The code I lifted from TGeoPolygon to help with the geometry of
76 // the modules, as well as to decide wether a hit is actually with
77 // in the real module shape. The point is, that the shape of the
78 // various ring modules are really polygons (much like the lid of a
79 // coffin), but it's segmented at constant radius. That is very
80 // hard to implement using GEANT 3.21 shapes, so instead the
81 // modules are implemented as TUBS (tube sections), and in the step
82 // procedure we do the test whether the track was inside the real
83 // shape of the module.
85 // * AliFMD1, AliFMD2, and AliFMD3
86 // These are specialisation of AliFMDSubDetector, that contains the
87 // particularities of each of the sub-detector system. It is
88 // envisioned that the classes should also define the support
89 // volumes and material for each of the detectors.
93 <img src="gif/AliFMDClass.gif">
97 The responsible person for this module is
98 <a href="mailto:Alla.Maevskaia@cern.ch">Alla Maevskaia</a>.
101 Many modifications by <a href="mailto:cholm@nbi.dk">Christian
102 Holm Christensen</a>.
108 #ifndef ROOT_TClonesArray
109 #include <TClonesArray.h>
111 #ifndef ROOT_TGeomtry
112 # include <TGeometry.h>
123 #ifndef ROOT_TVirtualMC
124 # include <TVirtualMC.h>
126 #ifndef ROOT_TBrowser
127 # include <TBrowser.h>
133 #ifndef ALIRUNDIGITIZER_H
134 # include "AliRunDigitizer.h"
137 # include "AliLoader.h"
149 # include "AliMagF.h"
154 #ifndef ALIFMDDIGIG_H
155 # include "AliFMDDigit.h"
158 # include "AliFMDHit.h"
160 #ifndef ALIFMDDIGITIZER_H
161 # include "AliFMDDigitizer.h"
164 # include "AliFMD1.h"
167 # include "AliFMD2.h"
170 # include "AliFMD3.h"
172 #ifndef ALIALTROBUFFER_H
173 # include "AliAltroBuffer.h"
176 //____________________________________________________________________
179 //____________________________________________________________________
188 // Default constructor for class AliFMD
190 AliDebug(0, "Default CTOR");
198 //____________________________________________________________________
199 AliFMD::AliFMD(const char *name, const char *title, bool detailed)
200 : AliDetector (name, title),
208 // Standard constructor for Forward Multiplicity Detector
210 AliDebug(0, "Standard CTOR");
212 // Initialise Hit array
214 gAlice->GetMCApp()->AddHitList(fHits);
216 // (S)Digits for the detectors disk
220 // CHC: What is this?
222 SetMarkerColor(kRed);
223 SetLineColor(kYellow);
226 // Create sub-volume managers
227 fInner = new AliFMDRing('I', detailed);
228 fOuter = new AliFMDRing('O', detailed);
229 fFMD1 = new AliFMD1();
230 fFMD2 = new AliFMD2();
231 fFMD3 = new AliFMD3();
233 // Specify parameters of sub-volume managers
234 fFMD1->SetInner(fInner);
237 fFMD2->SetInner(fInner);
238 fFMD2->SetOuter(fOuter);
240 fFMD3->SetInner(fInner);
241 fFMD3->SetOuter(fOuter);
248 fInner->SetLowR(4.3);
249 fInner->SetHighR(17.2);
250 fInner->SetWaferRadius(13.4/2);
251 fInner->SetTheta(36/2);
252 fInner->SetNStrips(512);
253 fInner->SetSiThickness(.03);
254 fInner->SetPrintboardThickness(.11);
255 fInner->SetBondingWidth(.5);
257 fOuter->SetLowR(15.6);
258 fOuter->SetHighR(28.0);
259 fOuter->SetWaferRadius(13.4/2);
260 fOuter->SetTheta(18/2);
261 fOuter->SetNStrips( 256);
262 fOuter->SetSiThickness(.03);
263 fOuter->SetPrintboardThickness(.1);
264 fOuter->SetBondingWidth(.5);
267 fFMD1->SetHoneycombThickness(1);
268 fFMD1->SetInnerZ(340.0);
270 fFMD2->SetHoneycombThickness(1);
271 fFMD2->SetInnerZ(83.4);
272 fFMD2->SetOuterZ(75.2);
274 fFMD3->SetHoneycombThickness(1);
275 fFMD3->SetInnerZ(-62.8);
276 fFMD3->SetOuterZ(-75.2);
279 //____________________________________________________________________
282 // Destructor for base class AliFMD
300 //====================================================================
302 // GEometry ANd Traking
304 //____________________________________________________________________
306 AliFMD::CreateGeometry()
309 // Create the geometry of Forward Multiplicity Detector version 0
311 // DebugGuard guard("AliFMD::CreateGeometry");
312 AliDebug(10, "Creating geometry");
320 par[0] = fLegRadius - .1;
322 par[2] = fLegLength / 2;
324 fShortLegId = gMC->Gsvolu(name.Data(),"TUBE",(*fIdtmed)[kPlasticId],par,3);
326 par[2] += fModuleSpacing / 2;
328 fLongLegId = gMC->Gsvolu(name.Data(),"TUBE",(*fIdtmed)[kPlasticId],par,3);
330 fInner->SetupGeometry((*fIdtmed)[kAirId],
333 fPrintboardRotationId,
334 fIdentityRotationId);
335 fOuter->SetupGeometry((*fIdtmed)[kAirId],
338 fPrintboardRotationId,
339 fIdentityRotationId);
341 fFMD1->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
342 fFMD2->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
343 fFMD3->SetupGeometry((*fIdtmed)[kAirId], (*fIdtmed)[kKaptionId]);
345 fFMD1->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);
346 fFMD2->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);
347 fFMD3->Geometry("ALIC", fPrintboardRotationId, fIdentityRotationId);
350 //____________________________________________________________________
351 void AliFMD::CreateMaterials()
353 // Register various materials and tracking mediums with the
356 // Currently defined materials and mediums are
358 // FMD Air Normal air
359 // FMD Si Active silicon of sensors
360 // FMD Carbon Normal carbon used in support, etc.
361 // FMD Kapton Carbon used in Honeycomb
362 // FMD PCB Printed circuit board material
363 // FMD Plastic Material for support legs
365 // Also defined are two rotation matricies.
367 // DebugGuard guard("AliFMD::CreateMaterials");
368 AliDebug(10, "Creating materials");
372 Double_t density = 0;
373 Double_t radiationLength = 0;
374 Double_t absorbtionLength = 999;
375 Int_t fieldType = gAlice->Field()->Integ(); // Field type
376 Double_t maxField = gAlice->Field()->Max(); // Field max.
377 Double_t maxBending = 0; // Max Angle
378 Double_t maxStepSize = 0.001; // Max step size
379 Double_t maxEnergyLoss = 1; // Max Delta E
380 Double_t precision = 0.001; // Precision
381 Double_t minStepSize = 0.001; // Minimum step size
386 density = fSiDensity;
387 radiationLength = 9.36;
393 AliMaterial(id, "FMD Si$", a, z, density, radiationLength, absorbtionLength);
394 AliMedium(kSiId, "FMD Si$",id,1,fieldType,maxField,maxBending,
395 maxStepSize,maxEnergyLoss,precision,minStepSize);
402 radiationLength = 18.8;
408 AliMaterial(id, "FMD Carbon$", a, z, density, radiationLength,
410 AliMedium(kCarbonId, "FMD Carbon$",id,0,fieldType,maxField,maxBending,
411 maxStepSize,maxEnergyLoss,precision,minStepSize);
415 Float_t as[] = { 12.0107, 14.0067, 15.9994,
416 1.00794, 28.0855, 107.8682 };
417 Float_t zs[] = { 6., 7., 8.,
419 Float_t ws[] = { 0.039730642, 0.001396798, 0.01169634,
420 0.004367771, 0.844665, 0.09814344903 };
427 AliMixture(id, "FMD Si Chip$", as, zs, density, 6, ws);
428 AliMedium(kSiChipId, "FMD Si Chip$", id, 0, fieldType, maxField,
429 maxBending, maxStepSize, maxEnergyLoss, precision, minStepSize);
435 Float_t as[] = { 1.00794, 12.0107, 14.010, 15.9994};
436 Float_t zs[] = { 1., 6., 7., 8.};
437 Float_t ws[] = { 0.026362, 0.69113, 0.07327, 0.209235};
444 AliMixture(id, "FMD Kaption$", as, zs, density, 4, ws);
445 AliMedium(kKaptionId, "FMD Kaption$",id,0,fieldType,maxField,maxBending,
446 maxStepSize,maxEnergyLoss,precision,minStepSize);
451 Float_t as[] = { 12.0107, 14.0067, 15.9994, 39.948 };
452 Float_t zs[] = { 6., 7., 8., 18. };
453 Float_t ws[] = { 0.000124, 0.755267, 0.231781, 0.012827 };
460 AliMixture(id, "FMD Air$", as, zs, density, 4, ws);
461 AliMedium(kAirId, "FMD Air$", id,0,fieldType,maxField,maxBending,
462 maxStepSize,maxEnergyLoss,precision,minStepSize);
467 Float_t zs[] = { 14., 20., 13., 12.,
471 Float_t as[] = { 28.0855, 40.078, 26.981538, 24.305,
472 10.811, 47.867, 22.98977, 39.0983,
473 55.845, 18.9984, 15.9994, 12.0107,
475 Float_t ws[] = { 0.15144894, 0.08147477, 0.04128158, 0.00904554,
476 0.01397570, 0.00287685, 0.00445114, 0.00498089,
477 0.00209828, 0.00420000, 0.36043788, 0.27529426,
478 0.01415852, 0.03427566};
485 AliMixture(id, "FMD PCB$", as, zs, density, 14, ws);
486 AliMedium(kPcbId, "FMD PCB$", id,1,fieldType,maxField,maxBending,
487 maxStepSize,maxEnergyLoss,precision,minStepSize);
492 Float_t as[] = { 1.01, 12.01 };
493 Float_t zs[] = { 1., 6. };
494 Float_t ws[] = { 1., 1. };
501 AliMixture(id, "FMD Plastic$", as, zs, density, -2, ws);
502 AliMedium(kPlasticId, "FMD Plastic$", id,0,fieldType,maxField,maxBending,
503 maxStepSize,maxEnergyLoss,precision,minStepSize);
505 AliMatrix(fPrintboardRotationId, 90, 90, 0, 90, 90, 0);
506 AliMatrix(fIdentityRotationId, 90, 0, 90, 90, 0, 0);
509 //____________________________________________________________________
514 // Initialis the FMD after it has been built
518 std::cout << "\n" << ClassName() << ": " << std::flush;
519 for (i = 0; i < 35; i++) std::cout << "*";
520 std::cout << " FMD_INIT ";
521 for (i = 0; i < 35; i++) std::cout << "*";
522 std::cout << "\n" << ClassName() << ": " << std::flush;
524 // Here the FMD initialisation code (if any!)
525 for (i = 0; i < 80; i++) std::cout << "*";
526 std::cout << std::endl;
532 //====================================================================
534 // Graphics and event display
536 //____________________________________________________________________
538 AliFMD::BuildGeometry()
541 // Build simple ROOT TNode geometry for event display
543 // Build a simplified geometry of the FMD used for event display
545 AliDebug(10, "Creating a simplified geometry");
547 TNode* top = gAlice->GetGeometry()->GetNode("alice");
549 fFMD1->SimpleGeometry(fNodes, top, GetLineColor(), 0);
550 fFMD2->SimpleGeometry(fNodes, top, GetLineColor(), 0);
551 fFMD3->SimpleGeometry(fNodes, top, GetLineColor(), 0);
554 //____________________________________________________________________
556 AliFMD::DrawDetector()
559 // Draw a shaded view of the Forward multiplicity detector version 0
561 // DebugGuard guard("AliFMD::DrawDetector");
562 AliDebug(10, "Draw detector");
564 //Set ALIC mother transparent
565 gMC->Gsatt("ALIC","SEEN",0);
567 //Set volumes visible
575 gMC->Gdopt("hide", "on");
576 gMC->Gdopt("shad", "on");
577 gMC->Gsatt("*", "fill", 7);
578 gMC->SetClipBox(".");
579 gMC->SetClipBox("*", 0, 1000, -1000, 1000, -1000, 1000);
581 gMC->Gdraw("alic", 40, 30, 0, 12, 12, .055, .055);
582 gMC->Gdhead(1111, "Forward Multiplicity Detector");
583 gMC->Gdman(16, 10, "MAN");
584 gMC->Gdopt("hide", "off");
587 //____________________________________________________________________
589 AliFMD::DistanceToPrimitive(Int_t, Int_t)
592 // Calculate the distance from the mouse to the FMD on the screen
598 //====================================================================
600 // Hit and Digit managment
602 //____________________________________________________________________
604 AliFMD::MakeBranch(Option_t * option)
606 // Create Tree branches for the FMD.
607 const Int_t kBufferSize = 16000;
608 TString branchname(GetName());
611 if (opt.Contains("H", TString::kIgnoreCase)) {
613 AliDetector::MakeBranch(option);
615 if (opt.Contains("D", TString::kIgnoreCase)) {
617 MakeBranchInTree(fLoader->TreeD(), branchname.Data(),
618 &fDigits, kBufferSize, 0);
620 if (opt.Contains("S", TString::kIgnoreCase)) {
622 MakeBranchInTree(fLoader->TreeS(), branchname.Data(),
623 &fSDigits, kBufferSize, 0);
627 //____________________________________________________________________
629 AliFMD::SetTreeAddress()
631 // Set branch address for the Hits and Digits Tree.
633 if (fLoader->TreeH()) HitsArray();
634 AliDetector::SetTreeAddress();
636 TTree *treeD = fLoader->TreeD();
639 TBranch* branch = treeD->GetBranch ("FMD");
640 if (branch) branch->SetAddress(&fDigits);
643 TTree *treeS = fLoader->TreeS();
646 TBranch* branch = treeS->GetBranch ("FMD");
647 if (branch) branch->SetAddress(&fSDigits);
653 //____________________________________________________________________
655 AliFMD::SetHitsAddressBranch(TBranch *b)
657 b->SetAddress(&fHits);
660 //____________________________________________________________________
662 AliFMD::AddHit(Int_t track, Int_t *vol, Float_t *hits)
664 // Add a hit to the hits tree
666 // The information of the two arrays are decoded as
670 // ivol[0] [UShort_t ] Detector #
671 // ivol[1] [Char_t ] Ring ID
672 // ivol[2] [UShort_t ] Sector #
673 // ivol[3] [UShort_t ] Strip #
674 // hits[0] [Float_t ] Track's X-coordinate at hit
675 // hits[1] [Float_t ] Track's Y-coordinate at hit
676 // hits[3] [Float_t ] Track's Z-coordinate at hit
677 // hits[4] [Float_t ] X-component of track's momentum
678 // hits[5] [Float_t ] Y-component of track's momentum
679 // hits[6] [Float_t ] Z-component of track's momentum
680 // hits[7] [Float_t ] Energy deposited by track
681 // hits[8] [Int_t ] Track's particle Id #
682 // hits[9] [Float_t ] Time when the track hit
684 UShort_t(vol[0]), // Detector #
685 Char_t(vol[1]), // Ring ID
686 UShort_t(vol[2]), // Sector #
687 UShort_t(vol[3]), // Strip #
694 hits[6], // Energy loss
695 Int_t(hits[7]), // PDG
699 //____________________________________________________________________
701 AliFMD::AddHit(Int_t track,
717 // Add a hit to the list
722 // detector Detector # (1, 2, or 3)
723 // ring Ring ID ('I' or 'O')
724 // sector Sector # (For inner/outer rings: 0-19/0-39)
725 // strip Strip # (For inner/outer rings: 0-511/0-255)
726 // x Track's X-coordinate at hit
727 // y Track's Y-coordinate at hit
728 // z Track's Z-coordinate at hit
729 // px X-component of track's momentum
730 // py Y-component of track's momentum
731 // pz Z-component of track's momentum
732 // edep Energy deposited by track
733 // pdg Track's particle Id #
734 // t Time when the track hit
736 TClonesArray& a = *(HitsArray());
737 // Search through the list of already registered hits, and see if we
738 // find a hit with the same parameters. If we do, then don't create
739 // a new hit, but rather update the energy deposited in the hit.
740 // This is done, so that a FLUKA based simulation will get the
741 // number of hits right, not just the enerrgy deposition.
742 for (Int_t i = 0; i < fNhits; i++) {
743 if (!a.At(i)) continue;
744 AliFMDHit* hit = static_cast<AliFMDHit*>(a.At(i));
745 if (hit->Detector() == detector
746 && hit->Ring() == ring
747 && hit->Sector() == sector
748 && hit->Strip() == strip
749 && hit->Track() == track) {
750 hit->SetEdep(hit->Edep() + edep);
754 // If hit wasn't already registered, do so know.
755 new (a[fNhits]) AliFMDHit(fIshunt, track, detector, ring, sector, strip,
756 x, y, z, px, py, pz, edep, pdg, t);
760 //____________________________________________________________________
762 AliFMD::AddDigit(Int_t* digits)
764 // Add a digit to the Digit tree
768 // digits[0] [UShort_t] Detector #
769 // digits[1] [Char_t] Ring ID
770 // digits[2] [UShort_t] Sector #
771 // digits[3] [UShort_t] Strip #
772 // digits[4] [UShort_t] ADC Count
773 // digits[5] [Short_t] ADC Count, -1 if not used
774 // digits[6] [Short_t] ADC Count, -1 if not used
776 AddDigit(UShort_t(digits[0]), // Detector #
777 Char_t(digits[1]), // Ring ID
778 UShort_t(digits[2]), // Sector #
779 UShort_t(digits[3]), // Strip #
780 UShort_t(digits[4]), // ADC Count1
781 Short_t(digits[5]), // ADC Count2
782 Short_t(digits[6])); // ADC Count3
785 //____________________________________________________________________
787 AliFMD::AddDigit(UShort_t detector,
795 // add a real digit - as coming from data
799 // detector Detector # (1, 2, or 3)
800 // ring Ring ID ('I' or 'O')
801 // sector Sector # (For inner/outer rings: 0-19/0-39)
802 // strip Strip # (For inner/outer rings: 0-511/0-255)
803 // count1 ADC count (a 10-bit word)
804 // count2 ADC count (a 10-bit word), or -1 if not used
805 // count3 ADC count (a 10-bit word), or -1 if not used
806 TClonesArray& a = *(DigitsArray());
809 AliFMDDigit(detector, ring, sector, strip, count1, count2, count3);
812 //____________________________________________________________________
814 AliFMD::AddSDigit(Int_t* digits)
816 // Add a digit to the SDigit tree
820 // digits[0] [UShort_t] Detector #
821 // digits[1] [Char_t] Ring ID
822 // digits[2] [UShort_t] Sector #
823 // digits[3] [UShort_t] Strip #
824 // digits[4] [Float_t] Total energy deposited
825 // digits[5] [UShort_t] ADC Count
826 // digits[6] [Short_t] ADC Count, -1 if not used
827 // digits[7] [Short_t] ADC Count, -1 if not used
829 AddSDigit(UShort_t(digits[0]), // Detector #
830 Char_t(digits[1]), // Ring ID
831 UShort_t(digits[2]), // Sector #
832 UShort_t(digits[3]), // Strip #
833 Float_t(digits[4]), // Edep
834 UShort_t(digits[5]), // ADC Count1
835 Short_t(digits[6]), // ADC Count2
836 Short_t(digits[7])); // ADC Count3
839 //____________________________________________________________________
841 AliFMD::AddSDigit(UShort_t detector,
850 // add a summable digit
854 // detector Detector # (1, 2, or 3)
855 // ring Ring ID ('I' or 'O')
856 // sector Sector # (For inner/outer rings: 0-19/0-39)
857 // strip Strip # (For inner/outer rings: 0-511/0-255)
858 // edep Total energy deposited
859 // count1 ADC count (a 10-bit word)
860 // count2 ADC count (a 10-bit word), or -1 if not used
861 // count3 ADC count (a 10-bit word), or -1 if not used
862 TClonesArray& a = *(SDigitsArray());
865 AliFMDSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
868 //____________________________________________________________________
870 AliFMD::ResetSDigits()
873 // Reset number of digits and the digits array for this detector
876 if (fSDigits) fSDigits->Clear();
880 //____________________________________________________________________
884 // Initialize hit array if not already, and return pointer to it.
886 fHits = new TClonesArray("AliFMDHit", 1000);
892 //____________________________________________________________________
894 AliFMD::DigitsArray()
896 // Initialize digit array if not already, and return pointer to it.
898 fDigits = new TClonesArray("AliFMDDigit", 1000);
904 //____________________________________________________________________
906 AliFMD::SDigitsArray()
908 // Initialize digit array if not already, and return pointer to it.
910 fSDigits = new TClonesArray("AliFMDSDigit", 1000);
916 //====================================================================
920 //____________________________________________________________________
922 AliFMD::Hits2Digits()
924 AliRunDigitizer* manager = new AliRunDigitizer(1, 1);
925 manager->SetInputStream(0, "galice.root");
926 manager->SetOutputFile("H2Dfile");
928 /* AliDigitizer* dig =*/ CreateDigitizer(manager);
932 //____________________________________________________________________
934 AliFMD::Hits2SDigits()
936 AliDigitizer* sdig = new AliFMDSDigitizer("galice.root");
941 //____________________________________________________________________
943 AliFMD::CreateDigitizer(AliRunDigitizer* manager) const
945 // Create a digitizer object
946 return new AliFMDDigitizer(manager);
949 //====================================================================
951 // Raw data simulation
953 //__________________________________________________________________
957 AliFMD* fmd = static_cast<AliFMD*>(gAlice->GetDetector(GetName()));
958 fLoader->LoadDigits();
959 TTree* digitTree = fLoader->TreeD();
961 Error("Digits2Raw", "no digit tree");
965 TClonesArray* digits = new TClonesArray("AliFMDDigit", 1000);
966 fmd->SetTreeAddress();
967 TBranch* digitBranch = digitTree->GetBranch(GetName());
969 Error("Digits2Raw", "no branch for %s", GetName());
972 digitBranch->SetAddress(&digits);
974 Int_t nEvents = Int_t(digitTree->GetEntries());
975 for (Int_t event = 0; event < nEvents; event++) {
977 digitTree->GetEvent(event);
979 Int_t nDigits = digits->GetEntries();
980 if (nDigits < 1) continue;
983 UShort_t prevDetector = 0;
984 Char_t prevRing = '\0';
985 UShort_t prevSector = 0;
986 // UShort_t prevStrip = 0;
988 // The first seen strip number for a channel
989 UShort_t startStrip = 0;
991 // Which channel number in the ALTRO channel we're at
994 // How many times the ALTRO Samples one VA1_ALICE channel
995 Int_t sampleRate = 1;
997 // A buffer to hold 1 ALTRO channel - Normally, one ALTRO channel
998 // holds 128 VA1_ALICE channels, sampled at a rate of `sampleRate'
999 TArrayI channel(128 * sampleRate);
1002 AliAltroBuffer* altro = 0;
1004 // Loop over the digits in the event. Note, that we assume the
1005 // the digits are in order in the branch. If they were not, we'd
1006 // have to cache all channels before we could write the data to
1007 // the ALTRO buffer, or we'd have to set up a map of the digits.
1008 for (Int_t i = 0; i < nDigits; i++) {
1010 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
1012 UShort_t det = digit->Detector();
1013 Char_t ring = digit->Ring();
1014 UShort_t sector = digit->Sector();
1015 UShort_t strip = digit->Strip();
1016 if (det != prevDetector) {
1017 AliDebug(10, Form("FMD: New DDL, was %d, now %d",
1018 kBaseDDL + prevDetector - 1,
1019 kBaseDDL + det - 1));
1020 // If an altro exists, delete the object, flushing the data to
1021 // disk, and closing the file.
1023 // When the first argument is false, we write the real
1025 AliDebug(10, Form("New altro: Write channel at %d Strip: %d "
1026 "Sector: %d Ring: %d",
1027 i, startStrip, prevSector, prevRing));
1028 // TPC to FMD translations
1031 // ----------+-----------
1036 altro->WriteChannel(Int_t(startStrip),
1038 Int_t((prevRing == 'I' ? 0 : 1)),
1039 channel.fN, channel.fArray, 0);
1041 altro->WriteDataHeader(kFALSE, kFALSE);
1047 // Need to open a new DDL!
1048 Int_t ddlId = kBaseDDL + det - 1;
1049 TString filename(Form("%s_%d.ddl", GetName(), ddlId));
1051 AliDebug(10, Form("New altro buffer with DDL file %s",
1053 AliDebug(10, Form("New altro at %d", i));
1054 // Create a new altro buffer - a `1' as the second argument
1055 // means `write mode'
1056 altro = new AliAltroBuffer(filename.Data(), 1);
1058 // Write a dummy (first argument is true) header to the DDL
1059 // file - later on, when we close the file, we write the real
1061 altro->WriteDataHeader(kTRUE, kFALSE);
1063 // Figure out the sample rate
1064 if (digit->Count2() > 0) sampleRate = 2;
1065 if (digit->Count3() > 0) sampleRate = 3;
1067 channel.Set(128 * sampleRate);
1070 prevSector = sector;
1073 else if (offset == 128
1074 || digit->Ring() != prevRing
1075 || digit->Sector() != prevSector) {
1076 // Force a new Altro channel
1077 AliDebug(10, Form("Flushing channel to disk because %s",
1078 (offset == 128 ? "channel is full" :
1079 (ring != prevRing ? "new ring up" :
1080 "new sector up"))));
1081 AliDebug(10, Form("New Channel: Write channel at %d Strip: %d "
1082 "Sector: %d Ring: %d",
1083 i, startStrip, prevSector, prevRing));
1084 altro->WriteChannel(Int_t(startStrip),
1086 Int_t((prevRing == 'I' ? 0 : 1)),
1087 channel.fN, channel.fArray, 0);
1088 // Reset and update channel variables
1093 prevSector = sector;
1096 // Store the counts of the ADC in the channel buffer
1097 channel[offset * sampleRate] = digit->Count1();
1099 channel[offset * sampleRate + 1] = digit->Count2();
1101 channel[offset * sampleRate + 2] = digit->Count3();
1104 // Finally, we need to close the final ALTRO buffer if it wasn't
1108 altro->WriteDataHeader(kFALSE, kFALSE);
1112 fLoader->UnloadDigits();
1115 //==================================================================
1117 // Various setter functions for the common paramters
1120 //__________________________________________________________________
1122 AliFMD::SetLegLength(Double_t length)
1124 // Set lenght of plastic legs that hold the hybrid (print board and
1125 // silicon sensor) onto the honeycomp support
1127 // DebugGuard guard("AliFMD::SetLegLength");
1128 AliDebug(10, "AliFMD::SetLegLength");
1129 fLegLength = length;
1130 fInner->SetLegLength(fLegLength);
1131 fOuter->SetLegLength(fLegLength);
1134 //__________________________________________________________________
1136 AliFMD::SetLegOffset(Double_t offset)
1138 // Set offset from edge of hybrid to plastic legs that hold the
1139 // hybrid (print board and silicon sensor) onto the honeycomp
1142 // DebugGuard guard("AliFMD::SetLegOffset");
1143 AliDebug(10, "AliFMD::SetLegOffset");
1144 fInner->SetLegOffset(offset);
1145 fOuter->SetLegOffset(offset);
1148 //__________________________________________________________________
1150 AliFMD::SetLegRadius(Double_t radius)
1152 // Set the diameter of the plastic legs that hold the hybrid (print
1153 // board and silicon sensor) onto the honeycomp support
1155 // DebugGuard guard("AliFMD::SetLegRadius");
1156 AliDebug(10, "AliFMD::SetLegRadius");
1157 fLegRadius = radius;
1158 fInner->SetLegRadius(fLegRadius);
1159 fOuter->SetLegRadius(fLegRadius);
1162 //__________________________________________________________________
1164 AliFMD::SetModuleSpacing(Double_t spacing)
1166 // Set the distance between the front and back sensor modules
1167 // (module staggering).
1169 // DebugGuard guard("AliFMD::SetModuleSpacing");
1170 AliDebug(10, "AliFMD::SetModuleSpacing");
1171 fModuleSpacing = spacing;
1172 fInner->SetModuleSpacing(fModuleSpacing);
1173 fOuter->SetModuleSpacing(fModuleSpacing);
1176 //====================================================================
1180 //__________________________________________________________________
1182 AliFMD::Browse(TBrowser* b)
1184 AliDebug(10, "AliFMD::Browse");
1185 AliDetector::Browse(b);
1186 if (fInner) b->Add(fInner, "Inner Ring");
1187 if (fOuter) b->Add(fOuter, "Outer Ring");
1188 if (fFMD1) b->Add(fFMD1, "FMD1 SubDetector");
1189 if (fFMD2) b->Add(fFMD2, "FMD2 SubDetector");
1190 if (fFMD3) b->Add(fFMD3, "FMD3 SubDetector");
1194 //********************************************************************
1198 //__________________________________________________________________
1202 //********************************************************************
1206 //__________________________________________________________________
1211 //_//____________________________________________________________________
1213 AliFMDv1::StepManager()
1216 // Called for every step in the Forward Multiplicity Detector
1218 // The procedure is as follows:
1220 // - IF NOT track is alive THEN RETURN ENDIF
1221 // - IF NOT particle is charged THEN RETURN ENDIF
1222 // - IF NOT volume name is "STRI" or "STRO" THEN RETURN ENDIF
1223 // - Get strip number (volume copy # minus 1)
1224 // - Get phi division number (mother volume copy #)
1225 // - Get module number (grand-mother volume copy #)
1226 // - section # = 2 * module # + phi division # - 1
1227 // - Get ring Id from volume name
1228 // - Get detector # from grand-grand-grand-mother volume name
1229 // - Get pointer to sub-detector object.
1230 // - Get track position
1231 // - IF track is entering volume AND track is inside real shape THEN
1232 // - Reset energy deposited
1233 // - Get track momentum
1234 // - Get particle ID #
1236 // - IF track is inside volume AND inside real shape THEN
1237 /// - Update energy deposited
1239 // - IF track is inside real shape AND (track is leaving volume,
1240 // or it died, or it is stopped THEN
1245 // DebugGuard guard("AliFMDv1::StepManager");
1246 AliDebug(10, "AliFMDv1::StepManager");
1249 // If the track is gone, return
1250 if (!gMC->IsTrackAlive()) return;
1252 // Only process charged particles
1253 if(TMath::Abs(gMC->TrackCharge()) <= 0) return;
1255 // Only do stuff is the track is in one of the strips.
1256 TString vol(gMC->CurrentVolName());
1257 if (!vol.Contains("STR")) return;
1260 // Get the strip number. Note, that GEANT numbers divisions from 1,
1261 // so we subtract one
1263 gMC->CurrentVolID(strip);
1266 // Get the phi division of the module
1267 Int_t phiDiv; // * The phi division number (1 or 2)
1268 gMC->CurrentVolOffID(1, phiDiv); // in the module
1270 // Active volume number - not used.
1272 // gMC->CurrentVolOffID(2, active);
1274 // Get the module number in the ring.
1276 gMC->CurrentVolOffID(3, module);
1278 // Ring copy number - the same as the detector number - not used
1279 // Int_t ringCopy; // * Ring copy number
1280 // gMC->CurrentVolOffID(4, ringCopy); // Same as detector number
1282 // Get the detector number from the path name
1283 Int_t detector = Int_t((gMC->CurrentVolOffName(5)[3]) - 48);
1285 // The sector number, calculated from module and phi division #
1286 Int_t sector = 2 * module + phiDiv - 1;
1288 // The ring ID is encoded in the volume name
1289 Char_t ring = vol[3];
1291 // Get a pointer to the sub detector structure
1292 AliFMDSubDetector* det = 0;
1294 case 1: det = fFMD1; break;
1295 case 2: det = fFMD2; break;
1296 case 3: det = fFMD3; break;
1300 // Get the current track position
1302 gMC->TrackPosition(v);
1303 // Check that the track is actually within the active area
1304 Bool_t isWithin = det->CheckHit(ring, module, v.X(), v.Y());
1305 Bool_t entering = gMC->IsTrackEntering() && isWithin;
1306 Bool_t inside = gMC->IsTrackInside() && isWithin;
1307 Bool_t out = (gMC->IsTrackExiting()
1308 || gMC->IsTrackDisappeared()
1309 || gMC->IsTrackStop()
1311 // Reset the energy deposition for this track, and update some of
1316 // Get production vertex and momentum of the track
1318 gMC->TrackMomentum(fCurrentP);
1319 fCurrentPdg = gMC->IdFromPDG(gMC->TrackPid());
1322 // fAnalyser->Update(detector, ring, isWithin, v.X(), v.Y());
1325 // If the track is inside, then update the energy deposition
1326 if (inside && fCurrentDeltaE >= 0)
1327 fCurrentDeltaE += 1000 * gMC->Edep();
1329 // The track exits the volume, or it disappeared in the volume, or
1330 // the track is stopped because it no longer fulfills the cuts
1331 // defined, then we create a hit.
1332 if (out && fCurrentDeltaE >= 0) {
1333 fCurrentDeltaE += 1000 * gMC->Edep();
1335 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber(),
1336 detector, ring, sector, strip,
1337 fCurrentV.X(), fCurrentV.Y(), fCurrentV.Z(),
1338 fCurrentP.X(), fCurrentP.Y(), fCurrentP.Z(),
1339 fCurrentDeltaE, fCurrentPdg, fCurrentV.T());
1340 fCurrentDeltaE = -1;
1343 //___________________________________________________________________