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