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