]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDHitDigitizer.cxx
add identifiers between MC and no MC data for histo titles
[u/mrichter/AliRoot.git] / FMD / AliFMDHitDigitizer.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: AliFMDHitDigitizer.cxx 28055 2008-08-18 00:33:20Z cholm $ */
16 /** @file    AliFMDHitDigitizer.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
26 // 
27 //  Digits consists of
28 //   - Detector #
29 //   - Ring ID                                             
30 //   - Sector #     
31 //   - Strip #
32 //   - ADC count in this channel                                  
33 //
34 // As the Digits and SDigits have so much in common, the classes
35 // AliFMDHitDigitizer and AliFMDSDigitizer are implemented via a base
36 // class AliFMDBaseDigitizer.
37 //
38 //                 +---------------------+
39 //                 | AliFMDBaseDigitizer |
40 //                 +---------------------+
41 //                           ^
42 //                           |
43 //                +----------+---------+
44 //                |                    |
45 //      +-----------------+     +------------------+
46 //      | AliFMDHitDigitizer |  | AliFMDSDigitizer |
47 //      +-----------------+     +------------------+
48 //
49 // These classes has several paramters: 
50 //
51 //     fPedestal
52 //     fPedestalWidth
53 //         (Only AliFMDHitDigitizer)
54 //         Mean and width of the pedestal.  The pedestal is simulated
55 //         by a Guassian, but derived classes my override MakePedestal
56 //         to simulate it differently (or pick it up from a database).
57 //
58 //     fVA1MipRange
59 //         The dymamic MIP range of the VA1_ALICE pre-amplifier chip 
60 //
61 //     fAltroChannelSize
62 //         The largest number plus one that can be stored in one
63 //         channel in one time step in the ALTRO ADC chip. 
64 //
65 //     fSampleRate
66 //         How many times the ALTRO ADC chip samples the VA1_ALICE
67 //         pre-amplifier signal.   The VA1_ALICE chip is read-out at
68 //         10MHz, while it's possible to drive the ALTRO chip at
69 //         25MHz.  That means, that the ALTRO chip can have time to
70 //         sample each VA1_ALICE signal up to 2 times.  Although it's
71 //         not certain this feature will be used in the production,
72 //         we'd like have the option, and so it should be reflected in
73 //         the code.
74 //
75 // These parameters are fetched from OCDB via the mananger AliFMDParameters.
76 //
77 // The shaping function of the VA1_ALICE is generally given by 
78 //
79 //      f(x) = A(1 - exp(-Bx))
80 //
81 // where A is the total charge collected in the pre-amp., and B is a
82 // paramter that depends on the shaping time of the VA1_ALICE circut.
83 // 
84 // When simulating the shaping function of the VA1_ALICe
85 // pre-amp. chip, we have to take into account, that the shaping
86 // function depends on the previous value of read from the pre-amp. 
87 //
88 // That results in the following algorithm:
89 //
90 //    last = 0;
91 //    FOR charge IN pre-amp. charge train DO 
92 //      IF last < charge THEN 
93 //        f(t) = (charge - last) * (1 - exp(-B * t)) + last
94 //      ELSE
95 //        f(t) = (last - charge) * exp(-B * t) + charge)
96 //      ENDIF
97 //      FOR i IN # samples DO 
98 //        adc_i = f(i / (# samples))
99 //      DONE
100 //      last = charge
101 //   DONE
102 //
103 // Here, 
104 //
105 //   pre-amp. charge train 
106 //       is a series of 128 charges read from the VA1_ALICE chip
107 //
108 //   # samples
109 //       is the number of times the ALTRO ADC samples each of the 128
110 //       charges from the pre-amp. 
111 //
112 // Where Q is the total charge collected by the VA1_ALICE
113 // pre-amplifier.   Q is then given by 
114 //
115 //           E S 
116 //      Q =  - -
117 //           e R
118 //
119 // where E is the total energy deposited in a silicon strip, R is the
120 // dynamic range of the VA1_ALICE pre-amp (fVA1MipRange), e is the
121 // energy deposited by a single MIP, and S ALTRO channel size in each
122 // time step (fAltroChannelSize).  
123 //
124 // The energy deposited per MIP is given by 
125 //
126 //      e = M * rho * w 
127 //
128 // where M is the universal number 1.664, rho is the density of
129 // silicon, and w is the depth of the silicon sensor. 
130 //
131 // The final ADC count is given by 
132 //
133 //      C' = C + P
134 //
135 // where P is the (randomized) pedestal (see MakePedestal)
136 //
137 // This class uses the class template AliFMDMap<Type> to make an
138 // internal cache of the energy deposted of the hits.  The class
139 // template is instantasized as 
140 //
141 //  typedef AliFMDMap<std::pair<Float_t, UShort_t> > AliFMDEdepMap;
142 //
143 // The first member of the values is the summed energy deposition in a
144 // given strip, while the second member of the values is the number of
145 // hits in a given strip.  Using the second member, it's possible to
146 // do some checks on just how many times a strip got hit, and what
147 // kind of error we get in our reconstructed hits.  Note, that this
148 // information is currently not written to the digits tree.  I think a
149 // QA (Quality Assurance) digit tree is better suited for that task.
150 // However, the information is there to be used in the future. 
151 //
152 //
153 // Latest changes by Christian Holm Christensen
154 //
155 //////////////////////////////////////////////////////////////////////////////
156
157 //      /1
158 //      |           A(-1 + B + exp(-B))
159 //      | f(x) dx = ------------------- = 1
160 //      |                    B
161 //      / 0
162 //
163 // and B is the a parameter defined by the shaping time (fShapingTime).  
164 //
165 // Solving the above equation, for A gives
166 //
167 //                 B
168 //      A = ----------------
169 //          -1 + B + exp(-B)
170 //
171 // So, if we define the function g: [0,1] -> [0:1] by 
172 //
173 //               / v
174 //               |              Bu + exp(-Bu) - Bv - exp(-Bv) 
175 //      g(u,v) = | f(x) dx = -A -----------------------------
176 //               |                            B
177 //               / u
178 //
179 // we can evaluate the ALTRO sample of the VA1_ALICE pre-amp between
180 // any two times (u, v), by 
181 //       
182 //
183 //                                B         Bu + exp(-Bu) - Bv - exp(-Bv)
184 //      C = Q g(u,v) = - Q ---------------- -----------------------------
185 //                         -1 + B + exp(-B)              B                  
186 //
187 //               Bu + exp(-Bu) - Bv - exp(-Bv) 
188 //        = -  Q -----------------------------
189 //                    -1 + B + exp(-B)
190 //
191
192 #include <TTree.h>              // ROOT_TTree
193 #include "AliFMDDebug.h"        // Better debug macros
194 #include "AliFMDHitDigitizer.h" // ALIFMDDIGITIZER_H
195 #include "AliFMD.h"             // ALIFMD_H
196 #include "AliFMDDigit.h"        // ALIFMDDIGIT_H
197 #include "AliFMDParameters.h"   // ALIFMDPARAMETERS_H
198 #include <AliRun.h>             // ALIRUN_H
199 #include <AliLoader.h>          // ALILOADER_H
200 #include <AliRunLoader.h>       // ALIRUNLOADER_H
201 #include <AliFMDHit.h>
202 #include <AliStack.h>
203 #include <TFile.h>
204 #include <TParticle.h>
205
206 //====================================================================
207 ClassImp(AliFMDHitDigitizer)    
208 #if 0
209 ;
210 #endif
211
212 //____________________________________________________________________
213 AliFMDHitDigitizer::AliFMDHitDigitizer(AliFMD* fmd, Output_t  output)
214   : AliFMDBaseDigitizer("FMD", (fOutput == kDigits ? 
215                                 "FMD Hit->Digit digitizer" :
216                                 "FMD Hit->SDigit digitizer")),
217     fOutput(output), 
218     fHoldTime(2e-6),
219     fStack(0)
220 {
221   fFMD = fmd;
222 }
223
224 //____________________________________________________________________
225 void
226 AliFMDHitDigitizer::Exec(Option_t* /*option*/)
227 {
228   // Run this digitizer 
229   // Get an inititialize parameter manager
230   AliFMDParameters::Instance()->Init();
231   if (AliLog::GetDebugLevel("FMD","") >= 10) 
232     AliFMDParameters::Instance()->Print("ALL");
233
234   // Get loader, and ask it to read in the hits 
235   AliLoader* loader = fFMD->GetLoader();
236   if (!loader) { 
237     AliError("Failed to get loader from detector object");
238     return;
239   }
240   loader->LoadHits("READ");
241   
242   // Get the run loader 
243   AliRunLoader* runLoader = loader->GetRunLoader();
244   if (!runLoader) {
245     AliError("Failed to get run loader from loader");
246     return;
247   }
248   
249   // Now loop over events
250   Int_t nEvents = runLoader->GetNumberOfEvents();
251   for (Int_t event = 0; event < nEvents; event++) { 
252     // Get the current event folder. 
253     TFolder* folder = loader->GetEventFolder();
254     if (!folder) { 
255       AliError("Failed to get event folder from loader");
256       return;
257     }
258
259     // Get the run-loader of this event. 
260     const char* loaderName = AliRunLoader::GetRunLoaderName();
261     AliRunLoader* thisLoader = 
262       static_cast<AliRunLoader*>(folder->FindObject(loaderName));
263     if (!thisLoader) { 
264       AliError(Form("Failed to get loader '%s' from event folder", loaderName));
265       return;
266     }
267     
268     // Read in the event
269     AliFMDDebug(1, ("Now digitizing (Hits->%s) event # %d", 
270                     (fOutput == kDigits ? "digits" : "sdigits"), event));
271     thisLoader->GetEvent(event);
272     
273     // Load kinematics to get primary information for SDigits
274     fStack = 0;
275     if (fOutput == kSDigits) {
276       if (thisLoader->LoadKinematics("READ")) {
277         AliError("Failed to get kinematics from event loader");
278         return;
279       }
280       AliFMDDebug(5, ("Loading stack of kinematics"));
281       fStack = thisLoader->Stack();
282     }
283
284     // Check that we have the hits 
285     if (!loader->TreeH() && loader->LoadHits()) {
286       AliError("Failed to load hits");
287       return;
288     }
289     TTree*   hitsTree = loader->TreeH();
290     TBranch* hitsBranch = hitsTree->GetBranch(fFMD->GetName());
291     if (!hitsBranch) { 
292       AliError("Failed to get hits branch in tree");
293       return;
294     }
295     // Check that we can make the output digits - This must come
296     // before AliFMD::SetBranchAddress
297     TTree* outTree = MakeOutputTree(loader);
298     if (!outTree) { 
299       AliError("Failed to get output tree");
300       return;
301     }
302     AliFMDDebug(5, ("Output tree name for %s is '%s'", 
303                     (fOutput == kDigits ? "digits" : "sdigits"),
304                     outTree->GetName()));
305     if (AliLog::GetDebugLevel("FMD","") >= 5) {
306       TFile* file = outTree->GetCurrentFile();
307       if (!file) {
308         AliWarning("Output tree has no file!");
309       }
310       else { 
311         AliFMDDebug(5, ("Output tree file %s content:", file->GetName()));
312         file->ls();
313       }
314     }
315
316     // Set-up the branch addresses 
317     fFMD->SetTreeAddress();
318     
319     // Now sum all contributions in cache 
320     SumContributions(hitsBranch);
321     loader->UnloadHits();
322
323     // And now digitize the hits 
324     DigitizeHits();
325     
326     // Write digits to tree
327     Int_t write = outTree->Fill();
328     AliFMDDebug(5, ("Wrote %d bytes to digit tree", write));
329
330     // Store the digits
331     StoreDigits(loader);
332
333   }  
334 }
335
336 //____________________________________________________________________
337 TTree*
338 AliFMDHitDigitizer::MakeOutputTree(AliLoader* loader)
339 {
340   if (fOutput == kDigits) 
341     return AliFMDBaseDigitizer::MakeOutputTree(loader);
342   
343   AliFMDDebug(5, ("Making sdigits tree"));
344   loader->LoadSDigits("UPDATE"); // RECREATE");
345   TTree* out = loader->TreeS();
346   if (!out) loader->MakeTree("S");
347   out = loader->TreeS(); 
348   if (out) { 
349     out->Reset();
350     fFMD->MakeBranch("S");
351   }
352   return out;
353 }
354
355
356 //____________________________________________________________________
357 void
358 AliFMDHitDigitizer::SumContributions(TBranch* hitsBranch) 
359 {
360   // Sum energy deposited contributions from each hit in a cache
361   // (fEdep).  
362   
363   // Clear array of deposited energies 
364   fEdep.Reset();
365   
366   // Get a list of hits from the FMD manager 
367   AliFMDDebug(5, ("Get array of FMD hits"));
368   TClonesArray *fmdHits = fFMD->Hits();
369   
370
371   // Get number of entries in the tree 
372   AliFMDDebug(5, ("Get # of tracks"));
373   Int_t ntracks  = Int_t(hitsBranch->GetEntries());
374   AliFMDDebug(5, ("We got %d tracks", ntracks));
375
376   Int_t read = 0;
377   // Loop over the tracks in the 
378   for (Int_t track = 0; track < ntracks; track++)  {
379     // Read in entry number `track' 
380     read += hitsBranch->GetEntry(track);
381
382     // Get the number of hits 
383     Int_t nhits = fmdHits->GetEntries ();
384     for (Int_t hit = 0; hit < nhits; hit++) {
385       // Get the hit number `hit'
386       AliFMDHit* fmdHit = 
387         static_cast<AliFMDHit*>(fmdHits->UncheckedAt(hit));
388
389       // Ignore hits that arrive too late
390       if (fmdHit->Time() > fHoldTime) continue;
391       
392
393       // Check if this is a primary particle
394       Bool_t isPrimary = kTRUE;
395       Int_t  trackno   = -1;
396       if (fStack) {
397         trackno = fmdHit->Track();
398         AliFMDDebug(10, ("Will get track # %d/%d from entry # %d", 
399                         trackno, fStack->GetNtrack(), track));
400         if (fStack->GetNtrack() < trackno) {
401           AliError(Form("Track number %d/%d out of bounds", 
402                         trackno, fStack->GetNtrack()));
403           continue;
404         }
405 #if 1
406         isPrimary = fStack->IsPhysicalPrimary(trackno);
407 #else // This is our hand-crafted code.  We use the ALICE definition
408         TParticle* part    = fStack->Particle(trackno);
409         isPrimary          = part->IsPrimary();
410         if (!isPrimary) { 
411           // Extended testing of mother status - this is for Pythia6.
412           Int_t      mother1   = part->GetFirstMother();
413           TParticle* mother    = fStack->Particle(mother1);
414           if (!mother || mother->GetStatusCode() > 1)
415             isPrimary = kTRUE;
416           AliFMDDebug(15,
417                       ("Track %d secondary, mother: %d - %s - status %d: %s", 
418                        trackno, mother1, 
419                        (mother ? "found"                 : "not found"), 
420                        (mother ? mother->GetStatusCode() : -1),
421                        (isPrimary ? "primary" : "secondary")));
422         }
423 #endif
424       }
425     
426       // Extract parameters 
427       AliFMDDebug(15,("Adding contribution %7.5f for FMD%d%c[%2d,%3d] "
428                       " for trackno %6d (%s)", 
429                       fmdHit->Edep(),
430                       fmdHit->Detector(), 
431                       fmdHit->Ring(),
432                       fmdHit->Sector(), 
433                       fmdHit->Strip(), 
434                       trackno, 
435                       (isPrimary ? "primary" : "secondary")));
436       AddContribution(fmdHit->Detector(),
437                       fmdHit->Ring(),
438                       fmdHit->Sector(),
439                       fmdHit->Strip(),
440                       fmdHit->Edep(), 
441                       isPrimary, 
442                       1, 
443                       &trackno);
444     }  // hit loop
445   } // track loop
446   AliFMDDebug(5, ("Size of cache: %d bytes, read %d bytes", 
447                   int(sizeof(fEdep)), read));
448 }
449
450
451 //____________________________________________________________________
452 UShort_t
453 AliFMDHitDigitizer::MakePedestal(UShort_t  detector, 
454                                  Char_t    ring, 
455                                  UShort_t  sector, 
456                                  UShort_t  strip) const 
457 {
458   // Make a pedestal 
459   if (fOutput == kSDigits) return 0;
460   return AliFMDBaseDigitizer::MakePedestal(detector, ring, sector, strip);
461 }
462
463
464
465 //____________________________________________________________________
466 void
467 AliFMDHitDigitizer::AddDigit(UShort_t        detector, 
468                              Char_t          ring,
469                              UShort_t        sector, 
470                              UShort_t        strip, 
471                              Float_t         edep, 
472                              UShort_t        count1, 
473                              Short_t         count2, 
474                              Short_t         count3,
475                              Short_t         count4, 
476                              UShort_t        ntotal,
477                              UShort_t        nprim, 
478                              const TArrayI&  refs) const
479 {
480   // Add a digit or summable digit
481   if (fOutput == kDigits) { 
482     AliFMDDebug(15,("Adding digit for FMD%d%c[%2d,%3d] = (%x,%x,%x,%x)",
483                     detector, ring, sector, strip, 
484                     count1, count2, count3, count4));
485     AliFMDBaseDigitizer::AddDigit(detector, ring, sector, strip, 0,
486                                   count1, count2, count3, count4, 
487                                   ntotal, nprim, refs);
488     return;
489   }
490   if (edep <= 0) { 
491     AliFMDDebug(15, ("Digit edep = %f <= 0 for FMD%d%c[%2d,%3d]", 
492                     edep, detector, ring, sector, strip));
493     return;
494   }
495   if (count1 == 0 && count2 <= 0 && count3 <= 0 && count4 <= 0) {
496     AliFMDDebug(15, ("Digit counts = (%x,%x,%x,%x) <= 0 for FMD%d%c[%2d,%3d]", 
497                     count1, count2, count3, count4, 
498                     detector, ring, sector, strip));
499     return;
500   }
501   AliFMDDebug(15, ("Adding sdigit for FMD%d%c[%2d,%3d] = "
502                    "(%x,%x,%x,%x) [%d/%d] %d",
503                    detector, ring, sector, strip, 
504                    count1, count2, count3, count4, nprim, ntotal, refs.fN));
505   fFMD->AddSDigitByFields(detector, ring, sector, strip, edep,
506                           count1, count2, count3, count4, 
507                           ntotal, nprim, fStoreTrackRefs ? refs.fArray : 0);
508 }
509
510 //____________________________________________________________________
511 void
512 AliFMDHitDigitizer::CheckDigit(AliFMDDigit*    digit,
513                                UShort_t        nhits,
514                                const TArrayI&  counts) 
515 {
516   // Check that digit is consistent
517   AliFMDParameters* param = AliFMDParameters::Instance();
518   UShort_t          det   = digit->Detector();
519   Char_t            ring  = digit->Ring();
520   UShort_t          sec   = digit->Sector();
521   UShort_t          str   = digit->Strip();
522   Float_t           mean  = param->GetPedestal(det,ring,sec,str);
523   Float_t           width = param->GetPedestalWidth(det,ring,sec,str);
524   UShort_t          range = param->GetVA1MipRange();
525   UShort_t          size  = param->GetAltroChannelSize();
526   Int_t             integral = counts[0];
527   if (counts[1] >= 0) integral += counts[1];
528   if (counts[2] >= 0) integral += counts[2];
529   if (counts[3] >= 0) integral += counts[3];
530   integral -= Int_t(mean + 2 * width);
531   if (integral < 0) integral = 0;
532   
533   Float_t convF = Float_t(range) / size;
534   Float_t mips  = integral * convF;
535   if (mips > Float_t(nhits) + .5 || mips < Float_t(nhits) - .5) 
536     Warning("CheckDigit", "Digit -> %4.2f MIPS != %d +/- .5 hits", 
537             mips, nhits);
538 }
539
540 //____________________________________________________________________
541 void
542 AliFMDHitDigitizer::StoreDigits(AliLoader* loader)
543 {
544   if (fOutput == kDigits) { 
545     AliFMDBaseDigitizer::StoreDigits(loader);
546     return;
547   }
548   AliFMDDebug(5, ("Storing %d sdigits",   fFMD->SDigits()->GetEntries()));
549   // Write the digits to disk 
550   loader->WriteSDigits("OVERWRITE");
551   loader->UnloadSDigits();
552   // Reset the digits in the AliFMD object 
553   fFMD->ResetSDigits();
554 }
555
556
557 //____________________________________________________________________
558 //
559 // EOF
560 // 
561
562
563
564