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