]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDHitDigitizer.cxx
Several fixes (the Make-Federico-Happy-Commit):
[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 <TFile.h>
203
204 //====================================================================
205 ClassImp(AliFMDHitDigitizer)    
206 #if 0
207 ;
208 #endif
209
210 //____________________________________________________________________
211 AliFMDHitDigitizer::AliFMDHitDigitizer(AliFMD* fmd, Output_t  output)
212   : AliFMDBaseDigitizer("FMD", (fOutput == kDigits ? 
213                                 "FMD Hit->Digit digitizer" :
214                                 "FMD Hit->SDigit digitizer")),
215     fOutput(output)
216 {
217   fFMD = fmd;
218 }
219
220 //____________________________________________________________________
221 void
222 AliFMDHitDigitizer::Exec(Option_t* /*option*/)
223 {
224   // Run this digitizer 
225   // Get an inititialize parameter manager
226   AliFMDParameters::Instance()->Init();
227   if (AliLog::GetDebugLevel("FMD","") >= 10) 
228     AliFMDParameters::Instance()->Print("ALL");
229
230   // Get loader, and ask it to read in the hits 
231   AliLoader* loader = fFMD->GetLoader();
232   if (!loader) { 
233     AliError("Failed to get loader from detector object");
234     return;
235   }
236   loader->LoadHits("READ");
237   
238   // Get the run loader 
239   AliRunLoader* runLoader = loader->GetRunLoader();
240   if (!runLoader) {
241     AliError("Failed to get run loader from loader");
242     return;
243   }
244   
245   // Now loop over events
246   Int_t nEvents = runLoader->GetNumberOfEvents();
247   for (Int_t event = 0; event < nEvents; event++) { 
248     // Get the current event folder. 
249     TFolder* folder = loader->GetEventFolder();
250     if (!folder) { 
251       AliError("Failed to get event folder from loader");
252       return;
253     }
254
255     // Get the run-loader of this event. 
256     const char* loaderName = AliRunLoader::GetRunLoaderName();
257     AliRunLoader* thisLoader = 
258       static_cast<AliRunLoader*>(folder->FindObject(loaderName));
259     if (!thisLoader) { 
260       AliError(Form("Failed to get loader '%s' from event folder", loaderName));
261       return;
262     }
263     
264     // Read in the event
265     AliFMDDebug(5, ("Now digitizing (Hits->%s) event # %d", 
266                     (fOutput == kDigits ? "digits" : "sdigits"), event));
267     thisLoader->GetEvent(event);
268     
269     // Check that we have the hits 
270     if (!loader->TreeH() && loader->LoadHits()) {
271       AliError("Failed to load hits");
272       return;
273     }
274     TTree*   hitsTree = loader->TreeH();
275     TBranch* hitsBranch = hitsTree->GetBranch(fFMD->GetName());
276     if (!hitsBranch) { 
277       AliError("Failed to get hits branch in tree");
278       return;
279     }
280     // Check that we can make the output digits - This must come
281     // before AliFMD::SetBranchAddress
282     TTree* outTree = MakeOutputTree(loader);
283     if (!outTree) { 
284       AliError("Failed to get output tree");
285       return;
286     }
287     AliFMDDebug(5, ("Output tree name for %s is '%s'", 
288                     (fOutput == kDigits ? "digits" : "sdigits"),
289                     outTree->GetName()));
290     if (AliLog::GetDebugLevel("FMD","") >= 5) {
291       TFile* file = outTree->GetCurrentFile();
292       if (!file) {
293         AliWarning("Output tree has no file!");
294       }
295       else { 
296         AliFMDDebug(5, ("Output tree file %s content:", file->GetName()));
297         file->ls();
298       }
299     }
300
301     // Set-up the branch addresses 
302     fFMD->SetTreeAddress();
303     
304     // Now sum all contributions in cache 
305     SumContributions(hitsBranch);
306     loader->UnloadHits();
307
308     // And now digitize the hits 
309     DigitizeHits();
310     
311     // Write digits to tree
312     Int_t write = outTree->Fill();
313     AliFMDDebug(1, ("Wrote %d bytes to digit tree", write));
314
315     // Store the digits
316     StoreDigits(loader);
317
318   }  
319 }
320
321 //____________________________________________________________________
322 TTree*
323 AliFMDHitDigitizer::MakeOutputTree(AliLoader* loader)
324 {
325   if (fOutput == kDigits) 
326     return AliFMDBaseDigitizer::MakeOutputTree(loader);
327   
328   AliFMDDebug(5, ("Making sdigits tree"));
329   loader->LoadSDigits("UPDATE"); // RECREATE");
330   TTree* out = loader->TreeS();
331   if (!out) loader->MakeTree("S");
332   out = loader->TreeS(); 
333   if (out) { 
334     out->Reset();
335     fFMD->MakeBranch("S");
336   }
337   return out;
338 }
339
340
341 //____________________________________________________________________
342 void
343 AliFMDHitDigitizer::SumContributions(TBranch* hitsBranch) 
344 {
345   // Sum energy deposited contributions from each hit in a cache
346   // (fEdep).  
347   
348   // Clear array of deposited energies 
349   fEdep.Reset();
350   
351   // Get a list of hits from the FMD manager 
352   AliFMDDebug(5, ("Get array of FMD hits"));
353   TClonesArray *fmdHits = fFMD->Hits();
354   
355
356   // Get number of entries in the tree 
357   AliFMDDebug(5, ("Get # of tracks"));
358   Int_t ntracks  = Int_t(hitsBranch->GetEntries());
359   AliFMDDebug(5, ("We got %d tracks", ntracks));
360
361   Int_t read = 0;
362   // Loop over the tracks in the 
363   for (Int_t track = 0; track < ntracks; track++)  {
364     // Read in entry number `track' 
365     read += hitsBranch->GetEntry(track);
366     
367     // Get the number of hits 
368     Int_t nhits = fmdHits->GetEntries ();
369     for (Int_t hit = 0; hit < nhits; hit++) {
370       // Get the hit number `hit'
371       AliFMDHit* fmdHit = 
372         static_cast<AliFMDHit*>(fmdHits->UncheckedAt(hit));
373       
374       // Extract parameters 
375       AddContribution(fmdHit->Detector(),
376                       fmdHit->Ring(),
377                       fmdHit->Sector(),
378                       fmdHit->Strip(),
379                       fmdHit->Edep());
380     }  // hit loop
381   } // track loop
382   AliFMDDebug(5, ("Size of cache: %d bytes, read %d bytes", 
383                    sizeof(fEdep), read));
384 }
385
386
387 //____________________________________________________________________
388 UShort_t
389 AliFMDHitDigitizer::MakePedestal(UShort_t  detector, 
390                                  Char_t    ring, 
391                                  UShort_t  sector, 
392                                  UShort_t  strip) const 
393 {
394   // Make a pedestal 
395   if (fOutput == kSDigits) return 0;
396   return AliFMDBaseDigitizer::MakePedestal(detector, ring, sector, strip);
397 }
398
399
400
401 //____________________________________________________________________
402 void
403 AliFMDHitDigitizer::AddDigit(UShort_t  detector, 
404                              Char_t    ring,
405                              UShort_t  sector, 
406                              UShort_t  strip, 
407                              Float_t   edep, 
408                              UShort_t  count1, 
409                              Short_t   count2, 
410                              Short_t   count3,
411                              Short_t   count4) const
412 {
413   // Add a digit or summable digit
414   if (fOutput == kDigits) { 
415     AliFMDBaseDigitizer::AddDigit(detector, ring, sector, strip, 0,
416                                   count1, count2, count3, count4);
417     return;
418   }
419   if (edep <= 0) return;
420   if (count1 == 0 && count2 <= 0 && count3 <= 0 && count4 <= 0) return;
421   fFMD->AddSDigitByFields(detector, ring, sector, strip, edep,
422                           count1, count2, count3, count4);
423 }
424
425 //____________________________________________________________________
426 void
427 AliFMDHitDigitizer::CheckDigit(AliFMDDigit*    digit,
428                                UShort_t        nhits,
429                                const TArrayI&  counts) 
430 {
431   // Check that digit is consistent
432   AliFMDParameters* param = AliFMDParameters::Instance();
433   UShort_t          det   = digit->Detector();
434   Char_t            ring  = digit->Ring();
435   UShort_t          sec   = digit->Sector();
436   UShort_t          str   = digit->Strip();
437   Float_t           mean  = param->GetPedestal(det,ring,sec,str);
438   Float_t           width = param->GetPedestalWidth(det,ring,sec,str);
439   UShort_t          range = param->GetVA1MipRange();
440   UShort_t          size  = param->GetAltroChannelSize();
441   Int_t             integral = counts[0];
442   if (counts[1] >= 0) integral += counts[1];
443   if (counts[2] >= 0) integral += counts[2];
444   if (counts[3] >= 0) integral += counts[3];
445   integral -= Int_t(mean + 2 * width);
446   if (integral < 0) integral = 0;
447   
448   Float_t convF = Float_t(range) / size;
449   Float_t mips  = integral * convF;
450   if (mips > Float_t(nhits) + .5 || mips < Float_t(nhits) - .5) 
451     Warning("CheckDigit", "Digit -> %4.2f MIPS != %d +/- .5 hits", 
452             mips, nhits);
453 }
454
455 //____________________________________________________________________
456 void
457 AliFMDHitDigitizer::StoreDigits(AliLoader* loader)
458 {
459   if (fOutput == kDigits) { 
460     AliFMDBaseDigitizer::StoreDigits(loader);
461     return;
462   }
463   AliFMDDebug(5, ("Storing %d sdigits",   fFMD->SDigits()->GetEntries()));
464   // Write the digits to disk 
465   loader->WriteSDigits("OVERWRITE");
466   loader->UnloadSDigits();
467   // Reset the digits in the AliFMD object 
468   fFMD->ResetSDigits();
469 }
470
471
472 //____________________________________________________________________
473 //
474 // EOF
475 // 
476
477
478
479