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() {
312 if(fTreeB) delete fTreeB;
313 if(fInZR) delete [] fInZR;
314 if(fInZI) delete [] fInZI;
315 if(fOutZR) delete [] fOutZR;
316 if(fOutZI) delete [] fOutZI;
318 //______________________________________________________________________
319 void AliITSsimulationSDD::SDigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
320 // create maps to build the lists of tracks for each summable digit
322 TObjArray *fHits = mod->GetHits();
323 Int_t nhits = fHits->GetEntriesFast();
329 AliITSpList *pList = new AliITSpList(2*fSegmentation->Npz(),
330 fScaleSize*fSegmentation->Npx());
332 // inputs to ListOfFiredCells.
333 TObjArray *alist = new TObjArray();
334 fHitMap1->SetArray(alist);
335 static TClonesArray *padr = 0;
336 if(!padr) padr = new TClonesArray("TVector",1000);
338 HitsToAnalogDigits(mod,alist,padr,pList);
347 fHitMap1->ClearMap();
348 fHitMap2->ClearMap();
350 //______________________________________________________________________
351 void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
352 // create maps to build the lists of tracks for each digit
354 TObjArray *fHits = mod->GetHits();
355 Int_t nhits = fHits->GetEntriesFast();
359 if (!nhits && fCheckNoise) {
362 fHitMap2->ClearMap();
364 } else if (!nhits) return;
366 AliITSpList *pList = new AliITSpList(2*fSegmentation->Npz(),
367 fScaleSize*fSegmentation->Npx());
369 // inputs to ListOfFiredCells.
370 TObjArray *alist = new TObjArray();
371 fHitMap1->SetArray(alist);
372 static TClonesArray *padr = 0;
373 if(!padr) padr = new TClonesArray("TVector",1000);
375 HitsToAnalogDigits(mod,alist,padr,pList);
384 fHitMap1->ClearMap();
385 fHitMap2->ClearMap();
387 //______________________________________________________________________
388 void AliITSsimulationSDD::SDigitsToDigits(AliITSpList *pList){
389 // Take Summable digits and create Digits.
391 // inputs to ListOfFiredCells.
392 TObjArray *alist = new TObjArray();
393 fHitMap1->SetArray(alist);
394 static TClonesArray *padr = 0;
395 if(!padr) padr = new TClonesArray("TVector",1000);
396 Int_t arg[6] = {0,0,0,0,0,0};
397 Double_t timeAmplitude;
400 // Fill maps from pList.
401 for(i=0;i<pList->GetMaxIndex();i++){
402 pList->GetMapIndex(i,arg[0],arg[1]);
403 for(j=0;j<pList->GetNEnteries();j++){
404 timeAmplitude = pList->GetTSignal(arg[0],arg[1],j);
405 if(timeAmplitude>0.0) continue;
406 arg[2] = pList->GetTrack(arg[0],arg[1],j);
407 arg[3] = pList->GetHit(arg[0],arg[1],j);
408 ListOfFiredCells(arg,timeAmplitude,alist,padr);
410 // Make sure map has full signal in it.
411 fHitMap2->SetHit(arg[0],arg[1],pList->GetSignal(arg[0],arg[1]));
420 fHitMap1->ClearMap();
421 fHitMap2->ClearMap();
423 //______________________________________________________________________
424 void AliITSsimulationSDD::FinishDigits(TObjArray *alist){
425 // introduce the electronics effects and do zero-suppression if required
426 Int_t nentries=alist->GetEntriesFast();
428 if(!nentries) return;
430 const char *kopt=fResponse->ZeroSuppOption();
431 ZeroSuppression(kopt);
433 //______________________________________________________________________
434 void AliITSsimulationSDD::HitsToAnalogDigits(AliITSmodule *mod,TObjArray *alst,
437 // create maps to build the lists of tracks for each digit
439 TObjArray *fHits = mod->GetHits();
440 Int_t nhits = fHits->GetEntriesFast();
441 Int_t arg[6] = {0,0,0,0,0,0};
443 Int_t nofAnodes = fNofMaps/2;
444 Float_t sddLength = fSegmentation->Dx();
445 Float_t sddWidth = fSegmentation->Dz();
446 Float_t anodePitch = fSegmentation->Dpz(dummy);
447 Float_t timeStep = fSegmentation->Dpx(dummy);
448 Float_t driftSpeed = fResponse->DriftSpeed();
449 Float_t maxadc = fResponse->MaxAdc();
450 Float_t topValue = fResponse->DynamicRange();
451 Float_t cHloss = fResponse->ChargeLoss();
452 Float_t norm = maxadc/topValue;
453 Float_t dfCoeff, s1; fResponse->DiffCoeff(dfCoeff,s1); // Signal 2d Shape
454 Double_t eVpairs = 3.6; // electron pair energy eV.
455 Float_t nsigma = fResponse->NSigmaIntegration(); //
456 Int_t nlookups = fResponse->GausNLookUp(); //
458 // Piergiorgio's part (apart for few variables which I made float
459 // when i thought that can be done
460 // Fill detector maps with GEANT hits
461 // loop over hits in the module
463 const Float_t kconv = 1.0e+6; // GeV->KeV
465 Int_t hitDetector; // detector number (lay,lad,hitDetector)
466 Int_t iWing; // which detector wing/side.
467 Int_t detector; // 2*(detector-1)+iWing
468 Int_t ii,kk,ka,kt; // loop indexs
469 Int_t ia,it,index; // sub-pixel integration indexies
470 Int_t iAnode; // anode number.
471 Int_t timeSample; // time buckett.
472 Int_t anodeWindow; // anode direction charge integration width
473 Int_t timeWindow; // time direction charge integration width
474 Int_t jamin,jamax; // anode charge integration window
475 Int_t jtmin,jtmax; // time charge integration window
476 Int_t ndiv; // Anode window division factor.
477 Int_t nsplit; // the number of splits in anode and time windows==1.
478 Int_t nOfSplits; // number of times track length is split into
479 Float_t nOfSplitsF; // Floating point version of nOfSplits.
480 Float_t kkF; // Floating point version of loop index kk.
481 Float_t pathInSDD; // Track length in SDD.
482 Float_t drPath; // average position of track in detector. in microns
483 Float_t drTime; // Drift time
484 Float_t nmul; // drift time window multiplication factor.
485 Float_t avDrft; // x position of path length segment in cm.
486 Float_t avAnode; // Anode for path length segment in Anode number (float)
487 Float_t xAnode; // Floating point anode number.
488 Float_t driftPath; // avDrft in microns.
489 Float_t width; // width of signal at anodes.
490 Double_t depEnergy; // Energy deposited in this GEANT step.
491 Double_t xL[3],dxL[3]; // local hit coordinates and diff.
492 Double_t sigA; // sigma of signal at anode.
493 Double_t sigT; // sigma in time/drift direction for track segment
494 Double_t aStep,aConst; // sub-pixel size and offset anode
495 Double_t tStep,tConst; // sub-pixel size and offset time
496 Double_t amplitude; // signal amplitude for track segment in nanoAmpere
497 Double_t chargeloss; // charge loss for track segment.
498 Double_t anodeAmplitude; // signal amplitude in anode direction
499 Double_t aExpo; // exponent of Gaussian anode direction
500 Double_t timeAmplitude; // signal amplitude in time direction
501 Double_t tExpo; // exponent of Gaussian time direction
502 // Double_t tof; // Time of flight in ns of this step.
504 for(ii=0; ii<nhits; ii++) {
505 if(!mod->LineSegmentL(ii,xL[0],dxL[0],xL[1],dxL[1],xL[2],dxL[2],
506 depEnergy,itrack)) continue;
508 hitDetector = mod->GetDet();
509 //tof = 1.E+09*(mod->GetHit(ii)->GetTOF()); // tof in ns.
510 //if(tof>sddLength/driftSpeed) continue; // hit happed too late.
512 // scale path to simulate a perpendicular track
513 // continue if the particle did not lose energy
514 // passing through detector
516 Warning("HitsToAnalogDigits",
517 "fTrack = %d hit=%d module=%d This particle has"
518 " passed without losing energy!",
519 itrack,ii,mod->GetIndex());
521 } // end if !depEnergy
523 pathInSDD = TMath::Sqrt(dxL[0]*dxL[0]+dxL[1]*dxL[1]+dxL[2]*dxL[2]);
525 if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
526 drPath = 10000.*(dxL[0]+2.*xL[0])*0.5;
527 if(drPath < 0) drPath = -drPath;
528 drPath = sddLength-drPath;
530 Warning("HitsToAnalogDigits",
531 "negative drift path drPath=%e sddLength=%e dxL[0]=%e "
533 drPath,sddLength,dxL[0],xL[0]);
535 } // end if drPath < 0
537 // Compute number of segments to brake step path into
538 drTime = drPath/driftSpeed; // Drift Time
539 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);// Sigma along the anodes
540 // calcuate the number of time the path length should be split into.
541 nOfSplits = (Int_t) (1. + 10000.*pathInSDD/sigA);
542 if(fFlag) nOfSplits = 1;
544 // loop over path segments, init. some variables.
545 depEnergy /= nOfSplits;
546 nOfSplitsF = (Float_t) nOfSplits;
547 for(kk=0;kk<nOfSplits;kk++) { // loop over path segments
548 kkF = (Float_t) kk + 0.5;
549 avDrft = xL[0]+dxL[0]*kkF/nOfSplitsF;
550 avAnode = xL[2]+dxL[2]*kkF/nOfSplitsF;
551 driftPath = 10000.*avDrft;
553 iWing = 2; // Assume wing is 2
554 if(driftPath < 0) { // if wing is not 2 it is 1.
556 driftPath = -driftPath;
557 } // end if driftPath < 0
558 driftPath = sddLength-driftPath;
559 detector = 2*(hitDetector-1) + iWing;
561 Warning("HitsToAnalogDigits","negative drift path "
562 "driftPath=%e sddLength=%e avDrft=%e dxL[0]=%e "
563 "xL[0]=%e",driftPath,sddLength,avDrft,dxL[0],xL[0]);
565 } // end if driftPath < 0
568 drTime = driftPath/driftSpeed; // drift time for segment.
569 timeSample = (Int_t) (fScaleSize*drTime/timeStep + 1);
570 // compute time Sample including tof information. The tof only
571 // effects the time of the signal is recoreded and not the
573 // timeSample = (Int_t) (fScaleSize*(drTime+tof)/timeStep + 1);
574 if(timeSample > fScaleSize*fMaxNofSamples) {
575 Warning("HItsToAnalogDigits","Wrong Time Sample: %e",
578 } // end if timeSample > fScaleSize*fMaxNoofSamples
581 xAnode = 10000.*(avAnode)/anodePitch + nofAnodes/2; // +1?
582 if(xAnode*anodePitch > sddWidth || xAnode*anodePitch < 0.)
583 Warning("HitsToAnalogDigits",
584 "Exceedubg sddWidth=%e Z = %e",
585 sddWidth,xAnode*anodePitch);
586 iAnode = (Int_t) (1.+xAnode); // xAnode?
587 if(iAnode < 1 || iAnode > nofAnodes) {
588 Warning("HitToAnalogDigits","Wrong iAnode: 1<%d>%d",
591 } // end if iAnode < 1 || iAnode > nofAnodes
593 // store straight away the particle position in the array
594 // of particles and take idhit=ii only when part is entering (this
595 // requires FillModules() in the macro for analysis) :
597 // Sigma along the anodes for track segment.
598 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);
599 sigT = sigA/driftSpeed;
600 // Peak amplitude in nanoAmpere
601 amplitude = fScaleSize*160.*depEnergy/
602 (timeStep*eVpairs*2.*acos(-1.)*sigT*sigA);
603 amplitude *= timeStep/25.; // WARNING!!!!! Amplitude scaling to
604 // account for clock variations
605 // (reference value: 40 MHz)
606 chargeloss = 1.-cHloss*driftPath/1000;
607 amplitude *= chargeloss;
608 width = 2.*nsigma/(nlookups-1);
616 } // end if drTime > 1200.
618 nsplit = 4; // hard-wired //nsplit=4;nsplit = (nsplit+1)/2*2;
619 // Sub-pixel size see computation of aExpo and tExpo.
620 aStep = anodePitch/(nsplit*fScaleSize*sigA);
621 aConst = xAnode*anodePitch/sigA;
622 tStep = timeStep/(nsplit*fScaleSize*sigT);
623 tConst = drTime/sigT;
624 // Define SDD window corresponding to the hit
625 anodeWindow = (Int_t)(fScaleSize*nsigma*sigA/anodePitch+1);
626 timeWindow = (Int_t) (fScaleSize*nsigma*sigT/timeStep+1.);
627 jamin = (iAnode - anodeWindow/ndiv - 1)*fScaleSize*nsplit +1;
628 jamax = (iAnode + anodeWindow/ndiv)*fScaleSize*nsplit;
629 if(jamin <= 0) jamin = 1;
630 if(jamax > fScaleSize*nofAnodes*nsplit)
631 jamax = fScaleSize*nofAnodes*nsplit;
632 // jtmin and jtmax are Hard-wired
633 jtmin = (Int_t)(timeSample-timeWindow*nmul-1)*nsplit+1;
634 jtmax = (Int_t)(timeSample+timeWindow*nmul)*nsplit;
635 if(jtmin <= 0) jtmin = 1;
636 if(jtmax > fScaleSize*fMaxNofSamples*nsplit)
637 jtmax = fScaleSize*fMaxNofSamples*nsplit;
638 // Spread the charge in the anode-time window
639 for(ka=jamin; ka <=jamax; ka++) {
640 ia = (ka-1)/(fScaleSize*nsplit) + 1;
642 Warning("HitsToAnalogDigits","ia < 1: ");
645 if(ia > nofAnodes) ia = nofAnodes;
646 aExpo = (aStep*(ka-0.5)-aConst);
647 if(TMath::Abs(aExpo) > nsigma) anodeAmplitude = 0.;
649 dummy = (Int_t) ((aExpo+nsigma)/width);
650 anodeAmplitude = amplitude*fResponse->GausLookUp(dummy);
651 } // end if TMath::Abs(aEspo) > nsigma
652 // index starts from 0
653 index = ((detector+1)%2)*nofAnodes+ia-1;
654 if(anodeAmplitude) for(kt=jtmin; kt<=jtmax; kt++) {
655 it = (kt-1)/nsplit+1; // it starts from 1
657 Warning("HitsToAnalogDigits","it < 1:");
660 if(it>fScaleSize*fMaxNofSamples)
661 it = fScaleSize*fMaxNofSamples;
662 tExpo = (tStep*(kt-0.5)-tConst);
663 if(TMath::Abs(tExpo) > nsigma) timeAmplitude = 0.;
665 dummy = (Int_t) ((tExpo+nsigma)/width);
666 timeAmplitude = anodeAmplitude*
667 fResponse->GausLookUp(dummy);
668 } // end if TMath::Abs(tExpo) > nsigma
669 // build the list of digits for this module
672 arg[2] = itrack; // track number
673 arg[3] = ii-1; // hit number.
674 timeAmplitude *= norm;
676 ListOfFiredCells(arg,timeAmplitude,alst,padr);
677 pList->AddSignal(index,it,itrack,ii-1,
678 mod->GetIndex(),timeAmplitude);
679 } // end if anodeAmplitude and loop over time in window
680 } // loop over anodes in window
681 } // end loop over "sub-hits"
682 } // end loop over hits
684 //______________________________________________________________________
685 void AliITSsimulationSDD::ListOfFiredCells(Int_t *arg,Double_t timeAmplitude,
686 TObjArray *alist,TClonesArray *padr){
687 // Returns the list of "fired" cells.
689 Int_t index = arg[0];
691 Int_t idtrack = arg[2];
692 Int_t idhit = arg[3];
693 Int_t counter = arg[4];
694 Int_t countadr = arg[5];
695 Double_t charge = timeAmplitude;
696 charge += fHitMap2->GetSignal(index,ik-1);
697 fHitMap2->SetHit(index, ik-1, charge);
700 Int_t it = (Int_t)((ik-1)/fScaleSize);
703 digits[2] = (Int_t)timeAmplitude;
705 if (idtrack >= 0) phys = (Float_t)timeAmplitude;
708 Double_t cellcharge = 0.;
709 AliITSTransientDigit* pdigit;
710 // build the list of fired cells and update the info
711 if (!fHitMap1->TestHit(index, it)) {
712 new((*padr)[countadr++]) TVector(3);
713 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
714 trinfo(0) = (Float_t)idtrack;
715 trinfo(1) = (Float_t)idhit;
716 trinfo(2) = (Float_t)timeAmplitude;
718 alist->AddAtAndExpand(new AliITSTransientDigit(phys,digits),counter);
719 fHitMap1->SetHit(index, it, counter);
721 pdigit=(AliITSTransientDigit*)alist->At(alist->GetLast());
723 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
724 trlist->Add(&trinfo);
726 pdigit = (AliITSTransientDigit*) fHitMap1->GetHit(index, it);
727 for(Int_t kk=0;kk<fScaleSize;kk++) {
728 cellcharge += fHitMap2->GetSignal(index,fScaleSize*it+kk);
731 (*pdigit).fSignal = (Int_t)cellcharge;
732 (*pdigit).fPhysics += phys;
733 // update list of tracks
734 TObjArray* trlist = (TObjArray*)pdigit->TrackList();
735 Int_t lastentry = trlist->GetLast();
736 TVector *ptrkp = (TVector*)trlist->At(lastentry);
737 TVector &trinfo = *ptrkp;
738 Int_t lasttrack = Int_t(trinfo(0));
739 Float_t lastcharge=(trinfo(2));
740 if (lasttrack==idtrack ) {
741 lastcharge += (Float_t)timeAmplitude;
742 trlist->RemoveAt(lastentry);
743 trinfo(0) = lasttrack;
745 trinfo(2) = lastcharge;
746 trlist->AddAt(&trinfo,lastentry);
748 new((*padr)[countadr++]) TVector(3);
749 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
750 trinfo(0) = (Float_t)idtrack;
751 trinfo(1) = (Float_t)idhit;
752 trinfo(2) = (Float_t)timeAmplitude;
753 trlist->Add(&trinfo);
754 } // end if lasttrack==idtrack
757 // check the track list - debugging
758 Int_t trk[20], htrk[20];
760 Int_t nptracks = trlist->GetEntriesFast();
763 for (tr=0;tr<nptracks;tr++) {
764 TVector *pptrkp = (TVector*)trlist->At(tr);
765 TVector &pptrk = *pptrkp;
766 trk[tr] = Int_t(pptrk(0));
767 htrk[tr] = Int_t(pptrk(1));
768 chtrk[tr] = (pptrk(2));
769 cout << "nptracks "<<nptracks << endl;
775 // update counter and countadr for next call.
779 //____________________________________________
781 void AliITSsimulationSDD::AddDigit(Int_t i, Int_t j, Int_t signal){
783 // tag with -1 signals coming from background tracks
784 // tag with -2 signals coming from pure electronic noise
786 Int_t digits[3], tracks[3], hits[3];
787 Float_t phys, charges[3];
789 Int_t trk[20], htrk[20];
792 Bool_t do10to8=fResponse->Do10to8();
794 if(do10to8) signal=Convert8to10(signal);
795 AliITSTransientDigit *obj = (AliITSTransientDigit*)fHitMap1->GetHit(i,j);
807 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
810 TObjArray* trlist=(TObjArray*)obj->TrackList();
811 Int_t nptracks=trlist->GetEntriesFast();
813 Warning("AddDigit","nptracks=%d > 20 nptracks set to 20",nptracks);
815 } // end if nptracks > 20
817 for (tr=0;tr<nptracks;tr++) {
818 TVector &pp =*((TVector*)trlist->At(tr));
819 trk[tr]=Int_t(pp(0));
820 htrk[tr]=Int_t(pp(1));
824 SortTracks(trk,chtrk,htrk,nptracks);
825 } // end if nptracks > 1
828 for (i=0; i<nptracks; i++) {
833 for (i=nptracks; i<3; i++) {
839 for (i=0; i<3; i++) {
844 } // end if/else nptracks < 3
846 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
848 } // end if/else !obj
850 //______________________________________________________________________
851 void AliITSsimulationSDD::SortTracks(Int_t *tracks,Float_t *charges,
852 Int_t *hits,Int_t ntr){
853 // Sort the list of tracks contributing to a given digit
854 // Only the 3 most significant tracks are acctually sorted
855 // Loop over signals, only 3 times
859 Int_t idx[3] = {-3,-3,-3};
860 Float_t jch[3] = {-3,-3,-3};
861 Int_t jtr[3] = {-3,-3,-3};
862 Int_t jhit[3] = {-3,-3,-3};
865 if (ntr<3) imax = ntr;
871 if((i == 1 && j == idx[i-1] )
872 ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
873 if(charges[j] > qmax) {
876 } // end if charges[j]>qmax
880 jch[i] = charges[jmax];
881 jtr[i] = tracks[jmax];
882 jhit[i] = hits[jmax];
895 } // end if jtr[i] == -3
898 //______________________________________________________________________
899 void AliITSsimulationSDD::ChargeToSignal() {
900 // add baseline, noise, electronics and ADC saturation effects
902 char opt1[20], opt2[20];
903 fResponse->ParamOptions(opt1,opt2);
904 char *read = strstr(opt1,"file");
905 Float_t baseline, noise;
908 static Bool_t readfile=kTRUE;
909 //read baseline and noise from file
910 if (readfile) ReadBaseline();
912 } else fResponse->GetNoiseParam(noise,baseline);
916 Float_t maxadc = fResponse->MaxAdc();
918 for (i=0;i<fNofMaps;i++) {
919 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
920 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
921 fInZR[k] = fHitMap2->GetSignal(i,k);
922 contrib = (baseline + noise*gRandom->Gaus());
925 for(k=0; k<fMaxNofSamples; k++) {
926 Double_t newcont = 0.;
927 Double_t maxcont = 0.;
928 for(kk=0;kk<fScaleSize;kk++) {
929 newcont = fInZR[fScaleSize*k+kk];
930 if(newcont > maxcont) maxcont = newcont;
933 if (newcont >= maxadc) newcont = maxadc -1;
934 if(newcont >= baseline){
935 Warning("","newcont=%d>=baseline=%d",newcont,baseline);
938 fHitMap2->SetHit(i,k,newcont);
940 } // end for i loop over anodes
944 for (i=0;i<fNofMaps;i++) {
945 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
946 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
947 fInZR[k] = fHitMap2->GetSignal(i,k);
948 contrib = (baseline + noise*gRandom->Gaus());
952 FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1);
953 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
954 Double_t rw = fElectronics->GetTraFunReal(k);
955 Double_t iw = fElectronics->GetTraFunImag(k);
956 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
957 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
959 FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1);
960 for(k=0; k<fMaxNofSamples; k++) {
961 Double_t newcont1 = 0.;
962 Double_t maxcont1 = 0.;
963 for(kk=0;kk<fScaleSize;kk++) {
964 newcont1 = fOutZR[fScaleSize*k+kk];
965 if(newcont1 > maxcont1) maxcont1 = newcont1;
968 if (newcont1 >= maxadc) newcont1 = maxadc -1;
969 fHitMap2->SetHit(i,k,newcont1);
971 } // end for i loop over anodes
974 //______________________________________________________________________
975 void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline,
977 // Returns the Baseline for a particular anode.
978 baseline = fBaseline[i];
981 //______________________________________________________________________
982 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl,
984 // Returns the compression alogirthm parameters
985 Int_t size = fD.GetSize();
987 db=fD[i]; tl=fT1[i]; th=fT2[i];
989 if (size <= 2 && i>=fNofMaps/2) {
990 db=fD[1]; tl=fT1[1]; th=fT2[1];
992 db=fD[0]; tl=fT1[0]; th=fT2[0];
993 } // end if size <=2 && i>=fNofMaps/2
996 //______________________________________________________________________
997 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){
998 // returns the compression alogirthm parameters
999 Int_t size = fD.GetSize();
1002 db=fD[i]; tl=fT1[i];
1004 if (size <= 2 && i>=fNofMaps/2) {
1005 db=fD[1]; tl=fT1[1];
1007 db=fD[0]; tl=fT1[0];
1008 } // end if size <=2 && i>=fNofMaps/2
1009 } // end if size > 2
1011 //______________________________________________________________________
1012 void AliITSsimulationSDD::SetCompressParam(){
1013 // Sets the compression alogirthm parameters
1016 fResponse->GiveCompressParam(cp);
1017 for (i=0; i<2; i++) {
1024 //______________________________________________________________________
1025 void AliITSsimulationSDD::ReadBaseline(){
1026 // read baseline and noise from file - either a .root file and in this
1027 // case data should be organised in a tree with one entry for each
1028 // module => reading should be done accordingly
1029 // or a classic file and do smth. like this:
1030 // Read baselines and noise for SDD
1034 char input[100], base[100], param[100];
1037 fResponse->Filenames(input,base,param);
1040 filtmp = gSystem->ExpandPathName(fFileName.Data());
1041 FILE *bline = fopen(filtmp,"r");
1045 while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) {
1047 Error("ReadBaseline","Anode number not in increasing order!",
1050 } // end if pos != na+1
1056 Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",filtmp);
1063 //______________________________________________________________________
1064 Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) const {
1065 // To the 10 to 8 bit lossive compression.
1066 // code from Davide C. and Albert W.
1068 if (signal < 128) return signal;
1069 if (signal < 256) return (128+((signal-128)>>1));
1070 if (signal < 512) return (192+((signal-256)>>3));
1071 if (signal < 1024) return (224+((signal-512)>>4));
1074 //______________________________________________________________________
1075 Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) const {
1076 // Undo the lossive 10 to 8 bit compression.
1077 // code from Davide C. and Albert W.
1078 if (signal < 0 || signal > 255) {
1079 Warning("Convert8to10","out of range signal=%d",signal);
1081 } // end if signal <0 || signal >255
1083 if (signal < 128) return signal;
1085 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
1086 else return (128+((signal-128)<<1)+1);
1087 } // end if signal < 192
1089 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
1090 else return (256+((signal-192)<<3)+4);
1091 } // end if signal < 224
1092 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
1093 return (512+((signal-224)<<4)+7);
1095 //______________________________________________________________________
1096 AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){
1097 //Return the correct map.
1099 return ((i==0)? fHitMap1 : fHitMap2);
1101 //______________________________________________________________________
1102 void AliITSsimulationSDD::ZeroSuppression(const char *option) {
1103 // perform the zero suppresion
1105 if (strstr(option,"2D")) {
1106 //Init2D(); // activate if param change module by module
1108 } else if (strstr(option,"1D")) {
1109 //Init1D(); // activate if param change module by module
1111 } else StoreAllDigits();
1113 //______________________________________________________________________
1114 void AliITSsimulationSDD::Init2D(){
1115 // read in and prepare arrays: fD, fT1, fT2
1116 // savemu[nanodes], savesigma[nanodes]
1117 // read baseline and noise from file - either a .root file and in this
1118 // case data should be organised in a tree with one entry for each
1119 // module => reading should be done accordingly
1120 // or a classic file and do smth. like this ( code from Davide C. and
1122 // Read 2D zero-suppression parameters for SDD
1124 if (!strstr(fParam,"file")) return;
1126 Int_t na,pos,tempTh;
1128 Float_t *savemu = new Float_t [fNofMaps];
1129 Float_t *savesigma = new Float_t [fNofMaps];
1130 char input[100],basel[100],par[100];
1132 Int_t minval = fResponse->MinVal();
1134 fResponse->Filenames(input,basel,par);
1137 filtmp = gSystem->ExpandPathName(fFileName.Data());
1138 FILE *param = fopen(filtmp,"r");
1142 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1144 Error("Init2D","Anode number not in increasing order!",filtmp);
1146 } // end if pos != na+1
1148 savesigma[na] = sigma;
1149 if ((2.*sigma) < mu) {
1150 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1153 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1154 if (tempTh < 0) tempTh=0;
1156 tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval;
1157 if (tempTh < 0) tempTh=0;
1162 Error("Init2D","THE FILE %s DOES NOT EXIST !",filtmp);
1169 delete [] savesigma;
1171 //______________________________________________________________________
1172 void AliITSsimulationSDD::Compress2D(){
1173 // simple ITS cluster finder -- online zero-suppression conditions
1176 Int_t minval = fResponse->MinVal();
1177 Bool_t write = fResponse->OutputOption();
1178 Bool_t do10to8 = fResponse->Do10to8();
1179 Int_t nz, nl, nh, low, i, j;
1181 for (i=0; i<fNofMaps; i++) {
1182 CompressionParam(i,db,tl,th);
1187 for (j=0; j<fMaxNofSamples; j++) {
1188 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1189 signal -= db; // if baseline eq. is done here
1190 if (signal <= 0) {nz++; continue;}
1191 if ((signal - tl) < minval) low++;
1192 if ((signal - th) >= minval) {
1195 FindCluster(i,j,signal,minval,cond);
1197 ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)){
1198 if(do10to8) signal = Convert10to8(signal);
1199 AddDigit(i,j,signal);
1200 } // end if cond&&j&&()
1201 } else if ((signal - tl) >= minval) nl++;
1202 } // end for j loop time samples
1203 if (write) TreeB()->Fill(nz,nl,nh,low,i+1);
1204 } //end for i loop anodes
1208 sprintf(hname,"TNtuple%d_%d",fModule,fEvent);
1209 TreeB()->Write(hname);
1214 //______________________________________________________________________
1215 void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal,
1216 Int_t minval,Bool_t &cond){
1217 // Find clusters according to the online 2D zero-suppression algorithm
1218 Bool_t do10to8 = fResponse->Do10to8();
1219 Bool_t high = kFALSE;
1221 fHitMap2->FlagHit(i,j);
1223 // check the online zero-suppression conditions
1225 const Int_t kMaxNeighbours = 4;
1228 Int_t xList[kMaxNeighbours], yList[kMaxNeighbours];
1229 fSegmentation->Neighbours(i,j,&nn,xList,yList);
1231 for (in=0; in<nn; in++) {
1234 if (fHitMap2->TestHit(ix,iy)==kUnused) {
1235 CompressionParam(ix,dbx,tlx,thx);
1236 Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy));
1237 qn -= dbx; // if baseline eq. is done here
1238 if ((qn-tlx) < minval) {
1239 fHitMap2->FlagHit(ix,iy);
1242 if ((qn - thx) >= minval) high=kTRUE;
1244 if(do10to8) signal = Convert10to8(signal);
1245 AddDigit(i,j,signal);
1247 if(do10to8) qns = Convert10to8(qn);
1249 if (!high) AddDigit(ix,iy,qns);
1251 if(!high) fHitMap2->FlagHit(ix,iy);
1252 } // end if qn-tlx < minval
1254 } // end for in loop over neighbours
1256 //______________________________________________________________________
1257 void AliITSsimulationSDD::Init1D(){
1258 // this is just a copy-paste of input taken from 2D algo
1259 // Torino people should give input
1260 // Read 1D zero-suppression parameters for SDD
1262 if (!strstr(fParam,"file")) return;
1264 Int_t na,pos,tempTh;
1266 Float_t *savemu = new Float_t [fNofMaps];
1267 Float_t *savesigma = new Float_t [fNofMaps];
1268 char input[100],basel[100],par[100];
1270 Int_t minval = fResponse->MinVal();
1272 fResponse->Filenames(input,basel,par);
1275 // set first the disable and tol param
1278 filtmp = gSystem->ExpandPathName(fFileName.Data());
1279 FILE *param = fopen(filtmp,"r");
1283 fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]);
1284 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1286 Error("Init1D","Anode number not in increasing order!",filtmp);
1288 } // end if pos != na+1
1290 savesigma[na]=sigma;
1291 if ((2.*sigma) < mu) {
1292 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1295 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1296 if (tempTh < 0) tempTh=0;
1301 Error("Init1D","THE FILE %s DOES NOT EXIST !",filtmp);
1308 delete [] savesigma;
1310 //______________________________________________________________________
1311 void AliITSsimulationSDD::Compress1D(){
1312 // 1D zero-suppression algorithm (from Gianluca A.)
1313 Int_t dis,tol,thres,decr,diff;
1314 UChar_t *str=fStream->Stream();
1316 Bool_t do10to8=fResponse->Do10to8();
1320 for (k=0; k<2; k++) {
1323 for (i=0; i<fNofMaps/2; i++) {
1324 Bool_t firstSignal=kTRUE;
1325 Int_t idx=i+k*fNofMaps/2;
1326 CompressionParam(idx,decr,thres);
1327 for (j=0; j<fMaxNofSamples; j++) {
1328 Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j));
1329 signal -= decr; // if baseline eq.
1330 if(do10to8) signal = Convert10to8(signal);
1331 if (signal <= thres) {
1335 // write diff in the buffer for HuffT
1336 str[counter]=(UChar_t)diff;
1339 } // end if signal <= thres
1341 if (diff > 127) diff=127;
1342 if (diff < -128) diff=-128;
1344 // tol has changed to 8 possible cases ? - one can write
1345 // this if(TMath::Abs(diff)<tol) ... else ...
1346 if(TMath::Abs(diff)<tol) diff=0;
1347 // or keep it as it was before
1349 if (tol==1 && (diff >= -2 && diff <= 1)) diff=0;
1350 if (tol==2 && (diff >= -4 && diff <= 3)) diff=0;
1351 if (tol==3 && (diff >= -16 && diff <= 15)) diff=0;
1353 AddDigit(idx,j,last+diff);
1355 AddDigit(idx,j,signal);
1356 } // end if singal < dis
1358 // write diff in the buffer used to compute Huffman tables
1359 if (firstSignal) str[counter]=(UChar_t)signal;
1360 else str[counter]=(UChar_t)diff;
1364 } // end for j loop time samples
1365 } // end for i loop anodes one half of detector
1369 fStream->CheckCount(counter);
1371 // open file and write out the stream of diff's
1372 static Bool_t open=kTRUE;
1373 static TFile *outFile;
1374 Bool_t write = fResponse->OutputOption();
1378 SetFileName("stream.root");
1379 cout<<"filename "<<fFileName<<endl;
1380 outFile=new TFile(fFileName,"recreate");
1381 cout<<"I have opened "<<fFileName<<" file "<<endl;
1388 fStream->ClearStream();
1390 // back to galice.root file
1392 TTree *fAli=gAlice->TreeK();
1395 if (fAli) file =fAli->GetCurrentFile();
1398 //______________________________________________________________________
1399 void AliITSsimulationSDD::StoreAllDigits(){
1400 // if non-zero-suppressed data
1401 Bool_t do10to8 = fResponse->Do10to8();
1402 Int_t i, j, digits[3];
1404 for (i=0; i<fNofMaps; i++) {
1405 for (j=0; j<fMaxNofSamples; j++) {
1406 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1407 if(do10to8) signal = Convert10to8(signal);
1408 if(do10to8) signal = Convert8to10(signal);
1412 fITS->AddRealDigit(1,digits);
1416 //______________________________________________________________________
1417 void AliITSsimulationSDD::CreateHistograms(Int_t scale){
1418 // Creates histograms of maps for debugging
1421 fHis=new TObjArray(fNofMaps);
1422 for (i=0;i<fNofMaps;i++) {
1423 TString sddName("sdd_");
1425 sprintf(candNum,"%d",i+1);
1426 sddName.Append(candNum);
1427 fHis->AddAt(new TH1F(sddName.Data(),"SDD maps",scale*fMaxNofSamples,
1428 0.,(Float_t) scale*fMaxNofSamples), i);
1431 //______________________________________________________________________
1432 void AliITSsimulationSDD::FillHistograms(){
1433 // fill 1D histograms from map
1437 for( Int_t i=0; i<fNofMaps; i++) {
1438 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
1439 Int_t nsamples = hist->GetNbinsX();
1440 for( Int_t j=0; j<nsamples; j++) {
1441 Double_t signal=fHitMap2->GetSignal(i,j);
1442 hist->Fill((Float_t)j,signal);
1446 //______________________________________________________________________
1447 void AliITSsimulationSDD::ResetHistograms(){
1448 // Reset histograms for this detector
1451 for (i=0;i<fNofMaps;i++ ) {
1452 if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset();
1455 //______________________________________________________________________
1456 TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
1457 // Fills a histogram from a give anode.
1459 if (!fHis) return 0;
1461 if(wing <=0 || wing > 2) {
1462 Warning("GetAnode","Wrong wing number: %d",wing);
1464 } // end if wing <=0 || wing >2
1465 if(anode <=0 || anode > fNofMaps/2) {
1466 Warning("GetAnode","Wrong anode number: %d",anode);
1468 } // end if ampde <=0 || andoe > fNofMaps/2
1470 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
1471 return (TH1F*)(fHis->At(index));
1473 //______________________________________________________________________
1474 void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
1475 // Writes the histograms to a file
1481 for(i=0; i<fNofMaps; i++) fHis->At(i)->Write(); //fAdcs[i]->Write();
1484 //______________________________________________________________________
1485 Float_t AliITSsimulationSDD::GetNoise() {
1486 // Returns the noise value
1487 //Bool_t do10to8=fResponse->Do10to8();
1488 //noise will always be in the liniar part of the signal
1490 Int_t threshold = fT1[0];
1491 char opt1[20], opt2[20];
1493 fResponse->ParamOptions(opt1,opt2);
1495 char *same = strstr(opt1,"same");
1496 Float_t noise,baseline;
1498 fResponse->GetNoiseParam(noise,baseline);
1500 static Bool_t readfile=kTRUE;
1501 //read baseline and noise from file
1502 if (readfile) ReadBaseline();
1506 TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2");
1507 if(c2) delete c2->GetPrimitive("noisehist");
1508 if(c2) delete c2->GetPrimitive("anode");
1509 else c2=new TCanvas("c2");
1511 c2->SetFillColor(0);
1513 TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold);
1514 TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0.,
1515 (float)fMaxNofSamples);
1517 for (i=0;i<fNofMaps;i++) {
1518 CompressionParam(i,decr,threshold);
1519 if (!same) GetAnodeBaseline(i,baseline,noise);
1521 for (k=0;k<fMaxNofSamples;k++) {
1522 Float_t signal=(Float_t)fHitMap2->GetSignal(i,k);
1523 //if (signal <= (float)threshold) noisehist->Fill(signal-baseline);
1524 if (signal <= (float)threshold) noisehist->Fill(signal);
1525 anode->Fill((float)k,signal);
1530 TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold);
1531 noisehist->Fit("gnoise","RQ");
1534 Float_t mnoise = gnoise->GetParameter(1);
1535 cout << "mnoise : " << mnoise << endl;
1536 Float_t rnoise = gnoise->GetParameter(2);
1537 cout << "rnoise : " << rnoise << endl;
1540 }//______________________________________________________________________
1541 void AliITSsimulationSDD::WriteSDigits(AliITSpList *pList){
1542 // Fills the Summable digits Tree
1544 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
1546 pList->GetMaxMapIndex(ni,nj);
1547 for(i=0;i<ni;i++)for(j=0;j<nj;j++){
1548 if(pList->GetSignalOnly(i,j)>0.5*fT1[0]){ // above small threshold.
1549 aliITS->AddSumDigit(*(pList->GetpListItem(i,j)));
1550 // cout << "pListSDD: " << *(pList->GetpListItem(i,j)) << endl;
1555 //______________________________________________________________________
1556 void AliITSsimulationSDD::Print() {
1557 // Print SDD simulation Parameters
1559 cout << "**************************************************" << endl;
1560 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1561 cout << "**************************************************" << endl;
1562 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1563 cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl;
1564 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1565 cout << "Number pf Anodes used: " << fNofMaps << endl;
1566 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1567 cout << "Scale size factor: " << fScaleSize << endl;
1568 cout << "**************************************************" << endl;