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