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