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"
40 #include "AliRunLoader.h"
41 #include "AliLoader.h"
43 #include "AliTRDfeeParam.h"
44 #include "AliTRDtrapConfig.h"
45 #include "AliTRDdigitsManager.h"
46 #include "AliTRDarrayADC.h"
47 #include "AliTRDarrayDictionary.h"
48 #include "AliTRDtrackletMCM.h"
49 #include "AliTRDmcmSim.h"
51 ClassImp(AliTRDmcmSim)
53 Bool_t AliTRDmcmSim::fgApplyCut = kTRUE;
54 Int_t AliTRDmcmSim::fgAddBaseline = 0;
56 const Int_t AliTRDmcmSim::fgkFormatIndex = std::ios_base::xalloc();
58 const Int_t AliTRDmcmSim::fgkNADC = AliTRDfeeParam::GetNadcMcm();
59 const UShort_t AliTRDmcmSim::fgkFPshifts[4] = {11, 14, 17, 21};
62 AliTRDmcmSim::AliTRDmcmSim() :
75 fTrklBranchName("mcmtrklbranch"),
88 // AliTRDmcmSim default constructor
89 // By default, nothing is initialized.
90 // It is necessary to issue Init before use.
92 for (Int_t iDict = 0; iDict < 3; iDict++)
101 AliTRDmcmSim::~AliTRDmcmSim()
104 // AliTRDmcmSim destructor
108 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
109 delete [] fADCR[iAdc];
110 delete [] fADCF[iAdc];
118 delete [] fGainCounterA;
119 delete [] fGainCounterB;
120 delete [] fTailAmplLong;
121 delete [] fTailAmplShort;
124 fTrackletArray->Delete();
125 delete fTrackletArray;
129 void AliTRDmcmSim::Init( Int_t det, Int_t robPos, Int_t mcmPos, Bool_t /* newEvent */ )
132 // Initialize the class with new MCM position information
133 // memory is allocated in the first initialization
137 fFeeParam = AliTRDfeeParam::Instance();
138 fTrapConfig = AliTRDtrapConfig::Instance();
144 fNTimeBin = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC13CPUA);
145 fRow = fFeeParam->GetPadRowFromMCM( fRobPos, fMcmPos );
148 fADCR = new Int_t *[fgkNADC];
149 fADCF = new Int_t *[fgkNADC];
150 fZSMap = new Int_t [fgkNADC];
151 fGainCounterA = new UInt_t[fgkNADC];
152 fGainCounterB = new UInt_t[fgkNADC];
153 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
154 fADCR[iAdc] = new Int_t[fNTimeBin];
155 fADCF[iAdc] = new Int_t[fNTimeBin];
159 fPedAcc = new UInt_t[fgkNADC]; // accumulator for pedestal filter
160 fTailAmplLong = new UShort_t[fgkNADC];
161 fTailAmplShort = new UShort_t[fgkNADC];
163 // tracklet calculation
164 fFitReg = new FitReg_t[fgkNADC];
165 fTrackletArray = new TClonesArray("AliTRDtrackletMCM", fgkMaxTracklets);
167 fMCMT = new UInt_t[fgkMaxTracklets];
170 fInitialized = kTRUE;
175 void AliTRDmcmSim::Reset()
177 // Resets the data values and internal filter registers
178 // by re-initialising them
180 if( !CheckInitialized() )
183 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
184 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
188 fZSMap[iAdc] = -1; // Default unread, low active bit mask
189 fGainCounterA[iAdc] = 0;
190 fGainCounterB[iAdc] = 0;
193 for(Int_t i = 0; i < fgkMaxTracklets; i++) {
197 for (Int_t iDict = 0; iDict < 3; iDict++)
200 FilterPedestalInit();
205 void AliTRDmcmSim::SetNTimebins(Int_t ntimebins)
207 // Reallocate memory if a change in the number of timebins
208 // is needed (should not be the case for real data)
210 if( !CheckInitialized() )
213 fNTimeBin = ntimebins;
214 for( Int_t iAdc = 0 ; iAdc < fgkNADC; iAdc++ ) {
217 fADCR[iAdc] = new Int_t[fNTimeBin];
218 fADCF[iAdc] = new Int_t[fNTimeBin];
222 Bool_t AliTRDmcmSim::LoadMCM(AliRunLoader* const runloader, Int_t det, Int_t rob, Int_t mcm)
224 // loads the ADC data as obtained from the digitsManager for the specified MCM.
225 // This method is meant for rare execution, e.g. in the visualization. When called
226 // frequently use SetData(...) instead.
231 AliError("No Runloader given");
235 AliLoader *trdLoader = runloader->GetLoader("TRDLoader");
237 AliError("Could not get TRDLoader");
241 Bool_t retval = kTRUE;
242 trdLoader->LoadDigits();
243 fDigitsManager = 0x0;
244 AliTRDdigitsManager *digMgr = new AliTRDdigitsManager();
245 digMgr->SetSDigits(0);
246 digMgr->CreateArrays();
247 digMgr->ReadDigits(trdLoader->TreeD());
248 AliTRDarrayADC *digits = (AliTRDarrayADC*) digMgr->GetDigits(det);
249 if (digits->HasData()) {
252 if (fNTimeBin != digits->GetNtime()) {
253 AliWarning(Form("Changing no. of timebins from %i to %i", fNTimeBin, digits->GetNtime()));
254 SetNTimebins(digits->GetNtime());
267 void AliTRDmcmSim::NoiseTest(Int_t nsamples, Int_t mean, Int_t sigma, Int_t inputGain, Int_t inputTail)
269 // This function can be used to test the filters.
270 // It feeds nsamples of ADC values with a gaussian distribution specified by mean and sigma.
271 // The filter chain implemented here consists of:
272 // Pedestal -> Gain -> Tail
273 // With inputGain and inputTail the input to the gain and tail filter, respectively,
274 // can be chosen where
276 // 1: pedestal output
278 // The input has to be chosen from a stage before.
279 // The filter behaviour is controlled by the TRAP parameters from AliTRDtrapConfig in the
280 // same way as in normal simulation.
281 // The functions produces four histograms with the values at the different stages.
283 if( !CheckInitialized() )
286 TString nameInputGain;
287 TString nameInputTail;
291 nameInputGain = "Noise";
295 nameInputGain = "Pedestal";
299 AliError("Undefined input to tail cancellation filter");
305 nameInputTail = "Noise";
309 nameInputTail = "Pedestal";
313 nameInputTail = "Gain";
317 AliError("Undefined input to tail cancellation filter");
321 TH1F *h = new TH1F("noise", "Gaussian Noise;sample;ADC count",
322 nsamples, 0, nsamples);
323 TH1F *hfp = new TH1F("ped", "Noise #rightarrow Pedestal filter;sample;ADC count", nsamples, 0, nsamples);
324 TH1F *hfg = new TH1F("gain",
325 (nameInputGain + "#rightarrow Gain;sample;ADC count").Data(),
326 nsamples, 0, nsamples);
327 TH1F *hft = new TH1F("tail",
328 (nameInputTail + "#rightarrow Tail;sample;ADC count").Data(),
329 nsamples, 0, nsamples);
331 hfp->SetStats(kFALSE);
332 hfg->SetStats(kFALSE);
333 hft->SetStats(kFALSE);
335 Int_t value; // ADC count with noise (10 bit)
336 Int_t valuep; // pedestal filter output (12 bit)
337 Int_t valueg; // gain filter output (12 bit)
338 Int_t valuet; // tail filter value (12 bit)
340 for (Int_t i = 0; i < nsamples; i++) {
341 value = (Int_t) gRandom->Gaus(mean, sigma); // generate noise with gaussian distribution
342 h->SetBinContent(i, value);
344 valuep = FilterPedestalNextSample(1, 0, ((Int_t) value) << 2);
347 valueg = FilterGainNextSample(1, ((Int_t) value) << 2);
349 valueg = FilterGainNextSample(1, valuep);
352 valuet = FilterTailNextSample(1, ((Int_t) value) << 2);
353 else if (inputTail == 1)
354 valuet = FilterTailNextSample(1, valuep);
356 valuet = FilterTailNextSample(1, valueg);
358 hfp->SetBinContent(i, valuep >> 2);
359 hfg->SetBinContent(i, valueg >> 2);
360 hft->SetBinContent(i, valuet >> 2);
363 TCanvas *c = new TCanvas;
375 Bool_t AliTRDmcmSim::CheckInitialized() const
378 // Check whether object is initialized
382 AliError(Form ("AliTRDmcmSim is not initialized but function other than Init() is called."));
387 void AliTRDmcmSim::Print(Option_t* const option) const
389 // Prints the data stored and/or calculated for this MCM.
390 // The output is controlled by option which can be a sequence of any of
391 // the following characters:
392 // R - prints raw ADC data
393 // F - prints filtered data
394 // H - prints detected hits
395 // T - prints found tracklets
396 // The later stages are only meaningful after the corresponding calculations
397 // have been performed.
399 if ( !CheckInitialized() )
402 printf("MCM %i on ROB %i in detector %i\n", fMcmPos, fRobPos, fDetector);
404 TString opt = option;
405 if (opt.Contains("R") || opt.Contains("F")) {
409 if (opt.Contains("H")) {
410 printf("Found %i hits:\n", fNHits);
411 for (Int_t iHit = 0; iHit < fNHits; iHit++) {
412 printf("Hit %3i in timebin %2i, ADC %2i has charge %3i and position %3i\n",
413 iHit, fHits[iHit].fTimebin, fHits[iHit].fChannel, fHits[iHit].fQtot, fHits[iHit].fYpos);
417 if (opt.Contains("T")) {
418 printf("Tracklets:\n");
419 for (Int_t iTrkl = 0; iTrkl < fTrackletArray->GetEntriesFast(); iTrkl++) {
420 printf("tracklet %i: 0x%08x\n", iTrkl, ((AliTRDtrackletMCM*) (*fTrackletArray)[iTrkl])->GetTrackletWord());
425 void AliTRDmcmSim::Draw(Option_t* const option)
427 // Plots the data stored in a 2-dim. timebin vs. ADC channel plot.
428 // The option selects what data is plotted and can be a sequence of
429 // the following characters:
430 // R - plot raw data (default)
431 // F - plot filtered data (meaningless if R is specified)
432 // In addition to the ADC values:
434 // T - plot tracklets
436 if( !CheckInitialized() )
439 TString opt = option;
441 TH2F *hist = new TH2F("mcmdata", Form("Data of MCM %i on ROB %i in detector %i", \
442 fMcmPos, fRobPos, fDetector), \
443 fgkNADC, -0.5, fgkNADC-.5, fNTimeBin, -.5, fNTimeBin-.5);
444 hist->GetXaxis()->SetTitle("ADC Channel");
445 hist->GetYaxis()->SetTitle("Timebin");
446 hist->SetStats(kFALSE);
448 if (opt.Contains("R")) {
449 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
450 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
451 hist->SetBinContent(iAdc+1, iTimeBin+1, fADCR[iAdc][iTimeBin] >> fgkAddDigits);
456 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
457 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
458 hist->SetBinContent(iAdc+1, iTimeBin+1, fADCF[iAdc][iTimeBin] >> fgkAddDigits);
464 if (opt.Contains("H")) {
465 TGraph *grHits = new TGraph();
466 for (Int_t iHit = 0; iHit < fNHits; iHit++) {
467 grHits->SetPoint(iHit,
468 fHits[iHit].fChannel + 1 + fHits[iHit].fYpos/256.,
469 fHits[iHit].fTimebin);
474 if (opt.Contains("T")) {
475 TLine *trklLines = new TLine[4];
476 for (Int_t iTrkl = 0; iTrkl < fTrackletArray->GetEntries(); iTrkl++) {
477 AliTRDtrackletMCM *trkl = (AliTRDtrackletMCM*) (*fTrackletArray)[iTrkl];
478 Float_t padWidth = 0.635 + 0.03 * (fDetector % 6);
479 Float_t offset = padWidth/256. * ((((((fRobPos & 0x1) << 2) + (fMcmPos & 0x3)) * 18) << 8) - ((18*4*2 - 18*2 - 3) << 7)); // revert adding offset in FitTracklet
480 Int_t ndrift = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrNdrift, fDetector, fRobPos, fMcmPos) >> 5;
481 Float_t slope = trkl->GetdY() * 140e-4 / ndrift;
483 Int_t t0 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS);
484 Int_t t1 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFE);
486 trklLines[iTrkl].SetX1((offset - (trkl->GetY() - slope * t0)) / padWidth); // ??? sign?
487 trklLines[iTrkl].SetY1(t0);
488 trklLines[iTrkl].SetX2((offset - (trkl->GetY() - slope * t1)) / padWidth); // ??? sign?
489 trklLines[iTrkl].SetY2(t1);
490 trklLines[iTrkl].SetLineColor(2);
491 trklLines[iTrkl].SetLineWidth(2);
492 printf("Tracklet %i: y = %f, dy = %f, offset = %f\n", iTrkl, trkl->GetY(), (trkl->GetdY() * 140e-4), offset);
493 trklLines[iTrkl].Draw();
498 void AliTRDmcmSim::SetData( Int_t adc, Int_t* const data )
501 // Store ADC data into array of raw data
504 if( !CheckInitialized() ) return;
506 if( adc < 0 || adc >= fgkNADC ) {
507 AliError(Form ("Error: ADC %i is out of range (0 .. %d).", adc, fgkNADC-1));
511 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
512 fADCR[adc][it] = (Int_t) (data[it]) << fgkAddDigits;
513 fADCF[adc][it] = (Int_t) (data[it]) << fgkAddDigits;
517 void AliTRDmcmSim::SetData( Int_t adc, Int_t it, Int_t data )
520 // Store ADC data into array of raw data
523 if( !CheckInitialized() ) return;
525 if( adc < 0 || adc >= fgkNADC ) {
526 AliError(Form ("Error: ADC %i is out of range (0 .. %d).", adc, fgkNADC-1));
530 fADCR[adc][it] = data << fgkAddDigits;
531 fADCF[adc][it] = data << fgkAddDigits;
534 void AliTRDmcmSim::SetData(AliTRDarrayADC* const adcArray, AliTRDdigitsManager * const digitsManager)
536 // Set the ADC data from an AliTRDarrayADC
538 if( !CheckInitialized() )
541 fDigitsManager = digitsManager;
542 if (fDigitsManager) {
543 for (Int_t iDict = 0; iDict < 3; iDict++) {
544 AliTRDarrayDictionary *newDict = (AliTRDarrayDictionary*) fDigitsManager->GetDictionary(fDetector, iDict);
545 if (fDict[iDict] != 0x0 && newDict != 0x0) {
547 if (fDict[iDict] == newDict)
550 fDict[iDict] = newDict;
552 if (fDict[iDict]->GetDim() == 0) {
553 AliError(Form("Dictionary %i of det. %i has dim. 0", fDetector, iDict));
556 fDict[iDict]->Expand();
559 fDict[iDict] = newDict;
561 fDict[iDict]->Expand();
566 if (fNTimeBin != adcArray->GetNtime())
567 SetNTimebins(adcArray->GetNtime());
569 Int_t offset = (fMcmPos % 4 + 1) * 21 + (fRobPos % 2) * 84 - 1;
571 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
572 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
573 Int_t value = adcArray->GetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin);
574 if (value < 0 || (offset - iAdc < 1) || (offset - iAdc > 165)) {
575 fADCR[iAdc][iTimeBin] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP) + (fgAddBaseline << fgkAddDigits);
576 fADCF[iAdc][iTimeBin] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP) + (fgAddBaseline << fgkAddDigits);
580 fADCR[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
581 fADCF[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
587 void AliTRDmcmSim::SetDataByPad(AliTRDarrayADC* const adcArray, AliTRDdigitsManager * const digitsManager)
589 // Set the ADC data from an AliTRDarrayADC
590 // (by pad, to be used during initial reading in simulation)
592 if( !CheckInitialized() )
595 fDigitsManager = digitsManager;
596 if (fDigitsManager) {
597 for (Int_t iDict = 0; iDict < 3; iDict++) {
598 AliTRDarrayDictionary *newDict = (AliTRDarrayDictionary*) fDigitsManager->GetDictionary(fDetector, iDict);
599 if (fDict[iDict] != 0x0 && newDict != 0x0) {
601 if (fDict[iDict] == newDict)
604 fDict[iDict] = newDict;
606 if (fDict[iDict]->GetDim() == 0) {
607 AliError(Form("Dictionary %i of det. %i has dim. 0", fDetector, iDict));
610 fDict[iDict]->Expand();
613 fDict[iDict] = newDict;
615 fDict[iDict]->Expand();
620 if (fNTimeBin != adcArray->GetNtime())
621 SetNTimebins(adcArray->GetNtime());
623 Int_t offset = (fMcmPos % 4 + 1) * 18 + (fRobPos % 2) * 72 + 1;
625 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
626 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
628 Int_t pad = offset - iAdc;
629 if (pad > -1 && pad < 144)
630 value = adcArray->GetData(GetRow(), offset - iAdc, iTimeBin);
631 // Int_t value = adcArray->GetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin);
632 if (value < 0 || (offset - iAdc < 1) || (offset - iAdc > 165)) {
633 fADCR[iAdc][iTimeBin] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP) + (fgAddBaseline << fgkAddDigits);
634 fADCF[iAdc][iTimeBin] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP) + (fgAddBaseline << fgkAddDigits);
638 fADCR[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
639 fADCF[iAdc][iTimeBin] = (value << fgkAddDigits) + (fgAddBaseline << fgkAddDigits);
645 void AliTRDmcmSim::SetDataPedestal( Int_t adc )
648 // Store ADC data into array of raw data
651 if( !CheckInitialized() )
654 if( adc < 0 || adc >= fgkNADC ) {
658 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
659 fADCR[adc][it] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP) + (fgAddBaseline << fgkAddDigits);
660 fADCF[adc][it] = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP) + (fgAddBaseline << fgkAddDigits);
664 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
666 // retrieve the MC hit information (not available in TRAP hardware)
668 if (index < 0 || index >= fNHits)
671 channel = fHits[index].fChannel;
672 timebin = fHits[index].fTimebin;
673 qtot = fHits[index].fQtot;
674 ypos = fHits[index].fYpos;
675 y = (Float_t) ((((((fRobPos & 0x1) << 2) + (fMcmPos & 0x3)) * 18) << 8) - ((18*4*2 - 18*2 - 1) << 7) -
676 (channel << 8) - ypos)
677 * (0.635 + 0.03 * (fDetector % 6))
679 label = fHits[index].fLabel[0];
684 Int_t AliTRDmcmSim::GetCol( Int_t adc )
687 // Return column id of the pad for the given ADC channel
690 if( !CheckInitialized() )
693 Int_t col = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adc);
694 if (col < 0 || col >= fFeeParam->GetNcol())
700 Int_t AliTRDmcmSim::ProduceRawStream( UInt_t *buf, Int_t bufSize, UInt_t iEv) const
703 // Produce raw data stream from this MCM and put in buf
704 // Returns number of words filled, or negative value
705 // with -1 * number of overflowed words
708 if( !CheckInitialized() )
712 UInt_t mcmHeader = 0;
714 Int_t nw = 0; // Number of written words
715 Int_t of = 0; // Number of overflowed words
716 Int_t rawVer = fFeeParam->GetRAWversion();
718 Int_t nActiveADC = 0; // number of activated ADC bits in a word
720 if( !CheckInitialized() )
723 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBSF) != 0) // store unfiltered data
728 // Produce ADC mask : nncc cccm mmmm mmmm mmmm mmmm mmmm 1100
729 // n : unused , c : ADC count, m : selected ADCs
731 (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC15CPUA) & (1 << 13))) { // check for zs flag in TRAP configuration
732 for( Int_t iAdc = 0 ; iAdc < fgkNADC ; iAdc++ ) {
733 if( ~fZSMap[iAdc] != 0 ) { // 0 means not suppressed
734 adcMask |= (1 << (iAdc+4) ); // last 4 digit reserved for 1100=0xc
735 nActiveADC++; // number of 1 in mmm....m
739 if ((nActiveADC == 0) &&
740 (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kC15CPUA) & (1 << 8))) // check for DEH flag in TRAP configuration
743 // assemble adc mask word
744 adcMask |= (1 << 30) | ( ( 0x3FFFFFFC ) & (~(nActiveADC) << 25) ) | 0xC; // nn = 01, ccccc are inverted, 0xc=1100
748 mcmHeader = (1<<31) | (fRobPos << 28) | (fMcmPos << 24) | ((iEv % 0x100000) << 4) | 0xC;
750 buf[nw++] = mcmHeader;
762 // Produce ADC data. 3 timebins are packed into one 32 bits word
763 // In this version, different ADC channel will NOT share the same word
765 UInt_t aa=0, a1=0, a2=0, a3=0;
767 for (Int_t iAdc = 0; iAdc < 21; iAdc++ ) {
768 if( rawVer>= 3 && ~fZSMap[iAdc] == 0 ) continue; // Zero Suppression, 0 means not suppressed
769 aa = !(iAdc & 1) + 2;
770 for (Int_t iT = 0; iT < fNTimeBin; iT+=3 ) {
771 a1 = ((iT ) < fNTimeBin ) ? adc[iAdc][iT ] >> fgkAddDigits : 0;
772 a2 = ((iT + 1) < fNTimeBin ) ? adc[iAdc][iT+1] >> fgkAddDigits : 0;
773 a3 = ((iT + 2) < fNTimeBin ) ? adc[iAdc][iT+2] >> fgkAddDigits : 0;
774 x = (a3 << 22) | (a2 << 12) | (a1 << 2) | aa;
784 if( of != 0 ) return -of; else return nw;
787 Int_t AliTRDmcmSim::ProduceTrackletStream( UInt_t *buf, Int_t bufSize )
790 // Produce tracklet data stream from this MCM and put in buf
791 // Returns number of words filled, or negative value
792 // with -1 * number of overflowed words
795 if( !CheckInitialized() )
798 Int_t nw = 0; // Number of written words
799 Int_t of = 0; // Number of overflowed words
801 // Produce tracklet data. A maximum of four 32 Bit words will be written per MCM
802 // fMCMT is filled continuously until no more tracklet words available
804 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
806 buf[nw++] = ((AliTRDtrackletMCM*) (*fTrackletArray)[iTracklet])->GetTrackletWord();
811 if( of != 0 ) return -of; else return nw;
814 void AliTRDmcmSim::Filter()
817 // Filter the raw ADC values. The active filter stages and their
818 // parameters are taken from AliTRDtrapConfig.
819 // The raw data is stored separate from the filtered data. Thus,
820 // it is possible to run the filters on a set of raw values
821 // sequentially for parameter tuning.
824 if( !CheckInitialized() )
827 // Apply filters sequentially. Bypass is handled by filters
828 // since counters and internal registers may be updated even
829 // if the filter is bypassed.
830 // The first filter takes the data from fADCR and
833 // Non-linearity filter not implemented.
837 // Crosstalk filter not implemented.
840 void AliTRDmcmSim::FilterPedestalInit(Int_t baseline)
842 // Initializes the pedestal filter assuming that the input has
843 // been constant for a long time (compared to the time constant).
845 UShort_t fptc = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPTC); // 0..3, 0 - fastest, 3 - slowest
847 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++)
848 fPedAcc[iAdc] = (baseline << 2) * (1 << fgkFPshifts[fptc]);
851 UShort_t AliTRDmcmSim::FilterPedestalNextSample(Int_t adc, Int_t timebin, UShort_t value)
853 // Returns the output of the pedestal filter given the input value.
854 // The output depends on the internal registers and, thus, the
855 // history of the filter.
857 UShort_t fpnp = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP); // 0..511 -> 0..127.75, pedestal at the output
858 UShort_t fptc = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPTC); // 0..3, 0 - fastest, 3 - slowest
859 UShort_t fpby = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPBY); // 0..1 bypass, active low
861 UShort_t accumulatorShifted;
865 inpAdd = value + fpnp;
867 accumulatorShifted = (fPedAcc[adc] >> fgkFPshifts[fptc]) & 0x3FF; // 10 bits
868 if (timebin == 0) // the accumulator is disabled in the drift time
870 correction = (value & 0x3FF) - accumulatorShifted;
871 fPedAcc[adc] = (fPedAcc[adc] + correction) & 0x7FFFFFFF; // 31 bits
877 if (inpAdd <= accumulatorShifted)
881 inpAdd = inpAdd - accumulatorShifted;
889 void AliTRDmcmSim::FilterPedestal()
892 // Apply pedestal filter
894 // As the first filter in the chain it reads data from fADCR
895 // and outputs to fADCF.
896 // It has only an effect if previous samples have been fed to
897 // find the pedestal. Currently, the simulation assumes that
898 // the input has been stable for a sufficiently long time.
900 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
901 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
902 fADCF[iAdc][iTimeBin] = FilterPedestalNextSample(iAdc, iTimeBin, fADCR[iAdc][iTimeBin]);
907 void AliTRDmcmSim::FilterGainInit()
909 // Initializes the gain filter. In this case, only threshold
910 // counters are reset.
912 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
913 // these are counters which in hardware continue
914 // until maximum or reset
915 fGainCounterA[iAdc] = 0;
916 fGainCounterB[iAdc] = 0;
920 UShort_t AliTRDmcmSim::FilterGainNextSample(Int_t adc, UShort_t value)
922 // Apply the gain filter to the given value.
923 // BEGIN_LATEX O_{i}(t) = #gamma_{i} * I_{i}(t) + a_{i} END_LATEX
924 // The output depends on the internal registers and, thus, the
925 // history of the filter.
927 UShort_t fgby = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFGBY); // bypass, active low
928 UShort_t fgf = fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGF0 + adc)); // 0x700 + (0 & 0x1ff);
929 UShort_t fga = fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGA0 + adc)); // 40;
930 UShort_t fgta = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFGTA); // 20;
931 UShort_t fgtb = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFGTB); // 2060;
933 UInt_t corr; // corrected value
936 corr = (value * fgf) >> 11;
937 corr = corr > 0xfff ? 0xfff : corr;
938 corr = AddUintClipping(corr, fga, 12);
940 // Update threshold counters
941 // not really useful as they are cleared with every new event
942 if (!((fGainCounterA[adc] == 0x3FFFFFF) || (fGainCounterB[adc] == 0x3FFFFFF)))
946 fGainCounterB[adc]++;
947 else if (corr >= fgta)
948 fGainCounterA[adc]++;
957 void AliTRDmcmSim::FilterGain()
959 // Read data from fADCF and apply gain filter.
961 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
962 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
963 fADCF[iAdc][iTimeBin] = FilterGainNextSample(iAdc, fADCF[iAdc][iTimeBin]);
968 void AliTRDmcmSim::FilterTailInit(Int_t baseline)
970 // Initializes the tail filter assuming that the input has
971 // been at the baseline value (configured by FTFP) for a
972 // sufficiently long time.
974 // exponents and weight calculated from configuration
975 UShort_t alphaLong = 0x3ff & fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTAL); // the weight of the long component
976 UShort_t lambdaLong = (1 << 10) | (1 << 9) | (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTLL) & 0x1FF); // the multiplier
977 UShort_t lambdaShort = (0 << 10) | (1 << 9) | (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTLS) & 0x1FF); // the multiplier
979 Float_t lambdaL = lambdaLong * 1.0 / (1 << 11);
980 Float_t lambdaS = lambdaShort * 1.0 / (1 << 11);
981 Float_t alphaL = alphaLong * 1.0 / (1 << 11);
983 qup = (1 - lambdaL) * (1 - lambdaS);
984 qdn = 1 - lambdaS * alphaL - lambdaL * (1 - alphaL);
985 Float_t kdc = qup/qdn;
991 baseline = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFPNP);
993 ql = lambdaL * (1 - lambdaS) * alphaL;
994 qs = lambdaS * (1 - lambdaL) * (1 - alphaL);
996 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
997 Int_t value = baseline & 0xFFF;
998 Int_t corr = (value * fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGF0 + iAdc))) >> 11;
999 corr = corr > 0xfff ? 0xfff : corr;
1000 corr = AddUintClipping(corr, fTrapConfig->GetTrapReg(AliTRDtrapConfig::TrapReg_t(AliTRDtrapConfig::kFGA0 + iAdc)), 12);
1002 kt = kdc * baseline;
1003 aout = baseline - (UShort_t) kt;
1005 fTailAmplLong[iAdc] = (UShort_t) (aout * ql / (ql + qs));
1006 fTailAmplShort[iAdc] = (UShort_t) (aout * qs / (ql + qs));
1010 UShort_t AliTRDmcmSim::FilterTailNextSample(Int_t adc, UShort_t value)
1012 // Returns the output of the tail filter for the given input value.
1013 // The output depends on the internal registers and, thus, the
1014 // history of the filter.
1016 // exponents and weight calculated from configuration
1017 UShort_t alphaLong = 0x3ff & fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTAL); // the weight of the long component
1018 UShort_t lambdaLong = (1 << 10) | (1 << 9) | (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTLL) & 0x1FF); // the multiplier of the long component
1019 UShort_t lambdaShort = (0 << 10) | (1 << 9) | (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTLS) & 0x1FF); // the multiplier of the short component
1021 // intermediate signals
1027 UShort_t inpVolt = value & 0xFFF; // 12 bits
1029 // add the present generator outputs
1030 aQ = AddUintClipping(fTailAmplLong[adc], fTailAmplShort[adc], 12);
1032 // calculate the difference between the input and the generated signal
1034 aDiff = inpVolt - aQ;
1038 // the inputs to the two generators, weighted
1039 alInpv = (aDiff * alphaLong) >> 11;
1041 // the new values of the registers, used next time
1043 tmp = AddUintClipping(fTailAmplLong[adc], alInpv, 12);
1044 tmp = (tmp * lambdaLong) >> 11;
1045 fTailAmplLong[adc] = tmp & 0xFFF;
1047 tmp = AddUintClipping(fTailAmplShort[adc], aDiff - alInpv, 12);
1048 tmp = (tmp * lambdaShort) >> 11;
1049 fTailAmplShort[adc] = tmp & 0xFFF;
1051 // the output of the filter
1052 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kFTBY) == 0) // bypass mode, active low
1058 void AliTRDmcmSim::FilterTail()
1060 // Apply tail cancellation filter to all data.
1062 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1063 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1064 fADCF[iAdc][iTimeBin] = FilterTailNextSample(iAdc, fADCF[iAdc][iTimeBin]);
1069 void AliTRDmcmSim::ZSMapping()
1072 // Zero Suppression Mapping implemented in TRAP chip
1073 // only implemented for up to 30 timebins
1075 // See detail TRAP manual "Data Indication" section:
1076 // http://www.kip.uni-heidelberg.de/ti/TRD/doc/trap/TRAP-UserManual.pdf
1079 if( !CheckInitialized() )
1082 Int_t eBIS = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBIS);
1083 Int_t eBIT = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBIT);
1084 Int_t eBIL = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBIL);
1085 Int_t eBIN = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBIN);
1087 Int_t **adc = fADCF;
1089 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++)
1092 for( Int_t it = 0 ; it < fNTimeBin ; it++ ) {
1093 Int_t iAdc; // current ADC channel
1098 Int_t supp; // suppression of the current channel (low active)
1100 // ----- first channel -----
1104 ac = adc[iAdc ][it]; // current
1105 an = adc[iAdc+1][it]; // next
1107 mask = ( ac >= ap && ac >= an ) ? 0 : 0x1; // peak center detection
1108 mask += ( ap + ac + an > eBIT ) ? 0 : 0x2; // cluster
1109 mask += ( ac > eBIS ) ? 0 : 0x4; // absolute large peak
1111 supp = (eBIL >> mask) & 1;
1113 fZSMap[iAdc] &= ~((1-supp) << it);
1114 if( eBIN == 0 ) { // neighbour sensitivity
1115 fZSMap[iAdc+1] &= ~((1-supp) << it);
1118 // ----- last channel -----
1121 ap = adc[iAdc-1][it]; // previous
1122 ac = adc[iAdc ][it]; // current
1125 mask = ( ac >= ap && ac >= an ) ? 0 : 0x1; // peak center detection
1126 mask += ( ap + ac + an > eBIT ) ? 0 : 0x2; // cluster
1127 mask += ( ac > eBIS ) ? 0 : 0x4; // absolute large peak
1129 supp = (eBIL >> mask) & 1;
1131 fZSMap[iAdc] &= ~((1-supp) << it);
1132 if( eBIN == 0 ) { // neighbour sensitivity
1133 fZSMap[iAdc-1] &= ~((1-supp) << it);
1136 // ----- middle channels -----
1137 for( iAdc = 1 ; iAdc < fgkNADC-1; iAdc++ ) {
1138 ap = adc[iAdc-1][it]; // previous
1139 ac = adc[iAdc ][it]; // current
1140 an = adc[iAdc+1][it]; // next
1142 mask = ( ac >= ap && ac >= an ) ? 0 : 0x1; // peak center detection
1143 mask += ( ap + ac + an > eBIT ) ? 0 : 0x2; // cluster
1144 mask += ( ac > eBIS ) ? 0 : 0x4; // absolute large peak
1146 supp = (eBIL >> mask) & 1;
1148 fZSMap[iAdc] &= ~((1-supp) << it);
1149 if( eBIN == 0 ) { // neighbour sensitivity
1150 fZSMap[iAdc-1] &= ~((1-supp) << it);
1151 fZSMap[iAdc+1] &= ~((1-supp) << it);
1158 void AliTRDmcmSim::AddHitToFitreg(Int_t adc, UShort_t timebin, UShort_t qtot, Short_t ypos, Int_t label[])
1160 // Add the given hit to the fit register which is lateron used for
1161 // the tracklet calculation.
1162 // In addition to the fit sums in the fit register MC information
1165 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0)) &&
1166 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE0)))
1167 fFitReg[adc].fQ0 += qtot;
1169 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS1)) &&
1170 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1)))
1171 fFitReg[adc].fQ1 += qtot;
1173 if ((timebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS) ) &&
1174 (timebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFE)))
1176 fFitReg[adc].fSumX += timebin;
1177 fFitReg[adc].fSumX2 += timebin*timebin;
1178 fFitReg[adc].fNhits++;
1179 fFitReg[adc].fSumY += ypos;
1180 fFitReg[adc].fSumY2 += ypos*ypos;
1181 fFitReg[adc].fSumXY += timebin*ypos;
1184 // register hits (MC info)
1185 fHits[fNHits].fChannel = adc;
1186 fHits[fNHits].fQtot = qtot;
1187 fHits[fNHits].fYpos = ypos;
1188 fHits[fNHits].fTimebin = timebin;
1189 fHits[fNHits].fLabel[0] = label[0];
1190 fHits[fNHits].fLabel[1] = label[1];
1191 fHits[fNHits].fLabel[2] = label[2];
1195 void AliTRDmcmSim::CalcFitreg()
1198 // Detect the hits and fill the fit registers.
1199 // Requires 12-bit data from fADCF which means Filter()
1200 // has to be called before even if all filters are bypassed.
1202 //??? to be clarified:
1203 UInt_t adcMask = 0xffffffff;
1205 UShort_t timebin, adcch, adcLeft, adcCentral, adcRight, hitQual, timebin1, timebin2, qtotTemp;
1206 Short_t ypos, fromLeft, fromRight, found;
1207 UShort_t qTotal[19+1]; // the last is dummy
1208 UShort_t marked[6], qMarked[6], worse1, worse2;
1210 timebin1 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFS);
1211 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0)
1213 timebin1 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0);
1214 timebin2 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFE);
1215 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1)
1217 timebin2 = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1);
1219 // reset the fit registers
1221 for (adcch = 0; adcch < fgkNADC-2; adcch++) // due to border channels
1223 fFitReg[adcch].fNhits = 0;
1224 fFitReg[adcch].fQ0 = 0;
1225 fFitReg[adcch].fQ1 = 0;
1226 fFitReg[adcch].fSumX = 0;
1227 fFitReg[adcch].fSumY = 0;
1228 fFitReg[adcch].fSumX2 = 0;
1229 fFitReg[adcch].fSumY2 = 0;
1230 fFitReg[adcch].fSumXY = 0;
1233 for (timebin = timebin1; timebin < timebin2; timebin++)
1235 // first find the hit candidates and store the total cluster charge in qTotal array
1236 // in case of not hit store 0 there.
1237 for (adcch = 0; adcch < fgkNADC-2; adcch++) {
1238 if ( ( (adcMask >> adcch) & 7) == 7) //??? all 3 channels are present in case of ZS
1240 adcLeft = fADCF[adcch ][timebin];
1241 adcCentral = fADCF[adcch+1][timebin];
1242 adcRight = fADCF[adcch+2][timebin];
1243 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPVBY) == 1)
1244 hitQual = ( (adcLeft * adcRight) <
1245 (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPVT) * adcCentral) );
1248 // The accumulated charge is with the pedestal!!!
1249 qtotTemp = adcLeft + adcCentral + adcRight;
1251 (qtotTemp >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPHT)) &&
1252 (adcLeft <= adcCentral) &&
1253 (adcCentral > adcRight) )
1254 qTotal[adcch] = qtotTemp;
1259 qTotal[adcch] = 0; //jkl
1260 if (qTotal[adcch] != 0)
1261 AliDebug(10,Form("ch %2d qTotal %5d",adcch, qTotal[adcch]));
1267 marked[4] = 19; // invalid channel
1268 marked[5] = 19; // invalid channel
1270 while ((adcch < 16) && (found < 3))
1272 if (qTotal[adcch] > 0)
1275 marked[2*found+1]=adcch;
1284 while ((adcch > 2) && (found < 3))
1286 if (qTotal[adcch] > 0)
1288 marked[2*found]=adcch;
1295 AliDebug(10,Form("Fromleft=%d, Fromright=%d",fromLeft, fromRight));
1296 // here mask the hit candidates in the middle, if any
1297 if ((fromLeft >= 0) && (fromRight >= 0) && (fromLeft < fromRight))
1298 for (adcch = fromLeft+1; adcch < fromRight; adcch++)
1302 for (adcch = 0; adcch < 19; adcch++)
1303 if (qTotal[adcch] > 0) found++;
1306 if (found > 4) // sorting like in the TRAP in case of 5 or 6 candidates!
1308 if (marked[4] == marked[5]) marked[5] = 19;
1309 for (found=0; found<6; found++)
1311 qMarked[found] = qTotal[marked[found]] >> 4;
1312 AliDebug(10,Form("ch_%d qTotal %d qTotals %d",marked[found],qTotal[marked[found]],qMarked[found]));
1315 Sort6To2Worst(marked[0], marked[3], marked[4], marked[1], marked[2], marked[5],
1323 // Now mask the two channels with the smallest charge
1327 AliDebug(10,Form("Kill ch %d\n",worse1));
1332 AliDebug(10,Form("Kill ch %d\n",worse2));
1336 for (adcch = 0; adcch < 19; adcch++) {
1337 if (qTotal[adcch] > 0) // the channel is marked for processing
1339 adcLeft = fADCF[adcch ][timebin];
1340 adcCentral = fADCF[adcch+1][timebin];
1341 adcRight = fADCF[adcch+2][timebin];
1342 // hit detected, in TRAP we have 4 units and a hit-selection, here we proceed all channels!
1343 // subtract the pedestal TPFP, clipping instead of wrapping
1345 Int_t regTPFP = fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPFP);
1346 AliDebug(10, Form("Hit found, time=%d, adcch=%d/%d/%d, adc values=%d/%d/%d, regTPFP=%d, TPHT=%d\n",
1347 timebin, adcch, adcch+1, adcch+2, adcLeft, adcCentral, adcRight, regTPFP,
1348 fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPHT)));
1350 if (adcLeft < regTPFP) adcLeft = 0; else adcLeft -= regTPFP;
1351 if (adcCentral < regTPFP) adcCentral = 0; else adcCentral -= regTPFP;
1352 if (adcRight < regTPFP) adcRight = 0; else adcRight -= regTPFP;
1354 // Calculate the center of gravity
1355 // checking for adcCentral != 0 (in case of "bad" configuration)
1356 if (adcCentral == 0)
1358 ypos = 128*(adcLeft - adcRight) / adcCentral;
1359 if (ypos < 0) ypos = -ypos;
1360 // make the correction using the position LUT
1361 ypos = ypos + fTrapConfig->GetTrapReg((AliTRDtrapConfig::TrapReg_t) (AliTRDtrapConfig::kTPL00 + (ypos & 0x7F)),
1362 fDetector, fRobPos, fMcmPos);
1363 if (adcLeft > adcRight) ypos = -ypos;
1365 // label calculation (up to 3)
1366 Int_t mcLabel[] = {-1, -1, -1};
1367 if (fDigitsManager) {
1368 const Int_t maxLabels = 9;
1369 Int_t label[maxLabels] = { 0 }; // up to 9 different labels possible
1370 Int_t count[maxLabels] = { 0 };
1373 padcol[0] = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adcch);
1374 padcol[1] = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adcch+1);
1375 padcol[2] = fFeeParam->GetPadColFromADC(fRobPos, fMcmPos, adcch+2);
1376 Int_t padrow = fFeeParam->GetPadRowFromMCM(fRobPos, fMcmPos);
1377 for (Int_t iDict = 0; iDict < 3; iDict++) {
1380 for (Int_t iPad = 0; iPad < 3; iPad++) {
1381 if (padcol[iPad] < 0)
1383 Int_t currLabel = fDict[iDict]->GetData(padrow, padcol[iPad], timebin);
1384 AliDebug(10, Form("Read label: %4i for det: %3i, row: %i, col: %i, tb: %i\n", currLabel, fDetector, padrow, padcol[iPad], timebin));
1385 for (Int_t iLabel = 0; iLabel < nLabels; iLabel++) {
1386 if (currLabel == label[iLabel]) {
1392 if (currLabel >= 0) {
1393 label[nLabels] = currLabel;
1399 Int_t index[2*maxLabels];
1400 TMath::Sort(maxLabels, count, index);
1401 for (Int_t i = 0; i < 3; i++) {
1402 if (count[index[i]] <= 0)
1404 mcLabel[i] = label[index[i]];
1408 // add the hit to the fitregister
1409 AddHitToFitreg(adcch, timebin, qTotal[adcch], ypos, mcLabel);
1414 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1415 if (fFitReg[iAdc].fNhits != 0) {
1416 AliDebug(2, Form("fitreg[%i]: nHits = %i, sumX = %i, sumY = %i, sumX2 = %i, sumY2 = %i, sumXY = %i", iAdc,
1417 fFitReg[iAdc].fNhits,
1418 fFitReg[iAdc].fSumX,
1419 fFitReg[iAdc].fSumY,
1420 fFitReg[iAdc].fSumX2,
1421 fFitReg[iAdc].fSumY2,
1422 fFitReg[iAdc].fSumXY
1428 void AliTRDmcmSim::TrackletSelection()
1430 // Select up to 4 tracklet candidates from the fit registers
1431 // and assign them to the CPUs.
1433 UShort_t adcIdx, i, j, ntracks, tmp;
1434 UShort_t trackletCand[18][2]; // store the adcch[0] and number of hits[1] for all tracklet candidates
1437 for (adcIdx = 0; adcIdx < 18; adcIdx++) // ADCs
1438 if ( (fFitReg[adcIdx].fNhits
1439 >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPCL)) &&
1440 (fFitReg[adcIdx].fNhits+fFitReg[adcIdx+1].fNhits
1441 >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPCT)))
1443 trackletCand[ntracks][0] = adcIdx;
1444 trackletCand[ntracks][1] = fFitReg[adcIdx].fNhits+fFitReg[adcIdx+1].fNhits;
1445 AliDebug(10,Form("%d %2d %4d\n", ntracks, trackletCand[ntracks][0], trackletCand[ntracks][1]));
1449 for (i=0; i<ntracks;i++)
1450 AliDebug(10,Form("%d %d %d\n",i,trackletCand[i][0], trackletCand[i][1]));
1454 // primitive sorting according to the number of hits
1455 for (j = 0; j < (ntracks-1); j++)
1457 for (i = j+1; i < ntracks; i++)
1459 if ( (trackletCand[j][1] < trackletCand[i][1]) ||
1460 ( (trackletCand[j][1] == trackletCand[i][1]) && (trackletCand[j][0] < trackletCand[i][0]) ) )
1463 tmp = trackletCand[j][1];
1464 trackletCand[j][1] = trackletCand[i][1];
1465 trackletCand[i][1] = tmp;
1466 tmp = trackletCand[j][0];
1467 trackletCand[j][0] = trackletCand[i][0];
1468 trackletCand[i][0] = tmp;
1472 ntracks = 4; // cut the rest, 4 is the max
1474 // else is not necessary to sort
1476 // now sort, so that the first tracklet going to CPU0 corresponds to the highest adc channel - as in the TRAP
1477 for (j = 0; j < (ntracks-1); j++)
1479 for (i = j+1; i < ntracks; i++)
1481 if (trackletCand[j][0] < trackletCand[i][0])
1484 tmp = trackletCand[j][1];
1485 trackletCand[j][1] = trackletCand[i][1];
1486 trackletCand[i][1] = tmp;
1487 tmp = trackletCand[j][0];
1488 trackletCand[j][0] = trackletCand[i][0];
1489 trackletCand[i][0] = tmp;
1493 for (i = 0; i < ntracks; i++) // CPUs with tracklets.
1494 fFitPtr[i] = trackletCand[i][0]; // pointer to the left channel with tracklet for CPU[i]
1495 for (i = ntracks; i < 4; i++) // CPUs without tracklets
1496 fFitPtr[i] = 31; // pointer to the left channel with tracklet for CPU[i] = 31 (invalid)
1497 AliDebug(10,Form("found %i tracklet candidates\n", ntracks));
1498 for (i = 0; i < 4; i++)
1499 AliDebug(10,Form("fitPtr[%i]: %i\n", i, fFitPtr[i]));
1502 void AliTRDmcmSim::FitTracklet()
1504 // Perform the actual tracklet fit based on the fit sums
1505 // which have been filled in the fit registers.
1507 // parameters in fitred.asm (fit program)
1508 Int_t decPlaces = 5;
1511 rndAdd = (1 << (decPlaces-1)) + 1;
1512 else if (decPlaces == 1)
1514 Int_t ndriftDp = 5; // decimal places for drift time
1515 Long64_t shift = ((Long64_t) 1 << 32);
1517 // calculated in fitred.asm
1518 Int_t padrow = ((fRobPos >> 1) << 2) | (fMcmPos >> 2);
1519 Int_t yoffs = (((((fRobPos & 0x1) << 2) + (fMcmPos & 0x3)) * 18) << 8) -
1520 ((18*4*2 - 18*2 - 1) << 7);
1521 yoffs = yoffs << decPlaces; // holds position of ADC channel 1
1522 Int_t layer = fDetector % 6;
1523 UInt_t scaleY = (UInt_t) ((0.635 + 0.03 * layer)/(256.0 * 160.0e-4) * shift);
1524 UInt_t scaleD = (UInt_t) ((0.635 + 0.03 * layer)/(256.0 * 140.0e-4) * shift);
1526 Int_t deflCorr = (Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCorr, fDetector, fRobPos, fMcmPos);
1527 Int_t ndrift = (Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrNdrift, fDetector, fRobPos, fMcmPos);
1529 // local variables for calculation
1530 Long64_t mult, temp, denom; //???
1531 UInt_t q0, q1, pid; // charges in the two windows and total charge
1532 UShort_t nHits; // number of hits
1533 Int_t slope, offset; // slope and offset of the tracklet
1534 Int_t sumX, sumY, sumXY, sumX2; // fit sums from fit registers
1535 Int_t sumY2; // not used in the current TRAP program, now used for error calculation (simulation only)
1536 Float_t fitError, fitSlope, fitOffset;
1537 FitReg_t *fit0, *fit1; // pointers to relevant fit registers
1539 // const uint32_t OneDivN[32] = { // 2**31/N : exactly like in the TRAP, the simple division here gives the same result!
1540 // 0x00000000, 0x80000000, 0x40000000, 0x2AAAAAA0, 0x20000000, 0x19999990, 0x15555550, 0x12492490,
1541 // 0x10000000, 0x0E38E380, 0x0CCCCCC0, 0x0BA2E8B0, 0x0AAAAAA0, 0x09D89D80, 0x09249240, 0x08888880,
1542 // 0x08000000, 0x07878780, 0x071C71C0, 0x06BCA1A0, 0x06666660, 0x06186180, 0x05D17450, 0x0590B210,
1543 // 0x05555550, 0x051EB850, 0x04EC4EC0, 0x04BDA120, 0x04924920, 0x0469EE50, 0x04444440, 0x04210840};
1545 for (Int_t cpu = 0; cpu < 4; cpu++) {
1546 if (fFitPtr[cpu] == 31)
1548 fMCMT[cpu] = 0x10001000; //??? AliTRDfeeParam::GetTrackletEndmarker();
1552 fit0 = &fFitReg[fFitPtr[cpu] ];
1553 fit1 = &fFitReg[fFitPtr[cpu]+1]; // next channel
1556 mult = mult << (32 + decPlaces);
1560 nHits = fit0->fNhits + fit1->fNhits; // number of hits
1561 sumX = fit0->fSumX + fit1->fSumX;
1562 sumX2 = fit0->fSumX2 + fit1->fSumX2;
1563 denom = ((Long64_t) nHits)*((Long64_t) sumX2) - ((Long64_t) sumX)*((Long64_t) sumX);
1565 mult = mult / denom; // exactly like in the TRAP program
1566 q0 = fit0->fQ0 + fit1->fQ0;
1567 q1 = fit0->fQ1 + fit1->fQ1;
1568 sumY = fit0->fSumY + fit1->fSumY + 256*fit1->fNhits;
1569 sumXY = fit0->fSumXY + fit1->fSumXY + 256*fit1->fSumX;
1570 sumY2 = fit0->fSumY2 + fit1->fSumY2 + 512*fit1->fSumY + 256*256*fit1->fNhits;
1572 slope = nHits*sumXY - sumX * sumY;
1573 offset = sumX2*sumY - sumX * sumXY;
1574 temp = mult * slope;
1575 slope = temp >> 32; // take the upper 32 bits
1577 temp = mult * offset;
1578 offset = temp >> 32; // take the upper 32 bits
1580 offset = offset + yoffs;
1581 AliDebug(10, Form("slope = %i, slope * ndrift = %i, deflCorr: %i",
1582 slope, slope * ndrift, deflCorr));
1583 slope = ((slope * ndrift) >> ndriftDp) + deflCorr;
1584 offset = offset - (fFitPtr[cpu] << (8 + decPlaces));
1587 temp = temp * scaleD;
1588 slope = (temp >> 32);
1590 temp = temp * scaleY;
1591 offset = (temp >> 32);
1593 // rounding, like in the TRAP
1594 slope = (slope + rndAdd) >> decPlaces;
1595 offset = (offset + rndAdd) >> decPlaces;
1597 AliDebug(5, Form("Det: %3i, ROB: %i, MCM: %2i: deflection: %i, min: %i, max: %i",
1598 fDetector, fRobPos, fMcmPos, slope,
1599 (Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCutStart + 2*fFitPtr[cpu], fDetector, fRobPos, fMcmPos),
1600 (Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCutStart + 1 + 2*fFitPtr[cpu], fDetector, fRobPos, fMcmPos)));
1602 AliDebug(5, Form("Fit sums: x = %i, X = %i, y = %i, Y = %i, Z = %i",
1603 sumX, sumX2, sumY, sumY2, sumXY));
1605 fitSlope = (Float_t) (nHits * sumXY - sumX * sumY) / (nHits * sumX2 - sumX*sumX);
1607 fitOffset = (Float_t) (sumX2 * sumY - sumX * sumXY) / (nHits * sumX2 - sumX*sumX);
1609 Float_t sx = (Float_t) sumX;
1610 Float_t sx2 = (Float_t) sumX2;
1611 Float_t sy = (Float_t) sumY;
1612 Float_t sy2 = (Float_t) sumY2;
1613 Float_t sxy = (Float_t) sumXY;
1614 fitError = sy2 - (sx2 * sy*sy - 2 * sx * sxy * sy + nHits * sxy*sxy) / (nHits * sx2 - sx*sx);
1615 //fitError = (Float_t) sumY2 - (Float_t) (sumY*sumY) / nHits - fitSlope * ((Float_t) (sumXY - sumX*sumY) / nHits);
1617 Bool_t rejected = kFALSE;
1618 // deflection range table from DMEM
1619 if ((slope < ((Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCutStart + 2*fFitPtr[cpu], fDetector, fRobPos, fMcmPos))) ||
1620 (slope > ((Int_t) fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrDeflCutStart + 1 + 2*fFitPtr[cpu], fDetector, fRobPos, fMcmPos))))
1623 if (rejected && GetApplyCut())
1625 fMCMT[cpu] = 0x10001000; //??? AliTRDfeeParam::GetTrackletEndmarker();
1629 if (slope > 63 || slope < -64) { // wrapping in TRAP!
1630 AliError(Form("Overflow in slope: %i, tracklet discarded!", slope));
1631 fMCMT[cpu] = 0x10001000;
1635 slope = slope & 0x7F; // 7 bit
1637 if (offset > 0xfff || offset < -0xfff)
1638 AliWarning("Overflow in offset");
1639 offset = offset & 0x1FFF; // 13 bit
1641 pid = GetPID(q0 >> fgkAddDigits, q1 >> fgkAddDigits); // divided by 4 because in simulation there are two additional decimal places
1644 AliWarning("Overflow in PID");
1645 pid = pid & 0xFF; // 8 bit, exactly like in the TRAP program
1647 // assemble and store the tracklet word
1648 fMCMT[cpu] = (pid << 24) | (padrow << 20) | (slope << 13) | offset;
1650 // calculate MC label
1651 Int_t mcLabel[] = { -1, -1, -1};
1654 if (fDigitsManager) {
1655 const Int_t maxLabels = 30;
1656 Int_t label[maxLabels] = {0}; // up to 30 different labels possible
1657 Int_t count[maxLabels] = {0};
1659 for (Int_t iHit = 0; iHit < fNHits; iHit++) {
1660 if ((fHits[iHit].fChannel - fFitPtr[cpu] < 0) ||
1661 (fHits[iHit].fChannel - fFitPtr[cpu] > 1))
1664 // counting contributing hits
1665 if (fHits[iHit].fTimebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS0) &&
1666 fHits[iHit].fTimebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE0))
1668 if (fHits[iHit].fTimebin >= fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQS1) &&
1669 fHits[iHit].fTimebin < fTrapConfig->GetTrapReg(AliTRDtrapConfig::kTPQE1))
1672 for (Int_t i = 0; i < 3; i++) {
1673 Int_t currLabel = fHits[iHit].fLabel[i];
1674 for (Int_t iLabel = 0; iLabel < nLabels; iLabel++) {
1675 if (currLabel == label[iLabel]) {
1681 if (currLabel >= 0 && nLabels < maxLabels) {
1682 label[nLabels] = currLabel;
1688 Int_t index[2*maxLabels];
1689 TMath::Sort(maxLabels, count, index);
1690 for (Int_t i = 0; i < 3; i++) {
1691 if (count[index[i]] <= 0)
1693 mcLabel[i] = label[index[i]];
1696 new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletMCM((UInt_t) fMCMT[cpu], fDetector*2 + fRobPos%2, fRobPos, fMcmPos);
1697 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetLabel(mcLabel);
1700 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetNHits(fit0->fNhits + fit1->fNhits);
1701 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetNHits0(nHits0);
1702 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetNHits1(nHits1);
1703 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetQ0(q0);
1704 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetQ1(q1);
1705 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetSlope(fitSlope);
1706 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetOffset(fitOffset);
1707 ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetError(TMath::Sqrt(TMath::Abs(fitError)/nHits));
1709 // // cluster information
1710 // Float_t *res = new Float_t[nHits];
1711 // Float_t *qtot = new Float_t[nHits];
1713 // for (Int_t iHit = 0; iHit < fNHits; iHit++) {
1714 // // check if hit contributes
1715 // if (fHits[iHit].fChannel == fFitPtr[cpu]) {
1716 // res[nCls] = fHits[iHit].fYpos - (fitSlope * fHits[iHit].fTimebin + fitOffset);
1717 // qtot[nCls] = fHits[iHit].fQtot;
1720 // else if (fHits[iHit].fChannel == fFitPtr[cpu] + 1) {
1721 // res[nCls] = fHits[iHit].fYpos + 256 - (fitSlope * fHits[iHit].fTimebin + fitOffset);
1722 // qtot[nCls] = fHits[iHit].fQtot;
1726 // ((AliTRDtrackletMCM*) (*fTrackletArray)[fTrackletArray->GetEntriesFast()-1])->SetClusters(res, qtot, nCls);
1731 AliError(Form("Strange fit error: %f from Sx: %i, Sy: %i, Sxy: %i, Sx2: %i, Sy2: %i, nHits: %i",
1732 fitError, sumX, sumY, sumXY, sumX2, sumY2, nHits));
1733 AliDebug(3, Form("fit slope: %f, offset: %f, error: %f",
1734 fitSlope, fitOffset, TMath::Sqrt(TMath::Abs(fitError)/nHits)));
1740 void AliTRDmcmSim::Tracklet()
1742 // Run the tracklet calculation by calling sequentially:
1743 // CalcFitreg(); TrackletSelection(); FitTracklet()
1744 // and store the tracklets
1746 if (!fInitialized) {
1747 AliError("Called uninitialized! Nothing done!");
1751 fTrackletArray->Delete();
1756 TrackletSelection();
1760 Bool_t AliTRDmcmSim::StoreTracklets()
1762 // store the found tracklets via the loader
1764 if (fTrackletArray->GetEntriesFast() == 0)
1767 AliRunLoader *rl = AliRunLoader::Instance();
1768 AliDataLoader *dl = 0x0;
1770 dl = rl->GetLoader("TRDLoader")->GetDataLoader("tracklets");
1772 AliError("Could not get the tracklets data loader!");
1776 TTree *trackletTree = dl->Tree();
1777 if (!trackletTree) {
1779 trackletTree = dl->Tree();
1782 AliTRDtrackletMCM *trkl = 0x0;
1783 TBranch *trkbranch = trackletTree->GetBranch(fTrklBranchName.Data());
1785 trkbranch = trackletTree->Branch(fTrklBranchName.Data(), "AliTRDtrackletMCM", &trkl, 32000);
1787 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
1788 trkl = ((AliTRDtrackletMCM*) (*fTrackletArray)[iTracklet]);
1789 trkbranch->SetAddress(&trkl);
1796 void AliTRDmcmSim::WriteData(AliTRDarrayADC *digits)
1798 // write back the processed data configured by EBSF
1799 // EBSF = 1: unfiltered data; EBSF = 0: filtered data
1800 // zero-suppressed valued are written as -1 to digits
1802 if( !CheckInitialized() )
1805 Int_t offset = (fMcmPos % 4 + 1) * 21 + (fRobPos % 2) * 84 - 1;
1807 if (fTrapConfig->GetTrapReg(AliTRDtrapConfig::kEBSF) != 0) // store unfiltered data
1809 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1810 if (~fZSMap[iAdc] == 0) {
1811 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1812 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, -1);
1815 else if (iAdc < 2 || iAdc == 20) {
1816 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1817 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, (fADCR[iAdc][iTimeBin] >> fgkAddDigits) - fgAddBaseline);
1823 for (Int_t iAdc = 0; iAdc < fgkNADC; iAdc++) {
1824 if (~fZSMap[iAdc] != 0) {
1825 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1826 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, (fADCF[iAdc][iTimeBin] >> fgkAddDigits) - fgAddBaseline);
1830 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
1831 digits->SetDataByAdcCol(GetRow(), offset - iAdc, iTimeBin, -1);
1839 // ******************************
1842 // Memory area for the LUT: 0xC100 to 0xC3FF
1844 // The addresses for the parameters (the order is optimized for maximum calculation speed in the MCMs):
1846 // 0xC029: nBins(sF)
1848 // 0xC02B: TableLength
1849 // Defined in AliTRDtrapConfig.h
1851 // The algorithm implemented in the TRAP program of the MCMs (Venelin Angelov)
1852 // 1) set the read pointer to the beginning of the Parameters in DMEM
1853 // 2) shift right the FitReg with the Q0 + (Q1 << 16) to get Q1
1854 // 3) read cor1 with rpointer++
1856 // 5) read nBins with rpointer++
1857 // 6) start nBins*cor1*Q1
1858 // 7) read cor0 with rpointer++
1859 // 8) swap hi-low parts in FitReg, now is Q1 + (Q0 << 16)
1860 // 9) shift right to get Q0
1861 // 10) start cor0*Q0
1862 // 11) read TableLength
1863 // 12) compare cor0*Q0 with nBins
1864 // 13) if >=, clip cor0*Q0 to nBins-1
1865 // 14) add cor0*Q0 to nBins*cor1*Q1
1866 // 15) compare the result with TableLength
1867 // 16) if >=, clip to TableLength-1
1868 // 17) read from the LUT 8 bits
1871 Int_t AliTRDmcmSim::GetPID(Int_t q0, Int_t q1)
1873 // return PID calculated from charges accumulated in two time windows
1878 UInt_t nBinsQ0 = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTnbins); // number of bins in q0 / 4 !!
1879 UInt_t pidTotalSize = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTLength);
1880 if(nBinsQ0==0 || pidTotalSize==0) // make sure we don't run into trouble if one of the values is not configured
1883 ULong_t corrQ0 = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTcor0, fDetector, fRobPos, fMcmPos);
1884 ULong_t corrQ1 = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTcor1, fDetector, fRobPos, fMcmPos);
1885 if(corrQ0==0 || corrQ1==0) // make sure we don't run into trouble if one of the values is not configured
1889 addrQ0 = (((addrQ0*q0)>>16)>>16); // because addrQ0 = (q0 * corrQ0) >> 32; does not work for unknown reasons
1891 if(addrQ0 >= nBinsQ0) { // check for overflow
1892 AliDebug(5,Form("Overflow in q0: %llu/4 is bigger then %u", addrQ0, nBinsQ0));
1893 addrQ0 = nBinsQ0 -1;
1897 addr = (((addr*q1)>>16)>>16);
1898 addr = addrQ0 + nBinsQ0*addr; // because addr = addrQ0 + nBinsQ0* (((corrQ1*q1)>>32); does not work
1900 if(addr >= pidTotalSize) {
1901 AliDebug(5,Form("Overflow in q1. Address %llu/4 is bigger then %u", addr, pidTotalSize));
1902 addr = pidTotalSize -1;
1905 // 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)
1907 UInt_t result = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTStart+(addr/4));
1908 return (result>>((addr%4)*8)) & 0xFF;
1913 // help functions, to be cleaned up
1915 UInt_t AliTRDmcmSim::AddUintClipping(UInt_t a, UInt_t b, UInt_t nbits) const
1918 // This function adds a and b (unsigned) and clips to
1919 // the specified number of bits.
1925 UInt_t maxv = (1 << nbits) - 1;;
1931 if ((sum < a) || (sum < b))
1937 void AliTRDmcmSim::Sort2(UShort_t idx1i, UShort_t idx2i, \
1938 UShort_t val1i, UShort_t val2i, \
1939 UShort_t * const idx1o, UShort_t * const idx2o, \
1940 UShort_t * const val1o, UShort_t * const val2o) const
1942 // sorting for tracklet selection
1960 void AliTRDmcmSim::Sort3(UShort_t idx1i, UShort_t idx2i, UShort_t idx3i, \
1961 UShort_t val1i, UShort_t val2i, UShort_t val3i, \
1962 UShort_t * const idx1o, UShort_t * const idx2o, UShort_t * const idx3o, \
1963 UShort_t * const val1o, UShort_t * const val2o, UShort_t * const val3o)
1965 // sorting for tracklet selection
1970 if (val1i > val2i) sel=4; else sel=0;
1971 if (val2i > val3i) sel=sel + 2;
1972 if (val3i > val1i) sel=sel + 1;
1975 case 6 : // 1 > 2 > 3 => 1 2 3
1976 case 0 : // 1 = 2 = 3 => 1 2 3 : in this case doesn't matter, but so is in hardware!
1985 case 4 : // 1 > 2, 2 <= 3, 3 <= 1 => 1 3 2
1994 case 2 : // 1 <= 2, 2 > 3, 3 <= 1 => 2 1 3
2003 case 3 : // 1 <= 2, 2 > 3, 3 > 1 => 2 3 1
2012 case 1 : // 1 <= 2, 2 <= 3, 3 > 1 => 3 2 1
2021 case 5 : // 1 > 2, 2 <= 3, 3 > 1 => 3 1 2
2030 default: // the rest should NEVER happen!
2031 AliError("ERROR in Sort3!!!\n");
2036 void AliTRDmcmSim::Sort6To4(UShort_t idx1i, UShort_t idx2i, UShort_t idx3i, UShort_t idx4i, UShort_t idx5i, UShort_t idx6i, \
2037 UShort_t val1i, UShort_t val2i, UShort_t val3i, UShort_t val4i, UShort_t val5i, UShort_t val6i, \
2038 UShort_t * const idx1o, UShort_t * const idx2o, UShort_t * const idx3o, UShort_t * const idx4o, \
2039 UShort_t * const val1o, UShort_t * const val2o, UShort_t * const val3o, UShort_t * const val4o)
2041 // sorting for tracklet selection
2043 UShort_t idx21s, idx22s, idx23s, dummy;
2044 UShort_t val21s, val22s, val23s;
2045 UShort_t idx23as, idx23bs;
2046 UShort_t val23as, val23bs;
2048 Sort3(idx1i, idx2i, idx3i, val1i, val2i, val3i,
2049 idx1o, &idx21s, &idx23as,
2050 val1o, &val21s, &val23as);
2052 Sort3(idx4i, idx5i, idx6i, val4i, val5i, val6i,
2053 idx2o, &idx22s, &idx23bs,
2054 val2o, &val22s, &val23bs);
2056 Sort2(idx23as, idx23bs, val23as, val23bs, &idx23s, &dummy, &val23s, &dummy);
2058 Sort3(idx21s, idx22s, idx23s, val21s, val22s, val23s,
2059 idx3o, idx4o, &dummy,
2060 val3o, val4o, &dummy);
2064 void AliTRDmcmSim::Sort6To2Worst(UShort_t idx1i, UShort_t idx2i, UShort_t idx3i, UShort_t idx4i, UShort_t idx5i, UShort_t idx6i, \
2065 UShort_t val1i, UShort_t val2i, UShort_t val3i, UShort_t val4i, UShort_t val5i, UShort_t val6i, \
2066 UShort_t * const idx5o, UShort_t * const idx6o)
2068 // sorting for tracklet selection
2070 UShort_t idx21s, idx22s, idx23s, dummy1, dummy2, dummy3, dummy4, dummy5;
2071 UShort_t val21s, val22s, val23s;
2072 UShort_t idx23as, idx23bs;
2073 UShort_t val23as, val23bs;
2075 Sort3(idx1i, idx2i, idx3i, val1i, val2i, val3i,
2076 &dummy1, &idx21s, &idx23as,
2077 &dummy2, &val21s, &val23as);
2079 Sort3(idx4i, idx5i, idx6i, val4i, val5i, val6i,
2080 &dummy1, &idx22s, &idx23bs,
2081 &dummy2, &val22s, &val23bs);
2083 Sort2(idx23as, idx23bs, val23as, val23bs, &idx23s, idx5o, &val23s, &dummy1);
2085 Sort3(idx21s, idx22s, idx23s, val21s, val22s, val23s,
2086 &dummy1, &dummy2, idx6o,
2087 &dummy3, &dummy4, &dummy5);
2091 // ----- I/O implementation -----
2093 ostream& AliTRDmcmSim::Text(ostream& os)
2095 // manipulator to activate output in text format (default)
2097 os.iword(fgkFormatIndex) = 0;
2101 ostream& AliTRDmcmSim::Cfdat(ostream& os)
2103 // manipulator to activate output in CFDAT format
2104 // to send to the FEE via SCSN
2106 os.iword(fgkFormatIndex) = 1;
2110 ostream& AliTRDmcmSim::Raw(ostream& os)
2112 // manipulator to activate output as raw data dump
2114 os.iword(fgkFormatIndex) = 2;
2118 ostream& operator<<(ostream& os, const AliTRDmcmSim& mcm)
2120 // output implementation
2122 // no output for non-initialized MCM
2123 if (!mcm.CheckInitialized())
2126 // ----- human-readable output -----
2127 if (os.iword(AliTRDmcmSim::fgkFormatIndex) == 0) {
2129 os << "MCM " << mcm.fMcmPos << " on ROB " << mcm.fRobPos <<
2130 " in detector " << mcm.fDetector << std::endl;
2132 os << "----- Unfiltered ADC data (10 bit) -----" << std::endl;
2134 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++)
2135 os << std::setw(5) << iChannel;
2137 for (Int_t iTimeBin = 0; iTimeBin < mcm.fNTimeBin; iTimeBin++) {
2138 os << "tb " << std::setw(2) << iTimeBin << ":";
2139 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++) {
2140 os << std::setw(5) << (mcm.fADCR[iChannel][iTimeBin] >> mcm.fgkAddDigits);
2145 os << "----- Filtered ADC data (10+2 bit) -----" << std::endl;
2147 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++)
2148 os << std::setw(4) << iChannel
2149 << ((~mcm.fZSMap[iChannel] != 0) ? "!" : " ");
2151 for (Int_t iTimeBin = 0; iTimeBin < mcm.fNTimeBin; iTimeBin++) {
2152 os << "tb " << std::setw(2) << iTimeBin << ":";
2153 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++) {
2154 os << std::setw(4) << (mcm.fADCF[iChannel][iTimeBin])
2155 << (((mcm.fZSMap[iChannel] & (1 << iTimeBin)) == 0) ? "!" : " ");
2161 // ----- CFDAT output -----
2162 else if(os.iword(AliTRDmcmSim::fgkFormatIndex) == 1) {
2164 Int_t addrOffset = 0x2000;
2165 Int_t addrStep = 0x80;
2167 for (Int_t iTimeBin = 0; iTimeBin < mcm.fNTimeBin; iTimeBin++) {
2168 for (Int_t iChannel = 0; iChannel < mcm.fgkNADC; iChannel++) {
2169 os << std::setw(5) << 10
2170 << std::setw(5) << addrOffset + iChannel * addrStep + iTimeBin
2171 << std::setw(5) << (mcm.fADCF[iChannel][iTimeBin])
2172 << std::setw(5) << dest << std::endl;
2178 // ----- raw data ouptut -----
2179 else if (os.iword(AliTRDmcmSim::fgkFormatIndex) == 2) {
2180 Int_t bufSize = 300;
2181 UInt_t *buf = new UInt_t[bufSize];
2183 Int_t bufLength = mcm.ProduceRawStream(&buf[0], bufSize);
2185 for (Int_t i = 0; i < bufLength; i++)
2186 std::cout << "0x" << std::hex << buf[i] << std::endl;
2192 os << "unknown format set" << std::endl;
2199 void AliTRDmcmSim::PrintFitRegXml(ostream& os) const
2201 // print fit registres in XML format
2203 bool tracklet=false;
2205 for (Int_t cpu = 0; cpu < 4; cpu++) {
2206 if(fFitPtr[cpu] != 31)
2210 if(tracklet==true) {
2211 os << "<nginject>" << std::endl;
2212 os << "<ack roc=\""<< fDetector << "\" cmndid=\"0\">" << std::endl;
2213 os << "<dmem-readout>" << std::endl;
2214 os << "<d det=\"" << fDetector << "\">" << std::endl;
2215 os << " <ro-board rob=\"" << fRobPos << "\">" << std::endl;
2216 os << " <m mcm=\"" << fMcmPos << "\">" << std::endl;
2218 for(int cpu=0; cpu<4; cpu++) {
2219 os << " <c cpu=\"" << cpu << "\">" << std::endl;
2220 if(fFitPtr[cpu] != 31) {
2221 for(int adcch=fFitPtr[cpu]; adcch<fFitPtr[cpu]+2; adcch++) {
2222 os << " <ch chnr=\"" << adcch << "\">"<< std::endl;
2223 os << " <hits>" << fFitReg[adcch].fNhits << "</hits>"<< std::endl;
2224 os << " <q0>" << fFitReg[adcch].fQ0/4 << "</q0>"<< std::endl; // divided by 4 because in simulation we have 2 additional decimal places
2225 os << " <q1>" << fFitReg[adcch].fQ1/4 << "</q1>"<< std::endl; // in the output
2226 os << " <sumx>" << fFitReg[adcch].fSumX << "</sumx>"<< std::endl;
2227 os << " <sumxsq>" << fFitReg[adcch].fSumX2 << "</sumxsq>"<< std::endl;
2228 os << " <sumy>" << fFitReg[adcch].fSumY << "</sumy>"<< std::endl;
2229 os << " <sumysq>" << fFitReg[adcch].fSumY2 << "</sumysq>"<< std::endl;
2230 os << " <sumxy>" << fFitReg[adcch].fSumXY << "</sumxy>"<< std::endl;
2231 os << " </ch>" << std::endl;
2234 os << " </c>" << std::endl;
2236 os << " </m>" << std::endl;
2237 os << " </ro-board>" << std::endl;
2238 os << "</d>" << std::endl;
2239 os << "</dmem-readout>" << std::endl;
2240 os << "</ack>" << std::endl;
2241 os << "</nginject>" << std::endl;
2246 void AliTRDmcmSim::PrintTrackletsXml(ostream& os) const
2248 // print tracklets in XML format
2250 os << "<nginject>" << std::endl;
2251 os << "<ack roc=\""<< fDetector << "\" cmndid=\"0\">" << std::endl;
2252 os << "<dmem-readout>" << std::endl;
2253 os << "<d det=\"" << fDetector << "\">" << std::endl;
2254 os << " <ro-board rob=\"" << fRobPos << "\">" << std::endl;
2255 os << " <m mcm=\"" << fMcmPos << "\">" << std::endl;
2257 Int_t pid, padrow, slope, offset;
2258 for(Int_t cpu=0; cpu<4; cpu++) {
2259 if(fMCMT[cpu] == 0x10001000) {
2266 pid = (fMCMT[cpu] & 0xFF000000) >> 24;
2267 padrow = (fMCMT[cpu] & 0xF00000 ) >> 20;
2268 slope = (fMCMT[cpu] & 0xFE000 ) >> 13;
2269 offset = (fMCMT[cpu] & 0x1FFF ) ;
2272 os << " <trk> <pid>" << pid << "</pid>" << " <padrow>" << padrow << "</padrow>"
2273 << " <slope>" << slope << "</slope>" << " <offset>" << offset << "</offset>" << "</trk>" << std::endl;
2276 os << " </m>" << std::endl;
2277 os << " </ro-board>" << std::endl;
2278 os << "</d>" << std::endl;
2279 os << "</dmem-readout>" << std::endl;
2280 os << "</ack>" << std::endl;
2281 os << "</nginject>" << std::endl;
2285 void AliTRDmcmSim::PrintAdcDatHuman(ostream& os) const
2287 // print ADC data in human-readable format
2289 os << "MCM " << fMcmPos << " on ROB " << fRobPos <<
2290 " in detector " << fDetector << std::endl;
2292 os << "----- Unfiltered ADC data (10 bit) -----" << std::endl;
2294 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++)
2295 os << std::setw(5) << iChannel;
2297 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
2298 os << "tb " << std::setw(2) << iTimeBin << ":";
2299 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2300 os << std::setw(5) << (fADCR[iChannel][iTimeBin] >> fgkAddDigits);
2305 os << "----- Filtered ADC data (10+2 bit) -----" << std::endl;
2307 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++)
2308 os << std::setw(4) << iChannel
2309 << ((~fZSMap[iChannel] != 0) ? "!" : " ");
2311 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
2312 os << "tb " << std::setw(2) << iTimeBin << ":";
2313 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2314 os << std::setw(4) << (fADCF[iChannel][iTimeBin])
2315 << (((fZSMap[iChannel] & (1 << iTimeBin)) == 0) ? "!" : " ");
2322 void AliTRDmcmSim::PrintAdcDatXml(ostream& os) const
2324 // print ADC data in XML format
2326 os << "<nginject>" << std::endl;
2327 os << "<ack roc=\""<< fDetector << "\" cmndid=\"0\">" << std::endl;
2328 os << "<dmem-readout>" << std::endl;
2329 os << "<d det=\"" << fDetector << "\">" << std::endl;
2330 os << " <ro-board rob=\"" << fRobPos << "\">" << std::endl;
2331 os << " <m mcm=\"" << fMcmPos << "\">" << std::endl;
2333 for(Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2334 os << " <ch chnr=\"" << iChannel << "\">" << std::endl;
2335 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
2336 os << "<tb>" << fADCF[iChannel][iTimeBin]/4 << "</tb>";
2338 os << " </ch>" << std::endl;
2341 os << " </m>" << std::endl;
2342 os << " </ro-board>" << std::endl;
2343 os << "</d>" << std::endl;
2344 os << "</dmem-readout>" << std::endl;
2345 os << "</ack>" << std::endl;
2346 os << "</nginject>" << std::endl;
2351 void AliTRDmcmSim::PrintAdcDatDatx(ostream& os, Bool_t broadcast) const
2353 // print ADC data in datx format (to send to FEE)
2355 fTrapConfig->PrintDatx(os, 2602, 1, 0, 127); // command to enable the ADC clock - necessary to write ADC values to MCM
2358 Int_t addrOffset = 0x2000;
2359 Int_t addrStep = 0x80;
2360 Int_t addrOffsetEBSIA = 0x20;
2362 for (Int_t iTimeBin = 0; iTimeBin < fNTimeBin; iTimeBin++) {
2363 for (Int_t iChannel = 0; iChannel < fgkNADC; iChannel++) {
2364 if(broadcast==kFALSE)
2365 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, (fADCF[iChannel][iTimeBin]/4), GetRobPos(), GetMcmPos());
2367 fTrapConfig->PrintDatx(os, addrOffset+iChannel*addrStep+addrOffsetEBSIA+iTimeBin, (fADCF[iChannel][iTimeBin]/4), 0, 127);
2374 void AliTRDmcmSim::PrintPidLutHuman()
2376 // print PID LUT in human readable format
2380 UInt_t addrEnd = AliTRDtrapConfig::fgkDmemAddrLUTStart + fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTLength)/4; // /4 because each addr contains 4 values
2381 UInt_t nBinsQ0 = fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTnbins);
2383 std::cout << "nBinsQ0: " << nBinsQ0 << std::endl;
2384 std::cout << "LUT table length: " << fTrapConfig->GetDmemUnsigned(AliTRDtrapConfig::fgkDmemAddrLUTLength) << std::endl;
2386 for(UInt_t addr=AliTRDtrapConfig::fgkDmemAddrLUTStart; addr< addrEnd; addr++) {
2387 result = fTrapConfig->GetDmemUnsigned(addr);
2388 std::cout << addr << " # x: " << ((addr-AliTRDtrapConfig::fgkDmemAddrLUTStart)%((nBinsQ0)/4))*4 << ", y: " <<(addr-AliTRDtrapConfig::fgkDmemAddrLUTStart)/(nBinsQ0/4)
2389 << " # " <<((result>>0)&0xFF)
2390 << " | " << ((result>>8)&0xFF)
2391 << " | " << ((result>>16)&0xFF)
2392 << " | " << ((result>>24)&0xFF) << std::endl;