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