]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDBaseDigitizer.cxx
Starting a collection of QA/Comparison macros
[u/mrichter/AliRoot.git] / FMD / AliFMDBaseDigitizer.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 /* $Id$ */
16 /** @file    AliFMDBaseDigitizer.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:38:26 2006
19     @brief   FMD Digitizers implementation
20     @ingroup FMD_sim
21 */
22 //////////////////////////////////////////////////////////////////////////////
23 //
24 //  This class contains the procedures simulation ADC  signal for the
25 //  Forward Multiplicity detector  : Hits->Digits and Hits->SDigits
26 // 
27 //  Digits consists of
28 //   - Detector #
29 //   - Ring ID                                             
30 //   - Sector #     
31 //   - Strip #
32 //   - ADC count in this channel                                  
33 //
34 //  Digits consists of
35 //   - Detector #
36 //   - Ring ID                                             
37 //   - Sector #     
38 //   - Strip #
39 //   - Total energy deposited in the strip
40 //   - ADC count in this channel                                  
41 //
42 // As the Digits and SDigits have so much in common, the classes
43 // AliFMDDigitizer and AliFMDSDigitizer are implemented via a base
44 // class AliFMDBaseDigitizer.
45 //
46 //                 +---------------------+
47 //                 | AliFMDBaseDigitizer |
48 //                 +---------------------+
49 //                           ^
50 //                           |
51 //                +----------+---------+
52 //                |                    |
53 //      +-----------------+     +------------------+
54 //      | AliFMDDigitizer |     | AliFMDSDigitizer |
55 //      +-----------------+     +------------------+
56 //
57 // These classes has several paramters: 
58 //
59 //     fPedestal
60 //     fPedestalWidth
61 //         (Only AliFMDDigitizer)
62 //         Mean and width of the pedestal.  The pedestal is simulated
63 //         by a Guassian, but derived classes my override MakePedestal
64 //         to simulate it differently (or pick it up from a database).
65 //
66 //     fVA1MipRange
67 //         The dymamic MIP range of the VA1_ALICE pre-amplifier chip 
68 //
69 //     fAltroChannelSize
70 //         The largest number plus one that can be stored in one
71 //         channel in one time step in the ALTRO ADC chip. 
72 //
73 //     fSampleRate
74 //         How many times the ALTRO ADC chip samples the VA1_ALICE
75 //         pre-amplifier signal.   The VA1_ALICE chip is read-out at
76 //         10MHz, while it's possible to drive the ALTRO chip at
77 //         25MHz.  That means, that the ALTRO chip can have time to
78 //         sample each VA1_ALICE signal up to 2 times.  Although it's
79 //         not certain this feature will be used in the production,
80 //         we'd like have the option, and so it should be reflected in
81 //         the code.
82 //
83 //
84 // The shaping function of the VA1_ALICE is generally given by 
85 //
86 //      f(x) = A(1 - exp(-Bx))
87 //
88 // where A is the total charge collected in the pre-amp., and B is a
89 // paramter that depends on the shaping time of the VA1_ALICE circut.
90 // 
91 // When simulating the shaping function of the VA1_ALICe
92 // pre-amp. chip, we have to take into account, that the shaping
93 // function depends on the previous value of read from the pre-amp. 
94 //
95 // That results in the following algorithm:
96 //
97 //    last = 0;
98 //    FOR charge IN pre-amp. charge train DO 
99 //      IF last < charge THEN 
100 //        f(t) = (charge - last) * (1 - exp(-B * t)) + last
101 //      ELSE
102 //        f(t) = (last - charge) * exp(-B * t) + charge)
103 //      ENDIF
104 //      FOR i IN # samples DO 
105 //        adc_i = f(i / (# samples))
106 //      DONE
107 //      last = charge
108 //   DONE
109 //
110 // Here, 
111 //
112 //   pre-amp. charge train 
113 //       is a series of 128 charges read from the VA1_ALICE chip
114 //
115 //   # samples
116 //       is the number of times the ALTRO ADC samples each of the 128
117 //       charges from the pre-amp. 
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 //      /1
165 //      |           A(-1 + B + exp(-B))
166 //      | f(x) dx = ------------------- = 1
167 //      |                    B
168 //      / 0
169 //
170 // and B is the a parameter defined by the shaping time (fShapingTime).  
171 //
172 // Solving the above equation, for A gives
173 //
174 //                 B
175 //      A = ----------------
176 //          -1 + B + exp(-B)
177 //
178 // So, if we define the function g: [0,1] -> [0:1] by 
179 //
180 //               / v
181 //               |              Bu + exp(-Bu) - Bv - exp(-Bv) 
182 //      g(u,v) = | f(x) dx = -A -----------------------------
183 //               |                            B
184 //               / u
185 //
186 // we can evaluate the ALTRO sample of the VA1_ALICE pre-amp between
187 // any two times (u, v), by 
188 //       
189 //
190 //                                B         Bu + exp(-Bu) - Bv - exp(-Bv)
191 //      C = Q g(u,v) = - Q ---------------- -----------------------------
192 //                         -1 + B + exp(-B)              B                  
193 //
194 //               Bu + exp(-Bu) - Bv - exp(-Bv) 
195 //        = -  Q -----------------------------
196 //                    -1 + B + exp(-B)
197 //
198
199 #include <TMath.h>
200 #include <TTree.h>              // ROOT_TTree
201 //#include <TRandom.h>          // ROOT_TRandom
202 // #include <AliLog.h>          // ALILOG_H
203 #include "AliFMDDebug.h" // Better debug macros
204 #include "AliFMDBaseDigitizer.h" // ALIFMDDIGITIZER_H
205 #include "AliFMD.h"             // ALIFMD_H
206 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
207 #include "AliFMDDetector.h"     // ALIFMDDETECTOR_H
208 #include "AliFMDRing.h"         // ALIFMDRING_H
209 #include "AliFMDHit.h"          // ALIFMDHIT_H
210 // #include "AliFMDDigit.h"     // ALIFMDDIGIT_H
211 #include "AliFMDParameters.h"   // ALIFMDPARAMETERS_H
212 // #include <AliDigitizationInput.h>    // ALIRUNDIGITIZER_H
213 //#include <AliRun.h>           // ALIRUN_H
214 #include <AliLoader.h>          // ALILOADER_H
215 #include <AliRun.h>             // ALILOADER_H
216 #include <AliRunLoader.h>       // ALIRUNLOADER_H
217 #include <TRandom.h>
218     
219 //====================================================================
220 ClassImp(AliFMDBaseDigitizer)
221 #if 0
222   ; // This is here to keep Emacs for indenting the next line
223 #endif
224
225 //____________________________________________________________________
226 AliFMDBaseDigitizer::AliFMDBaseDigitizer()  
227   : fFMD(0),
228     fRunLoader(0),
229     fEdep(AliFMDMap::kMaxDetectors, 
230           AliFMDMap::kMaxRings, 
231           AliFMDMap::kMaxSectors, 
232           AliFMDMap::kMaxStrips),
233     fShapingTime(6),
234     fStoreTrackRefs(kTRUE), 
235     fIgnoredLabels(0)
236 {
237   AliFMDDebug(1, ("Constructed"));
238   // Default ctor - don't use it
239 }
240
241 //____________________________________________________________________
242 AliFMDBaseDigitizer::AliFMDBaseDigitizer(AliDigitizationInput* digInput) 
243   : AliDigitizer(digInput, "AliFMDBaseDigitizer", "FMD Digitizer base class"), 
244     fFMD(0),
245     fRunLoader(0),
246     fEdep(0),        // nDet==0 means 51200 slots
247     fShapingTime(6),
248     fStoreTrackRefs(kTRUE), 
249     fIgnoredLabels(0)
250 {
251   // Normal CTOR
252   AliFMDDebug(1, ("Constructed"));
253   SetShapingTime();
254 }
255
256 //____________________________________________________________________
257 AliFMDBaseDigitizer::AliFMDBaseDigitizer(const Char_t* name, 
258                                          const Char_t* title) 
259   : AliDigitizer(name, title),
260     fFMD(0),
261     fRunLoader(0),
262     fEdep(0),        // nDet==0 means 51200 slots
263     fShapingTime(6),
264     fStoreTrackRefs(kTRUE), 
265     fIgnoredLabels(0)
266 {
267   // Normal CTOR
268   AliFMDDebug(1, (" Constructed"));
269   SetShapingTime();
270 }
271
272 //____________________________________________________________________
273 AliFMDBaseDigitizer::~AliFMDBaseDigitizer()
274 {
275   // Destructor
276 }
277
278 //____________________________________________________________________
279 AliFMDBaseDigitizer&
280 AliFMDBaseDigitizer::operator=(const AliFMDBaseDigitizer& o) 
281
282   // 
283   // Assignment operator
284   // 
285   // Return:
286   //    Reference to this object 
287   //
288   if (&o == this) return *this; 
289   AliDigitizer::operator=(o);
290   fRunLoader      = o.fRunLoader;
291   fEdep           = o.fEdep;
292   fShapingTime    = o.fShapingTime;
293   fStoreTrackRefs = o.fStoreTrackRefs;
294   fIgnoredLabels  = o.fIgnoredLabels;
295   return *this; 
296 }
297
298 //____________________________________________________________________
299 Bool_t 
300 AliFMDBaseDigitizer::Init()
301 {
302   // Initialization.   Get a pointer to the parameter manager, and
303   // initialize it.  
304   AliFMDParameters::Instance()->Init();
305   if (AliLog::GetDebugLevel("FMD","") >= 15) 
306     AliFMDParameters::Instance()->Print("");
307   return kTRUE;
308 }
309
310 //____________________________________________________________________
311 UShort_t
312 AliFMDBaseDigitizer::MakePedestal(UShort_t detector, 
313                                   Char_t   ring, 
314                                   UShort_t sector, 
315                                   UShort_t strip) const 
316
317   // Make a pedestal.  The pedestal value is drawn from a Gaussian
318   // distribution.  The mean of the distribution is the measured
319   // pedestal, and the width is the measured noise. 
320   AliFMDParameters* param =AliFMDParameters::Instance();
321   Float_t           mean  =param->GetPedestal(detector,ring,sector,strip);
322   Float_t           width =param->GetPedestalWidth(detector,ring,sector,strip);
323   return UShort_t(TMath::Max(gRandom->Gaus(mean, width), 0.));
324 }
325
326 //____________________________________________________________________
327 void
328 AliFMDBaseDigitizer::AddContribution(UShort_t detector, 
329                                      Char_t   ring, 
330                                      UShort_t sector, 
331                                      UShort_t strip, 
332                                      Float_t  edep, 
333                                      Bool_t   isPrimary,
334                                      Int_t    nTrack,
335                                      Int_t*   tracknos)
336 {
337   // Add edep contribution from (detector,ring,sector,strip) to cache
338   AliFMDParameters* param = AliFMDParameters::Instance();
339   AliFMDDebug(10, ("Adding contribution %7.5f for FMD%d%c[%2d,%3d] "
340                   " from %d tracks (%s)", 
341                   edep,
342                   detector, 
343                   ring,
344                   sector, 
345                   strip, 
346                   nTrack, 
347                   (isPrimary ? "primary" : "secondary")));
348   // Check if strip is `dead' 
349   if (param->IsDead(detector, ring, sector, strip)) { 
350     AliFMDDebug(5, ("FMD%d%c[%2d,%3d] is marked as dead", 
351                     detector, ring, sector, strip));
352     return;
353   }
354   // Check if strip is out-side read-out range 
355   // if (strip < minstrip || strip > maxstrip) {
356   //   AliFMDDebug(5, ("FMD%d%c[%2d,%3d] is outside range [%3d,%3d]", 
357   //                detector,ring,sector,strip,minstrip,maxstrip));
358   //   continue;
359   // }
360   
361   AliFMDEdepHitPair& entry = fEdep(detector, ring, sector, strip);
362
363   // Give warning in case of double sdigit 
364   if (entry.fEdep != 0)
365     AliFMDDebug(5, ("Double digit in FMD%d%c[%2d,%3d]", 
366                     detector, ring, sector, strip));
367       
368   // Sum energy deposition
369   Int_t oldN  =  entry.fN;
370   entry.fEdep += edep;
371   entry.fN    += nTrack;
372   if (isPrimary) entry.fNPrim += nTrack;
373   if (fStoreTrackRefs) { 
374     if (entry.fLabels.fN < entry.fN) {
375       AliFMDDebug(15, ("== New label array size %d, was %d, added %d", 
376                        entry.fN, entry.fLabels.fN, nTrack));
377       entry.fLabels.Set(entry.fN);
378     }
379     for (Int_t i = 0; i < nTrack; i++) {
380       AliFMDDebug(15, ("=> Setting track label # %d", oldN+i));
381       entry.fLabels[oldN + i] = tracknos[i];
382       AliFMDDebug(15, ("<= Setting track label # %d", oldN+i));
383     }
384   }
385   AliFMDDebug(15,("Adding contribution %f to FMD%d%c[%2d,%3d] (%f) track %d", 
386                   edep, detector, ring, sector, strip,
387                   entry.fEdep, (nTrack > 0 ? tracknos[0] : -1)));
388   
389 }
390
391 //____________________________________________________________________
392 void
393 AliFMDBaseDigitizer::DigitizeHits() const
394 {
395   // For the stored energy contributions in the cache (fEdep), convert
396   // the energy signal to ADC counts, and store the created digit in
397   // the digits array (AliFMD::fDigits)
398   //
399   AliFMDDebug(5, ("Will now digitize all the summed signals"));
400   fIgnoredLabels = 0;
401   AliFMDGeometry* geometry = AliFMDGeometry::Instance();
402   
403   TArrayI counts(4);
404   for (UShort_t detector=1; detector <= 3; detector++) {
405     AliFMDDebug(10, ("Processing hits in FMD%d", detector));
406     // Get pointer to subdetector 
407     AliFMDDetector* det = geometry->GetDetector(detector);
408     if (!det) continue;
409     for (UShort_t ringi = 0; ringi <= 1; ringi++) {
410       Char_t ring = ringi == 0 ? 'I' : 'O';
411       AliFMDDebug(10, (" Processing hits in FMD%d%c", detector,ring));
412       // Get pointer to Ring
413       AliFMDRing* r = det->GetRing(ring);
414       if (!r) continue;
415       
416       // Get number of sectors 
417       UShort_t nSectors = UShort_t(360. / r->GetTheta());
418       // Loop over the number of sectors 
419       for (UShort_t sector = 0; sector < nSectors; sector++) {
420         AliFMDDebug(10, ("  Processing hits in FMD%d%c[%2d]", 
421                         detector,ring,sector));
422         // Get number of strips 
423         UShort_t nStrips = r->GetNStrips();
424         // Loop over the stips 
425         Float_t last = 0;
426         for (UShort_t strip = 0; strip < nStrips; strip++) {
427           // Reset the counter array to the invalid value -1 
428           counts.Reset(-1);
429           // Reset the last `ADC' value when we've get to the end of a
430           // VA1_ALICE channel. 
431           if (strip % 128 == 0) last = 0;
432           
433           const AliFMDEdepHitPair& entry  = fEdep(detector,ring,sector,strip);
434           Float_t                  edep   = entry.fEdep;
435           UShort_t                 ntot   = entry.fN;
436           UShort_t                 nprim  = entry.fNPrim;
437           const TArrayI&           labels = entry.fLabels;
438           if (edep > 0)
439             AliFMDDebug(15, ("Edep = %f for FMD%d%c[%2d,%3d]", 
440                              edep, detector, ring, sector, strip));
441           ConvertToCount(edep, last, detector, ring, sector, strip, counts);
442           last = edep;
443           
444
445           // The following line was introduced - wrongly - by Peter
446           // Hristov.  It _will_ break the digitisation and the
447           // following reconstruction.  The behviour of the
448           // digitisation models exactly the front-end as it should
449           // (no matter what memory concuption it may entail).  The
450           // check should be on zero suppression, since that's what
451           // models the front-end - if zero suppression is turned on
452           // in the front-end, then we can suppress empty digits -
453           // otherwise we shoud never do that.  Note, that the line
454           // affects _both_ normal digitisation and digitisation for
455           // summable digits, since the condition is on the energy
456           // deposition and not on the actual number of counts.  If
457           // this line should go anywhere, it should be in the
458           // possible overloaded AliFMDSDigitizer::AddDigit - not
459           // here. 
460           // 
461           //   if (edep<=0) continue;
462           AddDigit(detector, ring, sector, strip, edep, 
463                    UShort_t(counts[0]), Short_t(counts[1]), 
464                    Short_t(counts[2]), Short_t(counts[3]), 
465                    ntot, nprim, labels);
466           AliFMDDebug(15, ("   Adding digit in FMD%d%c[%2d,%3d]=%d", 
467                            detector,ring,sector,strip,counts[0]));
468 #if 0
469           // This checks if the digit created will give the `right'
470           // number of particles when reconstructed, using a naiive
471           // approach.  It's here only as a quality check - nothing
472           // else. 
473           CheckDigit(digit, fEdep(detector, ring, sector, strip).fN,
474                      counts);
475 #endif
476         } // Strip
477       } // Sector 
478     } // Ring 
479   } // Detector 
480   if (fIgnoredLabels > 0) 
481     AliWarning(Form("%d track labels could not be associated with digits "
482                     "due to limited storage facilities in AliDigit", 
483                     fIgnoredLabels));
484 }
485
486 //____________________________________________________________________
487 void
488 AliFMDBaseDigitizer::ConvertToCount(Float_t   edep, 
489                                     Float_t   last,
490                                     UShort_t  detector, 
491                                     Char_t    ring, 
492                                     UShort_t  sector, 
493                                     UShort_t  strip,
494                                     TArrayI&  counts) const
495 {
496   // Convert the total energy deposited to a (set of) ADC count(s). 
497   // 
498   // This is done by 
499   // 
500   //               Energy_Deposited      ALTRO_Channel_Size
501   //    ADC = -------------------------- ------------------- + pedestal
502   //          Energy_Deposition_Of_1_MIP VA1_ALICE_MIP_Range
503   //
504   //               Energy_Deposited             fAltroChannelSize
505   //        = --------------------------------- ----------------- + pedestal 
506   //          1.664 * Si_Thickness * Si_Density   fVA1MipRange   
507   //          
508   // 
509   //        = Energy_Deposited * ConversionFactor + pedestal
510   // 
511   // However, this is modified by the response function of the
512   // VA1_ALICE pre-amp. chip in case we are doing oversampling of the
513   // VA1_ALICE output. 
514   // 
515   // In that case, we get N=fSampleRate values of the ADC, and the
516   // `EnergyDeposited' is a function of which sample where are
517   // calculating the ADC for 
518   // 
519   //     ADC_i = f(EnergyDeposited, i/N, Last) * ConversionFactor + pedestal 
520   // 
521   // where Last is the Energy deposited in the previous strip. 
522   // 
523   // Here, f is the shaping function of the VA1_ALICE.   This is given
524   // by 
525   //                       
526   //                    |   (E - l) * (1 - exp(-B * t) + l   if E > l
527   //       f(E, t, l) = <
528   //                    |   (l - E) * exp(-B * t) + E        otherwise
529   //                       
530   // 
531   //                  = E + (l - E) * ext(-B * t)
532   // 
533   AliFMDParameters* param = AliFMDParameters::Instance();
534   Float_t  convF          = (param->GetDACPerMIP() / param->GetEdepMip() *
535                              param->GetPulseGain(detector,ring,sector,strip));
536   Int_t    ped            = MakePedestal(detector,ring,sector,strip);
537   Int_t    maxAdc         = param->GetAltroChannelSize()-1;
538   if (maxAdc < 0) {
539     AliWarning(Form("Maximum ADC is %d < 0, forcing it to 1023", maxAdc));
540     maxAdc = 1023;
541   }
542   UShort_t rate           = param->GetSampleRate(detector,ring,sector,strip);
543   AliFMDDebug(15, ("Sample rate for FMD%d%c[%2d,%3d] = %d", 
544                    detector, ring, sector, strip, rate));
545   if (rate < 1 || rate > 4) {
546     AliWarning(Form("Invalid sample rate for for FMD%d%c[%2d,%3d] = %d", 
547                     detector, ring, sector, strip, rate));
548     rate = 1;
549   }
550
551   // In case we don't oversample, just return the end value. 
552   if (rate == 1) {
553     Float_t    a = edep * convF + ped;
554     if (a < 0) a = 0;
555     counts[0]    = UShort_t(TMath::Min(a, Float_t(maxAdc)));
556     AliFMDDebug(15, ("FMD%d%c[%2d,%3d]: converting ELoss %f to "
557                      "ADC %4d (%f,%d)",
558                      detector,ring,sector,strip,edep,counts[0],convF,ped));
559     return;
560   }
561
562   
563   // Create a pedestal 
564   Float_t b = fShapingTime;
565   for (Ssiz_t i = 0; i < rate;  i++) {
566     Float_t t  = Float_t(i) / rate + 1./rate;
567     Float_t s  = edep + (last - edep) * TMath::Exp(-b * t);
568     Float_t a  = Int_t(s * convF + ped);
569     if (a < 0) a = 0;
570     counts[i]  = UShort_t(TMath::Min(a, Float_t(maxAdc)));
571   }
572   AliFMDDebug(15, ("Converted edep = %f to ADC (%x,%x,%x,%x) "
573                    "[gain: %f=(%f/%f*%f), pedestal: %d, rate: %d]", 
574                    edep, counts[0], counts[1], counts[2], counts[3], 
575                    convF, param->GetDACPerMIP(),param->GetEdepMip(),
576                    param->GetPulseGain(detector,ring,sector,strip), 
577                    ped, rate));
578 }
579
580 //____________________________________________________________________
581 void
582 AliFMDBaseDigitizer::AddDigit(UShort_t        detector, 
583                               Char_t          ring,
584                               UShort_t        sector, 
585                               UShort_t        strip, 
586                               Float_t         /* edep */, 
587                               UShort_t        count1, 
588                               Short_t         count2, 
589                               Short_t         count3,
590                               Short_t         count4,
591                               UShort_t        ntot, 
592                               UShort_t        /* nprim */,
593                               const TArrayI&  refs) const
594 {
595   // Add a digit or summable digit
596   fFMD->AddDigitByFields(detector, ring, sector, strip, 
597                          count1, count2, count3, count4, 
598                          ntot, fStoreTrackRefs ? refs.fArray : 0);
599   if (fStoreTrackRefs && ntot > 3) fIgnoredLabels += ntot - 3;
600 }
601
602 //____________________________________________________________________
603 TTree*
604 AliFMDBaseDigitizer::MakeOutputTree(AliLoader* loader)
605 {
606   // Create output tree using loader.   If the passed loader differs
607   // from the currently set loader in the FMD object, reset the FMD
608   // loader to be the passed loader.   This is for the cases wher the
609   // output is different from the output. 
610   AliFMDDebug(5, ("Making digits tree"));
611   loader->LoadDigits("UPDATE"); // "RECREATE");
612   TTree* out = loader->TreeD();
613   if (!out) loader->MakeTree("D");
614   out = loader->TreeD(); 
615   if (out) { 
616     out->Reset();
617     if (loader != fFMD->GetLoader()) 
618       fFMD->SetLoader(loader);
619     fFMD->MakeBranch("D");
620   }
621   return out;
622 }
623
624 //____________________________________________________________________
625 void
626 AliFMDBaseDigitizer::StoreDigits(const AliLoader* loader)
627 {
628   // Write the digits to disk 
629   AliFMDDebug(5, ("Storing %d digits",   fFMD->Digits()->GetEntries()));
630   loader->WriteDigits("OVERWRITE");
631   loader->UnloadDigits();
632   // Reset the digits in the AliFMD object 
633   fFMD->ResetDigits();
634 }
635
636 //____________________________________________________________________
637 //
638 // EOF
639 // 
640
641
642
643