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 **************************************************************************/
23 #include <TStopwatch.h>
35 #include "AliITShit.h"
36 #include "AliITSdigit.h"
37 #include "AliITSmodule.h"
38 #include "AliITSpList.h"
39 #include "AliITSMapA1.h"
40 #include "AliITSMapA2.h"
41 #include "AliITSetfSDD.h"
42 #include "AliITSRawData.h"
43 #include "AliITSHuffman.h"
44 #include "AliITSsegmentation.h"
45 #include "AliITSresponse.h"
46 #include "AliITSsegmentationSDD.h"
47 #include "AliITSresponseSDD.h"
48 #include "AliITSsimulationSDD.h"
50 ClassImp(AliITSsimulationSDD)
51 ////////////////////////////////////////////////////////////////////////
53 // Written by Piergiorgio Cerello
56 // AliITSsimulationSDD is the simulation of SDDs.
60 <img src="picts/ITS/AliITShit_Class_Diagram.gif">
63 <font size=+2 color=red>
64 <p>This show the relasionships between the ITS hit class and the rest of Aliroot.
69 //______________________________________________________________________
70 Int_t power(Int_t b, Int_t e) {
71 // compute b to the e power, where both b and e are Int_ts.
74 for(i=0; i<e; i++) power *= b;
77 //______________________________________________________________________
78 void FastFourierTransform(AliITSetfSDD *alisddetf,Double_t *real,
79 Double_t *imag,Int_t direction) {
80 // Do a Fast Fourier Transform
82 Int_t samples = alisddetf->GetSamples();
83 Int_t l = (Int_t) ((log((Float_t) samples)/log(2.))+0.5);
86 Int_t m2 = samples/m1;
89 for(j=0; j<samples; j += m1) {
91 for(k=j; k<= j+m-1; k++) {
92 Double_t wsr = alisddetf->GetWeightReal(p);
93 Double_t wsi = alisddetf->GetWeightImag(p);
94 if(direction == -1) wsi = -wsi;
95 Double_t xr = *(real+k+m);
96 Double_t xi = *(imag+k+m);
97 *(real+k+m) = wsr*(*(real+k)-xr) - wsi*(*(imag+k)-xi);
98 *(imag+k+m) = wsr*(*(imag+k)-xi) + wsi*(*(real+k)-xr);
109 for(j=0; j<samples; j++) {
113 for(i1=1; i1<=l; i1++) {
116 p = p + p + j2 - j1 - j1;
119 Double_t xr = *(real+j);
120 Double_t xi = *(imag+j);
121 *(real+j) = *(real+p);
122 *(imag+j) = *(imag+p);
127 if(direction == -1) {
128 for(i=0; i<samples; i++) {
129 *(real+i) /= samples;
130 *(imag+i) /= samples;
132 } // end if direction == -1
135 //______________________________________________________________________
136 AliITSsimulationSDD::AliITSsimulationSDD(){
137 // Default constructor
161 SetPerpendTracksFlag();
165 //______________________________________________________________________
166 AliITSsimulationSDD::AliITSsimulationSDD(AliITSsimulationSDD &source){
167 // Copy constructor to satify Coding roules only.
169 if(this==&source) return;
170 Error("AliITSsimulationSSD","Not allowed to make a copy of "
171 "AliITSsimulationSDD Using default creater instead");
172 AliITSsimulationSDD();
174 //______________________________________________________________________
175 AliITSsimulationSDD& AliITSsimulationSDD::operator=(AliITSsimulationSDD &src){
176 // Assignment operator to satify Coding roules only.
178 if(this==&src) return *this;
179 Error("AliITSsimulationSSD","Not allowed to make a = with "
180 "AliITSsimulationSDD Using default creater instead");
183 //______________________________________________________________________
184 AliITSsimulationSDD::AliITSsimulationSDD(AliITSsegmentation *seg,
185 AliITSresponse *resp){
186 // Standard Constructor
205 Init((AliITSsegmentationSDD*)seg,(AliITSresponseSDD*)resp);
207 //______________________________________________________________________
208 void AliITSsimulationSDD::Init(AliITSsegmentationSDD *seg,
209 AliITSresponseSDD *resp){
210 // Standard Constructor
215 SetPerpendTracksFlag();
219 fHitMap2 = new AliITSMapA2(fSegmentation,fScaleSize,1);
220 fHitMap1 = new AliITSMapA1(fSegmentation);
222 fNofMaps = fSegmentation->Npz();
223 fMaxNofSamples = fSegmentation->Npx();
225 Float_t sddLength = fSegmentation->Dx();
226 Float_t sddWidth = fSegmentation->Dz();
229 Float_t anodePitch = fSegmentation->Dpz(dummy);
230 Double_t timeStep = (Double_t)fSegmentation->Dpx(dummy);
231 Float_t driftSpeed = fResponse->DriftSpeed();
233 if(anodePitch*(fNofMaps/2) > sddWidth) {
234 Warning("AliITSsimulationSDD",
235 "Too many anodes %d or too big pitch %f \n",
236 fNofMaps/2,anodePitch);
239 if(timeStep*fMaxNofSamples < sddLength/driftSpeed) {
240 Error("AliITSsimulationSDD",
241 "Time Interval > Allowed Time Interval: exit\n");
245 fElectronics = new AliITSetfSDD(timeStep/fScaleSize,
246 fResponse->Electronics());
248 char opt1[20], opt2[20];
249 fResponse->ParamOptions(opt1,opt2);
251 char *same = strstr(opt1,"same");
256 fNoise.Set(fNofMaps);
257 fBaseline.Set(fNofMaps);
260 const char *kopt=fResponse->ZeroSuppOption();
261 if (strstr(fParam,"file") ) {
264 if (strstr(kopt,"2D")) {
267 Init2D(); // desactivate if param change module by module
268 } else if(strstr(kopt,"1D")) {
271 Init1D(); // desactivate if param change module by module
279 } // end if else strstr
281 Bool_t write = fResponse->OutputOption();
282 if(write && strstr(kopt,"2D")) MakeTreeB();
284 // call here if baseline does not change by module
287 fITS = (AliITS*)gAlice->GetModule("ITS");
288 Int_t size = fNofMaps*fMaxNofSamples;
289 fStream = new AliITSInStream(size);
291 fInZR = new Double_t [fScaleSize*fMaxNofSamples];
292 fInZI = new Double_t [fScaleSize*fMaxNofSamples];
293 fOutZR = new Double_t [fScaleSize*fMaxNofSamples];
294 fOutZI = new Double_t [fScaleSize*fMaxNofSamples];
297 //______________________________________________________________________
298 AliITSsimulationSDD::~AliITSsimulationSDD() {
318 if(fTreeB) delete fTreeB;
319 if(fInZR) delete [] fInZR;
320 if(fInZI) delete [] fInZI;
321 if(fOutZR) delete [] fOutZR;
322 if(fOutZI) delete [] fOutZI;
324 //______________________________________________________________________
325 void AliITSsimulationSDD::SDigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
326 // create maps to build the lists of tracks for each summable digit
328 TObjArray *fHits = mod->GetHits();
329 Int_t nhits = fHits->GetEntriesFast();
335 AliITSpList *pList = new AliITSpList(2*fSegmentation->Npz(),
336 fScaleSize*fSegmentation->Npx());
338 // inputs to ListOfFiredCells.
339 TObjArray *alist = new TObjArray();
340 fHitMap1->SetArray(alist);
341 static TClonesArray *padr = 0;
342 if(!padr) padr = new TClonesArray("TVector",1000);
344 HitsToAnalogDigits(mod,alist,padr,pList);
352 fHitMap1->ClearMap();
353 fHitMap2->ClearMap();
355 //______________________________________________________________________
356 void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
357 // create maps to build the lists of tracks for each digit
359 TObjArray *fHits = mod->GetHits();
360 Int_t nhits = fHits->GetEntriesFast();
364 if (!nhits && fCheckNoise) {
367 fHitMap2->ClearMap();
369 } else if (!nhits) return;
371 AliITSpList *pList = new AliITSpList(2*fSegmentation->Npz(),
372 fScaleSize*fSegmentation->Npx());
374 // inputs to ListOfFiredCells.
375 TObjArray *alist = new TObjArray();
376 fHitMap1->SetArray(alist);
377 static TClonesArray *padr = 0;
378 if(!padr) padr = new TClonesArray("TVector",1000);
380 HitsToAnalogDigits(mod,alist,padr,pList);
388 fHitMap1->ClearMap();
389 fHitMap2->ClearMap();
391 //______________________________________________________________________
392 void AliITSsimulationSDD::SDigitsToDigits(AliITSpList *pList){
393 // Take Summable digits and create Digits.
395 // inputs to ListOfFiredCells.
396 TObjArray *alist = new TObjArray();
397 fHitMap1->SetArray(alist);
398 static TClonesArray *padr = 0;
399 if(!padr) padr = new TClonesArray("TVector",1000);
400 Int_t arg[6] = {0,0,0,0,0,0};
401 Double_t timeAmplitude;
404 // Fill maps from pList.
405 for(i=0;i<pList->GetMaxIndex();i++){
406 pList->GetMapIndex(i,arg[0],arg[1]);
407 for(j=0;j<pList->GetNEnteries();j++){
408 timeAmplitude = pList->GetTSignal(arg[0],arg[1],j);
409 if(timeAmplitude>0.0) continue;
410 arg[2] = pList->GetTrack(arg[0],arg[1],j);
411 arg[3] = pList->GetHit(arg[0],arg[1],j);
412 ListOfFiredCells(arg,timeAmplitude,alist,padr);
414 // Make sure map has full signal in it.
415 fHitMap2->SetHit(arg[0],arg[1],pList->GetSignal(arg[0],arg[1]));
424 fHitMap1->ClearMap();
425 fHitMap2->ClearMap();
427 //______________________________________________________________________
428 void AliITSsimulationSDD::FinishDigits(TObjArray *alist){
429 // introduce the electronics effects and do zero-suppression if required
430 Int_t nentries=alist->GetEntriesFast();
432 if(!nentries) return;
434 const char *kopt=fResponse->ZeroSuppOption();
435 ZeroSuppression(kopt);
437 //______________________________________________________________________
438 void AliITSsimulationSDD::HitsToAnalogDigits(AliITSmodule *mod,TObjArray *alst,
441 // create maps to build the lists of tracks for each digit
443 TObjArray *fHits = mod->GetHits();
444 Int_t nhits = fHits->GetEntriesFast();
445 Int_t arg[6] = {0,0,0,0,0,0};
447 Int_t nofAnodes = fNofMaps/2;
448 Float_t sddLength = fSegmentation->Dx();
449 Float_t sddWidth = fSegmentation->Dz();
450 Float_t anodePitch = fSegmentation->Dpz(dummy);
451 Float_t timeStep = fSegmentation->Dpx(dummy);
452 Float_t driftSpeed = fResponse->DriftSpeed();
453 Float_t maxadc = fResponse->MaxAdc();
454 Float_t topValue = fResponse->DynamicRange();
455 Float_t cHloss = fResponse->ChargeLoss();
456 Float_t norm = maxadc/topValue;
457 Float_t dfCoeff, s1; fResponse->DiffCoeff(dfCoeff,s1); // Signal 2d Shape
458 Double_t eVpairs = 3.6; // electron pair energy eV.
459 Float_t nsigma = fResponse->NSigmaIntegration(); //
460 Int_t nlookups = fResponse->GausNLookUp(); //
462 // Piergiorgio's part (apart for few variables which I made float
463 // when i thought that can be done
464 // Fill detector maps with GEANT hits
465 // loop over hits in the module
467 const Float_t kconv = 1.0e+6; // GeV->KeV
469 Int_t hitDetector; // detector number (lay,lad,hitDetector)
470 Int_t iWing; // which detector wing/side.
471 Int_t detector; // 2*(detector-1)+iWing
472 Int_t ii,kk,ka,kt; // loop indexs
473 Int_t ia,it,index; // sub-pixel integration indexies
474 Int_t iAnode; // anode number.
475 Int_t timeSample; // time buckett.
476 Int_t anodeWindow; // anode direction charge integration width
477 Int_t timeWindow; // time direction charge integration width
478 Int_t jamin,jamax; // anode charge integration window
479 Int_t jtmin,jtmax; // time charge integration window
480 Int_t ndiv; // Anode window division factor.
481 Int_t nsplit; // the number of splits in anode and time windows==1.
482 Int_t nOfSplits; // number of times track length is split into
483 Float_t nOfSplitsF; // Floating point version of nOfSplits.
484 Float_t kkF; // Floating point version of loop index kk.
485 Float_t pathInSDD; // Track length in SDD.
486 Float_t drPath; // average position of track in detector. in microns
487 Float_t drTime; // Drift time
488 Float_t nmul; // drift time window multiplication factor.
489 Float_t avDrft; // x position of path length segment in cm.
490 Float_t avAnode; // Anode for path length segment in Anode number (float)
491 Float_t xAnode; // Floating point anode number.
492 Float_t driftPath; // avDrft in microns.
493 Float_t width; // width of signal at anodes.
494 Double_t depEnergy; // Energy deposited in this GEANT step.
495 Double_t xL[3],dxL[3]; // local hit coordinates and diff.
496 Double_t sigA; // sigma of signal at anode.
497 Double_t sigT; // sigma in time/drift direction for track segment
498 Double_t aStep,aConst; // sub-pixel size and offset anode
499 Double_t tStep,tConst; // sub-pixel size and offset time
500 Double_t amplitude; // signal amplitude for track segment in nanoAmpere
501 Double_t chargeloss; // charge loss for track segment.
502 Double_t anodeAmplitude; // signal amplitude in anode direction
503 Double_t aExpo; // exponent of Gaussian anode direction
504 Double_t timeAmplitude; // signal amplitude in time direction
505 Double_t tExpo; // exponent of Gaussian time direction
506 // Double_t tof; // Time of flight in ns of this step.
508 for(ii=0; ii<nhits; ii++) {
509 if(!mod->LineSegmentL(ii,xL[0],dxL[0],xL[1],dxL[1],xL[2],dxL[2],
510 depEnergy,itrack)) continue;
512 hitDetector = mod->GetDet();
513 //tof = 1.E+09*(mod->GetHit(ii)->GetTOF()); // tof in ns.
514 //if(tof>sddLength/driftSpeed) continue; // hit happed too late.
516 // scale path to simulate a perpendicular track
517 // continue if the particle did not lose energy
518 // passing through detector
520 Warning("HitsToAnalogDigits",
521 "This particle has passed without losing energy!");
523 } // end if !depEnergy
525 pathInSDD = TMath::Sqrt(dxL[0]*dxL[0]+dxL[1]*dxL[1]+dxL[2]*dxL[2]);
527 if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
528 drPath = 10000.*(dxL[0]+2.*xL[0])*0.5;
529 if(drPath < 0) drPath = -drPath;
530 drPath = sddLength-drPath;
532 Warning("HitsToAnalogDigits","negative drift path %e",drPath);
534 } // end if drPath < 0
536 // Compute number of segments to brake step path into
537 drTime = drPath/driftSpeed; // Drift Time
538 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);// Sigma along the anodes
539 // calcuate the number of time the path length should be split into.
540 nOfSplits = (Int_t) (1. + 10000.*pathInSDD/sigA);
541 if(fFlag) nOfSplits = 1;
543 // loop over path segments, init. some variables.
544 depEnergy /= nOfSplits;
545 nOfSplitsF = (Float_t) nOfSplits;
546 for(kk=0;kk<nOfSplits;kk++) { // loop over path segments
547 kkF = (Float_t) kk + 0.5;
548 avDrft = xL[0]+dxL[0]*kkF/nOfSplitsF;
549 avAnode = xL[2]+dxL[2]*kkF/nOfSplitsF;
550 driftPath = 10000.*avDrft;
552 iWing = 2; // Assume wing is 2
553 if(driftPath < 0) { // if wing is not 2 it is 1.
555 driftPath = -driftPath;
556 } // end if driftPath < 0
557 driftPath = sddLength-driftPath;
558 detector = 2*(hitDetector-1) + iWing;
560 Warning("HitsToAnalogDigits","Warning: negative drift path %e",
563 } // end if driftPath < 0
566 drTime = driftPath/driftSpeed; // drift time for segment.
567 timeSample = (Int_t) (fScaleSize*drTime/timeStep + 1);
568 // compute time Sample including tof information. The tof only
569 // effects the time of the signal is recoreded and not the
571 // timeSample = (Int_t) (fScaleSize*(drTime+tof)/timeStep + 1);
572 if(timeSample > fScaleSize*fMaxNofSamples) {
573 Warning("HItsToAnalogDigits","Wrong Time Sample: %e",
576 } // end if timeSample > fScaleSize*fMaxNoofSamples
579 xAnode = 10000.*(avAnode)/anodePitch + nofAnodes/2; // +1?
580 if(xAnode*anodePitch > sddWidth || xAnode*anodePitch < 0.)
581 Warning("HitsToAnalogDigits","Z = %e",
583 iAnode = (Int_t) (1.+xAnode); // xAnode?
584 if(iAnode < 1 || iAnode > nofAnodes) {
585 Warning("HitToAnalogDigits","Wrong iAnode: %d",iAnode);
587 } // end if iAnode < 1 || iAnode > nofAnodes
589 // store straight away the particle position in the array
590 // of particles and take idhit=ii only when part is entering (this
591 // requires FillModules() in the macro for analysis) :
593 // Sigma along the anodes for track segment.
594 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);
595 sigT = sigA/driftSpeed;
596 // Peak amplitude in nanoAmpere
597 amplitude = fScaleSize*160.*depEnergy/
598 (timeStep*eVpairs*2.*acos(-1.)*sigT*sigA);
599 amplitude *= timeStep/25.; // WARNING!!!!! Amplitude scaling to
600 // account for clock variations
601 // (reference value: 40 MHz)
602 chargeloss = 1.-cHloss*driftPath/1000;
603 amplitude *= chargeloss;
604 width = 2.*nsigma/(nlookups-1);
612 } // end if drTime > 1200.
614 nsplit = 4; // hard-wired //nsplit=4;nsplit = (nsplit+1)/2*2;
615 // Sub-pixel size see computation of aExpo and tExpo.
616 aStep = anodePitch/(nsplit*fScaleSize*sigA);
617 aConst = xAnode*anodePitch/sigA;
618 tStep = timeStep/(nsplit*fScaleSize*sigT);
619 tConst = drTime/sigT;
620 // Define SDD window corresponding to the hit
621 anodeWindow = (Int_t)(fScaleSize*nsigma*sigA/anodePitch+1);
622 timeWindow = (Int_t) (fScaleSize*nsigma*sigT/timeStep+1.);
623 jamin = (iAnode - anodeWindow/ndiv - 1)*fScaleSize*nsplit +1;
624 jamax = (iAnode + anodeWindow/ndiv)*fScaleSize*nsplit;
625 if(jamin <= 0) jamin = 1;
626 if(jamax > fScaleSize*nofAnodes*nsplit)
627 jamax = fScaleSize*nofAnodes*nsplit;
628 // jtmin and jtmax are Hard-wired
629 jtmin = (Int_t)(timeSample-timeWindow*nmul-1)*nsplit+1;
630 jtmax = (Int_t)(timeSample+timeWindow*nmul)*nsplit;
631 if(jtmin <= 0) jtmin = 1;
632 if(jtmax > fScaleSize*fMaxNofSamples*nsplit)
633 jtmax = fScaleSize*fMaxNofSamples*nsplit;
634 // Spread the charge in the anode-time window
635 for(ka=jamin; ka <=jamax; ka++) {
636 ia = (ka-1)/(fScaleSize*nsplit) + 1;
638 Warning("HitsToAnalogDigits","ia < 1: ");
641 if(ia > nofAnodes) ia = nofAnodes;
642 aExpo = (aStep*(ka-0.5)-aConst);
643 if(TMath::Abs(aExpo) > nsigma) anodeAmplitude = 0.;
645 dummy = (Int_t) ((aExpo+nsigma)/width);
646 anodeAmplitude = amplitude*fResponse->GausLookUp(dummy);
647 } // end if TMath::Abs(aEspo) > nsigma
648 // index starts from 0
649 index = ((detector+1)%2)*nofAnodes+ia-1;
650 if(anodeAmplitude) for(kt=jtmin; kt<=jtmax; kt++) {
651 it = (kt-1)/nsplit+1; // it starts from 1
653 Warning("HitsToAnalogDigits","it < 1:");
656 if(it>fScaleSize*fMaxNofSamples)
657 it = fScaleSize*fMaxNofSamples;
658 tExpo = (tStep*(kt-0.5)-tConst);
659 if(TMath::Abs(tExpo) > nsigma) timeAmplitude = 0.;
661 dummy = (Int_t) ((tExpo+nsigma)/width);
662 timeAmplitude = anodeAmplitude*
663 fResponse->GausLookUp(dummy);
664 } // end if TMath::Abs(tExpo) > nsigma
665 // build the list of digits for this module
668 arg[2] = itrack; // track number
669 arg[3] = ii-1; // hit number.
670 timeAmplitude *= norm;
672 ListOfFiredCells(arg,timeAmplitude,alst,padr);
673 pList->AddSignal(index,it,itrack,ii-1,
674 mod->GetIndex(),timeAmplitude);
675 } // end if anodeAmplitude and loop over time in window
676 } // loop over anodes in window
677 } // end loop over "sub-hits"
678 } // end loop over hits
680 //______________________________________________________________________
681 void AliITSsimulationSDD::ListOfFiredCells(Int_t *arg,Double_t timeAmplitude,
682 TObjArray *alist,TClonesArray *padr){
683 // Returns the list of "fired" cells.
685 Int_t index = arg[0];
687 Int_t idtrack = arg[2];
688 Int_t idhit = arg[3];
689 Int_t counter = arg[4];
690 Int_t countadr = arg[5];
691 Double_t charge = timeAmplitude;
692 charge += fHitMap2->GetSignal(index,ik-1);
693 fHitMap2->SetHit(index, ik-1, charge);
696 Int_t it = (Int_t)((ik-1)/fScaleSize);
699 digits[2] = (Int_t)timeAmplitude;
701 if (idtrack >= 0) phys = (Float_t)timeAmplitude;
704 Double_t cellcharge = 0.;
705 AliITSTransientDigit* pdigit;
706 // build the list of fired cells and update the info
707 if (!fHitMap1->TestHit(index, it)) {
708 new((*padr)[countadr++]) TVector(3);
709 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
710 trinfo(0) = (Float_t)idtrack;
711 trinfo(1) = (Float_t)idhit;
712 trinfo(2) = (Float_t)timeAmplitude;
714 alist->AddAtAndExpand(new AliITSTransientDigit(phys,digits),counter);
715 fHitMap1->SetHit(index, it, counter);
717 pdigit=(AliITSTransientDigit*)alist->At(alist->GetLast());
719 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
720 trlist->Add(&trinfo);
722 pdigit = (AliITSTransientDigit*) fHitMap1->GetHit(index, it);
723 for(Int_t kk=0;kk<fScaleSize;kk++) {
724 cellcharge += fHitMap2->GetSignal(index,fScaleSize*it+kk);
727 (*pdigit).fSignal = (Int_t)cellcharge;
728 (*pdigit).fPhysics += phys;
729 // update list of tracks
730 TObjArray* trlist = (TObjArray*)pdigit->TrackList();
731 Int_t lastentry = trlist->GetLast();
732 TVector *ptrkp = (TVector*)trlist->At(lastentry);
733 TVector &trinfo = *ptrkp;
734 Int_t lasttrack = Int_t(trinfo(0));
735 Float_t lastcharge=(trinfo(2));
736 if (lasttrack==idtrack ) {
737 lastcharge += (Float_t)timeAmplitude;
738 trlist->RemoveAt(lastentry);
739 trinfo(0) = lasttrack;
741 trinfo(2) = lastcharge;
742 trlist->AddAt(&trinfo,lastentry);
744 new((*padr)[countadr++]) TVector(3);
745 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
746 trinfo(0) = (Float_t)idtrack;
747 trinfo(1) = (Float_t)idhit;
748 trinfo(2) = (Float_t)timeAmplitude;
749 trlist->Add(&trinfo);
750 } // end if lasttrack==idtrack
753 // check the track list - debugging
754 Int_t trk[20], htrk[20];
756 Int_t nptracks = trlist->GetEntriesFast();
759 for (tr=0;tr<nptracks;tr++) {
760 TVector *pptrkp = (TVector*)trlist->At(tr);
761 TVector &pptrk = *pptrkp;
762 trk[tr] = Int_t(pptrk(0));
763 htrk[tr] = Int_t(pptrk(1));
764 chtrk[tr] = (pptrk(2));
765 cout << "nptracks "<<nptracks << endl;
771 // update counter and countadr for next call.
775 //____________________________________________
777 void AliITSsimulationSDD::AddDigit(Int_t i, Int_t j, Int_t signal){
779 // tag with -1 signals coming from background tracks
780 // tag with -2 signals coming from pure electronic noise
782 Int_t digits[3], tracks[3], hits[3];
783 Float_t phys, charges[3];
785 Int_t trk[20], htrk[20];
788 Bool_t do10to8=fResponse->Do10to8();
790 if(do10to8) signal=Convert8to10(signal);
791 AliITSTransientDigit *obj = (AliITSTransientDigit*)fHitMap1->GetHit(i,j);
803 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
806 TObjArray* trlist=(TObjArray*)obj->TrackList();
807 Int_t nptracks=trlist->GetEntriesFast();
809 Warning("AddDigit","nptracks=%d > 20 nptracks set to 20",nptracks);
811 } // end if nptracks > 20
813 for (tr=0;tr<nptracks;tr++) {
814 TVector &pp =*((TVector*)trlist->At(tr));
815 trk[tr]=Int_t(pp(0));
816 htrk[tr]=Int_t(pp(1));
820 SortTracks(trk,chtrk,htrk,nptracks);
821 } // end if nptracks > 1
824 for (i=0; i<nptracks; i++) {
829 for (i=nptracks; i<3; i++) {
835 for (i=0; i<3; i++) {
840 } // end if/else nptracks < 3
842 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
844 } // end if/else !obj
846 //______________________________________________________________________
847 void AliITSsimulationSDD::SortTracks(Int_t *tracks,Float_t *charges,
848 Int_t *hits,Int_t ntr){
849 // Sort the list of tracks contributing to a given digit
850 // Only the 3 most significant tracks are acctually sorted
851 // Loop over signals, only 3 times
855 Int_t idx[3] = {-3,-3,-3};
856 Float_t jch[3] = {-3,-3,-3};
857 Int_t jtr[3] = {-3,-3,-3};
858 Int_t jhit[3] = {-3,-3,-3};
861 if (ntr<3) imax = ntr;
867 if((i == 1 && j == idx[i-1] )
868 ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
869 if(charges[j] > qmax) {
872 } // end if charges[j]>qmax
876 jch[i] = charges[jmax];
877 jtr[i] = tracks[jmax];
878 jhit[i] = hits[jmax];
891 } // end if jtr[i] == -3
894 //______________________________________________________________________
895 void AliITSsimulationSDD::ChargeToSignal() {
896 // add baseline, noise, electronics and ADC saturation effects
898 char opt1[20], opt2[20];
899 fResponse->ParamOptions(opt1,opt2);
900 char *read = strstr(opt1,"file");
901 Float_t baseline, noise;
904 static Bool_t readfile=kTRUE;
905 //read baseline and noise from file
906 if (readfile) ReadBaseline();
908 } else fResponse->GetNoiseParam(noise,baseline);
912 Float_t maxadc = fResponse->MaxAdc();
914 for (i=0;i<fNofMaps;i++) {
915 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
916 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
917 fInZR[k] = fHitMap2->GetSignal(i,k);
918 contrib = (baseline + noise*gRandom->Gaus());
921 for(k=0; k<fMaxNofSamples; k++) {
922 Double_t newcont = 0.;
923 Double_t maxcont = 0.;
924 for(kk=0;kk<fScaleSize;kk++) {
925 newcont = fInZR[fScaleSize*k+kk];
926 if(newcont > maxcont) maxcont = newcont;
929 if (newcont >= maxadc) newcont = maxadc -1;
930 if(newcont >= baseline){
931 Warning("","newcont=%d>=baseline=%d",newcont,baseline);
934 fHitMap2->SetHit(i,k,newcont);
936 } // end for i loop over anodes
940 for (i=0;i<fNofMaps;i++) {
941 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
942 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
943 fInZR[k] = fHitMap2->GetSignal(i,k);
944 contrib = (baseline + noise*gRandom->Gaus());
948 FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1);
949 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
950 Double_t rw = fElectronics->GetTraFunReal(k);
951 Double_t iw = fElectronics->GetTraFunImag(k);
952 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
953 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
955 FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1);
956 for(k=0; k<fMaxNofSamples; k++) {
957 Double_t newcont1 = 0.;
958 Double_t maxcont1 = 0.;
959 for(kk=0;kk<fScaleSize;kk++) {
960 newcont1 = fOutZR[fScaleSize*k+kk];
961 if(newcont1 > maxcont1) maxcont1 = newcont1;
964 if (newcont1 >= maxadc) newcont1 = maxadc -1;
965 fHitMap2->SetHit(i,k,newcont1);
967 } // end for i loop over anodes
970 //______________________________________________________________________
971 void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline,
973 // Returns the Baseline for a particular anode.
974 baseline = fBaseline[i];
977 //______________________________________________________________________
978 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl,
980 // Returns the compression alogirthm parameters
981 Int_t size = fD.GetSize();
983 db=fD[i]; tl=fT1[i]; th=fT2[i];
985 if (size <= 2 && i>=fNofMaps/2) {
986 db=fD[1]; tl=fT1[1]; th=fT2[1];
988 db=fD[0]; tl=fT1[0]; th=fT2[0];
989 } // end if size <=2 && i>=fNofMaps/2
992 //______________________________________________________________________
993 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){
994 // returns the compression alogirthm parameters
995 Int_t size = fD.GetSize();
1000 if (size <= 2 && i>=fNofMaps/2) {
1001 db=fD[1]; tl=fT1[1];
1003 db=fD[0]; tl=fT1[0];
1004 } // end if size <=2 && i>=fNofMaps/2
1005 } // end if size > 2
1007 //______________________________________________________________________
1008 void AliITSsimulationSDD::SetCompressParam(){
1009 // Sets the compression alogirthm parameters
1012 fResponse->GiveCompressParam(cp);
1013 for (i=0; i<2; i++) {
1020 //______________________________________________________________________
1021 void AliITSsimulationSDD::ReadBaseline(){
1022 // read baseline and noise from file - either a .root file and in this
1023 // case data should be organised in a tree with one entry for each
1024 // module => reading should be done accordingly
1025 // or a classic file and do smth. like this:
1026 // Read baselines and noise for SDD
1030 char input[100], base[100], param[100];
1033 fResponse->Filenames(input,base,param);
1036 filtmp = gSystem->ExpandPathName(fFileName.Data());
1037 FILE *bline = fopen(filtmp,"r");
1041 while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) {
1043 Error("ReadBaseline","Anode number not in increasing order!",
1046 } // end if pos != na+1
1052 Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",filtmp);
1059 //______________________________________________________________________
1060 Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) const {
1061 // To the 10 to 8 bit lossive compression.
1062 // code from Davide C. and Albert W.
1064 if (signal < 128) return signal;
1065 if (signal < 256) return (128+((signal-128)>>1));
1066 if (signal < 512) return (192+((signal-256)>>3));
1067 if (signal < 1024) return (224+((signal-512)>>4));
1070 //______________________________________________________________________
1071 Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) const {
1072 // Undo the lossive 10 to 8 bit compression.
1073 // code from Davide C. and Albert W.
1074 if (signal < 0 || signal > 255) {
1075 Warning("Convert8to10","out of range signal=%d",signal);
1077 } // end if signal <0 || signal >255
1079 if (signal < 128) return signal;
1081 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
1082 else return (128+((signal-128)<<1)+1);
1083 } // end if signal < 192
1085 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
1086 else return (256+((signal-192)<<3)+4);
1087 } // end if signal < 224
1088 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
1089 return (512+((signal-224)<<4)+7);
1091 //______________________________________________________________________
1092 AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){
1093 //Return the correct map.
1095 return ((i==0)? fHitMap1 : fHitMap2);
1097 //______________________________________________________________________
1098 void AliITSsimulationSDD::ZeroSuppression(const char *option) {
1099 // perform the zero suppresion
1101 if (strstr(option,"2D")) {
1102 //Init2D(); // activate if param change module by module
1104 } else if (strstr(option,"1D")) {
1105 //Init1D(); // activate if param change module by module
1107 } else StoreAllDigits();
1109 //______________________________________________________________________
1110 void AliITSsimulationSDD::Init2D(){
1111 // read in and prepare arrays: fD, fT1, fT2
1112 // savemu[nanodes], savesigma[nanodes]
1113 // read baseline and noise from file - either a .root file and in this
1114 // case data should be organised in a tree with one entry for each
1115 // module => reading should be done accordingly
1116 // or a classic file and do smth. like this ( code from Davide C. and
1118 // Read 2D zero-suppression parameters for SDD
1120 if (!strstr(fParam,"file")) return;
1122 Int_t na,pos,tempTh;
1124 Float_t *savemu = new Float_t [fNofMaps];
1125 Float_t *savesigma = new Float_t [fNofMaps];
1126 char input[100],basel[100],par[100];
1128 Int_t minval = fResponse->MinVal();
1130 fResponse->Filenames(input,basel,par);
1133 filtmp = gSystem->ExpandPathName(fFileName.Data());
1134 FILE *param = fopen(filtmp,"r");
1138 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1140 Error("Init2D","Anode number not in increasing order!",filtmp);
1142 } // end if pos != na+1
1144 savesigma[na] = sigma;
1145 if ((2.*sigma) < mu) {
1146 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1149 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1150 if (tempTh < 0) tempTh=0;
1152 tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval;
1153 if (tempTh < 0) tempTh=0;
1158 Error("Init2D","THE FILE %s DOES NOT EXIST !",filtmp);
1165 delete [] savesigma;
1167 //______________________________________________________________________
1168 void AliITSsimulationSDD::Compress2D(){
1169 // simple ITS cluster finder -- online zero-suppression conditions
1172 Int_t minval = fResponse->MinVal();
1173 Bool_t write = fResponse->OutputOption();
1174 Bool_t do10to8 = fResponse->Do10to8();
1175 Int_t nz, nl, nh, low, i, j;
1177 for (i=0; i<fNofMaps; i++) {
1178 CompressionParam(i,db,tl,th);
1183 for (j=0; j<fMaxNofSamples; j++) {
1184 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1185 signal -= db; // if baseline eq. is done here
1186 if (signal <= 0) {nz++; continue;}
1187 if ((signal - tl) < minval) low++;
1188 if ((signal - th) >= minval) {
1191 FindCluster(i,j,signal,minval,cond);
1193 ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)){
1194 if(do10to8) signal = Convert10to8(signal);
1195 AddDigit(i,j,signal);
1196 } // end if cond&&j&&()
1197 } else if ((signal - tl) >= minval) nl++;
1198 } // end for j loop time samples
1199 if (write) TreeB()->Fill(nz,nl,nh,low,i+1);
1200 } //end for i loop anodes
1204 sprintf(hname,"TNtuple%d_%d",fModule,fEvent);
1205 TreeB()->Write(hname);
1210 //______________________________________________________________________
1211 void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal,
1212 Int_t minval,Bool_t &cond){
1213 // Find clusters according to the online 2D zero-suppression algorithm
1214 Bool_t do10to8 = fResponse->Do10to8();
1215 Bool_t high = kFALSE;
1217 fHitMap2->FlagHit(i,j);
1219 // check the online zero-suppression conditions
1221 const Int_t kMaxNeighbours = 4;
1224 Int_t xList[kMaxNeighbours], yList[kMaxNeighbours];
1225 fSegmentation->Neighbours(i,j,&nn,xList,yList);
1227 for (in=0; in<nn; in++) {
1230 if (fHitMap2->TestHit(ix,iy)==kUnused) {
1231 CompressionParam(ix,dbx,tlx,thx);
1232 Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy));
1233 qn -= dbx; // if baseline eq. is done here
1234 if ((qn-tlx) < minval) {
1235 fHitMap2->FlagHit(ix,iy);
1238 if ((qn - thx) >= minval) high=kTRUE;
1240 if(do10to8) signal = Convert10to8(signal);
1241 AddDigit(i,j,signal);
1243 if(do10to8) qns = Convert10to8(qn);
1245 if (!high) AddDigit(ix,iy,qns);
1247 if(!high) fHitMap2->FlagHit(ix,iy);
1248 } // end if qn-tlx < minval
1250 } // end for in loop over neighbours
1252 //______________________________________________________________________
1253 void AliITSsimulationSDD::Init1D(){
1254 // this is just a copy-paste of input taken from 2D algo
1255 // Torino people should give input
1256 // Read 1D zero-suppression parameters for SDD
1258 if (!strstr(fParam,"file")) return;
1260 Int_t na,pos,tempTh;
1262 Float_t *savemu = new Float_t [fNofMaps];
1263 Float_t *savesigma = new Float_t [fNofMaps];
1264 char input[100],basel[100],par[100];
1266 Int_t minval = fResponse->MinVal();
1268 fResponse->Filenames(input,basel,par);
1271 // set first the disable and tol param
1274 filtmp = gSystem->ExpandPathName(fFileName.Data());
1275 FILE *param = fopen(filtmp,"r");
1279 fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]);
1280 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1282 Error("Init1D","Anode number not in increasing order!",filtmp);
1284 } // end if pos != na+1
1286 savesigma[na]=sigma;
1287 if ((2.*sigma) < mu) {
1288 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1291 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1292 if (tempTh < 0) tempTh=0;
1297 Error("Init1D","THE FILE %s DOES NOT EXIST !",filtmp);
1304 delete [] savesigma;
1306 //______________________________________________________________________
1307 void AliITSsimulationSDD::Compress1D(){
1308 // 1D zero-suppression algorithm (from Gianluca A.)
1309 Int_t dis,tol,thres,decr,diff;
1310 UChar_t *str=fStream->Stream();
1312 Bool_t do10to8=fResponse->Do10to8();
1316 for (k=0; k<2; k++) {
1319 for (i=0; i<fNofMaps/2; i++) {
1320 Bool_t firstSignal=kTRUE;
1321 Int_t idx=i+k*fNofMaps/2;
1322 CompressionParam(idx,decr,thres);
1323 for (j=0; j<fMaxNofSamples; j++) {
1324 Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j));
1325 signal -= decr; // if baseline eq.
1326 if(do10to8) signal = Convert10to8(signal);
1327 if (signal <= thres) {
1331 // write diff in the buffer for HuffT
1332 str[counter]=(UChar_t)diff;
1335 } // end if signal <= thres
1337 if (diff > 127) diff=127;
1338 if (diff < -128) diff=-128;
1340 // tol has changed to 8 possible cases ? - one can write
1341 // this if(TMath::Abs(diff)<tol) ... else ...
1342 if(TMath::Abs(diff)<tol) diff=0;
1343 // or keep it as it was before
1345 if (tol==1 && (diff >= -2 && diff <= 1)) diff=0;
1346 if (tol==2 && (diff >= -4 && diff <= 3)) diff=0;
1347 if (tol==3 && (diff >= -16 && diff <= 15)) diff=0;
1349 AddDigit(idx,j,last+diff);
1351 AddDigit(idx,j,signal);
1352 } // end if singal < dis
1354 // write diff in the buffer used to compute Huffman tables
1355 if (firstSignal) str[counter]=(UChar_t)signal;
1356 else str[counter]=(UChar_t)diff;
1360 } // end for j loop time samples
1361 } // end for i loop anodes one half of detector
1365 fStream->CheckCount(counter);
1367 // open file and write out the stream of diff's
1368 static Bool_t open=kTRUE;
1369 static TFile *outFile;
1370 Bool_t write = fResponse->OutputOption();
1374 SetFileName("stream.root");
1375 cout<<"filename "<<fFileName<<endl;
1376 outFile=new TFile(fFileName,"recreate");
1377 cout<<"I have opened "<<fFileName<<" file "<<endl;
1384 fStream->ClearStream();
1386 // back to galice.root file
1388 TTree *fAli=gAlice->TreeK();
1391 if (fAli) file =fAli->GetCurrentFile();
1394 //______________________________________________________________________
1395 void AliITSsimulationSDD::StoreAllDigits(){
1396 // if non-zero-suppressed data
1397 Bool_t do10to8 = fResponse->Do10to8();
1398 Int_t i, j, digits[3];
1400 for (i=0; i<fNofMaps; i++) {
1401 for (j=0; j<fMaxNofSamples; j++) {
1402 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1403 if(do10to8) signal = Convert10to8(signal);
1404 if(do10to8) signal = Convert8to10(signal);
1408 fITS->AddRealDigit(1,digits);
1412 //______________________________________________________________________
1413 void AliITSsimulationSDD::CreateHistograms(Int_t scale){
1414 // Creates histograms of maps for debugging
1417 fHis=new TObjArray(fNofMaps);
1418 for (i=0;i<fNofMaps;i++) {
1419 TString sddName("sdd_");
1421 sprintf(candNum,"%d",i+1);
1422 sddName.Append(candNum);
1423 fHis->AddAt(new TH1F(sddName.Data(),"SDD maps",scale*fMaxNofSamples,
1424 0.,(Float_t) scale*fMaxNofSamples), i);
1427 //______________________________________________________________________
1428 void AliITSsimulationSDD::FillHistograms(){
1429 // fill 1D histograms from map
1433 for( Int_t i=0; i<fNofMaps; i++) {
1434 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
1435 Int_t nsamples = hist->GetNbinsX();
1436 for( Int_t j=0; j<nsamples; j++) {
1437 Double_t signal=fHitMap2->GetSignal(i,j);
1438 hist->Fill((Float_t)j,signal);
1442 //______________________________________________________________________
1443 void AliITSsimulationSDD::ResetHistograms(){
1444 // Reset histograms for this detector
1447 for (i=0;i<fNofMaps;i++ ) {
1448 if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset();
1451 //______________________________________________________________________
1452 TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
1453 // Fills a histogram from a give anode.
1455 if (!fHis) return 0;
1457 if(wing <=0 || wing > 2) {
1458 Warning("GetAnode","Wrong wing number: %d",wing);
1460 } // end if wing <=0 || wing >2
1461 if(anode <=0 || anode > fNofMaps/2) {
1462 Warning("GetAnode","Wrong anode number: %d",anode);
1464 } // end if ampde <=0 || andoe > fNofMaps/2
1466 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
1467 return (TH1F*)(fHis->At(index));
1469 //______________________________________________________________________
1470 void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
1471 // Writes the histograms to a file
1477 for(i=0; i<fNofMaps; i++) fHis->At(i)->Write(); //fAdcs[i]->Write();
1480 //______________________________________________________________________
1481 Float_t AliITSsimulationSDD::GetNoise() {
1482 // Returns the noise value
1483 //Bool_t do10to8=fResponse->Do10to8();
1484 //noise will always be in the liniar part of the signal
1486 Int_t threshold = fT1[0];
1487 char opt1[20], opt2[20];
1489 fResponse->ParamOptions(opt1,opt2);
1491 char *same = strstr(opt1,"same");
1492 Float_t noise,baseline;
1494 fResponse->GetNoiseParam(noise,baseline);
1496 static Bool_t readfile=kTRUE;
1497 //read baseline and noise from file
1498 if (readfile) ReadBaseline();
1502 TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2");
1503 if(c2) delete c2->GetPrimitive("noisehist");
1504 if(c2) delete c2->GetPrimitive("anode");
1505 else c2=new TCanvas("c2");
1507 c2->SetFillColor(0);
1509 TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold);
1510 TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0.,
1511 (float)fMaxNofSamples);
1513 for (i=0;i<fNofMaps;i++) {
1514 CompressionParam(i,decr,threshold);
1515 if (!same) GetAnodeBaseline(i,baseline,noise);
1517 for (k=0;k<fMaxNofSamples;k++) {
1518 Float_t signal=(Float_t)fHitMap2->GetSignal(i,k);
1519 //if (signal <= (float)threshold) noisehist->Fill(signal-baseline);
1520 if (signal <= (float)threshold) noisehist->Fill(signal);
1521 anode->Fill((float)k,signal);
1526 TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold);
1527 noisehist->Fit("gnoise","RQ");
1530 Float_t mnoise = gnoise->GetParameter(1);
1531 cout << "mnoise : " << mnoise << endl;
1532 Float_t rnoise = gnoise->GetParameter(2);
1533 cout << "rnoise : " << rnoise << endl;
1536 }//______________________________________________________________________
1537 void AliITSsimulationSDD::WriteSDigits(AliITSpList *pList){
1538 // Fills the Summable digits Tree
1540 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
1542 pList->GetMaxMapIndex(ni,nj);
1543 for(i=0;i<ni;i++)for(j=0;j<nj;j++){
1544 if(pList->GetSignalOnly(i,j)>0.5*fT1[0]){ // above small threshold.
1545 aliITS->AddSumDigit(*(pList->GetpListItem(i,j)));
1546 // cout << "pListSDD: " << *(pList->GetpListItem(i,j)) << endl;
1551 //______________________________________________________________________
1552 void AliITSsimulationSDD::Print() {
1553 // Print SDD simulation Parameters
1555 cout << "**************************************************" << endl;
1556 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1557 cout << "**************************************************" << endl;
1558 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1559 cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl;
1560 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1561 cout << "Number pf Anodes used: " << fNofMaps << endl;
1562 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1563 cout << "Scale size factor: " << fScaleSize << endl;
1564 cout << "**************************************************" << endl;