]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDDigitizer.cxx
Fixed a few problems
[u/mrichter/AliRoot.git] / FMD / AliFMDDigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, 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 //  This class contains the procedures simulation ADC  signal for the
21 //  Forward Multiplicity detector  : Hits->Digits and Hits->SDigits
22 // 
23 //  Digits consists of
24 //   - Detector #
25 //   - Ring ID                                             
26 //   - Sector #     
27 //   - Strip #
28 //   - ADC count in this channel                                  
29 //
30 //  Digits consists of
31 //   - Detector #
32 //   - Ring ID                                             
33 //   - Sector #     
34 //   - Strip #
35 //   - Total energy deposited in the strip
36 //   - ADC count in this channel                                  
37 //
38 // As the Digits and SDigits have so much in common, the classes
39 // AliFMDDigitizer and AliFMDSDigitizer are implemented via a base
40 // class AliFMDBaseDigitizer.
41 //
42 //                 +---------------------+
43 //                 | AliFMDBaseDigitizer |
44 //                 +---------------------+
45 //                           ^
46 //                           |
47 //                +----------+---------+
48 //                |                    |
49 //      +-----------------+     +------------------+
50 //      | AliFMDDigitizer |     | AliFMDSDigitizer |
51 //      +-----------------+     +------------------+
52 //
53 // These classes has several paramters: 
54 //
55 //     fPedestal
56 //     fPedestalWidth
57 //         (Only AliFMDDigitizer)
58 //         Mean and width of the pedestal.  The pedestal is simulated
59 //         by a Guassian, but derived classes my override MakePedestal
60 //         to simulate it differently (or pick it up from a database).
61 //
62 //     fVA1MipRange
63 //         The dymamic MIP range of the VA1_ALICE pre-amplifier chip 
64 //
65 //     fAltroChannelSize
66 //         The largest number plus one that can be stored in one
67 //         channel in one time step in the ALTRO ADC chip. 
68 //
69 //     fSampleRate
70 //         How many times the ALTRO ADC chip samples the VA1_ALICE
71 //         pre-amplifier signal.   The VA1_ALICE chip is read-out at
72 //         10MHz, while it's possible to drive the ALTRO chip at
73 //         25MHz.  That means, that the ALTRO chip can have time to
74 //         sample each VA1_ALICE signal up to 2 times.  Although it's
75 //         not certain this feature will be used in the production,
76 //         we'd like have the option, and so it should be reflected in
77 //         the code.
78 //
79 // The shaping function of the VA1_ALICE is given by 
80 //
81 //      f(x) = A(1 - exp(-Bx))
82 //
83 // Where A is a normalization constant, tuned so that the integral 
84 //
85 //      /1
86 //      |           A(-1 + B + exp(-B))
87 //      | f(x) dx = ------------------- = 1
88 //      |                    B
89 //      / 0
90 //
91 // and B is the a parameter defined by the shaping time (fShapingTime).  
92 //
93 // Solving the above equation, for A gives
94 //
95 //                 B
96 //      A = ----------------
97 //          -1 + B + exp(-B)
98 //
99 // So, if we define the function g: [0,1] -> [0:1] by 
100 //
101 //               / v
102 //               |              Bu + exp(-Bu) - Bv - exp(-Bv) 
103 //      g(u,v) = | f(x) dx = -A -----------------------------
104 //               |                            B
105 //               / u
106 //
107 // we can evaluate the ALTRO sample of the VA1_ALICE pre-amp between
108 // any two times (u, v), by 
109 //       
110 //
111 //                                B         Bu + exp(-Bu) - Bv - exp(-Bv)
112 //      C = Q g(u,v) = - Q ---------------- -----------------------------
113 //                         -1 + B + exp(-B)              B                  
114 //
115 //               Bu + exp(-Bu) - Bv - exp(-Bv) 
116 //        = -  Q -----------------------------
117 //                    -1 + B + exp(-B)
118 //
119 // Where Q is the total charge collected by the VA1_ALICE
120 // pre-amplifier.   Q is then given by 
121 //
122 //           E S 
123 //      Q =  - -
124 //           e R
125 //
126 // where E is the total energy deposited in a silicon strip, R is the
127 // dynamic range of the VA1_ALICE pre-amp (fVA1MipRange), e is the
128 // energy deposited by a single MIP, and S ALTRO channel size in each
129 // time step (fAltroChannelSize).  
130 //
131 // The energy deposited per MIP is given by 
132 //
133 //      e = M * rho * w 
134 //
135 // where M is the universal number 1.664, rho is the density of
136 // silicon, and w is the depth of the silicon sensor. 
137 //
138 // The final ADC count is given by 
139 //
140 //      C' = C + P
141 //
142 // where P is the (randomized) pedestal (see MakePedestal)
143 //
144 // This class uses the class template AliFMDMap<Type> to make an
145 // internal cache of the energy deposted of the hits.  The class
146 // template is instantasized as 
147 //
148 //  typedef AliFMDMap<std::pair<Float_t, UShort_t> > AliFMDEdepMap;
149 //
150 // The first member of the values is the summed energy deposition in a
151 // given strip, while the second member of the values is the number of
152 // hits in a given strip.  Using the second member, it's possible to
153 // do some checks on just how many times a strip got hit, and what
154 // kind of error we get in our reconstructed hits.  Note, that this
155 // information is currently not written to the digits tree.  I think a
156 // QA (Quality Assurance) digit tree is better suited for that task.
157 // However, the information is there to be used in the future. 
158 //
159 //
160 // Latest changes by Christian Holm Christensen
161 //
162 //////////////////////////////////////////////////////////////////////////////
163
164 #ifndef ROOT_TTree
165 # include <TTree.h>
166 #endif
167 #ifndef ROOT_TRandom
168 # include <TRandom.h>
169 #endif
170 #ifndef ALILOG_H
171 # include "AliLog.h"
172 #endif
173 #ifndef ALIFMDDIGITIZER_H
174 # include "AliFMDDigitizer.h"
175 #endif
176 #ifndef ALIFMD_H
177 # include "AliFMD.h"
178 #endif
179 #ifndef ALIFMDHIT_H
180 # include "AliFMDHit.h"
181 #endif
182 #ifndef ALIFMDDIGIT_H
183 # include "AliFMDDigit.h"
184 #endif
185 #ifndef ALIFMDDIGIT_H
186 # include "AliFMDSDigit.h"
187 #endif
188 #ifndef ALIRUNDIGITIZER_H
189 # include "AliRunDigitizer.h"
190 #endif
191 #ifndef ALIRUN_H
192 # include "AliRun.h"
193 #endif
194 #ifndef ALILOADER_H
195 # include "AliLoader.h"
196 #endif
197 #ifndef ALIRUNLOADER_H
198 # include "AliRunLoader.h"
199 #endif
200     
201 //____________________________________________________________________
202 ClassImp(AliFMDEdepMap);
203
204 //====================================================================
205 ClassImp(AliFMDBaseDigitizer);
206
207 //____________________________________________________________________
208 AliFMDBaseDigitizer::AliFMDBaseDigitizer()  
209   : fRunLoader(0)
210 {
211   // Default ctor - don't use it
212 }
213
214 //____________________________________________________________________
215 AliFMDBaseDigitizer::AliFMDBaseDigitizer(AliRunDigitizer* manager) 
216   : AliDigitizer(manager, "AliFMDBaseDigitizer", "FMD Digitizer base class"), 
217     fRunLoader(0),
218     fEdep(kMaxDetectors, kMaxRings, kMaxSectors, kMaxStrips)
219 {
220   // Normal CTOR
221   AliDebug(1," processed");
222   SetVA1MipRange();
223   SetAltroChannelSize();
224   SetSampleRate();
225   SetShapingTime();
226 }
227
228 //____________________________________________________________________
229 AliFMDBaseDigitizer::AliFMDBaseDigitizer(const Char_t* name, 
230                                          const Char_t* title) 
231   : AliDigitizer(name, title),
232     fRunLoader(0),
233     fEdep(kMaxDetectors, kMaxRings, kMaxSectors, kMaxStrips)
234 {
235   // Normal CTOR
236   AliDebug(1," processed");
237   SetVA1MipRange();
238   SetAltroChannelSize();
239   SetSampleRate();
240   SetShapingTime();
241 }
242
243 //____________________________________________________________________
244 AliFMDBaseDigitizer::~AliFMDBaseDigitizer()
245 {
246   // Destructor
247 }
248
249 //____________________________________________________________________
250 Bool_t 
251 AliFMDBaseDigitizer::Init()
252 {
253   // Initialization
254   return kTRUE;
255 }
256  
257
258 //____________________________________________________________________
259 void
260 AliFMDBaseDigitizer::SumContributions(AliFMD* fmd) 
261 {
262   // Sum energy deposited contributions from each hit in a cache
263   // (fEdep).  
264   if (!fRunLoader) 
265     Fatal("SumContributions", "no run loader");
266   
267   // Clear array of deposited energies 
268   fEdep.Clear();
269   
270   // Get the FMD loader 
271   AliLoader* inFMD = fRunLoader->GetLoader("FMDLoader");
272   // And load the hits 
273   inFMD->LoadHits("READ");
274   
275   // Get the tree of hits 
276   TTree* hitsTree = inFMD->TreeH();
277   if (!hitsTree)  {
278     // Try again 
279     inFMD->LoadHits("READ");
280     hitsTree = inFMD->TreeH();
281   }
282   
283   // Get the FMD branch 
284   TBranch* hitsBranch = hitsTree->GetBranch("FMD");
285   if (hitsBranch) fmd->SetHitsAddressBranch(hitsBranch);
286   else            AliFatal("Branch FMD hit not found");
287   
288   // Get a list of hits from the FMD manager 
289   TClonesArray *fmdHits = fmd->Hits();
290   
291   // Get number of entries in the tree 
292   Int_t ntracks  = Int_t(hitsTree->GetEntries());
293   
294   // Loop over the tracks in the 
295   for (Int_t track = 0; track < ntracks; track++)  {
296     // Read in entry number `track' 
297     hitsBranch->GetEntry(track);
298     
299     // Get the number of hits 
300     Int_t nhits = fmdHits->GetEntries ();
301     for (Int_t hit = 0; hit < nhits; hit++) {
302       // Get the hit number `hit'
303       AliFMDHit* fmdHit = 
304         static_cast<AliFMDHit*>(fmdHits->UncheckedAt(hit));
305       
306       // Extract parameters 
307       UShort_t detector = fmdHit->Detector();
308       Char_t   ring     = fmdHit->Ring();
309       UShort_t sector   = fmdHit->Sector();
310       UShort_t strip    = fmdHit->Strip();
311       Float_t  edep     = fmdHit->Edep();
312       if (fEdep(detector, ring, sector, strip).first != 0)
313         AliDebug(1, Form("Double hit in %d%c(%d,%d)", 
314                          detector, ring, sector, strip));
315       
316       fEdep(detector, ring, sector, strip).first  += edep;
317       fEdep(detector, ring, sector, strip).second += 1;
318       // Add this to the energy deposited for this strip
319     }  // hit loop
320   } // track loop
321 }
322
323 //____________________________________________________________________
324 void
325 AliFMDBaseDigitizer::DigitizeHits(AliFMD* fmd) const
326 {
327   // For the stored energy contributions in the cache (fEdep), convert
328   // the energy signal to ADC counts, and store the created digit in
329   // the digits array (AliFMD::fDigits)
330   //
331   TArrayI counts(3);
332   for (UShort_t detector=1; detector <= 3; detector++) {
333     // Get pointer to subdetector 
334     AliFMDSubDetector* det = 0;
335     switch (detector) {
336     case 1: det = fmd->GetFMD1(); break;
337     case 2: det = fmd->GetFMD2(); break;
338     case 3: det = fmd->GetFMD3(); break;
339     }
340     if (!det) continue;
341     for (UShort_t ringi = 0; ringi <= 1; ringi++) {
342       // Get pointer to Ring
343       AliFMDRing* r = 0;
344       switch (ringi) {
345       case 0: if (det->GetInner()) r = det->GetInner(); break;
346       case 1: if (det->GetOuter()) r = det->GetOuter(); break;
347       }
348       if (!r) continue;
349       
350       // Get number of sectors 
351       UShort_t nSectors = UShort_t(360. / r->GetTheta());
352       // Loop over the number of sectors 
353       for (UShort_t sector = 0; sector < nSectors; sector++) {
354         // Get number of strips 
355         UShort_t nStrips = r->GetNStrips();
356         // Loop over the stips 
357         for (UShort_t strip = 0; strip < nStrips; strip++) {
358           counts.Reset(-1);
359           Float_t edep = fEdep(detector, r->GetId(), sector, strip).first;
360           ConvertToCount(edep, r->GetSiThickness(), fmd->GetSiDensity(), 
361                          counts);
362           AddDigit(fmd, detector, r->GetId(), sector, strip, 
363                    edep, UShort_t(counts[0]), 
364                    Short_t(counts[1]), Short_t(counts[2]));
365 #if 0
366           // This checks if the digit created will give the `right'
367           // number of particles when reconstructed, using a naiive
368           // approach.  It's here only as a quality check - nothing
369           // else. 
370           CheckDigit(fEdep(detector, r->GetId(), sector, strip).first,
371                      fEdep(detector, r->GetId(), sector, strip).second,
372                      counts);
373 #endif
374         } // Strip
375       } // Sector 
376     } // Ring 
377   } // Detector 
378 }
379
380 //____________________________________________________________________
381 Float_t
382 AliFMDBaseDigitizer::ShapeIntegral(Float_t u, Float_t v) const
383 {
384   // Calculates the integral 
385   // 
386   //      / v
387   //      |               Bu + exp(-Bu) - Bv - exp(-Bv) 
388   //      | f(x) dx = - A -----------------------------
389   //      |                             B
390   //      / u
391   // 
392   // of the shaping function of the VA1_ALICE between times u and v
393   //
394   //      f(x) = A(1 - exp(-Bx))
395   //
396   // where A is a normalization constant, tuned so that the integral 
397   //
398   //      /1
399   //      |                              B          
400   //      | f(x) dx = 1   =>  A = ----------------
401   //      |                       -1 + B + exp(-B)
402   //      / 0
403   //
404   // and B is the a parameter defined by the shaping time (fShapingTime).  
405   // 
406   // That is, the function return the value 
407   // 
408   //        Bu + exp(-Bu) - Bv - exp(-Bv) 
409   //      - -----------------------------
410   //               -1 + B + exp(-B)
411   // 
412   // u,v should lie in the interval [0,1], and u < v
413   if (u == 0 && v == 1) return 1;
414   Float_t B = fShapingTime;
415   
416   // Calculate the integral 
417   Float_t res = - ((B * u + TMath::Exp(-B * u) - B * v - TMath::Exp(-B * v)) /
418                  (-1 + B + TMath::Exp(-B)));
419   return res;
420 }
421
422 //____________________________________________________________________
423 void
424 AliFMDBaseDigitizer::ConvertToCount(Float_t   edep, 
425                                     Float_t   siThickness, 
426                                     Float_t   siDensity, 
427                                     TArrayI&  counts) const
428 {
429   // Put noise and make ADC signal
430   // This is calculated as the product 
431   // 
432   //   DeltaEmip * SiThickness * SiDensity / Number 
433   //
434   // Where 
435   //  
436   //   DeltaEmip     is the energy loss of a MIP 
437   //   SiThickness   is the thickness of the silicon 
438   //   SiDensity     is the Silicon density 
439   //   Number        is # of e^- per MIP
440   //
441   // Note: Need to check this is correct. 
442   // 
443   const Float_t mipI = 1.664 * siThickness * siDensity;
444   // const Float_t mipI = 1.664 * 0.04 * 2.33 / 22400; // = 6.923e-6;
445   
446   // Create a pedestal 
447   UShort_t ped = MakePedestal();
448   
449   Float_t convf = 1 / mipI * Float_t(fAltroChannelSize) / fVA1MipRange;
450   Int_t n = fSampleRate;
451   for (Ssiz_t i = 0; i < n;  i++) {
452     Float_t w = ShapeIntegral(Float_t(i)/n, Float_t(i+1)/n);
453     counts[i] = UShort_t(TMath::Min(w * edep * convf + ped, 
454                                     Float_t(fAltroChannelSize))); 
455   }
456 }
457
458
459 //====================================================================
460 ClassImp(AliFMDDigitizer);
461
462 //____________________________________________________________________
463 AliFMDDigitizer::AliFMDDigitizer()  
464   : AliFMDBaseDigitizer()
465 {
466   // Default ctor - don't use it
467 }
468
469 //____________________________________________________________________
470 AliFMDDigitizer::AliFMDDigitizer(AliRunDigitizer* manager) 
471   : AliFMDBaseDigitizer(manager)
472 {
473   // Normal CTOR
474   AliDebug(1," processed");
475   SetPedestal();
476 }
477
478 //____________________________________________________________________
479 void
480 AliFMDDigitizer::Exec(Option_t*) 
481 {
482   // Get the output manager 
483   TString outFolder(fManager->GetOutputFolderName());
484   AliRunLoader* out = 
485     AliRunLoader::GetRunLoader(outFolder.Data());
486   // Get the FMD output manager 
487   AliLoader* outFMD = out->GetLoader("FMDLoader");
488
489   // Get the input loader 
490   TString inFolder(fManager->GetInputFolderName(0));
491   fRunLoader = 
492     AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
493   if (!fRunLoader) {
494     AliError("Can not find Run Loader for input stream 0");
495     return;
496   }
497   // Get the AliRun object 
498   if (!fRunLoader->GetAliRun()) fRunLoader->LoadgAlice();
499
500   // Get the AliFMD object 
501   AliFMD* fmd = static_cast<AliFMD*>(fRunLoader->GetAliRun()->GetDetector("FMD"));
502   if (!fmd) {
503     AliError("Can not get FMD from gAlice");
504     return;
505   }
506
507   Int_t nFiles= fManager->GetNinputs();
508   for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
509     AliDebug(1,Form(" Digitizing event number %d",
510                     fManager->GetOutputEventNr()));
511     // Get the current loader 
512     fRunLoader = 
513       AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
514     if (!fRunLoader) Fatal("Exec", "no run loader");
515     // Cache contriutions 
516     SumContributions(fmd);
517   }
518   // Digitize the event 
519   DigitizeHits(fmd);
520
521   // Load digits from the tree 
522   outFMD->LoadDigits("update");
523   // Get the tree of digits 
524   TTree* digitTree = outFMD->TreeD();
525   if (!digitTree) {
526     outFMD->MakeTree("D");
527     digitTree = outFMD->TreeD();
528   }
529   digitTree->Reset();
530   // Make a branch in the tree 
531   TClonesArray* digits = fmd->Digits();
532   fmd->MakeBranchInTree(digitTree, fmd->GetName(), &(digits), 4000, 0);
533   // TBranch* digitBranch = digitTree->GetBranch(fmd->GetName());
534   // Fill the tree 
535   digitTree->Fill();
536   
537   // Write the digits to disk 
538   outFMD->WriteDigits("OVERWRITE");
539   outFMD->UnloadHits();
540   outFMD->UnloadDigits();
541
542   // Reset the digits in the AliFMD object 
543   fmd->ResetDigits();
544 }
545
546
547 //____________________________________________________________________
548 UShort_t
549 AliFMDDigitizer::MakePedestal() const 
550 {
551   return UShort_t(TMath::Max(gRandom->Gaus(fPedestal, fPedestalWidth), 0.));
552 }
553
554 //____________________________________________________________________
555 void
556 AliFMDDigitizer::AddDigit(AliFMD*  fmd,
557                           UShort_t detector, 
558                           Char_t   ring,
559                           UShort_t sector, 
560                           UShort_t strip, 
561                           Float_t  /* edep */, 
562                           UShort_t count1, 
563                           Short_t  count2, 
564                           Short_t  count3) const
565 {
566   fmd->AddDigit(detector, ring, sector, strip, count1, count2, count3);
567 }
568
569 //____________________________________________________________________
570 void
571 AliFMDDigitizer::CheckDigit(Float_t         /* edep */, 
572                             UShort_t        nhits,
573                             const TArrayI&  counts) 
574 {
575   Int_t integral = counts[0];
576   if (counts[1] >= 0) integral += counts[1];
577   if (counts[2] >= 0) integral += counts[2];
578   integral -= Int_t(fPedestal + 2 * fPedestalWidth);
579   if (integral < 0) integral = 0;
580   
581   Float_t convf = Float_t(fVA1MipRange) / fAltroChannelSize;
582   Float_t mips  = integral * convf;
583   if (mips > Float_t(nhits) + .5 || mips < Float_t(nhits) - .5) 
584     Warning("CheckDigit", "Digit -> %4.2f MIPS != %d +/- .5 hits", 
585             mips, nhits);
586 }
587
588 //====================================================================
589 ClassImp(AliFMDSDigitizer);
590
591 //____________________________________________________________________
592 AliFMDSDigitizer::AliFMDSDigitizer()  
593 {
594   // Default ctor - don't use it
595 }
596
597 //____________________________________________________________________
598 AliFMDSDigitizer::AliFMDSDigitizer(const Char_t* headerFile, 
599                                    const Char_t* /* sdigfile */)
600   : AliFMDBaseDigitizer("FMDSDigitizer", "FMD SDigitizer")
601 {
602   // Normal CTOR
603   AliDebug(1," processed");
604
605   fRunLoader = AliRunLoader::GetRunLoader(); // Open(headerFile);
606   if (!fRunLoader) 
607     Fatal("AliFMDSDigitizer", "cannot open session, header file '%s'",
608           headerFile);
609   AliLoader* loader = fRunLoader->GetLoader("FMDLoader");
610   if (!loader) 
611     Fatal("AliFMDSDigitizer", "cannot find FMD loader in specified event");
612
613   // Add task to tasks folder 
614   loader->PostSDigitizer(this);
615
616 }
617
618 //____________________________________________________________________
619 AliFMDSDigitizer::~AliFMDSDigitizer() 
620 {
621   AliLoader* loader = fRunLoader->GetLoader("FMDLoader");
622   loader->CleanSDigitizer();
623 }
624
625 //____________________________________________________________________
626 void
627 AliFMDSDigitizer::Exec(Option_t*) 
628 {
629   // Get the output manager 
630   if (!fRunLoader) {
631     Error("Exec", "Run loader is not set");
632     return;
633   }
634   if (!fRunLoader->GetAliRun()) fRunLoader->LoadgAlice();
635   if (!fRunLoader->TreeE())     fRunLoader->LoadHeader();
636   
637   AliLoader* fmdLoader = fRunLoader->GetLoader("FMDLoader");
638   if (!fmdLoader) Fatal("Exec", "no FMD loader");
639   
640   // Get the AliFMD object 
641   AliFMD* fmd = 
642     static_cast<AliFMD*>(fRunLoader->GetAliRun()->GetDetector("FMD"));
643   if (!fmd) {
644     AliError("Can not get FMD from gAlice");
645     return;
646   }
647
648   Int_t nEvents = Int_t(fRunLoader->TreeE()->GetEntries());
649   for (Int_t event = 0; event < nEvents; event++) {
650     AliDebug(1,Form(" Digitizing event number %d", event));
651     // Get the current loader 
652     fRunLoader->GetEvent(event);
653
654     if (!fmdLoader->TreeS()) fmdLoader->MakeTree("S");
655     // Make a branch
656     fmd->MakeBranch("S");
657     
658     // Cache contriutions 
659     SumContributions(fmd);
660
661     // Digitize the event 
662     DigitizeHits(fmd);
663
664     fmdLoader->TreeS()->Reset();
665     fmdLoader->TreeS()->Fill();
666     fmdLoader->WriteSDigits("OVERWRITE");
667   }
668 }
669
670 //____________________________________________________________________
671 void
672 AliFMDSDigitizer::AddDigit(AliFMD*  fmd,
673                            UShort_t detector, 
674                            Char_t   ring,
675                            UShort_t sector, 
676                            UShort_t strip, 
677                            Float_t  edep, 
678                            UShort_t count1, 
679                            Short_t  count2, 
680                            Short_t  count3) const
681 {
682   fmd->AddSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
683 }
684
685
686
687 //____________________________________________________________________
688 //
689 // EOF
690 // 
691
692
693
694