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