]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDmcmSim.cxx
Better starting value for estimate of covariance matrix (Maksym, Silvia)
[u/mrichter/AliRoot.git] / TRD / AliTRDmcmSim.cxx
CommitLineData
dfd03fc3 1/**************************************************************************
2 * Copyright(c) 1998-1999, 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
5eba8ada 16/* $Id$ */
17
dfd03fc3 18///////////////////////////////////////////////////////////////////////////////
19// //
20// TRD MCM (Multi Chip Module) simulator //
ce4786b9 21// which simulates the TRAP processing after the AD-conversion. //
22// The relevant parameters (i.e. configuration settings of the TRAP) //
23// are taken from AliTRDtrapConfig. //
dfd03fc3 24// //
25///////////////////////////////////////////////////////////////////////////////
26
ce4786b9 27#include <iostream>
28#include <iomanip>
ecf39416 29
ce4786b9 30#include "TCanvas.h"
31#include "TH1F.h"
32#include "TH2F.h"
33#include "TGraph.h"
34#include "TLine.h"
35#include "TRandom.h"
36#include "TClonesArray.h"
ce51199c 37#include "TMath.h"
c93255fe 38#include <TTree.h>
0c349049 39
dfd03fc3 40#include "AliLog.h"
b0a41e80 41#include "AliRunLoader.h"
42#include "AliLoader.h"
0c349049 43
dfd03fc3 44#include "AliTRDfeeParam.h"
b0a41e80 45#include "AliTRDtrapConfig.h"
4cc89512 46#include "AliTRDdigitsManager.h"
b65e5048 47#include "AliTRDarrayADC.h"
40bd6ee4 48#include "AliTRDarrayDictionary.h"
52c19022 49#include "AliTRDtrackletMCM.h"
b0a41e80 50#include "AliTRDmcmSim.h"
1d93b218 51
dfd03fc3 52ClassImp(AliTRDmcmSim)
53
40bd6ee4 54Bool_t AliTRDmcmSim::fgApplyCut = kTRUE;
ce4786b9 55Int_t AliTRDmcmSim::fgAddBaseline = 0;
56
5f006bd7 57const Int_t AliTRDmcmSim::fgkFormatIndex = std::ios_base::xalloc();
ce4786b9 58
59const Int_t AliTRDmcmSim::fgkNADC = AliTRDfeeParam::GetNadcMcm();
5f006bd7 60const UShort_t AliTRDmcmSim::fgkFPshifts[4] = {11, 14, 17, 21};
ce4786b9 61
62
5f006bd7 63AliTRDmcmSim::AliTRDmcmSim() :
ce4786b9 64 TObject(),
65 fInitialized(kFALSE),
66 fDetector(-1),
67 fRobPos(-1),
68 fMcmPos(-1),
69 fRow (-1),
70 fNTimeBin(-1),
71 fADCR(NULL),
72 fADCF(NULL),
73 fMCMT(NULL),
74 fTrackletArray(NULL),
75 fZSMap(NULL),
6b094867 76 fTrklBranchName("mcmtrklbranch"),
ce4786b9 77 fFeeParam(NULL),
78 fTrapConfig(NULL),
79 fDigitsManager(NULL),
80 fPedAcc(NULL),
81 fGainCounterA(NULL),
82 fGainCounterB(NULL),
83 fTailAmplLong(NULL),
84 fTailAmplShort(NULL),
85 fNHits(0),
86 fFitReg(NULL)
dfd03fc3 87{
88 //
b0a41e80 89 // AliTRDmcmSim default constructor
dfd03fc3 90 // By default, nothing is initialized.
91 // It is necessary to issue Init before use.
54d34aac 92
93 for (Int_t iDict = 0; iDict < 3; iDict++)
94 fDict[iDict] = 0x0;
95
96 fFitPtr[0] = 0;
97 fFitPtr[1] = 0;
98 fFitPtr[2] = 0;
99 fFitPtr[3] = 0;
dfd03fc3 100}
101
5f006bd7 102AliTRDmcmSim::~AliTRDmcmSim()
dfd03fc3 103{
104 //
105 // AliTRDmcmSim destructor
106 //
0c349049 107
b0a41e80 108 if(fInitialized) {
ce4786b9 109 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
110 delete [] fADCR[iAdc];
111 delete [] fADCF[iAdc];
dfd03fc3 112 }
16e077d0 113 delete [] fADCR;
114 delete [] fADCF;
ce4786b9 115 delete [] fZSMap;
1d93b218 116 delete [] fMCMT;
5f006bd7 117
b0a41e80 118 delete [] fPedAcc;
119 delete [] fGainCounterA;
120 delete [] fGainCounterB;
121 delete [] fTailAmplLong;
122 delete [] fTailAmplShort;
123 delete [] fFitReg;
5f006bd7 124
b0a41e80 125 fTrackletArray->Delete();
126 delete fTrackletArray;
1d93b218 127 }
dfd03fc3 128}
129
5f006bd7 130void AliTRDmcmSim::Init( Int_t det, Int_t robPos, Int_t mcmPos, Bool_t /* newEvent */ )
dfd03fc3 131{
0c349049 132 //
ce4786b9 133 // Initialize the class with new MCM position information
134 // memory is allocated in the first initialization
0c349049 135 //
5f006bd7 136
b0a41e80 137 if (!fInitialized) {
138 fFeeParam = AliTRDfeeParam::Instance();
139 fTrapConfig = AliTRDtrapConfig::Instance();
b0a41e80 140 }
141
142 fDetector = det;
0c349049 143 fRobPos = robPos;
144 fMcmPos = mcmPos;
dfd03fc3 145 fRow = fFeeParam->GetPadRowFromMCM( fRobPos, fMcmPos );
5f006bd7 146
b0a41e80 147 if (!fInitialized) {
ce4786b9 148 fADCR = new Int_t *[fgkNADC];
149 fADCF = new Int_t *[fgkNADC];
150 fZSMap = new Int_t [fgkNADC];
151 fGainCounterA = new UInt_t[fgkNADC];
152 fGainCounterB = new UInt_t[fgkNADC];
759042e7 153 fNTimeBin = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC13CPUA, fDetector, fRobPos, fMcmPos);
ce4786b9 154 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
155 fADCR[iAdc] = new Int_t[fNTimeBin];
156 fADCF[iAdc] = new Int_t[fNTimeBin];
dfd03fc3 157 }
5f006bd7 158
b0a41e80 159 // filter registers
ce4786b9 160 fPedAcc = new UInt_t[fgkNADC]; // accumulator for pedestal filter
161 fTailAmplLong = new UShort_t[fgkNADC];
162 fTailAmplShort = new UShort_t[fgkNADC];
5f006bd7 163
b0a41e80 164 // tracklet calculation
5f006bd7 165 fFitReg = new FitReg_t[fgkNADC];
ce4786b9 166 fTrackletArray = new TClonesArray("AliTRDtrackletMCM", fgkMaxTracklets);
5f006bd7 167
ce4786b9 168 fMCMT = new UInt_t[fgkMaxTracklets];
dfd03fc3 169 }
170
b0a41e80 171 fInitialized = kTRUE;
172
173 Reset();
174}
175
176void AliTRDmcmSim::Reset()
177{
178 // Resets the data values and internal filter registers
179 // by re-initialising them
180
5f006bd7 181 if( !CheckInitialized() )
ce4786b9 182 return;
5896bc23 183
ce4786b9 184 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
185 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
186 fADCR[iAdc][it] = 0;
187 fADCF[iAdc][it] = 0;
dfd03fc3 188 }
ce4786b9 189 fZSMap[iAdc] = -1; // Default unread, low active bit mask
190 fGainCounterA[iAdc] = 0;
191 fGainCounterB[iAdc] = 0;
dfd03fc3 192 }
5f006bd7 193
ce4786b9 194 for(Int_t i = 0; i < fgkMaxTracklets; i++) {
1d93b218 195 fMCMT[i] = 0;
196 }
ce4786b9 197
198 for (Int_t iDict = 0; iDict < 3; iDict++)
199 fDict[iDict] = 0x0;
5f006bd7 200
b0a41e80 201 FilterPedestalInit();
202 FilterGainInit();
ce4786b9 203 FilterTailInit();
b0a41e80 204}
1d93b218 205
5f006bd7 206void AliTRDmcmSim::SetNTimebins(Int_t ntimebins)
4ff7ed2b 207{
5f006bd7 208 // Reallocate memory if a change in the number of timebins
ce4786b9 209 // is needed (should not be the case for real data)
210
5f006bd7 211 if( !CheckInitialized() )
ce4786b9 212 return;
213
4ff7ed2b 214 fNTimeBin = ntimebins;
ce4786b9 215 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
759042e7 216 delete [] fADCR[iAdc];
217 delete [] fADCF[iAdc];
ce4786b9 218 fADCR[iAdc] = new Int_t[fNTimeBin];
219 fADCF[iAdc] = new Int_t[fNTimeBin];
4ff7ed2b 220 }
221}
222
5f006bd7 223Bool_t AliTRDmcmSim::LoadMCM(AliRunLoader* const runloader, Int_t det, Int_t rob, Int_t mcm)
b0a41e80 224{
ce4786b9 225 // loads the ADC data as obtained from the digitsManager for the specified MCM.
5f006bd7 226 // This method is meant for rare execution, e.g. in the visualization. When called
227 // frequently use SetData(...) instead.
b0a41e80 228
64e3d742 229 Init(det, rob, mcm);
b0a41e80 230
231 if (!runloader) {
232 AliError("No Runloader given");
233 return kFALSE;
234 }
235
236 AliLoader *trdLoader = runloader->GetLoader("TRDLoader");
237 if (!trdLoader) {
238 AliError("Could not get TRDLoader");
239 return kFALSE;
240 }
241
5eba8ada 242 Bool_t retval = kTRUE;
b0a41e80 243 trdLoader->LoadDigits();
40bd6ee4 244 fDigitsManager = 0x0;
b0a41e80 245 AliTRDdigitsManager *digMgr = new AliTRDdigitsManager();
246 digMgr->SetSDigits(0);
247 digMgr->CreateArrays();
248 digMgr->ReadDigits(trdLoader->TreeD());
249 AliTRDarrayADC *digits = (AliTRDarrayADC*) digMgr->GetDigits(det);
5eba8ada 250 if (digits->HasData()) {
251 digits->Expand();
252
5896bc23 253 if (fNTimeBin != digits->GetNtime()) {
ce4786b9 254 AliWarning(Form("Changing no. of timebins from %i to %i", fNTimeBin, digits->GetNtime()));
4ff7ed2b 255 SetNTimebins(digits->GetNtime());
5896bc23 256 }
4ff7ed2b 257
ce4786b9 258 SetData(digits);
b0a41e80 259 }
5f006bd7 260 else
5eba8ada 261 retval = kFALSE;
5f006bd7 262
b0a41e80 263 delete digMgr;
5f006bd7 264
4ff7ed2b 265 return retval;
b0a41e80 266}
267
268void AliTRDmcmSim::NoiseTest(Int_t nsamples, Int_t mean, Int_t sigma, Int_t inputGain, Int_t inputTail)
269{
5f006bd7 270 // This function can be used to test the filters.
b0a41e80 271 // It feeds nsamples of ADC values with a gaussian distribution specified by mean and sigma.
272 // The filter chain implemented here consists of:
273 // Pedestal -> Gain -> Tail
5f006bd7 274 // With inputGain and inputTail the input to the gain and tail filter, respectively,
275 // can be chosen where
b0a41e80 276 // 0: noise input
277 // 1: pedestal output
278 // 2: gain output
5f006bd7 279 // The input has to be chosen from a stage before.
280 // The filter behaviour is controlled by the TRAP parameters from AliTRDtrapConfig in the
b0a41e80 281 // same way as in normal simulation.
282 // The functions produces four histograms with the values at the different stages.
283
5f006bd7 284 if( !CheckInitialized() )
ce4786b9 285 return;
286
287 TString nameInputGain;
5f006bd7 288 TString nameInputTail;
ce4786b9 289
290 switch (inputGain) {
291 case 0:
292 nameInputGain = "Noise";
293 break;
294
295 case 1:
296 nameInputGain = "Pedestal";
297 break;
298
299 default:
300 AliError("Undefined input to tail cancellation filter");
301 return;
302 }
303
304 switch (inputTail) {
305 case 0:
306 nameInputTail = "Noise";
307 break;
308
309 case 1:
310 nameInputTail = "Pedestal";
311 break;
312
313 case 2:
314 nameInputTail = "Gain";
315 break;
316
317 default:
318 AliError("Undefined input to tail cancellation filter");
319 return;
320 }
321
b0a41e80 322 TH1F *h = new TH1F("noise", "Gaussian Noise;sample;ADC count",
323 nsamples, 0, nsamples);
ce4786b9 324 TH1F *hfp = new TH1F("ped", "Noise #rightarrow Pedestal filter;sample;ADC count", nsamples, 0, nsamples);
5f006bd7 325 TH1F *hfg = new TH1F("gain",
326 (nameInputGain + "#rightarrow Gain;sample;ADC count").Data(),
ce4786b9 327 nsamples, 0, nsamples);
5f006bd7 328 TH1F *hft = new TH1F("tail",
329 (nameInputTail + "#rightarrow Tail;sample;ADC count").Data(),
ce4786b9 330 nsamples, 0, nsamples);
b0a41e80 331 h->SetStats(kFALSE);
332 hfp->SetStats(kFALSE);
333 hfg->SetStats(kFALSE);
334 hft->SetStats(kFALSE);
5f006bd7 335
b0a41e80 336 Int_t value; // ADC count with noise (10 bit)
337 Int_t valuep; // pedestal filter output (12 bit)
338 Int_t valueg; // gain filter output (12 bit)
339 Int_t valuet; // tail filter value (12 bit)
5f006bd7 340
b0a41e80 341 for (Int_t i = 0; i < nsamples; i++) {
5f006bd7 342 value = (Int_t) gRandom->Gaus(mean, sigma); // generate noise with gaussian distribution
b0a41e80 343 h->SetBinContent(i, value);
344
345 valuep = FilterPedestalNextSample(1, 0, ((Int_t) value) << 2);
5f006bd7 346
b0a41e80 347 if (inputGain == 0)
348 valueg = FilterGainNextSample(1, ((Int_t) value) << 2);
5f006bd7 349 else
350 valueg = FilterGainNextSample(1, valuep);
351
b0a41e80 352 if (inputTail == 0)
353 valuet = FilterTailNextSample(1, ((Int_t) value) << 2);
354 else if (inputTail == 1)
5f006bd7 355 valuet = FilterTailNextSample(1, valuep);
b0a41e80 356 else
5f006bd7 357 valuet = FilterTailNextSample(1, valueg);
b0a41e80 358
359 hfp->SetBinContent(i, valuep >> 2);
360 hfg->SetBinContent(i, valueg >> 2);
361 hft->SetBinContent(i, valuet >> 2);
362 }
363
5f006bd7 364 TCanvas *c = new TCanvas;
b0a41e80 365 c->Divide(2,2);
366 c->cd(1);
367 h->Draw();
368 c->cd(2);
369 hfp->Draw();
370 c->cd(3);
371 hfg->Draw();
372 c->cd(4);
373 hft->Draw();
dfd03fc3 374}
375
ce4786b9 376Bool_t AliTRDmcmSim::CheckInitialized() const
ecf39416 377{
0c349049 378 //
379 // Check whether object is initialized
380 //
381
5f006bd7 382 if( ! fInitialized )
ce4786b9 383 AliError(Form ("AliTRDmcmSim is not initialized but function other than Init() is called."));
384
ecf39416 385 return fInitialized;
386}
387
ab9f7002 388void AliTRDmcmSim::Print(Option_t* const option) const
b0a41e80 389{
390 // Prints the data stored and/or calculated for this MCM.
5f006bd7 391 // The output is controlled by option which can be a sequence of any of
b0a41e80 392 // the following characters:
393 // R - prints raw ADC data
5f006bd7 394 // F - prints filtered data
b0a41e80 395 // H - prints detected hits
396 // T - prints found tracklets
5f006bd7 397 // The later stages are only meaningful after the corresponding calculations
b0a41e80 398 // have been performed.
399
5f006bd7 400 if ( !CheckInitialized() )
ce4786b9 401 return;
402
b0a41e80 403 printf("MCM %i on ROB %i in detector %i\n", fMcmPos, fRobPos, fDetector);
404
405 TString opt = option;
ce4786b9 406 if (opt.Contains("R") || opt.Contains("F")) {
407 std::cout << *this;
1d93b218 408 }
409
b0a41e80 410 if (opt.Contains("H")) {
411 printf("Found %i hits:\n", fNHits);
412 for (Int_t iHit = 0; iHit < fNHits; iHit++) {
413 printf("Hit %3i in timebin %2i, ADC %2i has charge %3i and position %3i\n",
ab9f7002 414 iHit, fHits[iHit].fTimebin, fHits[iHit].fChannel, fHits[iHit].fQtot, fHits[iHit].fYpos);
b0a41e80 415 }
1d93b218 416 }
1d93b218 417
b0a41e80 418 if (opt.Contains("T")) {
419 printf("Tracklets:\n");
420 for (Int_t iTrkl = 0; iTrkl < fTrackletArray->GetEntriesFast(); iTrkl++) {
421 printf("tracklet %i: 0x%08x\n", iTrkl, ((AliTRDtrackletMCM*) (*fTrackletArray)[iTrkl])->GetTrackletWord());
422 }
1d93b218 423 }
b0a41e80 424}
1d93b218 425
5f006bd7 426void AliTRDmcmSim::Draw(Option_t* const option)
b0a41e80 427{
428 // Plots the data stored in a 2-dim. timebin vs. ADC channel plot.
5f006bd7 429 // The option selects what data is plotted and can be a sequence of
b0a41e80 430 // the following characters:
431 // R - plot raw data (default)
432 // F - plot filtered data (meaningless if R is specified)
433 // In addition to the ADC values:
5f006bd7 434 // H - plot hits
b0a41e80 435 // T - plot tracklets
436
5f006bd7 437 if( !CheckInitialized() )
ce4786b9 438 return;
439
b0a41e80 440 TString opt = option;
441
442 TH2F *hist = new TH2F("mcmdata", Form("Data of MCM %i on ROB %i in detector %i", \
443 fMcmPos, fRobPos, fDetector), \
ce4786b9 444 fgkNADC, -0.5, fgkNADC-.5, fNTimeBin, -.5, fNTimeBin-.5);
b0a41e80 445 hist->GetXaxis()->SetTitle("ADC Channel");
446 hist->GetYaxis()->SetTitle("Timebin");
447 hist->SetStats(kFALSE);
448
449 if (opt.Contains("R")) {
450 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
ce4786b9 451 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
b0a41e80 452 hist->SetBinContent(iAdc+1, iTimeBin+1, fADCR[iAdc][iTimeBin] >> fgkAddDigits);
453 }
1d93b218 454 }
b0a41e80 455 }
456 else {
457 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
ce4786b9 458 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
b0a41e80 459 hist->SetBinContent(iAdc+1, iTimeBin+1, fADCF[iAdc][iTimeBin] >> fgkAddDigits);
460 }
1d93b218 461 }
1d93b218 462 }
b0a41e80 463 hist->Draw("colz");
1d93b218 464
b0a41e80 465 if (opt.Contains("H")) {
466 TGraph *grHits = new TGraph();
467 for (Int_t iHit = 0; iHit < fNHits; iHit++) {
5f006bd7 468 grHits->SetPoint(iHit,
469 fHits[iHit].fChannel + 1 + fHits[iHit].fYpos/256.,
ab9f7002 470 fHits[iHit].fTimebin);
b0a41e80 471 }
472 grHits->Draw("*");
473 }
1d93b218 474
b0a41e80 475 if (opt.Contains("T")) {
476 TLine *trklLines = new TLine[4];
64e3d742 477 for (Int_t iTrkl = 0; iTrkl < fTrackletArray->GetEntries(); iTrkl++) {
b0a41e80 478 AliTRDtrackletMCM *trkl = (AliTRDtrackletMCM*) (*fTrackletArray)[iTrkl];
ce4786b9 479 Float_t padWidth = 0.635 + 0.03 * (fDetector % 6);
480 Float_t offset = padWidth/256. * ((((((fRobPos & 0x1) << 2) + (fMcmPos & 0x3)) * 18) << 8) - ((18*4*2 - 18*2 - 3) << 7)); // revert adding offset in FitTracklet
5ac2e3b1 481 Int_t ndrift = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrNdrift, fDetector, fRobPos, fMcmPos) >> 5;
5f006bd7 482 Float_t slope = trkl->GetdY() * 140e-4 / ndrift;
ce4786b9 483
759042e7 484 Int_t t0 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS, fDetector, fRobPos, fMcmPos);
485 Int_t t1 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFE, fDetector, fRobPos, fMcmPos);
ce4786b9 486
487 trklLines[iTrkl].SetX1((offset - (trkl->GetY() - slope * t0)) / padWidth); // ??? sign?
488 trklLines[iTrkl].SetY1(t0);
489 trklLines[iTrkl].SetX2((offset - (trkl->GetY() - slope * t1)) / padWidth); // ??? sign?
490 trklLines[iTrkl].SetY2(t1);
b0a41e80 491 trklLines[iTrkl].SetLineColor(2);
492 trklLines[iTrkl].SetLineWidth(2);
493 printf("Tracklet %i: y = %f, dy = %f, offset = %f\n", iTrkl, trkl->GetY(), (trkl->GetdY() * 140e-4), offset);
494 trklLines[iTrkl].Draw();
495 }
496 }
1d93b218 497}
498
ce4786b9 499void AliTRDmcmSim::SetData( Int_t adc, Int_t* const data )
dfd03fc3 500{
0c349049 501 //
dfd03fc3 502 // Store ADC data into array of raw data
0c349049 503 //
dfd03fc3 504
ecf39416 505 if( !CheckInitialized() ) return;
dfd03fc3 506
ce4786b9 507 if( adc < 0 || adc >= fgkNADC ) {
508 AliError(Form ("Error: ADC %i is out of range (0 .. %d).", adc, fgkNADC-1));
dfd03fc3 509 return;
510 }
511
4ff7ed2b 512 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
ce4786b9 513 fADCR[adc][it] = (Int_t) (data[it]) << fgkAddDigits;
514 fADCF[adc][it] = (Int_t) (data[it]) << fgkAddDigits;
dfd03fc3 515 }
516}
517
ce4786b9 518void AliTRDmcmSim::SetData( Int_t adc, Int_t it, Int_t data )
dfd03fc3 519{
0c349049 520 //
dfd03fc3 521 // Store ADC data into array of raw data
0c349049 522 //
dfd03fc3 523
ecf39416 524 if( !CheckInitialized() ) return;
dfd03fc3 525
ce4786b9 526 if( adc < 0 || adc >= fgkNADC ) {
527 AliError(Form ("Error: ADC %i is out of range (0 .. %d).", adc, fgkNADC-1));
dfd03fc3 528 return;
529 }
530
ce4786b9 531 fADCR[adc][it] = data << fgkAddDigits;
532 fADCF[adc][it] = data << fgkAddDigits;
b0a41e80 533}
534
6b094867 535void AliTRDmcmSim::SetData(AliTRDarrayADC* const adcArray, AliTRDdigitsManager * const digitsManager)
b0a41e80 536{
ab9f7002 537 // Set the ADC data from an AliTRDarrayADC
538
5f006bd7 539 if( !CheckInitialized() )
b0a41e80 540 return;
b0a41e80 541
40bd6ee4 542 fDigitsManager = digitsManager;
ce4786b9 543 if (fDigitsManager) {
544 for (Int_t iDict = 0; iDict < 3; iDict++) {
545 AliTRDarrayDictionary *newDict = (AliTRDarrayDictionary*) fDigitsManager->GetDictionary(fDetector, iDict);
546 if (fDict[iDict] != 0x0 && newDict != 0x0) {
5f006bd7 547
ce4786b9 548 if (fDict[iDict] == newDict)
549 continue;
40bd6ee4 550
ce4786b9 551 fDict[iDict] = newDict;
5f006bd7 552 fDict[iDict]->Expand();
ce4786b9 553 }
554 else {
555 fDict[iDict] = newDict;
556 if (fDict[iDict])
557 fDict[iDict]->Expand();
27a030ab 558 }
5f006bd7 559
560 // If there is no data, set dictionary to zero to avoid crashes
27a030ab 561 if (fDict[iDict]->GetDim() == 0) {
562 AliError(Form("Dictionary %i of det. %i has dim. 0", fDetector, iDict));
563 fDict[iDict] = 0x0;
ce4786b9 564 }
565 }
5896bc23 566 }
4ff7ed2b 567
ce4786b9 568 if (fNTimeBin != adcArray->GetNtime())
569 SetNTimebins(adcArray->GetNtime());
5f006bd7 570
ce4786b9 571 Int_t offset = (fMcmPos % 4 + 1) * 21 + (fRobPos % 2) * 84 - 1;
b0a41e80 572
573 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
ce4786b9 574 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
575 Int_t value = adcArray->GetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin);
a10c6981 576 // treat 0 as suppressed,
577 // this is not correct but reported like that from arrayADC
578 if (value <= 0 || (offset - iAdc < 1) || (offset - iAdc > 165)) {
759042e7 579 fADCR[iAdc][iTimeBin] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP, fDetector, fRobPos, fMcmPos) + (fgAddBaseline << fgkAddDigits);
580 fADCF[iAdc][iTimeBin] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP, fDetector, fRobPos, fMcmPos) + (fgAddBaseline << fgkAddDigits);
b0a41e80 581 }
582 else {
ce4786b9 583 fZSMap[iAdc] = 0;
584 fADCR[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
585 fADCF[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
b0a41e80 586 }
587 }
588 }
dfd03fc3 589}
590
6b094867 591void AliTRDmcmSim::SetDataByPad(AliTRDarrayADC* const adcArray, AliTRDdigitsManager * const digitsManager)
ce51199c 592{
5f006bd7 593 // Set the ADC data from an AliTRDarrayADC
ce51199c 594 // (by pad, to be used during initial reading in simulation)
595
5f006bd7 596 if( !CheckInitialized() )
ce51199c 597 return;
598
599 fDigitsManager = digitsManager;
600 if (fDigitsManager) {
601 for (Int_t iDict = 0; iDict < 3; iDict++) {
602 AliTRDarrayDictionary *newDict = (AliTRDarrayDictionary*) fDigitsManager->GetDictionary(fDetector, iDict);
603 if (fDict[iDict] != 0x0 && newDict != 0x0) {
5f006bd7 604
ce51199c 605 if (fDict[iDict] == newDict)
606 continue;
607
608 fDict[iDict] = newDict;
5f006bd7 609 fDict[iDict]->Expand();
ce51199c 610 }
611 else {
612 fDict[iDict] = newDict;
613 if (fDict[iDict])
614 fDict[iDict]->Expand();
615 }
5f006bd7 616
617 // If there is no data, set dictionary to zero to avoid crashes
27a030ab 618 if (fDict[iDict]->GetDim() == 0) {
619 AliError(Form("Dictionary %i of det. %i has dim. 0", fDetector, iDict));
620 fDict[iDict] = 0x0;
621 }
ce51199c 622 }
623 }
624
625 if (fNTimeBin != adcArray->GetNtime())
626 SetNTimebins(adcArray->GetNtime());
5f006bd7 627
ce51199c 628 Int_t offset = (fMcmPos % 4 + 1) * 18 + (fRobPos % 2) * 72 + 1;
629
630 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
631 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
632 Int_t value = -1;
633 Int_t pad = offset - iAdc;
5f006bd7 634 if (pad > -1 && pad < 144)
ce51199c 635 value = adcArray->GetData(GetRow(), offset - iAdc, iTimeBin);
636 // Int_t value = adcArray->GetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin);
637 if (value < 0 || (offset - iAdc < 1) || (offset - iAdc > 165)) {
759042e7 638 fADCR[iAdc][iTimeBin] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP, fDetector, fRobPos, fMcmPos) + (fgAddBaseline << fgkAddDigits);
639 fADCF[iAdc][iTimeBin] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP, fDetector, fRobPos, fMcmPos) + (fgAddBaseline << fgkAddDigits);
ce51199c 640 }
641 else {
642 fZSMap[iAdc] = 0;
643 fADCR[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
644 fADCF[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
645 }
646 }
647 }
648}
649
ce4786b9 650void AliTRDmcmSim::SetDataPedestal( Int_t adc )
dfd03fc3 651{
0c349049 652 //
dfd03fc3 653 // Store ADC data into array of raw data
0c349049 654 //
dfd03fc3 655
5f006bd7 656 if( !CheckInitialized() )
ce4786b9 657 return;
dfd03fc3 658
ce4786b9 659 if( adc < 0 || adc >= fgkNADC ) {
dfd03fc3 660 return;
661 }
662
663 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
759042e7 664 fADCR[adc][it] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP, fDetector, fRobPos, fMcmPos) + (fgAddBaseline << fgkAddDigits);
665 fADCF[adc][it] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP, fDetector, fRobPos, fMcmPos) + (fgAddBaseline << fgkAddDigits);
dfd03fc3 666 }
667}
668
ce4786b9 669Bool_t AliTRDmcmSim::GetHit(Int_t index, Int_t &channel, Int_t &timebin, Int_t &qtot, Int_t &ypos, Float_t &y, Int_t &label) const
670{
671 // retrieve the MC hit information (not available in TRAP hardware)
672
673 if (index < 0 || index >= fNHits)
674 return kFALSE;
5f006bd7 675
ce4786b9 676 channel = fHits[index].fChannel;
677 timebin = fHits[index].fTimebin;
678 qtot = fHits[index].fQtot;
679 ypos = fHits[index].fYpos;
680 y = (Float_t) ((((((fRobPos & 0x1) << 2) + (fMcmPos & 0x3)) * 18) << 8) - ((18*4*2 - 18*2 - 1) << 7) -
5f006bd7 681 (channel << 8) - ypos)
ce4786b9 682 * (0.635 + 0.03 * (fDetector % 6))
683 / 256.0;
25b41f6f 684 label = fHits[index].fLabel[0];
ce4786b9 685
686 return kTRUE;
687}
688
689Int_t AliTRDmcmSim::GetCol( Int_t adc )
dfd03fc3 690{
0c349049 691 //
dfd03fc3 692 // Return column id of the pad for the given ADC channel
0c349049 693 //
694
5f006bd7 695 if( !CheckInitialized() )
f793c83d 696 return -1;
dfd03fc3 697
ce4786b9 698 Int_t col = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adc);
5f006bd7 699 if (col < 0 || col >= fFeeParam->GetNcol())
a6d08b7f 700 return -1;
5f006bd7 701 else
a6d08b7f 702 return col;
dfd03fc3 703}
704
ce4786b9 705Int_t AliTRDmcmSim::ProduceRawStream( UInt_t *buf, Int_t bufSize, UInt_t iEv) const
dfd03fc3 706{
0c349049 707 //
dfd03fc3 708 // Produce raw data stream from this MCM and put in buf
5f006bd7 709 // Returns number of words filled, or negative value
0c349049 710 // with -1 * number of overflowed words
711 //
dfd03fc3 712
5f006bd7 713 if( !CheckInitialized() )
ce4786b9 714 return 0;
715
dfd03fc3 716 UInt_t x;
7d619a80 717 UInt_t mcmHeader = 0;
718 UInt_t adcMask = 0;
dfd03fc3 719 Int_t nw = 0; // Number of written words
720 Int_t of = 0; // Number of overflowed words
721 Int_t rawVer = fFeeParam->GetRAWversion();
722 Int_t **adc;
b0a41e80 723 Int_t nActiveADC = 0; // number of activated ADC bits in a word
dfd03fc3 724
5f006bd7 725 if( !CheckInitialized() )
ce4786b9 726 return 0;
ecf39416 727
759042e7 728 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBSF, fDetector, fRobPos, fMcmPos) != 0) // store unfiltered data
dfd03fc3 729 adc = fADCR;
5f006bd7 730 else
dfd03fc3 731 adc = fADCF;
5f006bd7 732
b0a41e80 733 // Produce ADC mask : nncc cccm mmmm mmmm mmmm mmmm mmmm 1100
734 // n : unused , c : ADC count, m : selected ADCs
7d619a80 735 if( rawVer >= 3 &&
759042e7 736 (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC15CPUA, fDetector, fRobPos, fMcmPos) & (1 << 13))) { // check for zs flag in TRAP configuration
ce4786b9 737 for( Int_t iAdc = 0 ; iAdc < fgkNADC ; iAdc++ ) {
738 if( ~fZSMap[iAdc] != 0 ) { // 0 means not suppressed
7d619a80 739 adcMask |= (1 << (iAdc+4) ); // last 4 digit reserved for 1100=0xc
740 nActiveADC++; // number of 1 in mmm....m
dfd03fc3 741 }
742 }
b0a41e80 743
7d619a80 744 if ((nActiveADC == 0) &&
759042e7 745 (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC15CPUA, fDetector, fRobPos, fMcmPos) & (1 << 8))) // check for DEH flag in TRAP configuration
7d619a80 746 return 0;
747
748 // assemble adc mask word
749 adcMask |= (1 << 30) | ( ( 0x3FFFFFFC ) & (~(nActiveADC) << 25) ) | 0xC; // nn = 01, ccccc are inverted, 0xc=1100
750 }
751
752 // MCM header
753 mcmHeader = (1<<31) | (fRobPos << 28) | (fMcmPos << 24) | ((iEv % 0x100000) << 4) | 0xC;
754 if (nw < bufSize)
755 buf[nw++] = mcmHeader;
756 else
757 of++;
758
759 // ADC mask
760 if( adcMask != 0 ) {
761 if (nw < bufSize)
762 buf[nw++] = adcMask;
763 else
dfd03fc3 764 of++;
dfd03fc3 765 }
766
767 // Produce ADC data. 3 timebins are packed into one 32 bits word
768 // In this version, different ADC channel will NOT share the same word
769
770 UInt_t aa=0, a1=0, a2=0, a3=0;
771
772 for (Int_t iAdc = 0; iAdc < 21; iAdc++ ) {
ce4786b9 773 if( rawVer>= 3 && ~fZSMap[iAdc] == 0 ) continue; // Zero Suppression, 0 means not suppressed
dfd03fc3 774 aa = !(iAdc & 1) + 2;
775 for (Int_t iT = 0; iT < fNTimeBin; iT+=3 ) {
b0a41e80 776 a1 = ((iT ) < fNTimeBin ) ? adc[iAdc][iT ] >> fgkAddDigits : 0;
777 a2 = ((iT + 1) < fNTimeBin ) ? adc[iAdc][iT+1] >> fgkAddDigits : 0;
778 a3 = ((iT + 2) < fNTimeBin ) ? adc[iAdc][iT+2] >> fgkAddDigits : 0;
ecf39416 779 x = (a3 << 22) | (a2 << 12) | (a1 << 2) | aa;
ce4786b9 780 if (nw < bufSize) {
b0a41e80 781 buf[nw++] = x;
ecf39416 782 }
783 else {
b0a41e80 784 of++;
ecf39416 785 }
dfd03fc3 786 }
787 }
788
789 if( of != 0 ) return -of; else return nw;
790}
791
ce4786b9 792Int_t AliTRDmcmSim::ProduceTrackletStream( UInt_t *buf, Int_t bufSize )
987ba9a3 793{
794 //
b0a41e80 795 // Produce tracklet data stream from this MCM and put in buf
5f006bd7 796 // Returns number of words filled, or negative value
987ba9a3 797 // with -1 * number of overflowed words
798 //
799
5f006bd7 800 if( !CheckInitialized() )
ce4786b9 801 return 0;
802
987ba9a3 803 Int_t nw = 0; // Number of written words
804 Int_t of = 0; // Number of overflowed words
5f006bd7 805
806 // Produce tracklet data. A maximum of four 32 Bit words will be written per MCM
b0a41e80 807 // fMCMT is filled continuously until no more tracklet words available
987ba9a3 808
f793c83d 809 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
5f006bd7 810 if (nw < bufSize)
f793c83d 811 buf[nw++] = ((AliTRDtrackletMCM*) (*fTrackletArray)[iTracklet])->GetTrackletWord();
5f006bd7 812 else
f793c83d 813 of++;
987ba9a3 814 }
5f006bd7 815
b0a41e80 816 if( of != 0 ) return -of; else return nw;
817}
987ba9a3 818
b0a41e80 819void AliTRDmcmSim::Filter()
820{
821 //
822 // Filter the raw ADC values. The active filter stages and their
823 // parameters are taken from AliTRDtrapConfig.
5f006bd7 824 // The raw data is stored separate from the filtered data. Thus,
825 // it is possible to run the filters on a set of raw values
b0a41e80 826 // sequentially for parameter tuning.
827 //
987ba9a3 828
5f006bd7 829 if( !CheckInitialized() )
b0a41e80 830 return;
987ba9a3 831
b0a41e80 832 // Apply filters sequentially. Bypass is handled by filters
5f006bd7 833 // since counters and internal registers may be updated even
b0a41e80 834 // if the filter is bypassed.
5f006bd7 835 // The first filter takes the data from fADCR and
836 // outputs to fADCF.
837
b0a41e80 838 // Non-linearity filter not implemented.
839 FilterPedestal();
840 FilterGain();
841 FilterTail();
842 // Crosstalk filter not implemented.
843}
987ba9a3 844
5f006bd7 845void AliTRDmcmSim::FilterPedestalInit(Int_t baseline)
b0a41e80 846{
5f006bd7 847 // Initializes the pedestal filter assuming that the input has
b0a41e80 848 // been constant for a long time (compared to the time constant).
987ba9a3 849
759042e7 850 UShort_t fptc = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPTC, fDetector, fRobPos, fMcmPos); // 0..3, 0 - fastest, 3 - slowest
987ba9a3 851
ce4786b9 852 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++)
5f006bd7 853 fPedAcc[iAdc] = (baseline << 2) * (1 << fgkFPshifts[fptc]);
987ba9a3 854}
855
b0a41e80 856UShort_t AliTRDmcmSim::FilterPedestalNextSample(Int_t adc, Int_t timebin, UShort_t value)
1d93b218 857{
b0a41e80 858 // Returns the output of the pedestal filter given the input value.
5f006bd7 859 // The output depends on the internal registers and, thus, the
b0a41e80 860 // history of the filter.
1d93b218 861
759042e7 862 UShort_t fpnp = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP, fDetector, fRobPos, fMcmPos); // 0..511 -> 0..127.75, pedestal at the output
863 UShort_t fptc = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPTC, fDetector, fRobPos, fMcmPos); // 0..3, 0 - fastest, 3 - slowest
864 UShort_t fpby = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPBY, fDetector, fRobPos, fMcmPos); // 0..1 bypass, active low
1d93b218 865
b0a41e80 866 UShort_t accumulatorShifted;
867 Int_t correction;
868 UShort_t inpAdd;
5f006bd7 869
b0a41e80 870 inpAdd = value + fpnp;
1d93b218 871
ce4786b9 872 accumulatorShifted = (fPedAcc[adc] >> fgkFPshifts[fptc]) & 0x3FF; // 10 bits
b0a41e80 873 if (timebin == 0) // the accumulator is disabled in the drift time
874 {
875 correction = (value & 0x3FF) - accumulatorShifted;
876 fPedAcc[adc] = (fPedAcc[adc] + correction) & 0x7FFFFFFF; // 31 bits
1d93b218 877 }
5f006bd7 878
ce4786b9 879 if (fpby == 0)
880 return value;
881
b0a41e80 882 if (inpAdd <= accumulatorShifted)
883 return 0;
884 else
885 {
886 inpAdd = inpAdd - accumulatorShifted;
5f006bd7 887 if (inpAdd > 0xFFF)
b0a41e80 888 return 0xFFF;
5f006bd7 889 else
b0a41e80 890 return inpAdd;
891 }
1d93b218 892}
893
b0a41e80 894void AliTRDmcmSim::FilterPedestal()
dfd03fc3 895{
0c349049 896 //
b0a41e80 897 // Apply pedestal filter
0c349049 898 //
5f006bd7 899 // As the first filter in the chain it reads data from fADCR
900 // and outputs to fADCF.
901 // It has only an effect if previous samples have been fed to
902 // find the pedestal. Currently, the simulation assumes that
b0a41e80 903 // the input has been stable for a sufficiently long time.
dfd03fc3 904
b0a41e80 905 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
ce4786b9 906 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
b0a41e80 907 fADCF[iAdc][iTimeBin] = FilterPedestalNextSample(iAdc, iTimeBin, fADCR[iAdc][iTimeBin]);
dfd03fc3 908 }
909 }
b0a41e80 910}
911
912void AliTRDmcmSim::FilterGainInit()
913{
5f006bd7 914 // Initializes the gain filter. In this case, only threshold
b0a41e80 915 // counters are reset.
dfd03fc3 916
ce4786b9 917 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
5f006bd7 918 // these are counters which in hardware continue
b0a41e80 919 // until maximum or reset
920 fGainCounterA[iAdc] = 0;
921 fGainCounterB[iAdc] = 0;
922 }
dfd03fc3 923}
924
b0a41e80 925UShort_t AliTRDmcmSim::FilterGainNextSample(Int_t adc, UShort_t value)
dfd03fc3 926{
b0a41e80 927 // Apply the gain filter to the given value.
928 // BEGIN_LATEX O_{i}(t) = #gamma_{i} * I_{i}(t) + a_{i} END_LATEX
5f006bd7 929 // The output depends on the internal registers and, thus, the
b0a41e80 930 // history of the filter.
23200400 931
759042e7 932 UShort_t fgby = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFGBY, fDetector, fRobPos, fMcmPos); // bypass, active low
933 UShort_t fgf = fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGF0 + adc), fDetector, fRobPos, fMcmPos); // 0x700 + (0 & 0x1ff);
934 UShort_t fga = fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGA0 + adc), fDetector, fRobPos, fMcmPos); // 40;
935 UShort_t fgta = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFGTA, fDetector, fRobPos, fMcmPos); // 20;
936 UShort_t fgtb = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFGTB, fDetector, fRobPos, fMcmPos); // 2060;
dfd03fc3 937
ce4786b9 938 UInt_t corr; // corrected value
dfd03fc3 939
b0a41e80 940 value &= 0xFFF;
ce4786b9 941 corr = (value * fgf) >> 11;
942 corr = corr > 0xfff ? 0xfff : corr;
943 corr = AddUintClipping(corr, fga, 12);
b0a41e80 944
5f006bd7 945 // Update threshold counters
b0a41e80 946 // not really useful as they are cleared with every new event
ce4786b9 947 if (!((fGainCounterA[adc] == 0x3FFFFFF) || (fGainCounterB[adc] == 0x3FFFFFF)))
948 // stop when full
b0a41e80 949 {
5f006bd7 950 if (corr >= fgtb)
b0a41e80 951 fGainCounterB[adc]++;
5f006bd7 952 else if (corr >= fgta)
b0a41e80 953 fGainCounterA[adc]++;
dfd03fc3 954 }
b0a41e80 955
ce4786b9 956 if (fgby == 1)
5f006bd7 957 return corr;
ce4786b9 958 else
959 return value;
dfd03fc3 960}
961
dfd03fc3 962void AliTRDmcmSim::FilterGain()
963{
b0a41e80 964 // Read data from fADCF and apply gain filter.
0c349049 965
ce4786b9 966 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
b0a41e80 967 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
968 fADCF[iAdc][iTimeBin] = FilterGainNextSample(iAdc, fADCF[iAdc][iTimeBin]);
969 }
970 }
dfd03fc3 971}
972
b0a41e80 973void AliTRDmcmSim::FilterTailInit(Int_t baseline)
dfd03fc3 974{
5f006bd7 975 // Initializes the tail filter assuming that the input has
976 // been at the baseline value (configured by FTFP) for a
b0a41e80 977 // sufficiently long time.
978
979 // exponents and weight calculated from configuration
759042e7 980 UShort_t alphaLong = 0x3ff & fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTAL, fDetector, fRobPos, fMcmPos); // the weight of the long component
981 UShort_t lambdaLong = (1 << 10) | (1 << 9) | (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTLL, fDetector, fRobPos, fMcmPos) & 0x1FF); // the multiplier
982 UShort_t lambdaShort = (0 << 10) | (1 << 9) | (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTLS, fDetector, fRobPos, fMcmPos) & 0x1FF); // the multiplier
b0a41e80 983
984 Float_t lambdaL = lambdaLong * 1.0 / (1 << 11);
985 Float_t lambdaS = lambdaShort * 1.0 / (1 << 11);
986 Float_t alphaL = alphaLong * 1.0 / (1 << 11);
987 Float_t qup, qdn;
988 qup = (1 - lambdaL) * (1 - lambdaS);
989 qdn = 1 - lambdaS * alphaL - lambdaL * (1 - alphaL);
990 Float_t kdc = qup/qdn;
991
992 Float_t kt, ql, qs;
993 UShort_t aout;
ce4786b9 994
995 if (baseline < 0)
759042e7 996 baseline = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP, fDetector, fRobPos, fMcmPos);
5f006bd7 997
b0a41e80 998 ql = lambdaL * (1 - lambdaS) * alphaL;
999 qs = lambdaS * (1 - lambdaL) * (1 - alphaL);
1000
ce4786b9 1001 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1002 Int_t value = baseline & 0xFFF;
759042e7 1003 Int_t corr = (value * fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGF0 + iAdc), fDetector, fRobPos, fMcmPos)) >> 11;
ce4786b9 1004 corr = corr > 0xfff ? 0xfff : corr;
759042e7 1005 corr = AddUintClipping(corr, fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGA0 + iAdc), fDetector, fRobPos, fMcmPos), 12);
ce4786b9 1006
1007 kt = kdc * baseline;
1008 aout = baseline - (UShort_t) kt;
1009
b0a41e80 1010 fTailAmplLong[iAdc] = (UShort_t) (aout * ql / (ql + qs));
1011 fTailAmplShort[iAdc] = (UShort_t) (aout * qs / (ql + qs));
1012 }
1013}
dfd03fc3 1014
b0a41e80 1015UShort_t AliTRDmcmSim::FilterTailNextSample(Int_t adc, UShort_t value)
1016{
5f006bd7 1017 // Returns the output of the tail filter for the given input value.
1018 // The output depends on the internal registers and, thus, the
b0a41e80 1019 // history of the filter.
1020
1021 // exponents and weight calculated from configuration
759042e7 1022 UShort_t alphaLong = 0x3ff & fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTAL, fDetector, fRobPos, fMcmPos); // the weight of the long component
1023 UShort_t lambdaLong = (1 << 10) | (1 << 9) | (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTLL, fDetector, fRobPos, fMcmPos) & 0x1FF); // the multiplier of the long component
1024 UShort_t lambdaShort = (0 << 10) | (1 << 9) | (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTLS, fDetector, fRobPos, fMcmPos) & 0x1FF); // the multiplier of the short component
b0a41e80 1025
ce4786b9 1026 // intermediate signals
1027 UInt_t aDiff;
1028 UInt_t alInpv;
b0a41e80 1029 UShort_t aQ;
ce4786b9 1030 UInt_t tmp;
5f006bd7 1031
ab9f7002 1032 UShort_t inpVolt = value & 0xFFF; // 12 bits
5f006bd7 1033
ce4786b9 1034 // add the present generator outputs
1035 aQ = AddUintClipping(fTailAmplLong[adc], fTailAmplShort[adc], 12);
1036
1037 // calculate the difference between the input and the generated signal
5f006bd7 1038 if (inpVolt > aQ)
ce4786b9 1039 aDiff = inpVolt - aQ;
5f006bd7 1040 else
ce4786b9 1041 aDiff = 0;
5f006bd7 1042
ce4786b9 1043 // the inputs to the two generators, weighted
1044 alInpv = (aDiff * alphaLong) >> 11;
5f006bd7 1045
ce4786b9 1046 // the new values of the registers, used next time
1047 // long component
1048 tmp = AddUintClipping(fTailAmplLong[adc], alInpv, 12);
1049 tmp = (tmp * lambdaLong) >> 11;
1050 fTailAmplLong[adc] = tmp & 0xFFF;
1051 // short component
1052 tmp = AddUintClipping(fTailAmplShort[adc], aDiff - alInpv, 12);
1053 tmp = (tmp * lambdaShort) >> 11;
1054 fTailAmplShort[adc] = tmp & 0xFFF;
5f006bd7 1055
ce4786b9 1056 // the output of the filter
759042e7 1057 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTBY, fDetector, fRobPos, fMcmPos) == 0) // bypass mode, active low
b0a41e80 1058 return value;
1059 else
b0a41e80 1060 return aDiff;
b0a41e80 1061}
dfd03fc3 1062
b0a41e80 1063void AliTRDmcmSim::FilterTail()
1064{
5f006bd7 1065 // Apply tail cancellation filter to all data.
dfd03fc3 1066
b0a41e80 1067 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
ce4786b9 1068 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
b0a41e80 1069 fADCF[iAdc][iTimeBin] = FilterTailNextSample(iAdc, fADCF[iAdc][iTimeBin]);
dfd03fc3 1070 }
dfd03fc3 1071 }
dfd03fc3 1072}
1073
dfd03fc3 1074void AliTRDmcmSim::ZSMapping()
1075{
0c349049 1076 //
dfd03fc3 1077 // Zero Suppression Mapping implemented in TRAP chip
ce4786b9 1078 // only implemented for up to 30 timebins
dfd03fc3 1079 //
1080 // See detail TRAP manual "Data Indication" section:
1081 // http://www.kip.uni-heidelberg.de/ti/TRD/doc/trap/TRAP-UserManual.pdf
0c349049 1082 //
dfd03fc3 1083
5f006bd7 1084 if( !CheckInitialized() )
ce4786b9 1085 return;
1086
759042e7 1087 Int_t eBIS = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBIS, fDetector, fRobPos, fMcmPos);
1088 Int_t eBIT = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBIT, fDetector, fRobPos, fMcmPos);
1089 Int_t eBIL = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBIL, fDetector, fRobPos, fMcmPos);
1090 Int_t eBIN = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBIN, fDetector, fRobPos, fMcmPos);
ecf39416 1091
b0a41e80 1092 Int_t **adc = fADCF;
dfd03fc3 1093
5f006bd7 1094 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++)
ce4786b9 1095 fZSMap[iAdc] = -1;
b0a41e80 1096
1097 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
ce4786b9 1098 Int_t iAdc; // current ADC channel
1099 Int_t ap;
1100 Int_t ac;
1101 Int_t an;
1102 Int_t mask;
1103 Int_t supp; // suppression of the current channel (low active)
5f006bd7 1104
ce4786b9 1105 // ----- first channel -----
1106 iAdc = 0;
5f006bd7 1107
ce4786b9 1108 ap = 0; // previous
1109 ac = adc[iAdc ][it]; // current
1110 an = adc[iAdc+1][it]; // next
5f006bd7 1111
ce4786b9 1112 mask = ( ac >= ap && ac >= an ) ? 0 : 0x1; // peak center detection
1113 mask += ( ap + ac + an > eBIT ) ? 0 : 0x2; // cluster
1114 mask += ( ac > eBIS ) ? 0 : 0x4; // absolute large peak
5f006bd7 1115
ce4786b9 1116 supp = (eBIL >> mask) & 1;
5f006bd7 1117
ce4786b9 1118 fZSMap[iAdc] &= ~((1-supp) << it);
1119 if( eBIN == 0 ) { // neighbour sensitivity
1120 fZSMap[iAdc+1] &= ~((1-supp) << it);
dfd03fc3 1121 }
5f006bd7 1122
ce4786b9 1123 // ----- last channel -----
1124 iAdc = fgkNADC - 1;
5f006bd7 1125
ce4786b9 1126 ap = adc[iAdc-1][it]; // previous
1127 ac = adc[iAdc ][it]; // current
1128 an = 0; // next
5f006bd7 1129
ce4786b9 1130 mask = ( ac >= ap && ac >= an ) ? 0 : 0x1; // peak center detection
1131 mask += ( ap + ac + an > eBIT ) ? 0 : 0x2; // cluster
1132 mask += ( ac > eBIS ) ? 0 : 0x4; // absolute large peak
5f006bd7 1133
ce4786b9 1134 supp = (eBIL >> mask) & 1;
5f006bd7 1135
ce4786b9 1136 fZSMap[iAdc] &= ~((1-supp) << it);
1137 if( eBIN == 0 ) { // neighbour sensitivity
1138 fZSMap[iAdc-1] &= ~((1-supp) << it);
ecf39416 1139 }
5f006bd7 1140
ce4786b9 1141 // ----- middle channels -----
1142 for( iAdc = 1 ; iAdc < fgkNADC-1; iAdc++ ) {
1143 ap = adc[iAdc-1][it]; // previous
1144 ac = adc[iAdc ][it]; // current
1145 an = adc[iAdc+1][it]; // next
5f006bd7 1146
ce4786b9 1147 mask = ( ac >= ap && ac >= an ) ? 0 : 0x1; // peak center detection
1148 mask += ( ap + ac + an > eBIT ) ? 0 : 0x2; // cluster
1149 mask += ( ac > eBIS ) ? 0 : 0x4; // absolute large peak
5f006bd7 1150
ce4786b9 1151 supp = (eBIL >> mask) & 1;
5f006bd7 1152
ce4786b9 1153 fZSMap[iAdc] &= ~((1-supp) << it);
1154 if( eBIN == 0 ) { // neighbour sensitivity
1155 fZSMap[iAdc-1] &= ~((1-supp) << it);
1156 fZSMap[iAdc+1] &= ~((1-supp) << it);
ecf39416 1157 }
dfd03fc3 1158 }
ce4786b9 1159
dfd03fc3 1160 }
1161}
1162
25b41f6f 1163void AliTRDmcmSim::AddHitToFitreg(Int_t adc, UShort_t timebin, UShort_t qtot, Short_t ypos, Int_t label[])
dfd03fc3 1164{
5f006bd7 1165 // Add the given hit to the fit register which is lateron used for
1166 // the tracklet calculation.
1167 // In addition to the fit sums in the fit register MC information
b0a41e80 1168 // is stored.
1169
759042e7 1170 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0, fDetector, fRobPos, fMcmPos)) &&
1171 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE0, fDetector, fRobPos, fMcmPos)))
ab9f7002 1172 fFitReg[adc].fQ0 += qtot;
5f006bd7 1173
759042e7 1174 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS1, fDetector, fRobPos, fMcmPos)) &&
1175 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1, fDetector, fRobPos, fMcmPos)))
ab9f7002 1176 fFitReg[adc].fQ1 += qtot;
5f006bd7 1177
759042e7 1178 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS, fDetector, fRobPos, fMcmPos) ) &&
1179 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFE, fDetector, fRobPos, fMcmPos)))
b0a41e80 1180 {
ab9f7002 1181 fFitReg[adc].fSumX += timebin;
1182 fFitReg[adc].fSumX2 += timebin*timebin;
1183 fFitReg[adc].fNhits++;
1184 fFitReg[adc].fSumY += ypos;
1185 fFitReg[adc].fSumY2 += ypos*ypos;
1186 fFitReg[adc].fSumXY += timebin*ypos;
a10c6981 1187 AliDebug(10, Form("fitreg[%2i] in timebin %2i: X=%i, X2=%i, N=%i, Y=%i, Y2=%i, XY=%i, Q0=%i, Q1=%i",
1188 adc, timebin, fFitReg[adc].fSumX, fFitReg[adc].fSumX2, fFitReg[adc].fNhits,
1189 fFitReg[adc].fSumY, fFitReg[adc].fSumY2, fFitReg[adc].fSumXY, fFitReg[adc].fQ0, fFitReg[adc].fQ1));
b0a41e80 1190 }
1191
1192 // register hits (MC info)
ab9f7002 1193 fHits[fNHits].fChannel = adc;
1194 fHits[fNHits].fQtot = qtot;
1195 fHits[fNHits].fYpos = ypos;
1196 fHits[fNHits].fTimebin = timebin;
25b41f6f 1197 fHits[fNHits].fLabel[0] = label[0];
1198 fHits[fNHits].fLabel[1] = label[1];
1199 fHits[fNHits].fLabel[2] = label[2];
b0a41e80 1200 fNHits++;
1201}
dfd03fc3 1202
5f006bd7 1203void AliTRDmcmSim::CalcFitreg()
b0a41e80 1204{
1205 // Preprocessing.
1206 // Detect the hits and fill the fit registers.
5f006bd7 1207 // Requires 12-bit data from fADCF which means Filter()
b0a41e80 1208 // has to be called before even if all filters are bypassed.
1209
b0a41e80 1210 //??? to be clarified:
64e3d742 1211 UInt_t adcMask = 0xffffffff;
5f006bd7 1212
ab9f7002 1213 UShort_t timebin, adcch, adcLeft, adcCentral, adcRight, hitQual, timebin1, timebin2, qtotTemp;
b0a41e80 1214 Short_t ypos, fromLeft, fromRight, found;
5ac2e3b1 1215 UShort_t qTotal[19+1]; // the last is dummy
ab9f7002 1216 UShort_t marked[6], qMarked[6], worse1, worse2;
5f006bd7 1217
759042e7 1218 timebin1 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS, fDetector, fRobPos, fMcmPos);
1219 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0, fDetector, fRobPos, fMcmPos)
b0a41e80 1220 < timebin1)
759042e7 1221 timebin1 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0, fDetector, fRobPos, fMcmPos);
1222 timebin2 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFE, fDetector, fRobPos, fMcmPos);
1223 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1, fDetector, fRobPos, fMcmPos)
b0a41e80 1224 > timebin2)
759042e7 1225 timebin2 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1, fDetector, fRobPos, fMcmPos);
b0a41e80 1226
1227 // reset the fit registers
5f006bd7 1228 fNHits = 0;
ce4786b9 1229 for (adcch = 0; adcch < fgkNADC-2; adcch++) // due to border channels
b0a41e80 1230 {
ab9f7002 1231 fFitReg[adcch].fNhits = 0;
1232 fFitReg[adcch].fQ0 = 0;
1233 fFitReg[adcch].fQ1 = 0;
1234 fFitReg[adcch].fSumX = 0;
1235 fFitReg[adcch].fSumY = 0;
1236 fFitReg[adcch].fSumX2 = 0;
1237 fFitReg[adcch].fSumY2 = 0;
1238 fFitReg[adcch].fSumXY = 0;
b0a41e80 1239 }
5f006bd7 1240
b0a41e80 1241 for (timebin = timebin1; timebin < timebin2; timebin++)
1242 {
ab9f7002 1243 // first find the hit candidates and store the total cluster charge in qTotal array
b0a41e80 1244 // in case of not hit store 0 there.
ce4786b9 1245 for (adcch = 0; adcch < fgkNADC-2; adcch++) {
ab9f7002 1246 if ( ( (adcMask >> adcch) & 7) == 7) //??? all 3 channels are present in case of ZS
b0a41e80 1247 {
ab9f7002 1248 adcLeft = fADCF[adcch ][timebin];
1249 adcCentral = fADCF[adcch+1][timebin];
1250 adcRight = fADCF[adcch+2][timebin];
759042e7 1251 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPVBY, fDetector, fRobPos, fMcmPos) == 1)
5f006bd7 1252 hitQual = ( (adcLeft * adcRight) <
759042e7 1253 (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPVT, fDetector, fRobPos, fMcmPos) * adcCentral) );
5f006bd7 1254 else
ab9f7002 1255 hitQual = 1;
b0a41e80 1256 // The accumulated charge is with the pedestal!!!
ab9f7002 1257 qtotTemp = adcLeft + adcCentral + adcRight;
1258 if ( (hitQual) &&
759042e7 1259 (qtotTemp >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPHT, fDetector, fRobPos, fMcmPos)) &&
ab9f7002 1260 (adcLeft <= adcCentral) &&
1261 (adcCentral > adcRight) )
1262 qTotal[adcch] = qtotTemp;
b0a41e80 1263 else
ab9f7002 1264 qTotal[adcch] = 0;
b0a41e80 1265 }
1266 else
ab9f7002 1267 qTotal[adcch] = 0; //jkl
5f006bd7 1268 if (qTotal[adcch] != 0)
ce4786b9 1269 AliDebug(10,Form("ch %2d qTotal %5d",adcch, qTotal[adcch]));
b0a41e80 1270 }
dfd03fc3 1271
b0a41e80 1272 fromLeft = -1;
1273 adcch = 0;
1274 found = 0;
1275 marked[4] = 19; // invalid channel
1276 marked[5] = 19; // invalid channel
ab9f7002 1277 qTotal[19] = 0;
b0a41e80 1278 while ((adcch < 16) && (found < 3))
1279 {
ab9f7002 1280 if (qTotal[adcch] > 0)
b0a41e80 1281 {
1282 fromLeft = adcch;
1283 marked[2*found+1]=adcch;
1284 found++;
1285 }
1286 adcch++;
1287 }
5f006bd7 1288
b0a41e80 1289 fromRight = -1;
1290 adcch = 18;
1291 found = 0;
1292 while ((adcch > 2) && (found < 3))
1293 {
ab9f7002 1294 if (qTotal[adcch] > 0)
b0a41e80 1295 {
1296 marked[2*found]=adcch;
1297 found++;
1298 fromRight = adcch;
1299 }
1300 adcch--;
1301 }
dfd03fc3 1302
4ff7ed2b 1303 AliDebug(10,Form("Fromleft=%d, Fromright=%d",fromLeft, fromRight));
b0a41e80 1304 // here mask the hit candidates in the middle, if any
1305 if ((fromLeft >= 0) && (fromRight >= 0) && (fromLeft < fromRight))
1306 for (adcch = fromLeft+1; adcch < fromRight; adcch++)
ab9f7002 1307 qTotal[adcch] = 0;
5f006bd7 1308
b0a41e80 1309 found = 0;
1310 for (adcch = 0; adcch < 19; adcch++)
ab9f7002 1311 if (qTotal[adcch] > 0) found++;
b0a41e80 1312 // NOT READY
1313
1314 if (found > 4) // sorting like in the TRAP in case of 5 or 6 candidates!
1315 {
1316 if (marked[4] == marked[5]) marked[5] = 19;
1317 for (found=0; found<6; found++)
1318 {
ab9f7002 1319 qMarked[found] = qTotal[marked[found]] >> 4;
4ff7ed2b 1320 AliDebug(10,Form("ch_%d qTotal %d qTotals %d",marked[found],qTotal[marked[found]],qMarked[found]));
b0a41e80 1321 }
5f006bd7 1322
b0a41e80 1323 Sort6To2Worst(marked[0], marked[3], marked[4], marked[1], marked[2], marked[5],
ab9f7002 1324 qMarked[0],
1325 qMarked[3],
1326 qMarked[4],
1327 qMarked[1],
1328 qMarked[2],
1329 qMarked[5],
b0a41e80 1330 &worse1, &worse2);
1331 // Now mask the two channels with the smallest charge
1332 if (worse1 < 19)
1333 {
ab9f7002 1334 qTotal[worse1] = 0;
4ff7ed2b 1335 AliDebug(10,Form("Kill ch %d\n",worse1));
b0a41e80 1336 }
1337 if (worse2 < 19)
1338 {
ab9f7002 1339 qTotal[worse2] = 0;
4ff7ed2b 1340 AliDebug(10,Form("Kill ch %d\n",worse2));
b0a41e80 1341 }
1342 }
5f006bd7 1343
b0a41e80 1344 for (adcch = 0; adcch < 19; adcch++) {
ab9f7002 1345 if (qTotal[adcch] > 0) // the channel is marked for processing
b0a41e80 1346 {
ab9f7002 1347 adcLeft = fADCF[adcch ][timebin];
1348 adcCentral = fADCF[adcch+1][timebin];
1349 adcRight = fADCF[adcch+2][timebin];
b0a41e80 1350 // hit detected, in TRAP we have 4 units and a hit-selection, here we proceed all channels!
1351 // subtract the pedestal TPFP, clipping instead of wrapping
5f006bd7 1352
759042e7 1353 Int_t regTPFP = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP, fDetector, fRobPos, fMcmPos);
4ff7ed2b 1354 AliDebug(10, Form("Hit found, time=%d, adcch=%d/%d/%d, adc values=%d/%d/%d, regTPFP=%d, TPHT=%d\n",
5f006bd7 1355 timebin, adcch, adcch+1, adcch+2, adcLeft, adcCentral, adcRight, regTPFP,
759042e7 1356 fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPHT, fDetector, fRobPos, fMcmPos)));
b0a41e80 1357
ab9f7002 1358 if (adcLeft < regTPFP) adcLeft = 0; else adcLeft -= regTPFP;
1359 if (adcCentral < regTPFP) adcCentral = 0; else adcCentral -= regTPFP;
1360 if (adcRight < regTPFP) adcRight = 0; else adcRight -= regTPFP;
f793c83d 1361
b0a41e80 1362 // Calculate the center of gravity
f793c83d 1363 // checking for adcCentral != 0 (in case of "bad" configuration)
1364 if (adcCentral == 0)
1365 continue;
a10c6981 1366 ypos = 128*(adcRight - adcLeft) / adcCentral;
b0a41e80 1367 if (ypos < 0) ypos = -ypos;
ce4786b9 1368 // make the correction using the position LUT
25b41f6f 1369 ypos = ypos + fTrapConfig->GetTrapReg((AliTRDtrapConfig::TrapReg_t) (AliTRDtrapConfig::kTPL00 + (ypos & 0x7F)),
1370 fDetector, fRobPos, fMcmPos);
ab9f7002 1371 if (adcLeft > adcRight) ypos = -ypos;
40bd6ee4 1372
25b41f6f 1373 // label calculation (up to 3)
1374 Int_t mcLabel[] = {-1, -1, -1};
40bd6ee4 1375 if (fDigitsManager) {
25b41f6f 1376 const Int_t maxLabels = 9;
1377 Int_t label[maxLabels] = { 0 }; // up to 9 different labels possible
1378 Int_t count[maxLabels] = { 0 };
40bd6ee4 1379 Int_t nLabels = 0;
5f006bd7 1380 Int_t padcol[3];
40bd6ee4 1381 padcol[0] = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adcch);
1382 padcol[1] = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adcch+1);
1383 padcol[2] = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adcch+2);
1384 Int_t padrow = fFeeParam->GetPadRowFromMCM(fRobPos, fMcmPos);
1385 for (Int_t iDict = 0; iDict < 3; iDict++) {
ce4786b9 1386 if (!fDict[iDict])
40bd6ee4 1387 continue;
40bd6ee4 1388 for (Int_t iPad = 0; iPad < 3; iPad++) {
5f006bd7 1389 if (padcol[iPad] < 0)
40bd6ee4 1390 continue;
25b41f6f 1391 Int_t currLabel = fDict[iDict]->GetData(padrow, padcol[iPad], timebin);
4ff7ed2b 1392 AliDebug(10, Form("Read label: %4i for det: %3i, row: %i, col: %i, tb: %i\n", currLabel, fDetector, padrow, padcol[iPad], timebin));
40bd6ee4 1393 for (Int_t iLabel = 0; iLabel < nLabels; iLabel++) {
1394 if (currLabel == label[iLabel]) {
1395 count[iLabel]++;
ce51199c 1396 currLabel = -1;
40bd6ee4 1397 break;
1398 }
5f006bd7 1399 }
ce51199c 1400 if (currLabel >= 0) {
25b41f6f 1401 label[nLabels] = currLabel;
1402 count[nLabels] = 1;
1403 nLabels++;
40bd6ee4 1404 }
1405 }
1406 }
637666cd 1407 Int_t index[2*maxLabels];
25b41f6f 1408 TMath::Sort(maxLabels, count, index);
1409 for (Int_t i = 0; i < 3; i++) {
1410 if (count[index[i]] <= 0)
1411 break;
1412 mcLabel[i] = label[index[i]];
1413 }
40bd6ee4 1414 }
1415
1416 // add the hit to the fitregister
a10c6981 1417 AddHitToFitreg(adcch, timebin, qTotal[adcch] >> fgkAddDigits, ypos, mcLabel);
b0a41e80 1418 }
dfd03fc3 1419 }
1420 }
ce4786b9 1421
1422 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1423 if (fFitReg[iAdc].fNhits != 0) {
1424 AliDebug(2, Form("fitreg[%i]: nHits = %i, sumX = %i, sumY = %i, sumX2 = %i, sumY2 = %i, sumXY = %i", iAdc,
1425 fFitReg[iAdc].fNhits,
1426 fFitReg[iAdc].fSumX,
1427 fFitReg[iAdc].fSumY,
1428 fFitReg[iAdc].fSumX2,
1429 fFitReg[iAdc].fSumY2,
1430 fFitReg[iAdc].fSumXY
1431 ));
1432 }
1433 }
dfd03fc3 1434}
1435
5f006bd7 1436void AliTRDmcmSim::TrackletSelection()
dfd03fc3 1437{
5f006bd7 1438 // Select up to 4 tracklet candidates from the fit registers
b0a41e80 1439 // and assign them to the CPUs.
1440
ab9f7002 1441 UShort_t adcIdx, i, j, ntracks, tmp;
1442 UShort_t trackletCand[18][2]; // store the adcch[0] and number of hits[1] for all tracklet candidates
b0a41e80 1443
1444 ntracks = 0;
ab9f7002 1445 for (adcIdx = 0; adcIdx < 18; adcIdx++) // ADCs
5f006bd7 1446 if ( (fFitReg[adcIdx].fNhits
759042e7 1447 >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPCL, fDetector, fRobPos, fMcmPos)) &&
ab9f7002 1448 (fFitReg[adcIdx].fNhits+fFitReg[adcIdx+1].fNhits
759042e7 1449 >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPCT, fDetector, fRobPos, fMcmPos)))
b0a41e80 1450 {
ab9f7002 1451 trackletCand[ntracks][0] = adcIdx;
1452 trackletCand[ntracks][1] = fFitReg[adcIdx].fNhits+fFitReg[adcIdx+1].fNhits;
4ff7ed2b 1453 AliDebug(10,Form("%d %2d %4d\n", ntracks, trackletCand[ntracks][0], trackletCand[ntracks][1]));
b0a41e80 1454 ntracks++;
1455 };
1456
5f006bd7 1457 for (i=0; i<ntracks;i++)
4ff7ed2b 1458 AliDebug(10,Form("%d %d %d\n",i,trackletCand[i][0], trackletCand[i][1]));
b0a41e80 1459
1460 if (ntracks > 4)
1461 {
1462 // primitive sorting according to the number of hits
1463 for (j = 0; j < (ntracks-1); j++)
1464 {
1465 for (i = j+1; i < ntracks; i++)
1466 {
ab9f7002 1467 if ( (trackletCand[j][1] < trackletCand[i][1]) ||
1468 ( (trackletCand[j][1] == trackletCand[i][1]) && (trackletCand[j][0] < trackletCand[i][0]) ) )
b0a41e80 1469 {
1470 // swap j & i
ab9f7002 1471 tmp = trackletCand[j][1];
1472 trackletCand[j][1] = trackletCand[i][1];
1473 trackletCand[i][1] = tmp;
1474 tmp = trackletCand[j][0];
1475 trackletCand[j][0] = trackletCand[i][0];
1476 trackletCand[i][0] = tmp;
b0a41e80 1477 }
1478 }
1479 }
1480 ntracks = 4; // cut the rest, 4 is the max
dfd03fc3 1481 }
b0a41e80 1482 // else is not necessary to sort
5f006bd7 1483
b0a41e80 1484 // now sort, so that the first tracklet going to CPU0 corresponds to the highest adc channel - as in the TRAP
1485 for (j = 0; j < (ntracks-1); j++)
1486 {
1487 for (i = j+1; i < ntracks; i++)
1488 {
ab9f7002 1489 if (trackletCand[j][0] < trackletCand[i][0])
b0a41e80 1490 {
1491 // swap j & i
ab9f7002 1492 tmp = trackletCand[j][1];
1493 trackletCand[j][1] = trackletCand[i][1];
1494 trackletCand[i][1] = tmp;
1495 tmp = trackletCand[j][0];
1496 trackletCand[j][0] = trackletCand[i][0];
1497 trackletCand[i][0] = tmp;
b0a41e80 1498 }
dfd03fc3 1499 }
b0a41e80 1500 }
1501 for (i = 0; i < ntracks; i++) // CPUs with tracklets.
ab9f7002 1502 fFitPtr[i] = trackletCand[i][0]; // pointer to the left channel with tracklet for CPU[i]
b0a41e80 1503 for (i = ntracks; i < 4; i++) // CPUs without tracklets
1504 fFitPtr[i] = 31; // pointer to the left channel with tracklet for CPU[i] = 31 (invalid)
4ff7ed2b 1505 AliDebug(10,Form("found %i tracklet candidates\n", ntracks));
1506 for (i = 0; i < 4; i++)
1507 AliDebug(10,Form("fitPtr[%i]: %i\n", i, fFitPtr[i]));
b0a41e80 1508}
dfd03fc3 1509
b0a41e80 1510void AliTRDmcmSim::FitTracklet()
1511{
5f006bd7 1512 // Perform the actual tracklet fit based on the fit sums
1513 // which have been filled in the fit registers.
b0a41e80 1514
1515 // parameters in fitred.asm (fit program)
b0a41e80 1516 Int_t rndAdd = 0;
78c94f0b 1517 Int_t decPlaces = 5; // must be larger than 1 or change the following code
1518 // if (decPlaces > 1)
b0a41e80 1519 rndAdd = (1 << (decPlaces-1)) + 1;
78c94f0b 1520 // else if (decPlaces == 1)
1521 // rndAdd = 1;
1522
4ff7ed2b 1523 Int_t ndriftDp = 5; // decimal places for drift time
1524 Long64_t shift = ((Long64_t) 1 << 32);
1525
4ff7ed2b 1526 // calculated in fitred.asm
1527 Int_t padrow = ((fRobPos >> 1) << 2) | (fMcmPos >> 2);
5f006bd7 1528 Int_t yoffs = (((((fRobPos & 0x1) << 2) + (fMcmPos & 0x3)) * 18) << 8) -
4ff7ed2b 1529 ((18*4*2 - 18*2 - 1) << 7);
1530 yoffs = yoffs << decPlaces; // holds position of ADC channel 1
1531 Int_t layer = fDetector % 6;
1532 UInt_t scaleY = (UInt_t) ((0.635 + 0.03 * layer)/(256.0 * 160.0e-4) * shift);
1533 UInt_t scaleD = (UInt_t) ((0.635 + 0.03 * layer)/(256.0 * 140.0e-4) * shift);
4ff7ed2b 1534
8ea391e3 1535 Int_t deflCorr = (Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCorr, fDetector, fRobPos, fMcmPos);
5f006bd7 1536 Int_t ndrift = (Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrNdrift, fDetector, fRobPos, fMcmPos);
b0a41e80 1537
1538 // local variables for calculation
1539 Long64_t mult, temp, denom; //???
8ea391e3 1540 UInt_t q0, q1, pid; // charges in the two windows and total charge
b0a41e80 1541 UShort_t nHits; // number of hits
1542 Int_t slope, offset; // slope and offset of the tracklet
1543 Int_t sumX, sumY, sumXY, sumX2; // fit sums from fit registers
ce51199c 1544 Int_t sumY2; // not used in the current TRAP program, now used for error calculation (simulation only)
1545 Float_t fitError, fitSlope, fitOffset;
b0a41e80 1546 FitReg_t *fit0, *fit1; // pointers to relevant fit registers
5f006bd7 1547
b0a41e80 1548// const uint32_t OneDivN[32] = { // 2**31/N : exactly like in the TRAP, the simple division here gives the same result!
1549// 0x00000000, 0x80000000, 0x40000000, 0x2AAAAAA0, 0x20000000, 0x19999990, 0x15555550, 0x12492490,
1550// 0x10000000, 0x0E38E380, 0x0CCCCCC0, 0x0BA2E8B0, 0x0AAAAAA0, 0x09D89D80, 0x09249240, 0x08888880,
1551// 0x08000000, 0x07878780, 0x071C71C0, 0x06BCA1A0, 0x06666660, 0x06186180, 0x05D17450, 0x0590B210,
1552// 0x05555550, 0x051EB850, 0x04EC4EC0, 0x04BDA120, 0x04924920, 0x0469EE50, 0x04444440, 0x04210840};
1553
1554 for (Int_t cpu = 0; cpu < 4; cpu++) {
1555 if (fFitPtr[cpu] == 31)
1556 {
5f006bd7 1557 fMCMT[cpu] = 0x10001000; //??? AliTRDfeeParam::GetTrackletEndmarker();
dfd03fc3 1558 }
b0a41e80 1559 else
1560 {
1561 fit0 = &fFitReg[fFitPtr[cpu] ];
1562 fit1 = &fFitReg[fFitPtr[cpu]+1]; // next channel
1563
1564 mult = 1;
1565 mult = mult << (32 + decPlaces);
1566 mult = -mult;
1567
1568 // Merging
ab9f7002 1569 nHits = fit0->fNhits + fit1->fNhits; // number of hits
1570 sumX = fit0->fSumX + fit1->fSumX;
1571 sumX2 = fit0->fSumX2 + fit1->fSumX2;
f5821bdb 1572 denom = ((Long64_t) nHits)*((Long64_t) sumX2) - ((Long64_t) sumX)*((Long64_t) sumX);
b0a41e80 1573
1574 mult = mult / denom; // exactly like in the TRAP program
ab9f7002 1575 q0 = fit0->fQ0 + fit1->fQ0;
1576 q1 = fit0->fQ1 + fit1->fQ1;
1577 sumY = fit0->fSumY + fit1->fSumY + 256*fit1->fNhits;
1578 sumXY = fit0->fSumXY + fit1->fSumXY + 256*fit1->fSumX;
ce51199c 1579 sumY2 = fit0->fSumY2 + fit1->fSumY2 + 512*fit1->fSumY + 256*256*fit1->fNhits;
b0a41e80 1580
1581 slope = nHits*sumXY - sumX * sumY;
1582 offset = sumX2*sumY - sumX * sumXY;
1583 temp = mult * slope;
1584 slope = temp >> 32; // take the upper 32 bits
4ff7ed2b 1585 slope = -slope;
b0a41e80 1586 temp = mult * offset;
1587 offset = temp >> 32; // take the upper 32 bits
1588
4ff7ed2b 1589 offset = offset + yoffs;
5f006bd7 1590 AliDebug(10, Form("slope = %i, slope * ndrift = %i, deflCorr: %i",
ce4786b9 1591 slope, slope * ndrift, deflCorr));
1592 slope = ((slope * ndrift) >> ndriftDp) + deflCorr;
b0a41e80 1593 offset = offset - (fFitPtr[cpu] << (8 + decPlaces));
5f006bd7 1594
4ff7ed2b 1595 temp = slope;
1596 temp = temp * scaleD;
1597 slope = (temp >> 32);
4ff7ed2b 1598 temp = offset;
1599 temp = temp * scaleY;
1600 offset = (temp >> 32);
5f006bd7 1601
4ff7ed2b 1602 // rounding, like in the TRAP
1603 slope = (slope + rndAdd) >> decPlaces;
4ff7ed2b 1604 offset = (offset + rndAdd) >> decPlaces;
1605
5f006bd7 1606 AliDebug(5, Form("Det: %3i, ROB: %i, MCM: %2i: deflection: %i, min: %i, max: %i",
1607 fDetector, fRobPos, fMcmPos, slope,
1608 (Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCutStart + 2*fFitPtr[cpu], fDetector, fRobPos, fMcmPos),
8ea391e3 1609 (Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCutStart + 1 + 2*fFitPtr[cpu], fDetector, fRobPos, fMcmPos)));
ce4786b9 1610
5f006bd7 1611 AliDebug(5, Form("Fit sums: x = %i, X = %i, y = %i, Y = %i, Z = %i",
ce51199c 1612 sumX, sumX2, sumY, sumY2, sumXY));
1613
1614 fitSlope = (Float_t) (nHits * sumXY - sumX * sumY) / (nHits * sumX2 - sumX*sumX);
1615
1616 fitOffset = (Float_t) (sumX2 * sumY - sumX * sumXY) / (nHits * sumX2 - sumX*sumX);
1617
1618 Float_t sx = (Float_t) sumX;
1619 Float_t sx2 = (Float_t) sumX2;
1620 Float_t sy = (Float_t) sumY;
1621 Float_t sy2 = (Float_t) sumY2;
1622 Float_t sxy = (Float_t) sumXY;
1623 fitError = sy2 - (sx2 * sy*sy - 2 * sx * sxy * sy + nHits * sxy*sxy) / (nHits * sx2 - sx*sx);
1624 //fitError = (Float_t) sumY2 - (Float_t) (sumY*sumY) / nHits - fitSlope * ((Float_t) (sumXY - sumX*sumY) / nHits);
1625
40bd6ee4 1626 Bool_t rejected = kFALSE;
ce4786b9 1627 // deflection range table from DMEM
5f006bd7 1628 if ((slope < ((Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCutStart + 2*fFitPtr[cpu], fDetector, fRobPos, fMcmPos))) ||
8ea391e3 1629 (slope > ((Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCutStart + 1 + 2*fFitPtr[cpu], fDetector, fRobPos, fMcmPos))))
40bd6ee4 1630 rejected = kTRUE;
4ff7ed2b 1631
1632 if (rejected && GetApplyCut())
b0a41e80 1633 {
1634 fMCMT[cpu] = 0x10001000; //??? AliTRDfeeParam::GetTrackletEndmarker();
1635 }
1636 else
1637 {
4ff7ed2b 1638 if (slope > 63 || slope < -64) { // wrapping in TRAP!
27a030ab 1639 AliDebug(1,Form("Overflow in slope: %i, tracklet discarded!", slope));
40bd6ee4 1640 fMCMT[cpu] = 0x10001000;
1641 continue;
1642 }
b0a41e80 1643
4ff7ed2b 1644 slope = slope & 0x7F; // 7 bit
5f006bd7 1645
1646 if (offset > 0xfff || offset < -0xfff)
b0a41e80 1647 AliWarning("Overflow in offset");
1648 offset = offset & 0x1FFF; // 13 bit
1649
a10c6981 1650 pid = GetPID(q0, q1);
4ff7ed2b 1651
8ea391e3 1652 if (pid > 0xff)
1653 AliWarning("Overflow in PID");
1654 pid = pid & 0xFF; // 8 bit, exactly like in the TRAP program
5f006bd7 1655
b0a41e80 1656 // assemble and store the tracklet word
8ea391e3 1657 fMCMT[cpu] = (pid << 24) | (padrow << 20) | (slope << 13) | offset;
40bd6ee4 1658
1659 // calculate MC label
25b41f6f 1660 Int_t mcLabel[] = { -1, -1, -1};
4ff7ed2b 1661 Int_t nHits0 = 0;
1662 Int_t nHits1 = 0;
40bd6ee4 1663 if (fDigitsManager) {
25b41f6f 1664 const Int_t maxLabels = 30;
1665 Int_t label[maxLabels] = {0}; // up to 30 different labels possible
1666 Int_t count[maxLabels] = {0};
40bd6ee4 1667 Int_t nLabels = 0;
1668 for (Int_t iHit = 0; iHit < fNHits; iHit++) {
1669 if ((fHits[iHit].fChannel - fFitPtr[cpu] < 0) ||
1670 (fHits[iHit].fChannel - fFitPtr[cpu] > 1))
1671 continue;
4ff7ed2b 1672
1673 // counting contributing hits
759042e7 1674 if (fHits[iHit].fTimebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0, fDetector, fRobPos, fMcmPos) &&
1675 fHits[iHit].fTimebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE0, fDetector, fRobPos, fMcmPos))
4ff7ed2b 1676 nHits0++;
759042e7 1677 if (fHits[iHit].fTimebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS1, fDetector, fRobPos, fMcmPos) &&
1678 fHits[iHit].fTimebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1, fDetector, fRobPos, fMcmPos))
4ff7ed2b 1679 nHits1++;
1680
25b41f6f 1681 for (Int_t i = 0; i < 3; i++) {
1682 Int_t currLabel = fHits[iHit].fLabel[i];
1683 for (Int_t iLabel = 0; iLabel < nLabels; iLabel++) {
1684 if (currLabel == label[iLabel]) {
1685 count[iLabel]++;
1686 currLabel = -1;
1687 break;
1688 }
1689 }
1690 if (currLabel >= 0 && nLabels < maxLabels) {
1691 label[nLabels] = currLabel;
1692 count[nLabels]++;
1693 nLabels++;
1694 }
1695 }
1696 }
637666cd 1697 Int_t index[2*maxLabels];
25b41f6f 1698 TMath::Sort(maxLabels, count, index);
1699 for (Int_t i = 0; i < 3; i++) {
1700 if (count[index[i]] <= 0)
1701 break;
1702 mcLabel[i] = label[index[i]];
1703 }
40bd6ee4 1704 }
f793c83d 1705 new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletMCM((UInt_t) fMCMT[cpu], fDetector*2 + fRobPos%2, fRobPos, fMcmPos);
40bd6ee4 1706 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetLabel(mcLabel);
4ff7ed2b 1707
5f006bd7 1708
4ff7ed2b 1709 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetNHits(fit0->fNhits + fit1->fNhits);
1710 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetNHits0(nHits0);
1711 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetNHits1(nHits1);
a10c6981 1712 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetQ0(q0);
1713 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetQ1(q1);
ce51199c 1714 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetSlope(fitSlope);
1715 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetOffset(fitOffset);
1716 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetError(TMath::Sqrt(TMath::Abs(fitError)/nHits));
1717
1718// // cluster information
1719// Float_t *res = new Float_t[nHits];
1720// Float_t *qtot = new Float_t[nHits];
1721// Int_t nCls = 0;
1722// for (Int_t iHit = 0; iHit < fNHits; iHit++) {
1723// // check if hit contributes
1724// if (fHits[iHit].fChannel == fFitPtr[cpu]) {
1725// res[nCls] = fHits[iHit].fYpos - (fitSlope * fHits[iHit].fTimebin + fitOffset);
1726// qtot[nCls] = fHits[iHit].fQtot;
1727// nCls++;
1728// }
1729// else if (fHits[iHit].fChannel == fFitPtr[cpu] + 1) {
1730// res[nCls] = fHits[iHit].fYpos + 256 - (fitSlope * fHits[iHit].fTimebin + fitOffset);
1731// qtot[nCls] = fHits[iHit].fQtot;
1732// nCls++;
1733// }
1734// }
1735// ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetClusters(res, qtot, nCls);
1736// delete [] res;
1737// delete [] qtot;
1738
1739 if (fitError < 0)
1740 AliError(Form("Strange fit error: %f from Sx: %i, Sy: %i, Sxy: %i, Sx2: %i, Sy2: %i, nHits: %i",
1741 fitError, sumX, sumY, sumXY, sumX2, sumY2, nHits));
5f006bd7 1742 AliDebug(3, Form("fit slope: %f, offset: %f, error: %f",
ce51199c 1743 fitSlope, fitOffset, TMath::Sqrt(TMath::Abs(fitError)/nHits)));
b0a41e80 1744 }
dfd03fc3 1745 }
dfd03fc3 1746 }
1747}
1748
b0a41e80 1749void AliTRDmcmSim::Tracklet()
dfd03fc3 1750{
ab9f7002 1751 // Run the tracklet calculation by calling sequentially:
1752 // CalcFitreg(); TrackletSelection(); FitTracklet()
5f006bd7 1753 // and store the tracklets
ab9f7002 1754
b0a41e80 1755 if (!fInitialized) {
ab9f7002 1756 AliError("Called uninitialized! Nothing done!");
b0a41e80 1757 return;
dfd03fc3 1758 }
1759
b0a41e80 1760 fTrackletArray->Delete();
dfd03fc3 1761
b0a41e80 1762 CalcFitreg();
40bd6ee4 1763 if (fNHits == 0)
1764 return;
b0a41e80 1765 TrackletSelection();
1766 FitTracklet();
c8b1590d 1767}
1768
5f006bd7 1769Bool_t AliTRDmcmSim::StoreTracklets()
c8b1590d 1770{
36dc3337 1771 // store the found tracklets via the loader
1772
5f006bd7 1773 if (fTrackletArray->GetEntriesFast() == 0)
c8b1590d 1774 return kTRUE;
dfd03fc3 1775
b0a41e80 1776 AliRunLoader *rl = AliRunLoader::Instance();
1777 AliDataLoader *dl = 0x0;
1778 if (rl)
1779 dl = rl->GetLoader("TRDLoader")->GetDataLoader("tracklets");
1780 if (!dl) {
1781 AliError("Could not get the tracklets data loader!");
c8b1590d 1782 return kFALSE;
dfd03fc3 1783 }
b0a41e80 1784
c8b1590d 1785 TTree *trackletTree = dl->Tree();
1786 if (!trackletTree) {
1787 dl->MakeTree();
1788 trackletTree = dl->Tree();
1789 }
5f006bd7 1790
c8b1590d 1791 AliTRDtrackletMCM *trkl = 0x0;
6b094867 1792 TBranch *trkbranch = trackletTree->GetBranch(fTrklBranchName.Data());
c8b1590d 1793 if (!trkbranch)
6b094867 1794 trkbranch = trackletTree->Branch(fTrklBranchName.Data(), "AliTRDtrackletMCM", &trkl, 32000);
5f006bd7 1795
c8b1590d 1796 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
1797 trkl = ((AliTRDtrackletMCM*) (*fTrackletArray)[iTracklet]);
1798 trkbranch->SetAddress(&trkl);
c8b1590d 1799 trkbranch->Fill();
b0a41e80 1800 }
c8b1590d 1801
1802 return kTRUE;
dfd03fc3 1803}
1804
b0a41e80 1805void AliTRDmcmSim::WriteData(AliTRDarrayADC *digits)
dfd03fc3 1806{
b0a41e80 1807 // write back the processed data configured by EBSF
1808 // EBSF = 1: unfiltered data; EBSF = 0: filtered data
1809 // zero-suppressed valued are written as -1 to digits
dfd03fc3 1810
5f006bd7 1811 if( !CheckInitialized() )
b0a41e80 1812 return;
dfd03fc3 1813
ce4786b9 1814 Int_t offset = (fMcmPos % 4 + 1) * 21 + (fRobPos % 2) * 84 - 1;
dfd03fc3 1815
759042e7 1816 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBSF, fDetector, fRobPos, fMcmPos) != 0) // store unfiltered data
b0a41e80 1817 {
ce4786b9 1818 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1819 if (~fZSMap[iAdc] == 0) {
b0a41e80 1820 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
ce4786b9 1821 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, -1);
b0a41e80 1822 }
1823 }
ce51199c 1824 else if (iAdc < 2 || iAdc == 20) {
1825 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1826 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, (fADCR[iAdc][iTimeBin] >> fgkAddDigits) - fgAddBaseline);
1827 }
1828 }
b0a41e80 1829 }
1830 }
1831 else {
ce4786b9 1832 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1833 if (~fZSMap[iAdc] != 0) {
b0a41e80 1834 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
ce4786b9 1835 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, (fADCF[iAdc][iTimeBin] >> fgkAddDigits) - fgAddBaseline);
b0a41e80 1836 }
1837 }
1838 else {
1839 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
ce4786b9 1840 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, -1);
b0a41e80 1841 }
1842 }
1843 }
dfd03fc3 1844 }
b0a41e80 1845}
dfd03fc3 1846
8ea391e3 1847
1848// ******************************
1849// PID section
1850//
1851// Memory area for the LUT: 0xC100 to 0xC3FF
1852//
1853// The addresses for the parameters (the order is optimized for maximum calculation speed in the MCMs):
1854// 0xC028: cor1
1855// 0xC029: nBins(sF)
1856// 0xC02A: cor0
1857// 0xC02B: TableLength
1858// Defined in AliTRDtrapConfig.h
1859//
1860// The algorithm implemented in the TRAP program of the MCMs (Venelin Angelov)
1861// 1) set the read pointer to the beginning of the Parameters in DMEM
1862// 2) shift right the FitReg with the Q0 + (Q1 << 16) to get Q1
1863// 3) read cor1 with rpointer++
1864// 4) start cor1*Q1
1865// 5) read nBins with rpointer++
1866// 6) start nBins*cor1*Q1
1867// 7) read cor0 with rpointer++
1868// 8) swap hi-low parts in FitReg, now is Q1 + (Q0 << 16)
1869// 9) shift right to get Q0
1870// 10) start cor0*Q0
1871// 11) read TableLength
1872// 12) compare cor0*Q0 with nBins
1873// 13) if >=, clip cor0*Q0 to nBins-1
1874// 14) add cor0*Q0 to nBins*cor1*Q1
1875// 15) compare the result with TableLength
1876// 16) if >=, clip to TableLength-1
1877// 17) read from the LUT 8 bits
1878
1879
1880Int_t AliTRDmcmSim::GetPID(Int_t q0, Int_t q1)
1881{
6b094867 1882 // return PID calculated from charges accumulated in two time windows
1883
8ea391e3 1884 ULong64_t addrQ0;
1885 ULong64_t addr;
1886
1887 UInt_t nBinsQ0 = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTnbins); // number of bins in q0 / 4 !!
1888 UInt_t pidTotalSize = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTLength);
27a030ab 1889 if(nBinsQ0==0 || pidTotalSize==0) // make sure we don't run into trouble if the value for Q0 is not configured
1890 return 0; // Q1 not configured is ok for 1D LUT
8ea391e3 1891
1892 ULong_t corrQ0 = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTcor0, fDetector, fRobPos, fMcmPos);
1893 ULong_t corrQ1 = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTcor1, fDetector, fRobPos, fMcmPos);
27a030ab 1894 if(corrQ0==0) // make sure we don't run into trouble if one of the values is not configured
5ac2e3b1 1895 return 0;
8ea391e3 1896
1897 addrQ0 = corrQ0;
1898 addrQ0 = (((addrQ0*q0)>>16)>>16); // because addrQ0 = (q0 * corrQ0) >> 32; does not work for unknown reasons
8ea391e3 1899
1900 if(addrQ0 >= nBinsQ0) { // check for overflow
e3bd81f7 1901 AliDebug(5,Form("Overflow in q0: %llu/4 is bigger then %u", addrQ0, nBinsQ0));
8ea391e3 1902 addrQ0 = nBinsQ0 -1;
5f006bd7 1903 }
8ea391e3 1904
1905 addr = corrQ1;
1906 addr = (((addr*q1)>>16)>>16);
1907 addr = addrQ0 + nBinsQ0*addr; // because addr = addrQ0 + nBinsQ0* (((corrQ1*q1)>>32); does not work
8ea391e3 1908
1909 if(addr >= pidTotalSize) {
e3bd81f7 1910 AliDebug(5,Form("Overflow in q1. Address %llu/4 is bigger then %u", addr, pidTotalSize));
8ea391e3 1911 addr = pidTotalSize -1;
5f006bd7 1912 }
8ea391e3 1913
1914 // For a LUT with 11 input and 8 output bits, the first memory address is set to LUT[0] | (LUT[1] << 8) | (LUT[2] << 16) | (LUT[3] << 24)
1915 // and so on
1916 UInt_t result = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTStart+(addr/4));
1917 return (result>>((addr%4)*8)) & 0xFF;
1918}
1919
1920
1921
b0a41e80 1922// help functions, to be cleaned up
1923
ab9f7002 1924UInt_t AliTRDmcmSim::AddUintClipping(UInt_t a, UInt_t b, UInt_t nbits) const
b0a41e80 1925{
5f006bd7 1926 //
1927 // This function adds a and b (unsigned) and clips to
1928 // the specified number of bits.
1929 //
b0a41e80 1930
1931 UInt_t sum = a + b;
1932 if (nbits < 32)
1933 {
1934 UInt_t maxv = (1 << nbits) - 1;;
5f006bd7 1935 if (sum > maxv)
b0a41e80 1936 sum = maxv;
1937 }
1938 else
1939 {
5f006bd7 1940 if ((sum < a) || (sum < b))
b0a41e80 1941 sum = 0xFFFFFFFF;
1942 }
1943 return sum;
dfd03fc3 1944}
1945
982869bc 1946void AliTRDmcmSim::Sort2(UShort_t idx1i, UShort_t idx2i, \
1947 UShort_t val1i, UShort_t val2i, \
6b094867 1948 UShort_t * const idx1o, UShort_t * const idx2o, \
1949 UShort_t * const val1o, UShort_t * const val2o) const
dfd03fc3 1950{
ab9f7002 1951 // sorting for tracklet selection
dfd03fc3 1952
b0a41e80 1953 if (val1i > val2i)
1954 {
1955 *idx1o = idx1i;
1956 *idx2o = idx2i;
1957 *val1o = val1i;
1958 *val2o = val2i;
1959 }
1960 else
1961 {
1962 *idx1o = idx2i;
1963 *idx2o = idx1i;
1964 *val1o = val2i;
1965 *val2o = val1i;
1966 }
1967}
1968
982869bc 1969void AliTRDmcmSim::Sort3(UShort_t idx1i, UShort_t idx2i, UShort_t idx3i, \
1970 UShort_t val1i, UShort_t val2i, UShort_t val3i, \
6b094867 1971 UShort_t * const idx1o, UShort_t * const idx2o, UShort_t * const idx3o, \
1972 UShort_t * const val1o, UShort_t * const val2o, UShort_t * const val3o)
b0a41e80 1973{
ab9f7002 1974 // sorting for tracklet selection
1975
4ff7ed2b 1976 Int_t sel;
dfd03fc3 1977
dfd03fc3 1978
b0a41e80 1979 if (val1i > val2i) sel=4; else sel=0;
1980 if (val2i > val3i) sel=sel + 2;
1981 if (val3i > val1i) sel=sel + 1;
b0a41e80 1982 switch(sel)
1983 {
1984 case 6 : // 1 > 2 > 3 => 1 2 3
1985 case 0 : // 1 = 2 = 3 => 1 2 3 : in this case doesn't matter, but so is in hardware!
1986 *idx1o = idx1i;
1987 *idx2o = idx2i;
1988 *idx3o = idx3i;
1989 *val1o = val1i;
1990 *val2o = val2i;
1991 *val3o = val3i;
1992 break;
1993
1994 case 4 : // 1 > 2, 2 <= 3, 3 <= 1 => 1 3 2
1995 *idx1o = idx1i;
1996 *idx2o = idx3i;
1997 *idx3o = idx2i;
1998 *val1o = val1i;
1999 *val2o = val3i;
2000 *val3o = val2i;
2001 break;
2002
2003 case 2 : // 1 <= 2, 2 > 3, 3 <= 1 => 2 1 3
2004 *idx1o = idx2i;
2005 *idx2o = idx1i;
2006 *idx3o = idx3i;
2007 *val1o = val2i;
2008 *val2o = val1i;
2009 *val3o = val3i;
2010 break;
2011
2012 case 3 : // 1 <= 2, 2 > 3, 3 > 1 => 2 3 1
2013 *idx1o = idx2i;
2014 *idx2o = idx3i;
2015 *idx3o = idx1i;
2016 *val1o = val2i;
2017 *val2o = val3i;
2018 *val3o = val1i;
2019 break;
2020
2021 case 1 : // 1 <= 2, 2 <= 3, 3 > 1 => 3 2 1
2022 *idx1o = idx3i;
2023 *idx2o = idx2i;
2024 *idx3o = idx1i;
2025 *val1o = val3i;
2026 *val2o = val2i;
2027 *val3o = val1i;
2028 break;
2029
2030 case 5 : // 1 > 2, 2 <= 3, 3 > 1 => 3 1 2
2031 *idx1o = idx3i;
2032 *idx2o = idx1i;
2033 *idx3o = idx2i;
2034 *val1o = val3i;
2035 *val2o = val1i;
2036 *val3o = val2i;
2037 break;
2038
2039 default: // the rest should NEVER happen!
40bd6ee4 2040 AliError("ERROR in Sort3!!!\n");
b0a41e80 2041 break;
2042 }
b0a41e80 2043}
dfd03fc3 2044
982869bc 2045void AliTRDmcmSim::Sort6To4(UShort_t idx1i, UShort_t idx2i, UShort_t idx3i, UShort_t idx4i, UShort_t idx5i, UShort_t idx6i, \
2046 UShort_t val1i, UShort_t val2i, UShort_t val3i, UShort_t val4i, UShort_t val5i, UShort_t val6i, \
6b094867 2047 UShort_t * const idx1o, UShort_t * const idx2o, UShort_t * const idx3o, UShort_t * const idx4o, \
2048 UShort_t * const val1o, UShort_t * const val2o, UShort_t * const val3o, UShort_t * const val4o)
b0a41e80 2049{
ab9f7002 2050 // sorting for tracklet selection
dfd03fc3 2051
982869bc 2052 UShort_t idx21s, idx22s, idx23s, dummy;
2053 UShort_t val21s, val22s, val23s;
2054 UShort_t idx23as, idx23bs;
2055 UShort_t val23as, val23bs;
dfd03fc3 2056
b0a41e80 2057 Sort3(idx1i, idx2i, idx3i, val1i, val2i, val3i,
2058 idx1o, &idx21s, &idx23as,
2059 val1o, &val21s, &val23as);
dfd03fc3 2060
b0a41e80 2061 Sort3(idx4i, idx5i, idx6i, val4i, val5i, val6i,
2062 idx2o, &idx22s, &idx23bs,
2063 val2o, &val22s, &val23bs);
2064
2065 Sort2(idx23as, idx23bs, val23as, val23bs, &idx23s, &dummy, &val23s, &dummy);
2066
2067 Sort3(idx21s, idx22s, idx23s, val21s, val22s, val23s,
2068 idx3o, idx4o, &dummy,
2069 val3o, val4o, &dummy);
dfd03fc3 2070
dfd03fc3 2071}
2072
982869bc 2073void AliTRDmcmSim::Sort6To2Worst(UShort_t idx1i, UShort_t idx2i, UShort_t idx3i, UShort_t idx4i, UShort_t idx5i, UShort_t idx6i, \
2074 UShort_t val1i, UShort_t val2i, UShort_t val3i, UShort_t val4i, UShort_t val5i, UShort_t val6i, \
6b094867 2075 UShort_t * const idx5o, UShort_t * const idx6o)
b0a41e80 2076{
ab9f7002 2077 // sorting for tracklet selection
1d93b218 2078
982869bc 2079 UShort_t idx21s, idx22s, idx23s, dummy1, dummy2, dummy3, dummy4, dummy5;
2080 UShort_t val21s, val22s, val23s;
2081 UShort_t idx23as, idx23bs;
2082 UShort_t val23as, val23bs;
1d93b218 2083
b0a41e80 2084 Sort3(idx1i, idx2i, idx3i, val1i, val2i, val3i,
2085 &dummy1, &idx21s, &idx23as,
2086 &dummy2, &val21s, &val23as);
1d93b218 2087
b0a41e80 2088 Sort3(idx4i, idx5i, idx6i, val4i, val5i, val6i,
2089 &dummy1, &idx22s, &idx23bs,
2090 &dummy2, &val22s, &val23bs);
1d93b218 2091
b0a41e80 2092 Sort2(idx23as, idx23bs, val23as, val23bs, &idx23s, idx5o, &val23s, &dummy1);
b65e5048 2093
b0a41e80 2094 Sort3(idx21s, idx22s, idx23s, val21s, val22s, val23s,
2095 &dummy1, &dummy2, idx6o,
2096 &dummy3, &dummy4, &dummy5);
0d64b05f 2097}
f793c83d 2098
2099
ce4786b9 2100// ----- I/O implementation -----
2101
59f78ad5 2102ostream& AliTRDmcmSim::Text(ostream& os)
ce4786b9 2103{
2104 // manipulator to activate output in text format (default)
2105
2106 os.iword(fgkFormatIndex) = 0;
2107 return os;
2108}
2109
59f78ad5 2110ostream& AliTRDmcmSim::Cfdat(ostream& os)
ce4786b9 2111{
5f006bd7 2112 // manipulator to activate output in CFDAT format
ce4786b9 2113 // to send to the FEE via SCSN
2114
5f006bd7 2115 os.iword(fgkFormatIndex) = 1;
ce4786b9 2116 return os;
2117}
2118
59f78ad5 2119ostream& AliTRDmcmSim::Raw(ostream& os)
ce4786b9 2120{
2121 // manipulator to activate output as raw data dump
2122
2123 os.iword(fgkFormatIndex) = 2;
2124 return os;
2125}
2126
2127ostream& operator<<(ostream& os, const AliTRDmcmSim& mcm)
2128{
2129 // output implementation
5f006bd7 2130
ce4786b9 2131 // no output for non-initialized MCM
2132 if (!mcm.CheckInitialized())
2133 return os;
2134
2135 // ----- human-readable output -----
2136 if (os.iword(AliTRDmcmSim::fgkFormatIndex) == 0) {
5f006bd7 2137
2138 os << "MCM " << mcm.fMcmPos << " on ROB " << mcm.fRobPos <<
ce4786b9 2139 " in detector " << mcm.fDetector << std::endl;
5f006bd7 2140
ce4786b9 2141 os << "----- Unfiltered ADC data (10 bit) -----" << std::endl;
2142 os << "ch ";
5f006bd7 2143 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++)
ce4786b9 2144 os << std::setw(5) << iChannel;
2145 os << std::endl;
2146 for (Int_t iTimeBin = 0; iTimeBin < mcm.fNTimeBin; iTimeBin++) {
2147 os << "tb " << std::setw(2) << iTimeBin << ":";
2148 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++) {
2149 os << std::setw(5) << (mcm.fADCR[iChannel][iTimeBin] >> mcm.fgkAddDigits);
2150 }
2151 os << std::endl;
2152 }
5f006bd7 2153
ce4786b9 2154 os << "----- Filtered ADC data (10+2 bit) -----" << std::endl;
2155 os << "ch ";
5f006bd7 2156 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++)
ce4786b9 2157 os << std::setw(4) << iChannel
2158 << ((~mcm.fZSMap[iChannel] != 0) ? "!" : " ");
2159 os << std::endl;
2160 for (Int_t iTimeBin = 0; iTimeBin < mcm.fNTimeBin; iTimeBin++) {
2161 os << "tb " << std::setw(2) << iTimeBin << ":";
2162 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++) {
2163 os << std::setw(4) << (mcm.fADCF[iChannel][iTimeBin])
2164 << (((mcm.fZSMap[iChannel] & (1 << iTimeBin)) == 0) ? "!" : " ");
2165 }
2166 os << std::endl;
2167 }
2168 }
2169
2170 // ----- CFDAT output -----
2171 else if(os.iword(AliTRDmcmSim::fgkFormatIndex) == 1) {
2172 Int_t dest = 127;
2173 Int_t addrOffset = 0x2000;
2174 Int_t addrStep = 0x80;
5f006bd7 2175
ce4786b9 2176 for (Int_t iTimeBin = 0; iTimeBin < mcm.fNTimeBin; iTimeBin++) {
2177 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++) {
5f006bd7 2178 os << std::setw(5) << 10
2179 << std::setw(5) << addrOffset + iChannel * addrStep + iTimeBin
ce4786b9 2180 << std::setw(5) << (mcm.fADCF[iChannel][iTimeBin])
2181 << std::setw(5) << dest << std::endl;
2182 }
2183 os << std::endl;
2184 }
2185 }
2186
2187 // ----- raw data ouptut -----
2188 else if (os.iword(AliTRDmcmSim::fgkFormatIndex) == 2) {
2189 Int_t bufSize = 300;
2190 UInt_t *buf = new UInt_t[bufSize];
5f006bd7 2191
ce4786b9 2192 Int_t bufLength = mcm.ProduceRawStream(&buf[0], bufSize);
5f006bd7 2193
2194 for (Int_t i = 0; i < bufLength; i++)
78c94f0b 2195 std::cout << "0x" << std::hex << buf[i] << std::dec << std::endl;
5f006bd7 2196
ce4786b9 2197 delete [] buf;
2198 }
2199
2200 else {
2201 os << "unknown format set" << std::endl;
2202 }
2203
2204 return os;
2205}
8ea391e3 2206
2207
2208void AliTRDmcmSim::PrintFitRegXml(ostream& os) const
2209{
6b094867 2210 // print fit registres in XML format
2211
8ea391e3 2212 bool tracklet=false;
2213
2214 for (Int_t cpu = 0; cpu < 4; cpu++) {
2215 if(fFitPtr[cpu] != 31)
2216 tracklet=true;
2217 }
2218
2219 if(tracklet==true) {
2220 os << "<nginject>" << std::endl;
2221 os << "<ack roc=\""<< fDetector << "\" cmndid=\"0\">" << std::endl;
2222 os << "<dmem-readout>" << std::endl;
2223 os << "<d det=\"" << fDetector << "\">" << std::endl;
2224 os << " <ro-board rob=\"" << fRobPos << "\">" << std::endl;
2225 os << " <m mcm=\"" << fMcmPos << "\">" << std::endl;
5f006bd7 2226
8ea391e3 2227 for(int cpu=0; cpu<4; cpu++) {
2228 os << " <c cpu=\"" << cpu << "\">" << std::endl;
2229 if(fFitPtr[cpu] != 31) {
2230 for(int adcch=fFitPtr[cpu]; adcch<fFitPtr[cpu]+2; adcch++) {
5ac2e3b1 2231 os << " <ch chnr=\"" << adcch << "\">"<< std::endl;
8ea391e3 2232 os << " <hits>" << fFitReg[adcch].fNhits << "</hits>"<< std::endl;
2233 os << " <q0>" << fFitReg[adcch].fQ0/4 << "</q0>"<< std::endl; // divided by 4 because in simulation we have 2 additional decimal places
5f006bd7 2234 os << " <q1>" << fFitReg[adcch].fQ1/4 << "</q1>"<< std::endl; // in the output
8ea391e3 2235 os << " <sumx>" << fFitReg[adcch].fSumX << "</sumx>"<< std::endl;
2236 os << " <sumxsq>" << fFitReg[adcch].fSumX2 << "</sumxsq>"<< std::endl;
2237 os << " <sumy>" << fFitReg[adcch].fSumY << "</sumy>"<< std::endl;
2238 os << " <sumysq>" << fFitReg[adcch].fSumY2 << "</sumysq>"<< std::endl;
2239 os << " <sumxy>" << fFitReg[adcch].fSumXY << "</sumxy>"<< std::endl;
2240 os << " </ch>" << std::endl;
2241 }
2242 }
2243 os << " </c>" << std::endl;
2244 }
2245 os << " </m>" << std::endl;
2246 os << " </ro-board>" << std::endl;
2247 os << "</d>" << std::endl;
2248 os << "</dmem-readout>" << std::endl;
2249 os << "</ack>" << std::endl;
2250 os << "</nginject>" << std::endl;
2251 }
2252}
2253
2254
2255void AliTRDmcmSim::PrintTrackletsXml(ostream& os) const
2256{
6b094867 2257 // print tracklets in XML format
2258
8ea391e3 2259 os << "<nginject>" << std::endl;
2260 os << "<ack roc=\""<< fDetector << "\" cmndid=\"0\">" << std::endl;
2261 os << "<dmem-readout>" << std::endl;
2262 os << "<d det=\"" << fDetector << "\">" << std::endl;
2263 os << " <ro-board rob=\"" << fRobPos << "\">" << std::endl;
2264 os << " <m mcm=\"" << fMcmPos << "\">" << std::endl;
2265
2266 Int_t pid, padrow, slope, offset;
2267 for(Int_t cpu=0; cpu<4; cpu++) {
2268 if(fMCMT[cpu] == 0x10001000) {
2269 pid=-1;
2270 padrow=-1;
2271 slope=-1;
2272 offset=-1;
2273 }
2274 else {
2275 pid = (fMCMT[cpu] & 0xFF000000) >> 24;
2276 padrow = (fMCMT[cpu] & 0xF00000 ) >> 20;
2277 slope = (fMCMT[cpu] & 0xFE000 ) >> 13;
2278 offset = (fMCMT[cpu] & 0x1FFF ) ;
2279
2280 }
5f006bd7 2281 os << " <trk> <pid>" << pid << "</pid>" << " <padrow>" << padrow << "</padrow>"
5ac2e3b1 2282 << " <slope>" << slope << "</slope>" << " <offset>" << offset << "</offset>" << "</trk>" << std::endl;
8ea391e3 2283 }
2284
2285 os << " </m>" << std::endl;
2286 os << " </ro-board>" << std::endl;
2287 os << "</d>" << std::endl;
2288 os << "</dmem-readout>" << std::endl;
2289 os << "</ack>" << std::endl;
2290 os << "</nginject>" << std::endl;
2291}
2292
2293
2294void AliTRDmcmSim::PrintAdcDatHuman(ostream& os) const
2295{
6b094867 2296 // print ADC data in human-readable format
2297
5f006bd7 2298 os << "MCM " << fMcmPos << " on ROB " << fRobPos <<
8ea391e3 2299 " in detector " << fDetector << std::endl;
5f006bd7 2300
8ea391e3 2301 os << "----- Unfiltered ADC data (10 bit) -----" << std::endl;
2302 os << "ch ";
5f006bd7 2303 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++)
8ea391e3 2304 os << std::setw(5) << iChannel;
2305 os << std::endl;
2306 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
2307 os << "tb " << std::setw(2) << iTimeBin << ":";
2308 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2309 os << std::setw(5) << (fADCR[iChannel][iTimeBin] >> fgkAddDigits);
2310 }
2311 os << std::endl;
2312 }
5f006bd7 2313
8ea391e3 2314 os << "----- Filtered ADC data (10+2 bit) -----" << std::endl;
2315 os << "ch ";
5f006bd7 2316 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++)
8ea391e3 2317 os << std::setw(4) << iChannel
2318 << ((~fZSMap[iChannel] != 0) ? "!" : " ");
2319 os << std::endl;
2320 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
2321 os << "tb " << std::setw(2) << iTimeBin << ":";
2322 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2323 os << std::setw(4) << (fADCF[iChannel][iTimeBin])
2324 << (((fZSMap[iChannel] & (1 << iTimeBin)) == 0) ? "!" : " ");
2325 }
2326 os << std::endl;
2327 }
2328}
2329
2330
2331void AliTRDmcmSim::PrintAdcDatXml(ostream& os) const
2332{
5f006bd7 2333 // print ADC data in XML format
6b094867 2334
8ea391e3 2335 os << "<nginject>" << std::endl;
2336 os << "<ack roc=\""<< fDetector << "\" cmndid=\"0\">" << std::endl;
2337 os << "<dmem-readout>" << std::endl;
2338 os << "<d det=\"" << fDetector << "\">" << std::endl;
2339 os << " <ro-board rob=\"" << fRobPos << "\">" << std::endl;
2340 os << " <m mcm=\"" << fMcmPos << "\">" << std::endl;
2341
2342 for(Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2343 os << " <ch chnr=\"" << iChannel << "\">" << std::endl;
2344 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
2345 os << "<tb>" << fADCF[iChannel][iTimeBin]/4 << "</tb>";
2346 }
2347 os << " </ch>" << std::endl;
2348 }
2349
2350 os << " </m>" << std::endl;
2351 os << " </ro-board>" << std::endl;
2352 os << "</d>" << std::endl;
2353 os << "</dmem-readout>" << std::endl;
2354 os << "</ack>" << std::endl;
2355 os << "</nginject>" << std::endl;
2356}
2357
2358
2359
a10c6981 2360void AliTRDmcmSim::PrintAdcDatDatx(ostream& os, Bool_t broadcast, Int_t timeBinOffset) const
8ea391e3 2361{
6b094867 2362 // print ADC data in datx format (to send to FEE)
2363
8ea391e3 2364 fTrapConfig->PrintDatx(os, 2602, 1, 0, 127); // command to enable the ADC clock - necessary to write ADC values to MCM
2365 os << std::endl;
2366
2367 Int_t addrOffset = 0x2000;
2368 Int_t addrStep = 0x80;
2369 Int_t addrOffsetEBSIA = 0x20;
5f006bd7 2370
8ea391e3 2371 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
a10c6981 2372 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2373 if ((iTimeBin < timeBinOffset) || (iTimeBin >= fNTimeBin+timeBinOffset)) {
8ea391e3 2374 if(broadcast==kFALSE)
a10c6981 2375 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, 10, GetRobPos(), GetMcmPos());
8ea391e3 2376 else
a10c6981 2377 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, 10, 0, 127);
2378 }
2379 else {
2380 if(broadcast==kFALSE)
2381 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, (fADCF[iChannel][iTimeBin-timeBinOffset]/4), GetRobPos(), GetMcmPos());
2382 else
2383 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, (fADCF[iChannel][iTimeBin-timeBinOffset]/4), 0, 127);
2384 }
2385 }
2386 os << std::endl;
8ea391e3 2387 }
2388}
2389
2390
2391void AliTRDmcmSim::PrintPidLutHuman()
2392{
6b094867 2393 // print PID LUT in human readable format
2394
8ea391e3 2395 UInt_t result;
2396
2397 UInt_t addrEnd = AliTRDtrapConfig::fgkDmemAddrLUTStart + fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTLength)/4; // /4 because each addr contains 4 values
2398 UInt_t nBinsQ0 = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTnbins);
2399
2400 std::cout << "nBinsQ0: " << nBinsQ0 << std::endl;
5ac2e3b1 2401 std::cout << "LUT table length: " << fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTLength) << std::endl;
5f006bd7 2402
8ea391e3 2403 for(UInt_t addr=AliTRDtrapConfig::fgkDmemAddrLUTStart; addr< addrEnd; addr++) {
2404 result = fTrapConfig->GetDmemUnsigned(addr);
5ac2e3b1 2405 std::cout << addr << " # x: " << ((addr-AliTRDtrapConfig::fgkDmemAddrLUTStart)%((nBinsQ0)/4))*4 << ", y: " <<(addr-AliTRDtrapConfig::fgkDmemAddrLUTStart)/(nBinsQ0/4)
2406 << " # " <<((result>>0)&0xFF)
2407 << " | " << ((result>>8)&0xFF)
2408 << " | " << ((result>>16)&0xFF)
2409 << " | " << ((result>>24)&0xFF) << std::endl;
8ea391e3 2410 }
2411}