1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // TRD MCM (Multi Chip Module) simulator //
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. //
25 ///////////////////////////////////////////////////////////////////////////////
36 #include "TClonesArray.h"
41 #include "AliRunLoader.h"
42 #include "AliLoader.h"
44 #include "AliTRDfeeParam.h"
45 #include "AliTRDtrapConfigHandler.h"
46 #include "AliTRDtrapConfig.h"
47 #include "AliTRDdigitsManager.h"
48 #include "AliTRDarrayADC.h"
49 #include "AliTRDarrayDictionary.h"
50 #include "AliTRDtrackletMCM.h"
51 #include "AliTRDmcmSim.h"
53 ClassImp(AliTRDmcmSim)
55 Bool_t AliTRDmcmSim::fgApplyCut = kTRUE;
56 Int_t AliTRDmcmSim::fgAddBaseline = 0;
58 const Int_t AliTRDmcmSim::fgkFormatIndex = std::ios_base::xalloc();
60 const Int_t AliTRDmcmSim::fgkNADC = AliTRDfeeParam::GetNadcMcm();
61 const UShort_t AliTRDmcmSim::fgkFPshifts[4] = {11, 14, 17, 21};
64 AliTRDmcmSim::AliTRDmcmSim() :
77 fTrklBranchName("mcmtrklbranch"),
90 // AliTRDmcmSim default constructor
91 // By default, nothing is initialized.
92 // It is necessary to issue Init before use.
94 for (Int_t iDict = 0; iDict < 3; iDict++)
103 AliTRDmcmSim::~AliTRDmcmSim()
106 // AliTRDmcmSim destructor
110 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
111 delete [] fADCR[iAdc];
112 delete [] fADCF[iAdc];
120 delete [] fGainCounterA;
121 delete [] fGainCounterB;
122 delete [] fTailAmplLong;
123 delete [] fTailAmplShort;
126 fTrackletArray->Delete();
127 delete fTrackletArray;
131 void AliTRDmcmSim::Init( Int_t det, Int_t robPos, Int_t mcmPos, Bool_t /* newEvent */ )
134 // Initialize the class with new MCM position information
135 // memory is allocated in the first initialization
139 fFeeParam = AliTRDfeeParam::Instance();
140 fTrapConfig = AliTRDtrapConfigHandler::GetTrapConfig();
146 fRow = fFeeParam->GetPadRowFromMCM( fRobPos, fMcmPos );
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];
154 fNTimeBin = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC13CPUA, fDetector, fRobPos, fMcmPos);
155 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
156 fADCR[iAdc] = new Int_t[fNTimeBin];
157 fADCF[iAdc] = new Int_t[fNTimeBin];
161 fPedAcc = new UInt_t[fgkNADC]; // accumulator for pedestal filter
162 fTailAmplLong = new UShort_t[fgkNADC];
163 fTailAmplShort = new UShort_t[fgkNADC];
165 // tracklet calculation
166 fFitReg = new FitReg_t[fgkNADC];
167 fTrackletArray = new TClonesArray("AliTRDtrackletMCM", fgkMaxTracklets);
169 fMCMT = new UInt_t[fgkMaxTracklets];
172 fInitialized = kTRUE;
177 void AliTRDmcmSim::Reset()
179 // Resets the data values and internal filter registers
180 // by re-initialising them
182 if( !CheckInitialized() )
185 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
186 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
190 fZSMap[iAdc] = -1; // Default unread, low active bit mask
191 fGainCounterA[iAdc] = 0;
192 fGainCounterB[iAdc] = 0;
195 for(Int_t i = 0; i < fgkMaxTracklets; i++) {
199 for (Int_t iDict = 0; iDict < 3; iDict++)
202 FilterPedestalInit();
207 void AliTRDmcmSim::SetNTimebins(Int_t ntimebins)
209 // Reallocate memory if a change in the number of timebins
210 // is needed (should not be the case for real data)
212 if( !CheckInitialized() )
215 fNTimeBin = ntimebins;
216 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
217 delete [] fADCR[iAdc];
218 delete [] fADCF[iAdc];
219 fADCR[iAdc] = new Int_t[fNTimeBin];
220 fADCF[iAdc] = new Int_t[fNTimeBin];
224 Bool_t AliTRDmcmSim::LoadMCM(AliRunLoader* const runloader, Int_t det, Int_t rob, Int_t mcm)
226 // loads the ADC data as obtained from the digitsManager for the specified MCM.
227 // This method is meant for rare execution, e.g. in the visualization. When called
228 // frequently use SetData(...) instead.
233 AliError("No Runloader given");
237 AliLoader *trdLoader = runloader->GetLoader("TRDLoader");
239 AliError("Could not get TRDLoader");
243 Bool_t retval = kTRUE;
244 trdLoader->LoadDigits();
245 fDigitsManager = 0x0;
246 AliTRDdigitsManager *digMgr = new AliTRDdigitsManager();
247 digMgr->SetSDigits(0);
248 digMgr->CreateArrays();
249 digMgr->ReadDigits(trdLoader->TreeD());
250 AliTRDarrayADC *digits = (AliTRDarrayADC*) digMgr->GetDigits(det);
251 if (digits->HasData()) {
254 if (fNTimeBin != digits->GetNtime()) {
255 AliWarning(Form("Changing no. of timebins from %i to %i", fNTimeBin, digits->GetNtime()));
256 SetNTimebins(digits->GetNtime());
269 void AliTRDmcmSim::NoiseTest(Int_t nsamples, Int_t mean, Int_t sigma, Int_t inputGain, Int_t inputTail)
271 // This function can be used to test the filters.
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
275 // With inputGain and inputTail the input to the gain and tail filter, respectively,
276 // can be chosen where
278 // 1: pedestal output
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
282 // same way as in normal simulation.
283 // The functions produces four histograms with the values at the different stages.
285 if( !CheckInitialized() )
288 TString nameInputGain;
289 TString nameInputTail;
293 nameInputGain = "Noise";
297 nameInputGain = "Pedestal";
301 AliError("Undefined input to tail cancellation filter");
307 nameInputTail = "Noise";
311 nameInputTail = "Pedestal";
315 nameInputTail = "Gain";
319 AliError("Undefined input to tail cancellation filter");
323 TH1F *h = new TH1F("noise", "Gaussian Noise;sample;ADC count",
324 nsamples, 0, nsamples);
325 TH1F *hfp = new TH1F("ped", "Noise #rightarrow Pedestal filter;sample;ADC count", nsamples, 0, nsamples);
326 TH1F *hfg = new TH1F("gain",
327 (nameInputGain + "#rightarrow Gain;sample;ADC count").Data(),
328 nsamples, 0, nsamples);
329 TH1F *hft = new TH1F("tail",
330 (nameInputTail + "#rightarrow Tail;sample;ADC count").Data(),
331 nsamples, 0, nsamples);
333 hfp->SetStats(kFALSE);
334 hfg->SetStats(kFALSE);
335 hft->SetStats(kFALSE);
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)
342 for (Int_t i = 0; i < nsamples; i++) {
343 value = (Int_t) gRandom->Gaus(mean, sigma); // generate noise with gaussian distribution
344 h->SetBinContent(i, value);
346 valuep = FilterPedestalNextSample(1, 0, ((Int_t) value) << 2);
349 valueg = FilterGainNextSample(1, ((Int_t) value) << 2);
351 valueg = FilterGainNextSample(1, valuep);
354 valuet = FilterTailNextSample(1, ((Int_t) value) << 2);
355 else if (inputTail == 1)
356 valuet = FilterTailNextSample(1, valuep);
358 valuet = FilterTailNextSample(1, valueg);
360 hfp->SetBinContent(i, valuep >> 2);
361 hfg->SetBinContent(i, valueg >> 2);
362 hft->SetBinContent(i, valuet >> 2);
365 TCanvas *c = new TCanvas;
377 Bool_t AliTRDmcmSim::CheckInitialized() const
380 // Check whether object is initialized
384 AliError(Form ("AliTRDmcmSim is not initialized but function other than Init() is called."));
389 void AliTRDmcmSim::Print(Option_t* const option) const
391 // Prints the data stored and/or calculated for this MCM.
392 // The output is controlled by option which can be a sequence of any of
393 // the following characters:
394 // R - prints raw ADC data
395 // F - prints filtered data
396 // H - prints detected hits
397 // T - prints found tracklets
398 // The later stages are only meaningful after the corresponding calculations
399 // have been performed.
401 if ( !CheckInitialized() )
404 printf("MCM %i on ROB %i in detector %i\n", fMcmPos, fRobPos, fDetector);
406 TString opt = option;
407 if (opt.Contains("R") || opt.Contains("F")) {
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",
415 iHit, fHits[iHit].fTimebin, fHits[iHit].fChannel, fHits[iHit].fQtot, fHits[iHit].fYpos);
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());
427 void AliTRDmcmSim::Draw(Option_t* const option)
429 // Plots the data stored in a 2-dim. timebin vs. ADC channel plot.
430 // The option selects what data is plotted and can be a sequence of
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:
436 // T - plot tracklets
438 if( !CheckInitialized() )
441 TString opt = option;
443 TH2F *hist = new TH2F("mcmdata", Form("Data of MCM %i on ROB %i in detector %i", \
444 fMcmPos, fRobPos, fDetector), \
445 fgkNADC, -0.5, fgkNADC-.5, fNTimeBin, -.5, fNTimeBin-.5);
446 hist->GetXaxis()->SetTitle("ADC Channel");
447 hist->GetYaxis()->SetTitle("Timebin");
448 hist->SetStats(kFALSE);
450 if (opt.Contains("R")) {
451 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
452 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
453 hist->SetBinContent(iAdc+1, iTimeBin+1, fADCR[iAdc][iTimeBin] >> fgkAddDigits);
458 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
459 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
460 hist->SetBinContent(iAdc+1, iTimeBin+1, fADCF[iAdc][iTimeBin] >> fgkAddDigits);
466 if (opt.Contains("H")) {
467 TGraph *grHits = new TGraph();
468 for (Int_t iHit = 0; iHit < fNHits; iHit++) {
469 grHits->SetPoint(iHit,
470 fHits[iHit].fChannel + 1 + fHits[iHit].fYpos/256.,
471 fHits[iHit].fTimebin);
476 if (opt.Contains("T")) {
477 TLine *trklLines = new TLine[4];
478 for (Int_t iTrkl = 0; iTrkl < fTrackletArray->GetEntries(); iTrkl++) {
479 AliTRDtrackletMCM *trkl = (AliTRDtrackletMCM*) (*fTrackletArray)[iTrkl];
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
482 Int_t ndrift = fTrapConfig->GetDmemUnsigned(fgkDmemAddrNdrift, fDetector, fRobPos, fMcmPos) >> 5;
485 slope = trkl->GetdY() * 140e-4 / ndrift;
487 Int_t t0 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS, fDetector, fRobPos, fMcmPos);
488 Int_t t1 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFE, fDetector, fRobPos, fMcmPos);
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);
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();
502 void AliTRDmcmSim::SetData( Int_t adc, const Int_t* const data )
505 // Store ADC data into array of raw data
508 if( !CheckInitialized() ) return;
510 if( adc < 0 || adc >= fgkNADC ) {
511 AliError(Form ("Error: ADC %i is out of range (0 .. %d).", adc, fgkNADC-1));
515 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
516 fADCR[adc][it] = (Int_t) (data[it]) << fgkAddDigits;
517 fADCF[adc][it] = (Int_t) (data[it]) << fgkAddDigits;
521 void AliTRDmcmSim::SetData( Int_t adc, Int_t it, Int_t data )
524 // Store ADC data into array of raw data
527 if( !CheckInitialized() ) return;
529 if( adc < 0 || adc >= fgkNADC ) {
530 AliError(Form ("Error: ADC %i is out of range (0 .. %d).", adc, fgkNADC-1));
534 fADCR[adc][it] = data << fgkAddDigits;
535 fADCF[adc][it] = data << fgkAddDigits;
538 void AliTRDmcmSim::SetData(AliTRDarrayADC* const adcArray, AliTRDdigitsManager * const digitsManager)
540 // Set the ADC data from an AliTRDarrayADC
542 if( !CheckInitialized() )
545 fDigitsManager = digitsManager;
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) {
551 if (fDict[iDict] == newDict)
554 fDict[iDict] = newDict;
555 if(fDict[iDict]->GetDim() != 0)
556 fDict[iDict]->Expand();
559 fDict[iDict] = newDict;
560 if (fDict[iDict] && (fDict[iDict]->GetDim() != 0) )
561 fDict[iDict]->Expand();
564 // If there is no data, set dictionary to zero to avoid crashes
565 if (fDict[iDict]->GetDim() == 0) {
566 // AliError(Form("Dictionary %i of det. %i has dim. 0", iDict, fDetector));
572 if (fNTimeBin != adcArray->GetNtime())
573 SetNTimebins(adcArray->GetNtime());
575 Int_t offset = (fMcmPos % 4 + 1) * 21 + (fRobPos % 2) * 84 - 1;
577 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
578 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
579 Int_t value = adcArray->GetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin);
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)) {
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);
588 fADCR[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
589 fADCF[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
595 void AliTRDmcmSim::SetDataByPad(const AliTRDarrayADC* const adcArray, AliTRDdigitsManager * const digitsManager)
597 // Set the ADC data from an AliTRDarrayADC
598 // (by pad, to be used during initial reading in simulation)
600 if( !CheckInitialized() )
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) {
609 if (fDict[iDict] == newDict)
612 fDict[iDict] = newDict;
613 fDict[iDict]->Expand();
616 fDict[iDict] = newDict;
618 fDict[iDict]->Expand();
621 // If there is no data, set dictionary to zero to avoid crashes
622 if (fDict[iDict]->GetDim() == 0) {
623 AliError(Form("Dictionary %i of det. %i has dim. 0", iDict, fDetector));
629 if (fNTimeBin != adcArray->GetNtime())
630 SetNTimebins(adcArray->GetNtime());
632 Int_t offset = (fMcmPos % 4 + 1) * 18 + (fRobPos % 2) * 72 + 1;
634 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
635 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
637 Int_t pad = offset - iAdc;
638 if (pad > -1 && pad < 144)
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)) {
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);
647 fADCR[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
648 fADCF[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
654 void AliTRDmcmSim::SetDataPedestal( Int_t adc )
657 // Store ADC data into array of raw data
660 if( !CheckInitialized() )
663 if( adc < 0 || adc >= fgkNADC ) {
667 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
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);
673 Bool_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
675 // retrieve the MC hit information (not available in TRAP hardware)
677 if (index < 0 || index >= fNHits)
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) -
685 (channel << 8) - ypos)
686 * (0.635 + 0.03 * (fDetector % 6))
688 label = fHits[index].fLabel[0];
693 Int_t AliTRDmcmSim::GetCol( Int_t adc )
696 // Return column id of the pad for the given ADC channel
699 if( !CheckInitialized() )
702 Int_t col = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adc);
703 if (col < 0 || col >= fFeeParam->GetNcol())
709 Int_t AliTRDmcmSim::ProduceRawStream( UInt_t *buf, Int_t bufSize, UInt_t iEv) const
712 // Produce raw data stream from this MCM and put in buf
713 // Returns number of words filled, or negative value
714 // with -1 * number of overflowed words
717 if( !CheckInitialized() )
721 UInt_t mcmHeader = 0;
723 Int_t nw = 0; // Number of written words
724 Int_t of = 0; // Number of overflowed words
725 Int_t rawVer = fFeeParam->GetRAWversion();
727 Int_t nActiveADC = 0; // number of activated ADC bits in a word
729 if( !CheckInitialized() )
732 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBSF, fDetector, fRobPos, fMcmPos) != 0) // store unfiltered data
737 // Produce ADC mask : nncc cccm mmmm mmmm mmmm mmmm mmmm 1100
738 // n : unused , c : ADC count, m : selected ADCs
740 (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC15CPUA, fDetector, fRobPos, fMcmPos) & (1 << 13))) { // check for zs flag in TRAP configuration
741 for( Int_t iAdc = 0 ; iAdc < fgkNADC ; iAdc++ ) {
742 if( ~fZSMap[iAdc] != 0 ) { // 0 means not suppressed
743 adcMask |= (1 << (iAdc+4) ); // last 4 digit reserved for 1100=0xc
744 nActiveADC++; // number of 1 in mmm....m
748 if ((nActiveADC == 0) &&
749 (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC15CPUA, fDetector, fRobPos, fMcmPos) & (1 << 8))) // check for DEH flag in TRAP configuration
752 // assemble adc mask word
753 adcMask |= (1 << 30) | ( ( 0x3FFFFFFC ) & (~(nActiveADC) << 25) ) | 0xC; // nn = 01, ccccc are inverted, 0xc=1100
757 mcmHeader = (1<<31) | (fRobPos << 28) | (fMcmPos << 24) | ((iEv % 0x100000) << 4) | 0xC;
759 buf[nw++] = mcmHeader;
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
774 UInt_t aa=0, a1=0, a2=0, a3=0;
776 for (Int_t iAdc = 0; iAdc < 21; iAdc++ ) {
777 if( rawVer>= 3 && ~fZSMap[iAdc] == 0 ) continue; // Zero Suppression, 0 means not suppressed
778 aa = !(iAdc & 1) + 2;
779 for (Int_t iT = 0; iT < fNTimeBin; iT+=3 ) {
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;
783 x = (a3 << 22) | (a2 << 12) | (a1 << 2) | aa;
793 if( of != 0 ) return -of; else return nw;
796 Int_t AliTRDmcmSim::ProduceTrackletStream( UInt_t *buf, Int_t bufSize )
799 // Produce tracklet data stream from this MCM and put in buf
800 // Returns number of words filled, or negative value
801 // with -1 * number of overflowed words
804 if( !CheckInitialized() )
807 Int_t nw = 0; // Number of written words
808 Int_t of = 0; // Number of overflowed words
810 // Produce tracklet data. A maximum of four 32 Bit words will be written per MCM
811 // fMCMT is filled continuously until no more tracklet words available
813 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
815 buf[nw++] = ((AliTRDtrackletMCM*) (*fTrackletArray)[iTracklet])->GetTrackletWord();
820 if( of != 0 ) return -of; else return nw;
823 void AliTRDmcmSim::Filter()
826 // Filter the raw ADC values. The active filter stages and their
827 // parameters are taken from AliTRDtrapConfig.
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
830 // sequentially for parameter tuning.
833 if( !CheckInitialized() )
836 // Apply filters sequentially. Bypass is handled by filters
837 // since counters and internal registers may be updated even
838 // if the filter is bypassed.
839 // The first filter takes the data from fADCR and
842 // Non-linearity filter not implemented.
846 // Crosstalk filter not implemented.
849 void AliTRDmcmSim::FilterPedestalInit(Int_t baseline)
851 // Initializes the pedestal filter assuming that the input has
852 // been constant for a long time (compared to the time constant).
854 UShort_t fptc = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPTC, fDetector, fRobPos, fMcmPos); // 0..3, 0 - fastest, 3 - slowest
856 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++)
857 fPedAcc[iAdc] = (baseline << 2) * (1 << fgkFPshifts[fptc]);
860 UShort_t AliTRDmcmSim::FilterPedestalNextSample(Int_t adc, Int_t timebin, UShort_t value)
862 // Returns the output of the pedestal filter given the input value.
863 // The output depends on the internal registers and, thus, the
864 // history of the filter.
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
870 UShort_t accumulatorShifted;
874 inpAdd = value + fpnp;
876 accumulatorShifted = (fPedAcc[adc] >> fgkFPshifts[fptc]) & 0x3FF; // 10 bits
877 if (timebin == 0) // the accumulator is disabled in the drift time
879 correction = (value & 0x3FF) - accumulatorShifted;
880 fPedAcc[adc] = (fPedAcc[adc] + correction) & 0x7FFFFFFF; // 31 bits
886 if (inpAdd <= accumulatorShifted)
890 inpAdd = inpAdd - accumulatorShifted;
898 void AliTRDmcmSim::FilterPedestal()
901 // Apply pedestal filter
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
907 // the input has been stable for a sufficiently long time.
909 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
910 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
911 fADCF[iAdc][iTimeBin] = FilterPedestalNextSample(iAdc, iTimeBin, fADCR[iAdc][iTimeBin]);
916 void AliTRDmcmSim::FilterGainInit()
918 // Initializes the gain filter. In this case, only threshold
919 // counters are reset.
921 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
922 // these are counters which in hardware continue
923 // until maximum or reset
924 fGainCounterA[iAdc] = 0;
925 fGainCounterB[iAdc] = 0;
929 UShort_t AliTRDmcmSim::FilterGainNextSample(Int_t adc, UShort_t value)
931 // Apply the gain filter to the given value.
932 // BEGIN_LATEX O_{i}(t) = #gamma_{i} * I_{i}(t) + a_{i} END_LATEX
933 // The output depends on the internal registers and, thus, the
934 // history of the filter.
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;
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)
946 UInt_t corr; // corrected value
949 corr = (value * fgfExtended) >> 11;
950 corr = corr > 0xfff ? 0xfff : corr;
951 corr = AddUintClipping(corr, fga, 12);
953 // Update threshold counters
954 // not really useful as they are cleared with every new event
955 if (!((fGainCounterA[adc] == 0x3FFFFFF) || (fGainCounterB[adc] == 0x3FFFFFF)))
959 fGainCounterB[adc]++;
960 else if (corr >= fgta)
961 fGainCounterA[adc]++;
970 void AliTRDmcmSim::FilterGain()
972 // Read data from fADCF and apply gain filter.
974 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
975 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
976 fADCF[iAdc][iTimeBin] = FilterGainNextSample(iAdc, fADCF[iAdc][iTimeBin]);
981 void AliTRDmcmSim::FilterTailInit(Int_t baseline)
983 // Initializes the tail filter assuming that the input has
984 // been at the baseline value (configured by FTFP) for a
985 // sufficiently long time.
987 // exponents and weight calculated from configuration
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
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);
996 qup = (1 - lambdaL) * (1 - lambdaS);
997 qdn = 1 - lambdaS * alphaL - lambdaL * (1 - alphaL);
998 Float_t kdc = qup/qdn;
1004 baseline = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP, fDetector, fRobPos, fMcmPos);
1006 ql = lambdaL * (1 - lambdaS) * alphaL;
1007 qs = lambdaS * (1 - lambdaL) * (1 - alphaL);
1009 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1010 Int_t value = baseline & 0xFFF;
1011 Int_t corr = (value * fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGF0 + iAdc), fDetector, fRobPos, fMcmPos)) >> 11;
1012 corr = corr > 0xfff ? 0xfff : corr;
1013 corr = AddUintClipping(corr, fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGA0 + iAdc), fDetector, fRobPos, fMcmPos), 12);
1015 kt = kdc * baseline;
1016 aout = baseline - (UShort_t) kt;
1018 fTailAmplLong[iAdc] = (UShort_t) (aout * ql / (ql + qs));
1019 fTailAmplShort[iAdc] = (UShort_t) (aout * qs / (ql + qs));
1023 UShort_t AliTRDmcmSim::FilterTailNextSample(Int_t adc, UShort_t value)
1025 // Returns the output of the tail filter for the given input value.
1026 // The output depends on the internal registers and, thus, the
1027 // history of the filter.
1029 // exponents and weight calculated from configuration
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
1034 // intermediate signals
1040 UShort_t inpVolt = value & 0xFFF; // 12 bits
1042 // add the present generator outputs
1043 aQ = AddUintClipping(fTailAmplLong[adc], fTailAmplShort[adc], 12);
1045 // calculate the difference between the input and the generated signal
1047 aDiff = inpVolt - aQ;
1051 // the inputs to the two generators, weighted
1052 alInpv = (aDiff * alphaLong) >> 11;
1054 // the new values of the registers, used next time
1056 tmp = AddUintClipping(fTailAmplLong[adc], alInpv, 12);
1057 tmp = (tmp * lambdaLong) >> 11;
1058 fTailAmplLong[adc] = tmp & 0xFFF;
1060 tmp = AddUintClipping(fTailAmplShort[adc], aDiff - alInpv, 12);
1061 tmp = (tmp * lambdaShort) >> 11;
1062 fTailAmplShort[adc] = tmp & 0xFFF;
1064 // the output of the filter
1065 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTBY, fDetector, fRobPos, fMcmPos) == 0) // bypass mode, active low
1071 void AliTRDmcmSim::FilterTail()
1073 // Apply tail cancellation filter to all data.
1075 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1076 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1077 fADCF[iAdc][iTimeBin] = FilterTailNextSample(iAdc, fADCF[iAdc][iTimeBin]);
1082 void AliTRDmcmSim::ZSMapping()
1085 // Zero Suppression Mapping implemented in TRAP chip
1086 // only implemented for up to 30 timebins
1088 // See detail TRAP manual "Data Indication" section:
1089 // http://www.kip.uni-heidelberg.de/ti/TRD/doc/trap/TRAP-UserManual.pdf
1092 if( !CheckInitialized() )
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);
1100 Int_t **adc = fADCF;
1102 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++)
1105 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
1106 Int_t iAdc; // current ADC channel
1111 Int_t supp; // suppression of the current channel (low active)
1113 // ----- first channel -----
1117 ac = adc[iAdc ][it]; // current
1118 an = adc[iAdc+1][it]; // next
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
1124 supp = (eBIL >> mask) & 1;
1126 fZSMap[iAdc] &= ~((1-supp) << it);
1127 if( eBIN == 0 ) { // neighbour sensitivity
1128 fZSMap[iAdc+1] &= ~((1-supp) << it);
1131 // ----- last channel -----
1134 ap = adc[iAdc-1][it]; // previous
1135 ac = adc[iAdc ][it]; // current
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
1142 supp = (eBIL >> mask) & 1;
1144 fZSMap[iAdc] &= ~((1-supp) << it);
1145 if( eBIN == 0 ) { // neighbour sensitivity
1146 fZSMap[iAdc-1] &= ~((1-supp) << it);
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
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
1159 supp = (eBIL >> mask) & 1;
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);
1171 void AliTRDmcmSim::AddHitToFitreg(Int_t adc, UShort_t timebin, UShort_t qtot, Short_t ypos, Int_t label[])
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
1178 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0, fDetector, fRobPos, fMcmPos)) &&
1179 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE0, fDetector, fRobPos, fMcmPos)))
1180 fFitReg[adc].fQ0 += qtot;
1182 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS1, fDetector, fRobPos, fMcmPos)) &&
1183 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1, fDetector, fRobPos, fMcmPos)))
1184 fFitReg[adc].fQ1 += qtot;
1186 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS, fDetector, fRobPos, fMcmPos) ) &&
1187 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFE, fDetector, fRobPos, fMcmPos)))
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;
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));
1200 // register hits (MC info)
1201 fHits[fNHits].fChannel = adc;
1202 fHits[fNHits].fQtot = qtot;
1203 fHits[fNHits].fYpos = ypos;
1204 fHits[fNHits].fTimebin = timebin;
1205 fHits[fNHits].fLabel[0] = label[0];
1206 fHits[fNHits].fLabel[1] = label[1];
1207 fHits[fNHits].fLabel[2] = label[2];
1211 void AliTRDmcmSim::CalcFitreg()
1214 // Detect the hits and fill the fit registers.
1215 // Requires 12-bit data from fADCF which means Filter()
1216 // has to be called before even if all filters are bypassed.
1218 //??? to be clarified:
1219 UInt_t adcMask = 0xffffffff;
1222 Int_t adcLeft, adcCentral, adcRight;
1223 UShort_t timebin, adcch, timebin1, timebin2, qtotTemp;
1224 Short_t ypos, fromLeft, fromRight, found;
1225 UShort_t qTotal[19+1]; // the last is dummy
1226 UShort_t marked[6], qMarked[6], worse1, worse2;
1228 timebin1 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS, fDetector, fRobPos, fMcmPos);
1229 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0, fDetector, fRobPos, fMcmPos)
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)
1235 timebin2 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1, fDetector, fRobPos, fMcmPos);
1237 // reset the fit registers
1239 for (adcch = 0; adcch < fgkNADC-2; adcch++) // due to border channels
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;
1251 for (timebin = timebin1; timebin < timebin2; timebin++)
1253 // first find the hit candidates and store the total cluster charge in qTotal array
1254 // in case of not hit store 0 there.
1255 for (adcch = 0; adcch < fgkNADC-2; adcch++) {
1256 if ( ( (adcMask >> adcch) & 7) == 7) //??? all 3 channels are present in case of ZS
1258 adcLeft = fADCF[adcch ][timebin];
1259 adcCentral = fADCF[adcch+1][timebin];
1260 adcRight = fADCF[adcch+2][timebin];
1262 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPVBY, fDetector, fRobPos, fMcmPos) == 0) {
1263 // bypass the cluster verification
1267 hitQual = ( (adcLeft * adcRight) <
1268 ((fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPVT, fDetector, fRobPos, fMcmPos) * adcCentral*adcCentral) >> 10) );
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));
1276 // The accumulated charge is with the pedestal!!!
1277 qtotTemp = adcLeft + adcCentral + adcRight;
1279 (qtotTemp >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPHT, fDetector, fRobPos, fMcmPos)) &&
1280 (adcLeft <= adcCentral) &&
1281 (adcCentral > adcRight) )
1282 qTotal[adcch] = qtotTemp;
1287 qTotal[adcch] = 0; //jkl
1288 if (qTotal[adcch] != 0)
1289 AliDebug(10,Form("ch %2d qTotal %5d",adcch, qTotal[adcch]));
1295 marked[4] = 19; // invalid channel
1296 marked[5] = 19; // invalid channel
1298 while ((adcch < 16) && (found < 3))
1300 if (qTotal[adcch] > 0)
1303 marked[2*found+1]=adcch;
1312 while ((adcch > 2) && (found < 3))
1314 if (qTotal[adcch] > 0)
1316 marked[2*found]=adcch;
1323 AliDebug(10,Form("Fromleft=%d, Fromright=%d",fromLeft, fromRight));
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++)
1330 for (adcch = 0; adcch < 19; adcch++)
1331 if (qTotal[adcch] > 0) found++;
1334 if (found > 4) // sorting like in the TRAP in case of 5 or 6 candidates!
1336 if (marked[4] == marked[5]) marked[5] = 19;
1337 for (found=0; found<6; found++)
1339 qMarked[found] = qTotal[marked[found]] >> 4;
1340 AliDebug(10,Form("ch_%d qTotal %d qTotals %d",marked[found],qTotal[marked[found]],qMarked[found]));
1343 Sort6To2Worst(marked[0], marked[3], marked[4], marked[1], marked[2], marked[5],
1351 // Now mask the two channels with the smallest charge
1355 AliDebug(10,Form("Kill ch %d\n",worse1));
1360 AliDebug(10,Form("Kill ch %d\n",worse2));
1364 for (adcch = 0; adcch < 19; adcch++) {
1365 if (qTotal[adcch] > 0) // the channel is marked for processing
1367 adcLeft = fADCF[adcch ][timebin];
1368 adcCentral = fADCF[adcch+1][timebin];
1369 adcRight = fADCF[adcch+2][timebin];
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
1373 Int_t regTPFP = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP, fDetector, fRobPos, fMcmPos);
1374 AliDebug(10, Form("Hit found, time=%d, adcch=%d/%d/%d, adc values=%d/%d/%d, regTPFP=%d, TPHT=%d\n",
1375 timebin, adcch, adcch+1, adcch+2, adcLeft, adcCentral, adcRight, regTPFP,
1376 fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPHT, fDetector, fRobPos, fMcmPos)));
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;
1382 // Calculate the center of gravity
1383 // checking for adcCentral != 0 (in case of "bad" configuration)
1384 if (adcCentral == 0)
1386 ypos = 128*(adcRight - adcLeft) / adcCentral;
1387 if (ypos < 0) ypos = -ypos;
1388 // make the correction using the position LUT
1389 ypos = ypos + fTrapConfig->GetTrapReg((AliTRDtrapConfig::TrapReg_t) (AliTRDtrapConfig::kTPL00 + (ypos & 0x7F)),
1390 fDetector, fRobPos, fMcmPos);
1391 if (adcLeft > adcRight) ypos = -ypos;
1393 // label calculation (up to 3)
1394 Int_t mcLabel[] = {-1, -1, -1};
1395 if (fDigitsManager) {
1396 const Int_t maxLabels = 9;
1397 Int_t label[maxLabels] = { 0 }; // up to 9 different labels possible
1398 Int_t count[maxLabels] = { 0 };
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++) {
1408 for (Int_t iPad = 0; iPad < 3; iPad++) {
1409 if (padcol[iPad] < 0)
1411 Int_t currLabel = fDict[iDict]->GetData(padrow, padcol[iPad], timebin);
1412 AliDebug(10, Form("Read label: %4i for det: %3i, row: %i, col: %i, tb: %i\n", currLabel, fDetector, padrow, padcol[iPad], timebin));
1413 for (Int_t iLabel = 0; iLabel < nLabels; iLabel++) {
1414 if (currLabel == label[iLabel]) {
1420 if (currLabel >= 0) {
1421 label[nLabels] = currLabel;
1427 Int_t index[2*maxLabels];
1428 TMath::Sort(maxLabels, count, index);
1429 for (Int_t i = 0; i < 3; i++) {
1430 if (count[index[i]] <= 0)
1432 mcLabel[i] = label[index[i]];
1436 // add the hit to the fitregister
1437 AddHitToFitreg(adcch, timebin, qTotal[adcch] >> fgkAddDigits, ypos, mcLabel);
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
1456 void AliTRDmcmSim::TrackletSelection()
1458 // Select up to 4 tracklet candidates from the fit registers
1459 // and assign them to the CPUs.
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
1465 for (adcIdx = 0; adcIdx < 18; adcIdx++) // ADCs
1466 if ( (fFitReg[adcIdx].fNhits
1467 >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPCL, fDetector, fRobPos, fMcmPos)) &&
1468 (fFitReg[adcIdx].fNhits+fFitReg[adcIdx+1].fNhits
1469 >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPCT, fDetector, fRobPos, fMcmPos)))
1471 trackletCand[ntracks][0] = adcIdx;
1472 trackletCand[ntracks][1] = fFitReg[adcIdx].fNhits+fFitReg[adcIdx+1].fNhits;
1473 AliDebug(10,Form("%d %2d %4d\n", ntracks, trackletCand[ntracks][0], trackletCand[ntracks][1]));
1477 for (i=0; i<ntracks;i++)
1478 AliDebug(10,Form("%d %d %d\n",i,trackletCand[i][0], trackletCand[i][1]));
1482 // primitive sorting according to the number of hits
1483 for (j = 0; j < (ntracks-1); j++)
1485 for (i = j+1; i < ntracks; i++)
1487 if ( (trackletCand[j][1] < trackletCand[i][1]) ||
1488 ( (trackletCand[j][1] == trackletCand[i][1]) && (trackletCand[j][0] < trackletCand[i][0]) ) )
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;
1500 ntracks = 4; // cut the rest, 4 is the max
1502 // else is not necessary to sort
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++)
1507 for (i = j+1; i < ntracks; i++)
1509 if (trackletCand[j][0] < trackletCand[i][0])
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;
1521 for (i = 0; i < ntracks; i++) // CPUs with tracklets.
1522 fFitPtr[i] = trackletCand[i][0]; // pointer to the left channel with tracklet for CPU[i]
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)
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]));
1530 void AliTRDmcmSim::FitTracklet()
1532 // Perform the actual tracklet fit based on the fit sums
1533 // which have been filled in the fit registers.
1535 // parameters in fitred.asm (fit program)
1537 Int_t decPlaces = 5; // must be larger than 1 or change the following code
1538 // if (decPlaces > 1)
1539 rndAdd = (1 << (decPlaces-1)) + 1;
1540 // else if (decPlaces == 1)
1543 Int_t ndriftDp = 5; // decimal places for drift time
1544 Long64_t shift = ((Long64_t) 1 << 32);
1546 // calculated in fitred.asm
1547 Int_t padrow = ((fRobPos >> 1) << 2) | (fMcmPos >> 2);
1548 Int_t yoffs = (((((fRobPos & 0x1) << 2) + (fMcmPos & 0x3)) * 18) << 8) -
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);
1555 Int_t deflCorr = (Int_t) fTrapConfig->GetDmemUnsigned(fgkDmemAddrDeflCorr, fDetector, fRobPos, fMcmPos);
1556 Int_t ndrift = (Int_t) fTrapConfig->GetDmemUnsigned(fgkDmemAddrNdrift, fDetector, fRobPos, fMcmPos);
1558 // local variables for calculation
1559 Long64_t mult, temp, denom; //???
1560 UInt_t q0, q1, pid; // charges in the two windows and total charge
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
1564 Int_t sumY2; // not used in the current TRAP program, now used for error calculation (simulation only)
1565 Float_t fitError, fitSlope, fitOffset;
1566 FitReg_t *fit0, *fit1; // pointers to relevant fit registers
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};
1574 for (Int_t cpu = 0; cpu < 4; cpu++) {
1575 if (fFitPtr[cpu] == 31)
1577 fMCMT[cpu] = 0x10001000; //??? AliTRDfeeParam::GetTrackletEndmarker();
1581 fit0 = &fFitReg[fFitPtr[cpu] ];
1582 fit1 = &fFitReg[fFitPtr[cpu]+1]; // next channel
1585 mult = mult << (32 + decPlaces);
1589 nHits = fit0->fNhits + fit1->fNhits; // number of hits
1590 sumX = fit0->fSumX + fit1->fSumX;
1591 sumX2 = fit0->fSumX2 + fit1->fSumX2;
1592 denom = ((Long64_t) nHits)*((Long64_t) sumX2) - ((Long64_t) sumX)*((Long64_t) sumX);
1594 mult = mult / denom; // exactly like in the TRAP program
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;
1599 sumY2 = fit0->fSumY2 + fit1->fSumY2 + 512*fit1->fSumY + 256*256*fit1->fNhits;
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
1606 temp = mult * offset;
1607 offset = temp >> 32; // take the upper 32 bits
1609 offset = offset + yoffs;
1610 AliDebug(10, Form("slope = %i, slope * ndrift = %i, deflCorr: %i",
1611 slope, slope * ndrift, deflCorr));
1612 slope = ((slope * ndrift) >> ndriftDp) + deflCorr;
1613 offset = offset - (fFitPtr[cpu] << (8 + decPlaces));
1616 temp = temp * scaleD;
1617 slope = (temp >> 32);
1619 temp = temp * scaleY;
1620 offset = (temp >> 32);
1622 // rounding, like in the TRAP
1623 slope = (slope + rndAdd) >> decPlaces;
1624 offset = (offset + rndAdd) >> decPlaces;
1626 AliDebug(5, Form("Det: %3i, ROB: %i, MCM: %2i: deflection: %i, min: %i, max: %i",
1627 fDetector, fRobPos, fMcmPos, slope,
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)));
1631 AliDebug(5, Form("Fit sums: x = %i, X = %i, y = %i, Y = %i, Z = %i",
1632 sumX, sumX2, sumY, sumY2, sumXY));
1634 fitSlope = (Float_t) (nHits * sumXY - sumX * sumY) / (nHits * sumX2 - sumX*sumX);
1636 fitOffset = (Float_t) (sumX2 * sumY - sumX * sumXY) / (nHits * sumX2 - sumX*sumX);
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);
1646 Bool_t rejected = kFALSE;
1647 // deflection range table from DMEM
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))))
1652 if (rejected && GetApplyCut())
1654 fMCMT[cpu] = 0x10001000; //??? AliTRDfeeParam::GetTrackletEndmarker();
1658 if (slope > 63 || slope < -64) { // wrapping in TRAP!
1659 AliDebug(1,Form("Overflow in slope: %i, tracklet discarded!", slope));
1660 fMCMT[cpu] = 0x10001000;
1664 slope = slope & 0x7F; // 7 bit
1666 if (offset > 0xfff || offset < -0xfff)
1667 AliWarning("Overflow in offset");
1668 offset = offset & 0x1FFF; // 13 bit
1670 pid = GetPID(q0, q1);
1673 AliWarning("Overflow in PID");
1674 pid = pid & 0xFF; // 8 bit, exactly like in the TRAP program
1676 // assemble and store the tracklet word
1677 fMCMT[cpu] = (pid << 24) | (padrow << 20) | (slope << 13) | offset;
1679 // calculate MC label
1680 Int_t mcLabel[] = { -1, -1, -1};
1683 if (fDigitsManager) {
1684 const Int_t maxLabels = 30;
1685 Int_t label[maxLabels] = {0}; // up to 30 different labels possible
1686 Int_t count[maxLabels] = {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))
1693 // counting contributing hits
1694 if (fHits[iHit].fTimebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0, fDetector, fRobPos, fMcmPos) &&
1695 fHits[iHit].fTimebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE0, fDetector, fRobPos, fMcmPos))
1697 if (fHits[iHit].fTimebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS1, fDetector, fRobPos, fMcmPos) &&
1698 fHits[iHit].fTimebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1, fDetector, fRobPos, fMcmPos))
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]) {
1710 if (currLabel >= 0 && nLabels < maxLabels) {
1711 label[nLabels] = currLabel;
1717 Int_t index[2*maxLabels];
1718 TMath::Sort(maxLabels, count, index);
1719 for (Int_t i = 0; i < 3; i++) {
1720 if (count[index[i]] <= 0)
1722 mcLabel[i] = label[index[i]];
1725 new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletMCM((UInt_t) fMCMT[cpu], fDetector*2 + fRobPos%2, fRobPos, fMcmPos);
1726 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetLabel(mcLabel);
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);
1732 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetQ0(q0);
1733 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetQ1(q1);
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));
1738 // // cluster information
1739 // Float_t *res = new Float_t[nHits];
1740 // Float_t *qtot = new Float_t[nHits];
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;
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;
1755 // ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetClusters(res, qtot, nCls);
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));
1762 AliDebug(3, Form("fit slope: %f, offset: %f, error: %f",
1763 fitSlope, fitOffset, TMath::Sqrt(TMath::Abs(fitError)/nHits)));
1769 void AliTRDmcmSim::Tracklet()
1771 // Run the tracklet calculation by calling sequentially:
1772 // CalcFitreg(); TrackletSelection(); FitTracklet()
1773 // and store the tracklets
1775 if (!fInitialized) {
1776 AliError("Called uninitialized! Nothing done!");
1780 fTrackletArray->Delete();
1785 TrackletSelection();
1789 Bool_t AliTRDmcmSim::StoreTracklets()
1791 // store the found tracklets via the loader
1793 if (fTrackletArray->GetEntriesFast() == 0)
1796 AliRunLoader *rl = AliRunLoader::Instance();
1797 AliDataLoader *dl = 0x0;
1799 dl = rl->GetLoader("TRDLoader")->GetDataLoader("tracklets");
1801 AliError("Could not get the tracklets data loader!");
1805 TTree *trackletTree = dl->Tree();
1806 if (!trackletTree) {
1808 trackletTree = dl->Tree();
1811 AliTRDtrackletMCM *trkl = 0x0;
1812 TBranch *trkbranch = trackletTree->GetBranch(fTrklBranchName.Data());
1814 trkbranch = trackletTree->Branch(fTrklBranchName.Data(), "AliTRDtrackletMCM", &trkl, 32000);
1816 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
1817 trkl = ((AliTRDtrackletMCM*) (*fTrackletArray)[iTracklet]);
1818 trkbranch->SetAddress(&trkl);
1825 void AliTRDmcmSim::WriteData(AliTRDarrayADC *digits)
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
1831 if( !CheckInitialized() )
1834 Int_t offset = (fMcmPos % 4 + 1) * 21 + (fRobPos % 2) * 84 - 1;
1836 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBSF, fDetector, fRobPos, fMcmPos) != 0) // store unfiltered data
1838 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1839 if (~fZSMap[iAdc] == 0) {
1840 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1841 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, -1);
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);
1852 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1853 if (~fZSMap[iAdc] != 0) {
1854 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1855 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, (fADCF[iAdc][iTimeBin] >> fgkAddDigits) - fgAddBaseline);
1859 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1860 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, -1);
1868 // ******************************
1871 // Memory area for the LUT: 0xC100 to 0xC3FF
1873 // The addresses for the parameters (the order is optimized for maximum calculation speed in the MCMs):
1875 // 0xC029: nBins(sF)
1877 // 0xC02B: TableLength
1878 // Defined in AliTRDtrapConfig.h
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++
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
1900 Int_t AliTRDmcmSim::GetPID(Int_t q0, Int_t q1)
1902 // return PID calculated from charges accumulated in two time windows
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);
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
1912 ULong_t corrQ0 = fTrapConfig->GetDmemUnsigned(fgkDmemAddrLUTcor0, fDetector, fRobPos, fMcmPos);
1913 ULong_t corrQ1 = fTrapConfig->GetDmemUnsigned(fgkDmemAddrLUTcor1, fDetector, fRobPos, fMcmPos);
1914 if(corrQ0==0) // make sure we don't run into trouble if one of the values is not configured
1918 addrQ0 = (((addrQ0*q0)>>16)>>16); // because addrQ0 = (q0 * corrQ0) >> 32; does not work for unknown reasons
1920 if(addrQ0 >= nBinsQ0) { // check for overflow
1921 AliDebug(5,Form("Overflow in q0: %llu/4 is bigger then %u", addrQ0, nBinsQ0));
1922 addrQ0 = nBinsQ0 -1;
1926 addr = (((addr*q1)>>16)>>16);
1927 addr = addrQ0 + nBinsQ0*addr; // because addr = addrQ0 + nBinsQ0* (((corrQ1*q1)>>32); does not work
1929 if(addr >= pidTotalSize) {
1930 AliDebug(5,Form("Overflow in q1. Address %llu/4 is bigger then %u", addr, pidTotalSize));
1931 addr = pidTotalSize -1;
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)
1936 UInt_t result = fTrapConfig->GetDmemUnsigned(fgkDmemAddrLUTStart+(addr/4), fDetector, fRobPos, fMcmPos);
1937 return (result>>((addr%4)*8)) & 0xFF;
1942 // help functions, to be cleaned up
1944 UInt_t AliTRDmcmSim::AddUintClipping(UInt_t a, UInt_t b, UInt_t nbits) const
1947 // This function adds a and b (unsigned) and clips to
1948 // the specified number of bits.
1954 UInt_t maxv = (1 << nbits) - 1;;
1960 if ((sum < a) || (sum < b))
1966 void AliTRDmcmSim::Sort2(UShort_t idx1i, UShort_t idx2i, \
1967 UShort_t val1i, UShort_t val2i, \
1968 UShort_t * const idx1o, UShort_t * const idx2o, \
1969 UShort_t * const val1o, UShort_t * const val2o) const
1971 // sorting for tracklet selection
1989 void AliTRDmcmSim::Sort3(UShort_t idx1i, UShort_t idx2i, UShort_t idx3i, \
1990 UShort_t val1i, UShort_t val2i, UShort_t val3i, \
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)
1994 // sorting for tracklet selection
1999 if (val1i > val2i) sel=4; else sel=0;
2000 if (val2i > val3i) sel=sel + 2;
2001 if (val3i > val1i) sel=sel + 1;
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!
2014 case 4 : // 1 > 2, 2 <= 3, 3 <= 1 => 1 3 2
2023 case 2 : // 1 <= 2, 2 > 3, 3 <= 1 => 2 1 3
2032 case 3 : // 1 <= 2, 2 > 3, 3 > 1 => 2 3 1
2041 case 1 : // 1 <= 2, 2 <= 3, 3 > 1 => 3 2 1
2050 case 5 : // 1 > 2, 2 <= 3, 3 > 1 => 3 1 2
2059 default: // the rest should NEVER happen!
2060 AliError("ERROR in Sort3!!!\n");
2065 void 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, \
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)
2070 // sorting for tracklet selection
2072 UShort_t idx21s, idx22s, idx23s, dummy;
2073 UShort_t val21s, val22s, val23s;
2074 UShort_t idx23as, idx23bs;
2075 UShort_t val23as, val23bs;
2077 Sort3(idx1i, idx2i, idx3i, val1i, val2i, val3i,
2078 idx1o, &idx21s, &idx23as,
2079 val1o, &val21s, &val23as);
2081 Sort3(idx4i, idx5i, idx6i, val4i, val5i, val6i,
2082 idx2o, &idx22s, &idx23bs,
2083 val2o, &val22s, &val23bs);
2085 Sort2(idx23as, idx23bs, val23as, val23bs, &idx23s, &dummy, &val23s, &dummy);
2087 Sort3(idx21s, idx22s, idx23s, val21s, val22s, val23s,
2088 idx3o, idx4o, &dummy,
2089 val3o, val4o, &dummy);
2093 void 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, \
2095 UShort_t * const idx5o, UShort_t * const idx6o)
2097 // sorting for tracklet selection
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;
2104 Sort3(idx1i, idx2i, idx3i, val1i, val2i, val3i,
2105 &dummy1, &idx21s, &idx23as,
2106 &dummy2, &val21s, &val23as);
2108 Sort3(idx4i, idx5i, idx6i, val4i, val5i, val6i,
2109 &dummy1, &idx22s, &idx23bs,
2110 &dummy2, &val22s, &val23bs);
2112 Sort2(idx23as, idx23bs, val23as, val23bs, &idx23s, idx5o, &val23s, &dummy1);
2114 Sort3(idx21s, idx22s, idx23s, val21s, val22s, val23s,
2115 &dummy1, &dummy2, idx6o,
2116 &dummy3, &dummy4, &dummy5);
2120 // ----- I/O implementation -----
2122 ostream& AliTRDmcmSim::Text(ostream& os)
2124 // manipulator to activate output in text format (default)
2126 os.iword(fgkFormatIndex) = 0;
2130 ostream& AliTRDmcmSim::Cfdat(ostream& os)
2132 // manipulator to activate output in CFDAT format
2133 // to send to the FEE via SCSN
2135 os.iword(fgkFormatIndex) = 1;
2139 ostream& AliTRDmcmSim::Raw(ostream& os)
2141 // manipulator to activate output as raw data dump
2143 os.iword(fgkFormatIndex) = 2;
2147 ostream& operator<<(ostream& os, const AliTRDmcmSim& mcm)
2149 // output implementation
2151 // no output for non-initialized MCM
2152 if (!mcm.CheckInitialized())
2155 // ----- human-readable output -----
2156 if (os.iword(AliTRDmcmSim::fgkFormatIndex) == 0) {
2158 os << "MCM " << mcm.fMcmPos << " on ROB " << mcm.fRobPos <<
2159 " in detector " << mcm.fDetector << std::endl;
2161 os << "----- Unfiltered ADC data (10 bit) -----" << std::endl;
2163 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++)
2164 os << std::setw(5) << iChannel;
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);
2174 os << "----- Filtered ADC data (10+2 bit) -----" << std::endl;
2176 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++)
2177 os << std::setw(4) << iChannel
2178 << ((~mcm.fZSMap[iChannel] != 0) ? "!" : " ");
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) ? "!" : " ");
2190 // ----- CFDAT output -----
2191 else if(os.iword(AliTRDmcmSim::fgkFormatIndex) == 1) {
2193 Int_t addrOffset = 0x2000;
2194 Int_t addrStep = 0x80;
2196 for (Int_t iTimeBin = 0; iTimeBin < mcm.fNTimeBin; iTimeBin++) {
2197 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++) {
2198 os << std::setw(5) << 10
2199 << std::setw(5) << addrOffset + iChannel * addrStep + iTimeBin
2200 << std::setw(5) << (mcm.fADCF[iChannel][iTimeBin])
2201 << std::setw(5) << dest << std::endl;
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];
2212 Int_t bufLength = mcm.ProduceRawStream(&buf[0], bufSize);
2214 for (Int_t i = 0; i < bufLength; i++)
2215 std::cout << "0x" << std::hex << buf[i] << std::dec << std::endl;
2221 os << "unknown format set" << std::endl;
2228 void AliTRDmcmSim::PrintFitRegXml(ostream& os) const
2230 // print fit registres in XML format
2232 bool tracklet=false;
2234 for (Int_t cpu = 0; cpu < 4; cpu++) {
2235 if(fFitPtr[cpu] != 31)
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;
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++) {
2251 os << " <ch chnr=\"" << adcch << "\">"<< std::endl;
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
2254 os << " <q1>" << fFitReg[adcch].fQ1/4 << "</q1>"<< std::endl; // in the output
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;
2263 os << " </c>" << std::endl;
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;
2275 void AliTRDmcmSim::PrintTrackletsXml(ostream& os) const
2277 // print tracklets in XML format
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;
2286 Int_t pid, padrow, slope, offset;
2287 for(Int_t cpu=0; cpu<4; cpu++) {
2288 if(fMCMT[cpu] == 0x10001000) {
2295 pid = (fMCMT[cpu] & 0xFF000000) >> 24;
2296 padrow = (fMCMT[cpu] & 0xF00000 ) >> 20;
2297 slope = (fMCMT[cpu] & 0xFE000 ) >> 13;
2298 offset = (fMCMT[cpu] & 0x1FFF ) ;
2301 os << " <trk> <pid>" << pid << "</pid>" << " <padrow>" << padrow << "</padrow>"
2302 << " <slope>" << slope << "</slope>" << " <offset>" << offset << "</offset>" << "</trk>" << std::endl;
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;
2314 void AliTRDmcmSim::PrintAdcDatHuman(ostream& os) const
2316 // print ADC data in human-readable format
2318 os << "MCM " << fMcmPos << " on ROB " << fRobPos <<
2319 " in detector " << fDetector << std::endl;
2321 os << "----- Unfiltered ADC data (10 bit) -----" << std::endl;
2323 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++)
2324 os << std::setw(5) << iChannel;
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);
2334 os << "----- Filtered ADC data (10+2 bit) -----" << std::endl;
2336 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++)
2337 os << std::setw(4) << iChannel
2338 << ((~fZSMap[iChannel] != 0) ? "!" : " ");
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) ? "!" : " ");
2351 void AliTRDmcmSim::PrintAdcDatXml(ostream& os) const
2353 // print ADC data in XML format
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;
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>";
2367 os << " </ch>" << std::endl;
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;
2380 void AliTRDmcmSim::PrintAdcDatDatx(ostream& os, Bool_t broadcast, Int_t timeBinOffset) const
2382 // print ADC data in datx format (to send to FEE)
2384 fTrapConfig->PrintDatx(os, 2602, 1, 0, 127); // command to enable the ADC clock - necessary to write ADC values to MCM
2387 Int_t addrOffset = 0x2000;
2388 Int_t addrStep = 0x80;
2389 Int_t addrOffsetEBSIA = 0x20;
2391 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
2392 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2393 if ((iTimeBin < timeBinOffset) || (iTimeBin >= fNTimeBin+timeBinOffset)) {
2394 if(broadcast==kFALSE)
2395 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, 10, GetRobPos(), GetMcmPos());
2397 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, 10, 0, 127);
2400 if(broadcast==kFALSE)
2401 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, (fADCF[iChannel][iTimeBin-timeBinOffset]/4), GetRobPos(), GetMcmPos());
2403 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, (fADCF[iChannel][iTimeBin-timeBinOffset]/4), 0, 127);
2411 void AliTRDmcmSim::PrintPidLutHuman()
2413 // print PID LUT in human readable format
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);
2420 std::cout << "nBinsQ0: " << nBinsQ0 << std::endl;
2421 std::cout << "LUT table length: " << fTrapConfig->GetDmemUnsigned(fgkDmemAddrLUTLength, fDetector, fRobPos, fMcmPos) << std::endl;
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)
2427 << " # " <<((result>>0)&0xFF)
2428 << " | " << ((result>>8)&0xFF)
2429 << " | " << ((result>>16)&0xFF)
2430 << " | " << ((result>>24)&0xFF) << std::endl;
2436 Bool_t AliTRDmcmSim::ReadPackedConfig(AliTRDtrapConfig *cfg, Int_t hc, UInt_t *data, Int_t size)
2438 // Read the packed configuration from the passed memory block
2440 // To be used to retrieve the TRAP configuration from the
2441 // configuration as sent in the raw data.
2443 AliDebugClass(1, "Reading packed configuration");
2449 Int_t step, bwidth, nwords, exitFlag, bitcnt;
2452 UInt_t dat, msk, header, dataHi;
2454 while (idx < size && *data != 0x00000000) {
2456 Int_t rob = (*data >> 28) & 0x7;
2457 Int_t mcm = (*data >> 24) & 0xf;
2459 AliDebugClass(1, Form("Config of det. %3i MCM %i:%02i (0x%08x)", det, rob, mcm, *data));
2462 while (idx < size && *data != 0x00000000) {
2468 AliDebugClass(5, Form("read: 0x%08x", header));
2470 if (header & 0x01) // single data
2472 dat = (header >> 2) & 0xFFFF; // 16 bit data
2473 caddr = (header >> 18) & 0x3FFF; // 14 bit address
2475 if (caddr != 0x1FFF) // temp!!! because the end marker was wrong
2477 if (header & 0x02) // check if > 16 bits
2480 AliDebugClass(5, Form("read: 0x%08x", dataHi));
2483 err += ((dataHi ^ (dat | 1)) & 0xFFFF) != 0;
2484 dat = (dataHi & 0xFFFF0000) | dat;
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));
2491 AliDebugClass(5, Form("(single-write): no more data, missing end marker\n"));
2497 AliDebugClass(5, Form("(single-write): address 0x%04x => old endmarker?\n", caddr));
2502 else // block of data
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);
2521 msk = (1 << bwidth) - 1;
2530 AliDebugClass(5, Form("read 0x%08x", header));
2533 err += (header & 1);
2534 header = header >> 1;
2535 bitcnt = 31 - bwidth;
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));
2542 header = header >> bwidth;
2545 AliDebugClass(5, Form("(block-write): no end marker! %d words read\n", idx));
2556 AliDebugClass(5, Form("read 0x%08x", header));
2560 err += (header & 1);
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));
2569 AliDebugClass(5, Form("no end marker! %d words read", idx));
2575 default: return err;
2580 AliDebugClass(5, Form("no end marker! %d words read", idx));
2581 return -err; // only if the max length of the block reached!