1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
27 #include <TStopwatch.h>
39 #include "AliITShit.h"
40 #include "AliITSdigit.h"
41 #include "AliITSmodule.h"
42 #include "AliITSpList.h"
43 #include "AliITSMapA1.h"
44 #include "AliITSMapA2.h"
45 #include "AliITSetfSDD.h"
46 #include "AliITSRawData.h"
47 #include "AliITSHuffman.h"
48 #include "AliITSgeom.h"
49 #include "AliITSsegmentation.h"
50 #include "AliITSresponse.h"
51 #include "AliITSsegmentationSDD.h"
52 #include "AliITSresponseSDD.h"
53 #include "AliITSsimulationSDD.h"
55 ClassImp(AliITSsimulationSDD)
56 ////////////////////////////////////////////////////////////////////////
58 // Written by Piergiorgio Cerello
61 // AliITSsimulationSDD is the simulation of SDDs.
65 <img src="picts/ITS/AliITShit_Class_Diagram.gif">
68 <font size=+2 color=red>
69 <p>This show the relasionships between the ITS hit class and the rest of Aliroot.
74 //______________________________________________________________________
75 Int_t power(Int_t b, Int_t e) {
76 // compute b to the e power, where both b and e are Int_ts.
79 for(i=0; i<e; i++) power *= b;
82 //______________________________________________________________________
83 void FastFourierTransform(AliITSetfSDD *alisddetf,Double_t *real,
84 Double_t *imag,Int_t direction) {
85 // Do a Fast Fourier Transform
87 Int_t samples = alisddetf->GetSamples();
88 Int_t l = (Int_t) ((log((Float_t) samples)/log(2.))+0.5);
91 Int_t m2 = samples/m1;
94 for(j=0; j<samples; j += m1) {
96 for(k=j; k<= j+m-1; k++) {
97 Double_t wsr = alisddetf->GetWeightReal(p);
98 Double_t wsi = alisddetf->GetWeightImag(p);
99 if(direction == -1) wsi = -wsi;
100 Double_t xr = *(real+k+m);
101 Double_t xi = *(imag+k+m);
102 *(real+k+m) = wsr*(*(real+k)-xr) - wsi*(*(imag+k)-xi);
103 *(imag+k+m) = wsr*(*(imag+k)-xi) + wsi*(*(real+k)-xr);
114 for(j=0; j<samples; j++) {
118 for(i1=1; i1<=l; i1++) {
121 p = p + p + j2 - j1 - j1;
124 Double_t xr = *(real+j);
125 Double_t xi = *(imag+j);
126 *(real+j) = *(real+p);
127 *(imag+j) = *(imag+p);
132 if(direction == -1) {
133 for(i=0; i<samples; i++) {
134 *(real+i) /= samples;
135 *(imag+i) /= samples;
137 } // end if direction == -1
140 //______________________________________________________________________
141 AliITSsimulationSDD::AliITSsimulationSDD(){
142 // Default constructor
162 SetPerpendTracksFlag();
167 //______________________________________________________________________
168 AliITSsimulationSDD::AliITSsimulationSDD(AliITSsimulationSDD &source){
169 // Copy constructor to satify Coding roules only.
171 if(this==&source) return;
172 Error("AliITSsimulationSSD","Not allowed to make a copy of "
173 "AliITSsimulationSDD Using default creater instead");
174 AliITSsimulationSDD();
176 //______________________________________________________________________
177 AliITSsimulationSDD& AliITSsimulationSDD::operator=(AliITSsimulationSDD &src){
178 // Assignment operator to satify Coding roules only.
180 if(this==&src) return *this;
181 Error("AliITSsimulationSSD","Not allowed to make a = with "
182 "AliITSsimulationSDD Using default creater instead");
185 //______________________________________________________________________
186 AliITSsimulationSDD::AliITSsimulationSDD(AliITSsegmentation *seg,
187 AliITSresponse *resp){
188 // Standard Constructor
208 Init((AliITSsegmentationSDD*)seg,(AliITSresponseSDD*)resp);
210 //______________________________________________________________________
211 void AliITSsimulationSDD::Init(AliITSsegmentationSDD *seg,
212 AliITSresponseSDD *resp){
213 // Standard Constructor
218 SetPerpendTracksFlag();
223 fpList = new AliITSpList( fSegmentation->Npz(),
224 fScaleSize*fSegmentation->Npx() );
225 fHitSigMap2 = new AliITSMapA2(fSegmentation,fScaleSize,1);
226 fHitNoiMap2 = new AliITSMapA2(fSegmentation,fScaleSize,1);
227 fHitMap2 = fHitSigMap2;
229 fNofMaps = fSegmentation->Npz();
230 fMaxNofSamples = fSegmentation->Npx();
232 Float_t sddLength = fSegmentation->Dx();
233 Float_t sddWidth = fSegmentation->Dz();
236 Float_t anodePitch = fSegmentation->Dpz(dummy);
237 Double_t timeStep = (Double_t)fSegmentation->Dpx(dummy);
238 Float_t driftSpeed = fResponse->DriftSpeed();
240 if(anodePitch*(fNofMaps/2) > sddWidth) {
241 Warning("AliITSsimulationSDD",
242 "Too many anodes %d or too big pitch %f \n",
243 fNofMaps/2,anodePitch);
246 if(timeStep*fMaxNofSamples < sddLength/driftSpeed) {
247 Error("AliITSsimulationSDD",
248 "Time Interval > Allowed Time Interval: exit\n");
252 fElectronics = new AliITSetfSDD(timeStep/fScaleSize,
253 fResponse->Electronics());
255 char opt1[20], opt2[20];
256 fResponse->ParamOptions(opt1,opt2);
258 char *same = strstr(opt1,"same");
263 fNoise.Set(fNofMaps);
264 fBaseline.Set(fNofMaps);
267 const char *kopt=fResponse->ZeroSuppOption();
268 if (strstr(fParam,"file") ) {
271 if (strstr(kopt,"2D")) {
274 Init2D(); // desactivate if param change module by module
275 } else if(strstr(kopt,"1D")) {
278 Init1D(); // desactivate if param change module by module
286 } // end if else strstr
288 Bool_t write = fResponse->OutputOption();
289 if(write && strstr(kopt,"2D")) MakeTreeB();
291 // call here if baseline does not change by module
294 fITS = (AliITS*)gAlice->GetModule("ITS");
295 Int_t size = fNofMaps*fMaxNofSamples;
296 fStream = new AliITSInStream(size);
298 fInZR = new Double_t [fScaleSize*fMaxNofSamples];
299 fInZI = new Double_t [fScaleSize*fMaxNofSamples];
300 fOutZR = new Double_t [fScaleSize*fMaxNofSamples];
301 fOutZI = new Double_t [fScaleSize*fMaxNofSamples];
304 //______________________________________________________________________
305 AliITSsimulationSDD::~AliITSsimulationSDD() {
320 if(fTreeB) delete fTreeB;
321 if(fInZR) delete [] fInZR;
322 if(fInZI) delete [] fInZI;
323 if(fOutZR) delete [] fOutZR;
324 if(fOutZI) delete [] fOutZI;
326 //______________________________________________________________________
327 void AliITSsimulationSDD::InitSimulationModule( Int_t module, Int_t event ) {
328 // create maps to build the lists of tracks for each summable digit
333 //______________________________________________________________________
334 void AliITSsimulationSDD::ClearMaps() {
337 fHitSigMap2->ClearMap();
338 fHitNoiMap2->ClearMap();
340 //______________________________________________________________________
341 void AliITSsimulationSDD::SDigitiseModule( AliITSmodule *mod, Int_t md, Int_t ev){
342 // digitize module using the "slow" detector simulator creating
345 TObjArray *fHits = mod->GetHits();
346 Int_t nhits = fHits->GetEntriesFast();
349 InitSimulationModule( md, ev );
350 HitsToAnalogDigits( mod );
351 ChargeToSignal( kFALSE ); // - Process signal without add noise
352 fHitMap2 = fHitNoiMap2; // - Swap to noise map
353 ChargeToSignal( kTRUE ); // - Process only noise
354 fHitMap2 = fHitSigMap2; // - Return to signal map
358 //______________________________________________________________________
359 Bool_t AliITSsimulationSDD::AddSDigitsToModule( TClonesArray *pItemArray, Int_t mask ) {
360 // Add Summable digits to module maps.
361 Int_t nItems = pItemArray->GetEntries();
362 Double_t maxadc = fResponse->MaxAdc();
365 // cout << "Adding "<< nItems <<" SDigits to module " << fModule << endl;
366 for( Int_t i=0; i<nItems; i++ ) {
367 AliITSpListItem * pItem = (AliITSpListItem *)(pItemArray->At( i ));
368 if( pItem->GetModule() != fModule ) {
369 Error( "AliITSsimulationSDD",
370 "Error reading, SDigits module %d != current module %d: exit\n",
371 pItem->GetModule(), fModule );
375 if(pItem->GetSignal()>0.0 ) sig = kTRUE;
376 fpList->AddItemTo( mask, pItem ); // Add SignalAfterElect + noise
377 AliITSpListItem * pItem2 = fpList->GetpListItem( pItem->GetIndex() );
378 Double_t sigAE = pItem2->GetSignalAfterElect();
379 if( sigAE >= maxadc ) sigAE = maxadc-1; // avoid overflow signal
382 fpList->GetMapIndex( pItem->GetIndex(), ia, it );
383 fHitMap2->SetHit( ia, it, sigAE );
387 //______________________________________________________________________
388 void AliITSsimulationSDD::FinishSDigitiseModule() {
389 // digitize module using the "slow" detector simulator from
390 // the sum of summable digits.
394 //______________________________________________________________________
395 void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
396 // create maps to build the lists of tracks for each digit
398 TObjArray *fHits = mod->GetHits();
399 Int_t nhits = fHits->GetEntriesFast();
401 InitSimulationModule( md, ev );
403 if( !nhits && fCheckNoise ) {
404 ChargeToSignal( kTRUE ); // process noise
411 HitsToAnalogDigits( mod );
412 ChargeToSignal( kTRUE ); // process signal + noise
414 for( Int_t i=0; i<fNofMaps; i++ ) {
415 for( Int_t j=0; j<fMaxNofSamples; j++ ) {
416 Int_t jdx = j*fScaleSize;
417 Int_t index = fpList->GetHitIndex( i, j );
418 AliITSpListItem pItemTmp2( fModule, index, 0. );
419 // put the fScaleSize analog digits in only one
420 for( Int_t ik=0; ik<fScaleSize; ik++ ) {
421 AliITSpListItem *pItemTmp = fpList->GetpListItem( i, jdx+ik );
422 if( pItemTmp == 0 ) continue;
423 pItemTmp2.Add( pItemTmp );
425 fpList->DeleteHit( i, j );
426 fpList->AddItemTo( 0, &pItemTmp2 );
433 //______________________________________________________________________
434 void AliITSsimulationSDD::FinishDigits() {
435 // introduce the electronics effects and do zero-suppression if required
438 if( fCrosstalkFlag ) ApplyCrosstalk();
440 const char *kopt = fResponse->ZeroSuppOption();
441 ZeroSuppression( kopt );
443 //______________________________________________________________________
444 void AliITSsimulationSDD::HitsToAnalogDigits( AliITSmodule *mod ) {
445 // create maps to build the lists of tracks for each digit
447 TObjArray *fHits = mod->GetHits();
448 Int_t nhits = fHits->GetEntriesFast();
449 // Int_t arg[6] = {0,0,0,0,0,0};
451 Int_t nofAnodes = fNofMaps/2;
452 Float_t sddLength = fSegmentation->Dx();
453 Float_t sddWidth = fSegmentation->Dz();
454 Float_t anodePitch = fSegmentation->Dpz(dummy);
455 Float_t timeStep = fSegmentation->Dpx(dummy);
456 Float_t driftSpeed = fResponse->DriftSpeed();
457 Float_t maxadc = fResponse->MaxAdc();
458 Float_t topValue = fResponse->DynamicRange();
459 Float_t cHloss = fResponse->ChargeLoss();
460 Float_t norm = maxadc/topValue;
461 Float_t dfCoeff, s1; fResponse->DiffCoeff(dfCoeff,s1); // Signal 2d Shape
462 Double_t eVpairs = 3.6; // electron pair energy eV.
463 Float_t nsigma = fResponse->NSigmaIntegration(); //
464 Int_t nlookups = fResponse->GausNLookUp(); //
465 Float_t jitter = ((AliITSresponseSDD*)fResponse)->JitterError(); //
467 // Piergiorgio's part (apart for few variables which I made float
468 // when i thought that can be done
469 // Fill detector maps with GEANT hits
470 // loop over hits in the module
472 const Float_t kconv = 1.0e+6; // GeV->KeV
474 Int_t hitDetector; // detector number (lay,lad,hitDetector)
475 Int_t iWing; // which detector wing/side.
476 Int_t detector; // 2*(detector-1)+iWing
477 Int_t ii,kk,ka,kt; // loop indexs
478 Int_t ia,it,index; // sub-pixel integration indexies
479 Int_t iAnode; // anode number.
480 Int_t timeSample; // time buckett.
481 Int_t anodeWindow; // anode direction charge integration width
482 Int_t timeWindow; // time direction charge integration width
483 Int_t jamin,jamax; // anode charge integration window
484 Int_t jtmin,jtmax; // time charge integration window
485 Int_t ndiv; // Anode window division factor.
486 Int_t nsplit; // the number of splits in anode and time windows==1.
487 Int_t nOfSplits; // number of times track length is split into
488 Float_t nOfSplitsF; // Floating point version of nOfSplits.
489 Float_t kkF; // Floating point version of loop index kk.
490 Float_t pathInSDD; // Track length in SDD.
491 Float_t drPath; // average position of track in detector. in microns
492 Float_t drTime; // Drift time
493 Float_t nmul; // drift time window multiplication factor.
494 Float_t avDrft; // x position of path length segment in cm.
495 Float_t avAnode; // Anode for path length segment in Anode number (float)
496 Float_t xAnode; // Floating point anode number.
497 Float_t driftPath; // avDrft in microns.
498 Float_t width; // width of signal at anodes.
499 Double_t depEnergy; // Energy deposited in this GEANT step.
500 Double_t xL[3],dxL[3]; // local hit coordinates and diff.
501 Double_t sigA; // sigma of signal at anode.
502 Double_t sigT; // sigma in time/drift direction for track segment
503 Double_t aStep,aConst; // sub-pixel size and offset anode
504 Double_t tStep,tConst; // sub-pixel size and offset time
505 Double_t amplitude; // signal amplitude for track segment in nanoAmpere
506 Double_t chargeloss; // charge loss for track segment.
507 Double_t anodeAmplitude; // signal amplitude in anode direction
508 Double_t aExpo; // exponent of Gaussian anode direction
509 Double_t timeAmplitude; // signal amplitude in time direction
510 Double_t tExpo; // exponent of Gaussian time direction
511 // Double_t tof; // Time of flight in ns of this step.
513 for(ii=0; ii<nhits; ii++) {
514 if(!mod->LineSegmentL(ii,xL[0],dxL[0],xL[1],dxL[1],xL[2],dxL[2],
515 depEnergy,itrack)) continue;
516 xL[0] += 0.0001*gRandom->Gaus( 0, jitter ); //
518 hitDetector = mod->GetDet();
519 //tof = 1.E+09*(mod->GetHit(ii)->GetTOF()); // tof in ns.
520 //if(tof>sddLength/driftSpeed) continue; // hit happed too late.
522 // scale path to simulate a perpendicular track
523 // continue if the particle did not lose energy
524 // passing through detector
526 Warning("HitsToAnalogDigits",
527 "fTrack = %d hit=%d module=%d This particle has"
528 " passed without losing energy!",
529 itrack,ii,mod->GetIndex());
531 } // end if !depEnergy
533 pathInSDD = TMath::Sqrt(dxL[0]*dxL[0]+dxL[1]*dxL[1]+dxL[2]*dxL[2]);
535 if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
536 drPath = 10000.*(dxL[0]+2.*xL[0])*0.5;
537 if(drPath < 0) drPath = -drPath;
538 drPath = sddLength-drPath;
540 Warning("HitsToAnalogDigits",
541 "negative drift path drPath=%e sddLength=%e dxL[0]=%e "
543 drPath,sddLength,dxL[0],xL[0]);
545 } // end if drPath < 0
547 // Compute number of segments to brake step path into
548 drTime = drPath/driftSpeed; // Drift Time
549 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);// Sigma along the anodes
550 // calcuate the number of time the path length should be split into.
551 nOfSplits = (Int_t) (1. + 10000.*pathInSDD/sigA);
552 if(fFlag) nOfSplits = 1;
554 // loop over path segments, init. some variables.
555 depEnergy /= nOfSplits;
556 nOfSplitsF = (Float_t) nOfSplits;
557 for(kk=0;kk<nOfSplits;kk++) { // loop over path segments
558 kkF = (Float_t) kk + 0.5;
559 avDrft = xL[0]+dxL[0]*kkF/nOfSplitsF;
560 avAnode = xL[2]+dxL[2]*kkF/nOfSplitsF;
561 driftPath = 10000.*avDrft;
563 iWing = 2; // Assume wing is 2
564 if(driftPath < 0) { // if wing is not 2 it is 1.
566 driftPath = -driftPath;
567 } // end if driftPath < 0
568 driftPath = sddLength-driftPath;
569 detector = 2*(hitDetector-1) + iWing;
571 Warning("HitsToAnalogDigits","negative drift path "
572 "driftPath=%e sddLength=%e avDrft=%e dxL[0]=%e "
573 "xL[0]=%e",driftPath,sddLength,avDrft,dxL[0],xL[0]);
575 } // end if driftPath < 0
578 drTime = driftPath/driftSpeed; // drift time for segment.
579 timeSample = (Int_t) (fScaleSize*drTime/timeStep + 1);
580 // compute time Sample including tof information. The tof only
581 // effects the time of the signal is recoreded and not the
583 // timeSample = (Int_t) (fScaleSize*(drTime+tof)/timeStep + 1);
584 if(timeSample > fScaleSize*fMaxNofSamples) {
585 Warning("HitsToAnalogDigits","Wrong Time Sample: %e",
588 } // end if timeSample > fScaleSize*fMaxNoofSamples
591 xAnode = 10000.*(avAnode)/anodePitch + nofAnodes/2; // +1?
592 if(xAnode*anodePitch > sddWidth || xAnode*anodePitch < 0.)
593 Warning("HitsToAnalogDigits",
594 "Exceedubg sddWidth=%e Z = %e",
595 sddWidth,xAnode*anodePitch);
596 iAnode = (Int_t) (1.+xAnode); // xAnode?
597 if(iAnode < 1 || iAnode > nofAnodes) {
598 Warning("HitToAnalogDigits","Wrong iAnode: 1<%d>%d",
601 } // end if iAnode < 1 || iAnode > nofAnodes
603 // store straight away the particle position in the array
604 // of particles and take idhit=ii only when part is entering (this
605 // requires FillModules() in the macro for analysis) :
607 // Sigma along the anodes for track segment.
608 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);
609 sigT = sigA/driftSpeed;
610 // Peak amplitude in nanoAmpere
611 amplitude = fScaleSize*160.*depEnergy/
612 (timeStep*eVpairs*2.*acos(-1.)*sigT*sigA);
613 amplitude *= timeStep/25.; // WARNING!!!!! Amplitude scaling to
614 // account for clock variations
615 // (reference value: 40 MHz)
616 chargeloss = 1.-cHloss*driftPath/1000;
617 amplitude *= chargeloss;
618 width = 2.*nsigma/(nlookups-1);
626 } // end if drTime > 1200.
628 nsplit = 4; // hard-wired //nsplit=4;nsplit = (nsplit+1)/2*2;
629 // Sub-pixel size see computation of aExpo and tExpo.
630 aStep = anodePitch/(nsplit*fScaleSize*sigA);
631 aConst = xAnode*anodePitch/sigA;
632 tStep = timeStep/(nsplit*fScaleSize*sigT);
633 tConst = drTime/sigT;
634 // Define SDD window corresponding to the hit
635 anodeWindow = (Int_t)(fScaleSize*nsigma*sigA/anodePitch+1);
636 timeWindow = (Int_t) (fScaleSize*nsigma*sigT/timeStep+1.);
637 jamin = (iAnode - anodeWindow/ndiv - 1)*fScaleSize*nsplit +1;
638 jamax = (iAnode + anodeWindow/ndiv)*fScaleSize*nsplit;
639 if(jamin <= 0) jamin = 1;
640 if(jamax > fScaleSize*nofAnodes*nsplit)
641 jamax = fScaleSize*nofAnodes*nsplit;
642 // jtmin and jtmax are Hard-wired
643 jtmin = (Int_t)(timeSample-timeWindow*nmul-1)*nsplit+1;
644 jtmax = (Int_t)(timeSample+timeWindow*nmul)*nsplit;
645 if(jtmin <= 0) jtmin = 1;
646 if(jtmax > fScaleSize*fMaxNofSamples*nsplit)
647 jtmax = fScaleSize*fMaxNofSamples*nsplit;
648 // Spread the charge in the anode-time window
649 for(ka=jamin; ka <=jamax; ka++) {
650 ia = (ka-1)/(fScaleSize*nsplit) + 1;
652 Warning("HitsToAnalogDigits","ia < 1: ");
655 if(ia > nofAnodes) ia = nofAnodes;
656 aExpo = (aStep*(ka-0.5)-aConst);
657 if(TMath::Abs(aExpo) > nsigma) anodeAmplitude = 0.;
659 dummy = (Int_t) ((aExpo+nsigma)/width);
660 anodeAmplitude = amplitude*fResponse->GausLookUp(dummy);
661 } // end if TMath::Abs(aEspo) > nsigma
662 // index starts from 0
663 index = ((detector+1)%2)*nofAnodes+ia-1;
664 if(anodeAmplitude) for(kt=jtmin; kt<=jtmax; kt++) {
665 it = (kt-1)/nsplit+1; // it starts from 1
667 Warning("HitsToAnalogDigits","it < 1:");
670 if(it>fScaleSize*fMaxNofSamples)
671 it = fScaleSize*fMaxNofSamples;
672 tExpo = (tStep*(kt-0.5)-tConst);
673 if(TMath::Abs(tExpo) > nsigma) timeAmplitude = 0.;
675 dummy = (Int_t) ((tExpo+nsigma)/width);
676 timeAmplitude = anodeAmplitude*
677 fResponse->GausLookUp(dummy);
678 } // end if TMath::Abs(tExpo) > nsigma
679 // build the list of Sdigits for this module
682 // arg[2] = itrack; // track number
683 // arg[3] = ii-1; // hit number.
684 timeAmplitude *= norm;
686 // ListOfFiredCells(arg,timeAmplitude,alst,padr);
687 Double_t charge = timeAmplitude;
688 charge += fHitMap2->GetSignal(index,it-1);
689 fHitMap2->SetHit(index, it-1, charge);
690 fpList->AddSignal(index,it-1,itrack,ii-1,
691 mod->GetIndex(),timeAmplitude);
692 } // end if anodeAmplitude and loop over time in window
693 } // loop over anodes in window
694 } // end loop over "sub-hits"
695 } // end loop over hits
699 //______________________________________________________________________
700 void AliITSsimulationSDD::ListOfFiredCells(Int_t *arg,Double_t timeAmplitude,
701 TObjArray *alist,TClonesArray *padr){
702 // Returns the list of "fired" cells.
704 Int_t index = arg[0];
706 Int_t idtrack = arg[2];
707 Int_t idhit = arg[3];
708 Int_t counter = arg[4];
709 Int_t countadr = arg[5];
710 Double_t charge = timeAmplitude;
711 charge += fHitMap2->GetSignal(index,ik-1);
712 fHitMap2->SetHit(index, ik-1, charge);
715 Int_t it = (Int_t)((ik-1)/fScaleSize);
718 digits[2] = (Int_t)timeAmplitude;
720 if (idtrack >= 0) phys = (Float_t)timeAmplitude;
723 Double_t cellcharge = 0.;
724 AliITSTransientDigit* pdigit;
725 // build the list of fired cells and update the info
726 if (!fHitMap1->TestHit(index, it)) {
727 new((*padr)[countadr++]) TVector(3);
728 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
729 trinfo(0) = (Float_t)idtrack;
730 trinfo(1) = (Float_t)idhit;
731 trinfo(2) = (Float_t)timeAmplitude;
733 alist->AddAtAndExpand(new AliITSTransientDigit(phys,digits),counter);
734 fHitMap1->SetHit(index, it, counter);
736 pdigit=(AliITSTransientDigit*)alist->At(alist->GetLast());
738 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
739 trlist->Add(&trinfo);
741 pdigit = (AliITSTransientDigit*) fHitMap1->GetHit(index, it);
742 for(Int_t kk=0;kk<fScaleSize;kk++) {
743 cellcharge += fHitMap2->GetSignal(index,fScaleSize*it+kk);
746 (*pdigit).fSignal = (Int_t)cellcharge;
747 (*pdigit).fPhysics += phys;
748 // update list of tracks
749 TObjArray* trlist = (TObjArray*)pdigit->TrackList();
750 Int_t lastentry = trlist->GetLast();
751 TVector *ptrkp = (TVector*)trlist->At(lastentry);
752 TVector &trinfo = *ptrkp;
753 Int_t lasttrack = Int_t(trinfo(0));
754 Float_t lastcharge=(trinfo(2));
755 if (lasttrack==idtrack ) {
756 lastcharge += (Float_t)timeAmplitude;
757 trlist->RemoveAt(lastentry);
758 trinfo(0) = lasttrack;
760 trinfo(2) = lastcharge;
761 trlist->AddAt(&trinfo,lastentry);
763 new((*padr)[countadr++]) TVector(3);
764 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
765 trinfo(0) = (Float_t)idtrack;
766 trinfo(1) = (Float_t)idhit;
767 trinfo(2) = (Float_t)timeAmplitude;
768 trlist->Add(&trinfo);
769 } // end if lasttrack==idtrack
772 // check the track list - debugging
773 Int_t trk[20], htrk[20];
775 Int_t nptracks = trlist->GetEntriesFast();
778 for (tr=0;tr<nptracks;tr++) {
779 TVector *pptrkp = (TVector*)trlist->At(tr);
780 TVector &pptrk = *pptrkp;
781 trk[tr] = Int_t(pptrk(0));
782 htrk[tr] = Int_t(pptrk(1));
783 chtrk[tr] = (pptrk(2));
784 cout << "nptracks "<<nptracks << endl;
790 // update counter and countadr for next call.
796 //____________________________________________
797 void AliITSsimulationSDD::AddDigit( Int_t i, Int_t j, Int_t signal ) {
799 Int_t digits[3], tracks[3], hits[3];
800 Float_t phys, charges[3];
802 if( fResponse->Do10to8() ) signal = Convert8to10( signal );
807 AliITSpListItem *pItem = fpList->GetpListItem( i, j );
810 for( Int_t l=0; l<3; l++ ) {
816 Int_t idtrack = pItem->GetTrack( 0 );
817 if( idtrack >= 0 ) phys = pItem->GetSignal();
820 for( Int_t l=0; l<3; l++ ) {
821 tracks[l] = pItem->GetTrack( l );
822 hits[l] = pItem->GetHit( l );
823 charges[l] = pItem->GetSignal( l );
827 fITS->AddSimDigit( 1, phys, digits, tracks, hits, charges );
831 //____________________________________________
832 void AliITSsimulationSDD::AddDigit(Int_t i, Int_t j, Int_t signal){
834 // tag with -1 signals coming from background tracks
835 // tag with -2 signals coming from pure electronic noise
837 Int_t digits[3], tracks[3], hits[3];
838 Float_t phys, charges[3];
840 Int_t trk[20], htrk[20];
843 Bool_t do10to8=fResponse->Do10to8();
845 if(do10to8) signal=Convert8to10(signal);
846 AliITSTransientDigit *obj = (AliITSTransientDigit*)fHitMap1->GetHit(i,j);
858 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
861 TObjArray* trlist=(TObjArray*)obj->TrackList();
862 Int_t nptracks=trlist->GetEntriesFast();
864 Warning("AddDigit","nptracks=%d > 20 nptracks set to 20",nptracks);
866 } // end if nptracks > 20
868 for (tr=0;tr<nptracks;tr++) {
869 TVector &pp =*((TVector*)trlist->At(tr));
870 trk[tr]=Int_t(pp(0));
871 htrk[tr]=Int_t(pp(1));
875 SortTracks(trk,chtrk,htrk,nptracks);
876 } // end if nptracks > 1
879 for (i=0; i<nptracks; i++) {
884 for (i=nptracks; i<3; i++) {
890 for (i=0; i<3; i++) {
895 } // end if/else nptracks < 3
897 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
899 } // end if/else !obj
903 //______________________________________________________________________
904 void AliITSsimulationSDD::SortTracks(Int_t *tracks,Float_t *charges,
905 Int_t *hits,Int_t ntr){
906 // Sort the list of tracks contributing to a given digit
907 // Only the 3 most significant tracks are acctually sorted
908 // Loop over signals, only 3 times
912 Int_t idx[3] = {-3,-3,-3};
913 Float_t jch[3] = {-3,-3,-3};
914 Int_t jtr[3] = {-3,-3,-3};
915 Int_t jhit[3] = {-3,-3,-3};
918 if (ntr<3) imax = ntr;
924 if((i == 1 && j == idx[i-1] )
925 ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
926 if(charges[j] > qmax) {
929 } // end if charges[j]>qmax
933 jch[i] = charges[jmax];
934 jtr[i] = tracks[jmax];
935 jhit[i] = hits[jmax];
948 } // end if jtr[i] == -3
952 //______________________________________________________________________
953 void AliITSsimulationSDD::ChargeToSignal(Bool_t bAddNoise) {
954 // add baseline, noise, electronics and ADC saturation effects
956 char opt1[20], opt2[20];
957 fResponse->ParamOptions(opt1,opt2);
958 char *read = strstr(opt1,"file");
959 Float_t baseline, noise;
962 static Bool_t readfile=kTRUE;
963 //read baseline and noise from file
964 if (readfile) ReadBaseline();
966 } else fResponse->GetNoiseParam(noise,baseline);
970 Float_t maxadc = fResponse->MaxAdc();
972 for (i=0;i<fNofMaps;i++) {
973 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
974 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
975 fInZR[k] = fHitMap2->GetSignal(i,k);
977 contrib = (baseline + noise*gRandom->Gaus());
981 for(k=0; k<fMaxNofSamples; k++) {
982 Double_t newcont = 0.;
983 Double_t maxcont = 0.;
984 for(kk=0;kk<fScaleSize;kk++) {
985 newcont = fInZR[fScaleSize*k+kk];
986 if(newcont > maxcont) maxcont = newcont;
989 if (newcont >= maxadc) newcont = maxadc -1;
990 if(newcont >= baseline){
991 Warning("","newcont=%d>=baseline=%d",newcont,baseline);
994 fHitMap2->SetHit(i,k,newcont);
996 } // end for i loop over anodes
1000 for (i=0;i<fNofMaps;i++) {
1001 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
1002 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
1003 fInZR[k] = fHitMap2->GetSignal(i,k);
1005 contrib = (baseline + noise*gRandom->Gaus());
1006 fInZR[k] += contrib;
1010 FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1);
1011 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
1012 Double_t rw = fElectronics->GetTraFunReal(k);
1013 Double_t iw = fElectronics->GetTraFunImag(k);
1014 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
1015 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
1017 FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1);
1018 for(k=0; k<fMaxNofSamples; k++) {
1019 Double_t newcont1 = 0.;
1020 Double_t maxcont1 = 0.;
1021 for(kk=0;kk<fScaleSize;kk++) {
1022 newcont1 = fOutZR[fScaleSize*k+kk];
1023 if(newcont1 > maxcont1) maxcont1 = newcont1;
1025 newcont1 = maxcont1;
1026 if (newcont1 >= maxadc) newcont1 = maxadc -1;
1027 fHitMap2->SetHit(i,k,newcont1);
1029 } // end for i loop over anodes
1032 //____________________________________________________________________
1033 void AliITSsimulationSDD::ApplyDeadChannels() {
1034 // Set dead channel signal to zero
1035 AliITSresponseSDD * response = (AliITSresponseSDD *)fResponse;
1038 if( response->GetDeadModules() == 0 &&
1039 response->GetDeadChips() == 0 &&
1040 response->GetDeadChannels() == 0 )
1043 static AliITS *iTS = (AliITS*)gAlice->GetModule( "ITS" );
1045 Int_t fMaxNofSamples = fSegmentation->Npx();
1046 AliITSgeom *geom = iTS->GetITSgeom();
1047 Int_t firstSDDMod = geom->GetStartDet( 1 );
1049 for( Int_t j=0; j<2; j++ ) {
1050 Int_t mod = (fModule-firstSDDMod)*2 + j;
1051 for( Int_t u=0; u<response->Chips(); u++ )
1052 for( Int_t v=0; v<response->Channels(); v++ ) {
1053 Float_t Gain = response->Gain( mod, u, v );
1054 for( Int_t k=0; k<fMaxNofSamples; k++ ) {
1055 Int_t i = j*response->Chips()*response->Channels() +
1056 u*response->Channels() +
1058 Double_t signal = Gain * fHitMap2->GetSignal( i, k );
1059 fHitMap2->SetHit( i, k, signal ); ///
1064 //______________________________________________________________________
1065 void AliITSsimulationSDD::ApplyCrosstalk() {
1066 // function add the crosstalk effect to signal
1067 // temporal function, should be checked...!!!
1069 Int_t fNofMaps = fSegmentation->Npz();
1070 Int_t fMaxNofSamples = fSegmentation->Npx();
1072 // create and inizialice crosstalk map
1073 Float_t* ctk = new Float_t[fNofMaps*fMaxNofSamples+1];
1075 Error( "ApplyCrosstalk", "no memory for temporal map: exit \n" );
1078 memset( ctk, 0, sizeof(Float_t)*(fNofMaps*fMaxNofSamples+1) );
1080 Float_t noise, baseline;
1081 fResponse->GetNoiseParam( noise, baseline );
1083 for( Int_t z=0; z<fNofMaps; z++ ) {
1089 for( Int_t l=0; l<fMaxNofSamples; l++ ) {
1090 Float_t fadc = (Float_t)fHitMap2->GetSignal( z, l );
1091 if( fadc > baseline ) {
1092 if( on == kFALSE && l<fMaxNofSamples-4 ) {
1093 Float_t fadc1 = (Float_t)fHitMap2->GetSignal( z, l+1 );
1094 if( fadc1 < fadc ) continue;
1101 else { // end fadc > baseline
1105 // make smooth derivative
1106 Float_t* dev = new Float_t[fMaxNofSamples+1];
1107 memset( dev, 0, sizeof(Float_t)*(fMaxNofSamples+1) );
1109 Error( "ApplyCrosstalk",
1110 "no memory for temporal array: exit \n" );
1113 for( Int_t i=tstart; i<tstop; i++ ) {
1114 if( i > 2 && i < fMaxNofSamples-2 )
1115 dev[i] = -0.2*fHitMap2->GetSignal( z,i-2 )
1116 -0.1*fHitMap2->GetSignal( z,i-1 )
1117 +0.1*fHitMap2->GetSignal( z,i+1 )
1118 +0.2*fHitMap2->GetSignal( z,i+2 );
1121 // add crosstalk contribution to neibourg anodes
1122 for( Int_t i=tstart; i<tstop; i++ ) {
1123 Int_t anode = z - 1;
1124 Int_t i1 = (Int_t)((i-tstart)*.61+tstart+0.5); //
1125 Float_t ctktmp = -dev[i1] * 0.25;
1127 ctk[anode*fMaxNofSamples+i] += ctktmp;
1130 if( anode < fNofMaps ) {
1131 ctk[anode*fMaxNofSamples+i] += ctktmp;
1136 } // if( nTsteps > 2 )
1138 } // if( on == kTRUE )
1143 for( Int_t a=0; a<fNofMaps; a++ )
1144 for( Int_t t=0; t<fMaxNofSamples; t++ ) {
1145 Float_t signal = fHitMap2->GetSignal( a, t ) + ctk[a*fMaxNofSamples+t];
1146 fHitMap2->SetHit( a, t, signal );
1151 //______________________________________________________________________
1152 void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline,
1154 // Returns the Baseline for a particular anode.
1155 baseline = fBaseline[i];
1158 //______________________________________________________________________
1159 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl,
1161 // Returns the compression alogirthm parameters
1162 Int_t size = fD.GetSize();
1164 db=fD[i]; tl=fT1[i]; th=fT2[i];
1166 if (size <= 2 && i>=fNofMaps/2) {
1167 db=fD[1]; tl=fT1[1]; th=fT2[1];
1169 db=fD[0]; tl=fT1[0]; th=fT2[0];
1170 } // end if size <=2 && i>=fNofMaps/2
1173 //______________________________________________________________________
1174 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){
1175 // returns the compression alogirthm parameters
1176 Int_t size = fD.GetSize();
1179 db=fD[i]; tl=fT1[i];
1181 if (size <= 2 && i>=fNofMaps/2) {
1182 db=fD[1]; tl=fT1[1];
1184 db=fD[0]; tl=fT1[0];
1185 } // end if size <=2 && i>=fNofMaps/2
1186 } // end if size > 2
1188 //______________________________________________________________________
1189 void AliITSsimulationSDD::SetCompressParam(){
1190 // Sets the compression alogirthm parameters
1193 fResponse->GiveCompressParam(cp);
1194 for (i=0; i<2; i++) {
1201 //______________________________________________________________________
1202 void AliITSsimulationSDD::ReadBaseline(){
1203 // read baseline and noise from file - either a .root file and in this
1204 // case data should be organised in a tree with one entry for each
1205 // module => reading should be done accordingly
1206 // or a classic file and do smth. like this:
1207 // Read baselines and noise for SDD
1211 char input[100], base[100], param[100];
1214 fResponse->Filenames(input,base,param);
1217 filtmp = gSystem->ExpandPathName(fFileName.Data());
1218 FILE *bline = fopen(filtmp,"r");
1222 while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) {
1224 Error("ReadBaseline","Anode number not in increasing order!",
1227 } // end if pos != na+1
1233 Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",filtmp);
1240 //______________________________________________________________________
1241 Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) const {
1242 // To the 10 to 8 bit lossive compression.
1243 // code from Davide C. and Albert W.
1245 if (signal < 128) return signal;
1246 if (signal < 256) return (128+((signal-128)>>1));
1247 if (signal < 512) return (192+((signal-256)>>3));
1248 if (signal < 1024) return (224+((signal-512)>>4));
1251 //______________________________________________________________________
1252 Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) const {
1253 // Undo the lossive 10 to 8 bit compression.
1254 // code from Davide C. and Albert W.
1255 if (signal < 0 || signal > 255) {
1256 Warning("Convert8to10","out of range signal=%d",signal);
1258 } // end if signal <0 || signal >255
1260 if (signal < 128) return signal;
1262 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
1263 else return (128+((signal-128)<<1)+1);
1264 } // end if signal < 192
1266 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
1267 else return (256+((signal-192)<<3)+4);
1268 } // end if signal < 224
1269 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
1270 return (512+((signal-224)<<4)+8);
1274 //______________________________________________________________________
1275 AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){
1276 //Return the correct map.
1278 return ((i==0)? fHitMap1 : fHitMap2);
1281 //______________________________________________________________________
1282 void AliITSsimulationSDD::ZeroSuppression(const char *option) {
1283 // perform the zero suppresion
1285 if (strstr(option,"2D")) {
1286 //Init2D(); // activate if param change module by module
1288 } else if (strstr(option,"1D")) {
1289 //Init1D(); // activate if param change module by module
1291 } else StoreAllDigits();
1293 //______________________________________________________________________
1294 void AliITSsimulationSDD::Init2D(){
1295 // read in and prepare arrays: fD, fT1, fT2
1296 // savemu[nanodes], savesigma[nanodes]
1297 // read baseline and noise from file - either a .root file and in this
1298 // case data should be organised in a tree with one entry for each
1299 // module => reading should be done accordingly
1300 // or a classic file and do smth. like this ( code from Davide C. and
1302 // Read 2D zero-suppression parameters for SDD
1304 if (!strstr(fParam,"file")) return;
1306 Int_t na,pos,tempTh;
1308 Float_t *savemu = new Float_t [fNofMaps];
1309 Float_t *savesigma = new Float_t [fNofMaps];
1310 char input[100],basel[100],par[100];
1312 Int_t minval = fResponse->MinVal();
1314 fResponse->Filenames(input,basel,par);
1317 filtmp = gSystem->ExpandPathName(fFileName.Data());
1318 FILE *param = fopen(filtmp,"r");
1322 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1324 Error("Init2D","Anode number not in increasing order!",filtmp);
1326 } // end if pos != na+1
1328 savesigma[na] = sigma;
1329 if ((2.*sigma) < mu) {
1330 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1333 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1334 if (tempTh < 0) tempTh=0;
1336 tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval;
1337 if (tempTh < 0) tempTh=0;
1342 Error("Init2D","THE FILE %s DOES NOT EXIST !",filtmp);
1349 delete [] savesigma;
1351 //______________________________________________________________________
1352 void AliITSsimulationSDD::Compress2D(){
1353 // simple ITS cluster finder -- online zero-suppression conditions
1356 Int_t minval = fResponse->MinVal();
1357 Bool_t write = fResponse->OutputOption();
1358 Bool_t do10to8 = fResponse->Do10to8();
1359 Int_t nz, nl, nh, low, i, j;
1361 for (i=0; i<fNofMaps; i++) {
1362 CompressionParam(i,db,tl,th);
1367 for (j=0; j<fMaxNofSamples; j++) {
1368 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1369 signal -= db; // if baseline eq. is done here
1370 if (signal <= 0) {nz++; continue;}
1371 if ((signal - tl) < minval) low++;
1372 if ((signal - th) >= minval) {
1375 FindCluster(i,j,signal,minval,cond);
1377 ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)){
1378 if(do10to8) signal = Convert10to8(signal);
1379 AddDigit(i,j,signal);
1380 } // end if cond&&j&&()
1381 } else if ((signal - tl) >= minval) nl++;
1382 } // end for j loop time samples
1383 if (write) TreeB()->Fill(nz,nl,nh,low,i+1);
1384 } //end for i loop anodes
1388 sprintf(hname,"TNtuple%d_%d",fModule,fEvent);
1389 TreeB()->Write(hname);
1394 //______________________________________________________________________
1395 void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal,
1396 Int_t minval,Bool_t &cond){
1397 // Find clusters according to the online 2D zero-suppression algorithm
1398 Bool_t do10to8 = fResponse->Do10to8();
1399 Bool_t high = kFALSE;
1401 fHitMap2->FlagHit(i,j);
1403 // check the online zero-suppression conditions
1405 const Int_t kMaxNeighbours = 4;
1408 Int_t xList[kMaxNeighbours], yList[kMaxNeighbours];
1409 fSegmentation->Neighbours(i,j,&nn,xList,yList);
1411 for (in=0; in<nn; in++) {
1414 if (fHitMap2->TestHit(ix,iy)==kUnused) {
1415 CompressionParam(ix,dbx,tlx,thx);
1416 Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy));
1417 qn -= dbx; // if baseline eq. is done here
1418 if ((qn-tlx) < minval) {
1419 fHitMap2->FlagHit(ix,iy);
1422 if ((qn - thx) >= minval) high=kTRUE;
1424 if(do10to8) signal = Convert10to8(signal);
1425 AddDigit(i,j,signal);
1427 if(do10to8) qns = Convert10to8(qn);
1429 if (!high) AddDigit(ix,iy,qns);
1431 if(!high) fHitMap2->FlagHit(ix,iy);
1432 } // end if qn-tlx < minval
1434 } // end for in loop over neighbours
1436 //______________________________________________________________________
1437 void AliITSsimulationSDD::Init1D(){
1438 // this is just a copy-paste of input taken from 2D algo
1439 // Torino people should give input
1440 // Read 1D zero-suppression parameters for SDD
1442 if (!strstr(fParam,"file")) return;
1444 Int_t na,pos,tempTh;
1446 Float_t *savemu = new Float_t [fNofMaps];
1447 Float_t *savesigma = new Float_t [fNofMaps];
1448 char input[100],basel[100],par[100];
1450 Int_t minval = fResponse->MinVal();
1452 fResponse->Filenames(input,basel,par);
1455 // set first the disable and tol param
1458 filtmp = gSystem->ExpandPathName(fFileName.Data());
1459 FILE *param = fopen(filtmp,"r");
1463 fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]);
1464 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1466 Error("Init1D","Anode number not in increasing order!",filtmp);
1468 } // end if pos != na+1
1470 savesigma[na]=sigma;
1471 if ((2.*sigma) < mu) {
1472 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1475 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1476 if (tempTh < 0) tempTh=0;
1481 Error("Init1D","THE FILE %s DOES NOT EXIST !",filtmp);
1488 delete [] savesigma;
1490 //______________________________________________________________________
1491 void AliITSsimulationSDD::Compress1D(){
1492 // 1D zero-suppression algorithm (from Gianluca A.)
1493 Int_t dis,tol,thres,decr,diff;
1494 UChar_t *str=fStream->Stream();
1496 Bool_t do10to8=fResponse->Do10to8();
1500 for (k=0; k<2; k++) {
1503 for (i=0; i<fNofMaps/2; i++) {
1504 Bool_t firstSignal=kTRUE;
1505 Int_t idx=i+k*fNofMaps/2;
1506 CompressionParam(idx,decr,thres);
1507 for (j=0; j<fMaxNofSamples; j++) {
1508 Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j));
1509 signal -= decr; // if baseline eq.
1510 if(do10to8) signal = Convert10to8(signal);
1511 if (signal <= thres) {
1515 // write diff in the buffer for HuffT
1516 str[counter]=(UChar_t)diff;
1519 } // end if signal <= thres
1521 if (diff > 127) diff=127;
1522 if (diff < -128) diff=-128;
1524 // tol has changed to 8 possible cases ? - one can write
1525 // this if(TMath::Abs(diff)<tol) ... else ...
1526 if(TMath::Abs(diff)<tol) diff=0;
1527 // or keep it as it was before
1529 if (tol==1 && (diff >= -2 && diff <= 1)) diff=0;
1530 if (tol==2 && (diff >= -4 && diff <= 3)) diff=0;
1531 if (tol==3 && (diff >= -16 && diff <= 15)) diff=0;
1533 AddDigit(idx,j,last+diff);
1535 AddDigit(idx,j,signal);
1536 } // end if singal < dis
1538 // write diff in the buffer used to compute Huffman tables
1539 if (firstSignal) str[counter]=(UChar_t)signal;
1540 else str[counter]=(UChar_t)diff;
1544 } // end for j loop time samples
1545 } // end for i loop anodes one half of detector
1549 fStream->CheckCount(counter);
1551 // open file and write out the stream of diff's
1552 static Bool_t open=kTRUE;
1553 static TFile *outFile;
1554 Bool_t write = fResponse->OutputOption();
1555 TDirectory *savedir = gDirectory;
1559 SetFileName("stream.root");
1560 cout<<"filename "<<fFileName<<endl;
1561 outFile=new TFile(fFileName,"recreate");
1562 cout<<"I have opened "<<fFileName<<" file "<<endl;
1569 fStream->ClearStream();
1571 // back to galice.root file
1572 if(savedir) savedir->cd();
1574 //______________________________________________________________________
1575 void AliITSsimulationSDD::StoreAllDigits(){
1576 // if non-zero-suppressed data
1577 Bool_t do10to8 = fResponse->Do10to8();
1578 Int_t i, j, digits[3];
1580 for (i=0; i<fNofMaps; i++) {
1581 for (j=0; j<fMaxNofSamples; j++) {
1582 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1583 if(do10to8) signal = Convert10to8(signal);
1584 if(do10to8) signal = Convert8to10(signal);
1588 fITS->AddRealDigit(1,digits);
1592 //______________________________________________________________________
1593 void AliITSsimulationSDD::CreateHistograms(Int_t scale){
1594 // Creates histograms of maps for debugging
1597 fHis=new TObjArray(fNofMaps);
1598 for (i=0;i<fNofMaps;i++) {
1599 TString sddName("sdd_");
1601 sprintf(candNum,"%d",i+1);
1602 sddName.Append(candNum);
1603 fHis->AddAt(new TH1F(sddName.Data(),"SDD maps",scale*fMaxNofSamples,
1604 0.,(Float_t) scale*fMaxNofSamples), i);
1607 //______________________________________________________________________
1608 void AliITSsimulationSDD::FillHistograms(){
1609 // fill 1D histograms from map
1613 for( Int_t i=0; i<fNofMaps; i++) {
1614 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
1615 Int_t nsamples = hist->GetNbinsX();
1616 for( Int_t j=0; j<nsamples; j++) {
1617 Double_t signal=fHitMap2->GetSignal(i,j);
1618 hist->Fill((Float_t)j,signal);
1622 //______________________________________________________________________
1623 void AliITSsimulationSDD::ResetHistograms(){
1624 // Reset histograms for this detector
1627 for (i=0;i<fNofMaps;i++ ) {
1628 if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset();
1631 //______________________________________________________________________
1632 TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
1633 // Fills a histogram from a give anode.
1635 if (!fHis) return 0;
1637 if(wing <=0 || wing > 2) {
1638 Warning("GetAnode","Wrong wing number: %d",wing);
1640 } // end if wing <=0 || wing >2
1641 if(anode <=0 || anode > fNofMaps/2) {
1642 Warning("GetAnode","Wrong anode number: %d",anode);
1644 } // end if ampde <=0 || andoe > fNofMaps/2
1646 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
1647 return (TH1F*)(fHis->At(index));
1649 //______________________________________________________________________
1650 void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
1651 // Writes the histograms to a file
1657 for(i=0; i<fNofMaps; i++) fHis->At(i)->Write(); //fAdcs[i]->Write();
1660 //______________________________________________________________________
1661 Float_t AliITSsimulationSDD::GetNoise() {
1662 // Returns the noise value
1663 //Bool_t do10to8=fResponse->Do10to8();
1664 //noise will always be in the liniar part of the signal
1666 Int_t threshold = fT1[0];
1667 char opt1[20], opt2[20];
1669 fResponse->ParamOptions(opt1,opt2);
1671 char *same = strstr(opt1,"same");
1672 Float_t noise,baseline;
1674 fResponse->GetNoiseParam(noise,baseline);
1676 static Bool_t readfile=kTRUE;
1677 //read baseline and noise from file
1678 if (readfile) ReadBaseline();
1682 TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2");
1683 if(c2) delete c2->GetPrimitive("noisehist");
1684 if(c2) delete c2->GetPrimitive("anode");
1685 else c2=new TCanvas("c2");
1687 c2->SetFillColor(0);
1689 TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold);
1690 TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0.,
1691 (float)fMaxNofSamples);
1693 for (i=0;i<fNofMaps;i++) {
1694 CompressionParam(i,decr,threshold);
1695 if (!same) GetAnodeBaseline(i,baseline,noise);
1697 for (k=0;k<fMaxNofSamples;k++) {
1698 Float_t signal=(Float_t)fHitMap2->GetSignal(i,k);
1699 //if (signal <= (float)threshold) noisehist->Fill(signal-baseline);
1700 if (signal <= (float)threshold) noisehist->Fill(signal);
1701 anode->Fill((float)k,signal);
1706 TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold);
1707 noisehist->Fit("gnoise","RQ");
1710 Float_t mnoise = gnoise->GetParameter(1);
1711 cout << "mnoise : " << mnoise << endl;
1712 Float_t rnoise = gnoise->GetParameter(2);
1713 cout << "rnoise : " << rnoise << endl;
1717 //______________________________________________________________________
1718 void AliITSsimulationSDD::WriteSDigits(){
1719 // Fills the Summable digits Tree
1720 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
1722 for( Int_t i=0; i<fNofMaps; i++ ) {
1723 for( Int_t j=0; j<fMaxNofSamples; j++ ) {
1724 Double_t sig = fHitMap2->GetSignal( i, j );
1726 Int_t jdx = j*fScaleSize;
1727 Int_t index = fpList->GetHitIndex( i, j );
1728 AliITSpListItem pItemTmp2( fModule, index, 0. );
1729 // put the fScaleSize analog digits in only one
1730 for( Int_t ik=0; ik<fScaleSize; ik++ ) {
1731 AliITSpListItem *pItemTmp = fpList->GetpListItem( i, jdx+ik );
1732 if( pItemTmp == 0 ) continue;
1733 pItemTmp2.Add( pItemTmp );
1735 pItemTmp2.AddSignalAfterElect( fModule, index, sig );
1736 pItemTmp2.AddNoise( fModule, index, fHitNoiMap2->GetSignal( i, j ) );
1737 aliITS->AddSumDigit( pItemTmp2 );
1738 } // end if (sig > 0.2)
1743 //______________________________________________________________________
1744 void AliITSsimulationSDD::Print() {
1745 // Print SDD simulation Parameters
1747 cout << "**************************************************" << endl;
1748 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1749 cout << "**************************************************" << endl;
1750 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1751 cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl;
1752 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1753 cout << "Number pf Anodes used: " << fNofMaps << endl;
1754 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1755 cout << "Scale size factor: " << fScaleSize << endl;
1756 cout << "**************************************************" << endl;