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