]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDDigitizer.cxx
Calibration object splitted in: pedestal + E calib + reco parameters
[u/mrichter/AliRoot.git] / FMD / AliFMDDigitizer.cxx
CommitLineData
4347b38f 1/**************************************************************************
2 * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
66d2ede1 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 **************************************************************************/
4347b38f 15/* $Id$ */
c2fc1258 16/** @file AliFMDDigitizer.cxx
17 @author Christian Holm Christensen <cholm@nbi.dk>
18 @date Mon Mar 27 12:38:26 2006
19 @brief FMD Digitizers implementation
02a27b50 20 @ingroup FMD_sim
c2fc1258 21*/
4347b38f 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//
4347b38f 83//
e802be3e 84// The shaping function of the VA1_ALICE is generally given by
4347b38f 85//
e802be3e 86// f(x) = A(1 - exp(-Bx))
4347b38f 87//
e802be3e 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.
4347b38f 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
e802be3e 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
56b1929b 199#include <TTree.h> // ROOT_TTree
e939a978 200#include <TRandom.h> // ROOT_TRandom
f95a63c4 201// #include <AliLog.h> // ALILOG_H
202#include "AliFMDDebug.h" // Better debug macros
e802be3e 203#include "AliFMDDigitizer.h" // ALIFMDDIGITIZER_H
204#include "AliFMD.h" // ALIFMD_H
6169f936 205// #include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
206// #include "AliFMDDetector.h" // ALIFMDDETECTOR_H
207// #include "AliFMDRing.h" // ALIFMDRING_H
208// #include "AliFMDHit.h" // ALIFMDHIT_H
e802be3e 209#include "AliFMDDigit.h" // ALIFMDDIGIT_H
8f6ee336 210#include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
56b1929b 211#include <AliRunDigitizer.h> // ALIRUNDIGITIZER_H
212#include <AliRun.h> // ALIRUN_H
213#include <AliLoader.h> // ALILOADER_H
214#include <AliRunLoader.h> // ALIRUNLOADER_H
4347b38f 215
4347b38f 216//====================================================================
925e6570 217ClassImp(AliFMDDigitizer)
66d2ede1 218
4347b38f 219//____________________________________________________________________
220AliFMDDigitizer::AliFMDDigitizer()
221 : AliFMDBaseDigitizer()
222{
223 // Default ctor - don't use it
224}
66d2ede1 225
4347b38f 226//____________________________________________________________________
227AliFMDDigitizer::AliFMDDigitizer(AliRunDigitizer* manager)
228 : AliFMDBaseDigitizer(manager)
229{
230 // Normal CTOR
f95a63c4 231 AliFMDDebug(1, (" processed"));
4347b38f 232}
66d2ede1 233
4347b38f 234//____________________________________________________________________
235void
236AliFMDDigitizer::Exec(Option_t*)
237{
238 // Get the output manager
239 TString outFolder(fManager->GetOutputFolderName());
240 AliRunLoader* out =
241 AliRunLoader::GetRunLoader(outFolder.Data());
242 // Get the FMD output manager
243 AliLoader* outFMD = out->GetLoader("FMDLoader");
244
245 // Get the input loader
246 TString inFolder(fManager->GetInputFolderName(0));
247 fRunLoader =
a3537838 248 AliRunLoader::GetRunLoader(inFolder.Data());
4347b38f 249 if (!fRunLoader) {
250 AliError("Can not find Run Loader for input stream 0");
251 return;
252 }
253 // Get the AliRun object
254 if (!fRunLoader->GetAliRun()) fRunLoader->LoadgAlice();
255
256 // Get the AliFMD object
8f6ee336 257 AliFMD* fmd =
258 static_cast<AliFMD*>(fRunLoader->GetAliRun()->GetDetector("FMD"));
4347b38f 259 if (!fmd) {
260 AliError("Can not get FMD from gAlice");
261 return;
262 }
263
264 Int_t nFiles= fManager->GetNinputs();
265 for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
f95a63c4 266 AliFMDDebug(1, (" Digitizing event number %d",
4347b38f 267 fManager->GetOutputEventNr()));
268 // Get the current loader
269 fRunLoader =
270 AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
271 if (!fRunLoader) Fatal("Exec", "no run loader");
272 // Cache contriutions
273 SumContributions(fmd);
274 }
275 // Digitize the event
276 DigitizeHits(fmd);
277
278 // Load digits from the tree
279 outFMD->LoadDigits("update");
280 // Get the tree of digits
281 TTree* digitTree = outFMD->TreeD();
282 if (!digitTree) {
283 outFMD->MakeTree("D");
284 digitTree = outFMD->TreeD();
285 }
286 digitTree->Reset();
287 // Make a branch in the tree
288 TClonesArray* digits = fmd->Digits();
289 fmd->MakeBranchInTree(digitTree, fmd->GetName(), &(digits), 4000, 0);
290 // TBranch* digitBranch = digitTree->GetBranch(fmd->GetName());
291 // Fill the tree
a1f80595 292 Int_t write = 0;
293 write = digitTree->Fill();
f95a63c4 294 AliFMDDebug(1, ("Wrote %d bytes to digit tree", write));
66d2ede1 295
4347b38f 296 // Write the digits to disk
297 outFMD->WriteDigits("OVERWRITE");
298 outFMD->UnloadHits();
299 outFMD->UnloadDigits();
4110645f 300
4347b38f 301 // Reset the digits in the AliFMD object
302 fmd->ResetDigits();
303}
3d44ce66 304
88cb7938 305
4347b38f 306//____________________________________________________________________
307UShort_t
8f6ee336 308AliFMDDigitizer::MakePedestal(UShort_t detector,
309 Char_t ring,
310 UShort_t sector,
311 UShort_t strip) const
4347b38f 312{
088f8e79 313 // Make a pedestal
8f6ee336 314 AliFMDParameters* param =AliFMDParameters::Instance();
315 Float_t mean =param->GetPedestal(detector,ring,sector,strip);
316 Float_t width =param->GetPedestalWidth(detector,ring,sector,strip);
317 return UShort_t(TMath::Max(gRandom->Gaus(mean, width), 0.));
4347b38f 318}
88cb7938 319
4347b38f 320//____________________________________________________________________
321void
322AliFMDDigitizer::AddDigit(AliFMD* fmd,
323 UShort_t detector,
324 Char_t ring,
325 UShort_t sector,
326 UShort_t strip,
327 Float_t /* edep */,
328 UShort_t count1,
329 Short_t count2,
330 Short_t count3) const
331{
088f8e79 332 // Add a digit
69b696b9 333 fmd->AddDigitByFields(detector, ring, sector, strip, count1, count2, count3);
4347b38f 334}
3d44ce66 335
4347b38f 336//____________________________________________________________________
337void
8f6ee336 338AliFMDDigitizer::CheckDigit(AliFMDDigit* digit,
4347b38f 339 UShort_t nhits,
340 const TArrayI& counts)
341{
088f8e79 342 // Check that digit is consistent
8f6ee336 343 AliFMDParameters* param = AliFMDParameters::Instance();
344 UShort_t det = digit->Detector();
345 Char_t ring = digit->Ring();
346 UShort_t sec = digit->Sector();
347 UShort_t str = digit->Strip();
348 Float_t mean = param->GetPedestal(det,ring,sec,str);
349 Float_t width = param->GetPedestalWidth(det,ring,sec,str);
350 UShort_t range = param->GetVA1MipRange();
351 UShort_t size = param->GetAltroChannelSize();
352 Int_t integral = counts[0];
4347b38f 353 if (counts[1] >= 0) integral += counts[1];
354 if (counts[2] >= 0) integral += counts[2];
8f6ee336 355 integral -= Int_t(mean + 2 * width);
4347b38f 356 if (integral < 0) integral = 0;
357
8f6ee336 358 Float_t convF = Float_t(range) / size;
088f8e79 359 Float_t mips = integral * convF;
4347b38f 360 if (mips > Float_t(nhits) + .5 || mips < Float_t(nhits) - .5)
361 Warning("CheckDigit", "Digit -> %4.2f MIPS != %d +/- .5 hits",
362 mips, nhits);
363}
66d2ede1 364
4347b38f 365//____________________________________________________________________
366//
367// EOF
368//
88cb7938 369
370
371
66d2ede1 372