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