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