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);
353 fHitMap1->ClearMap();
354 fHitMap2->ClearMap();
356 //______________________________________________________________________
357 void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
358 // create maps to build the lists of tracks for each digit
360 TObjArray *fHits = mod->GetHits();
361 Int_t nhits = fHits->GetEntriesFast();
365 if (!nhits && fCheckNoise) {
368 fHitMap2->ClearMap();
370 } else if (!nhits) return;
372 AliITSpList *pList = new AliITSpList(2*fSegmentation->Npz(),
373 fScaleSize*fSegmentation->Npx());
375 // inputs to ListOfFiredCells.
376 TObjArray *alist = new TObjArray();
377 fHitMap1->SetArray(alist);
378 static TClonesArray *padr = 0;
379 if(!padr) padr = new TClonesArray("TVector",1000);
381 HitsToAnalogDigits(mod,alist,padr,pList);
390 fHitMap1->ClearMap();
391 fHitMap2->ClearMap();
393 //______________________________________________________________________
394 void AliITSsimulationSDD::SDigitsToDigits(AliITSpList *pList){
395 // Take Summable digits and create Digits.
397 // inputs to ListOfFiredCells.
398 TObjArray *alist = new TObjArray();
399 fHitMap1->SetArray(alist);
400 static TClonesArray *padr = 0;
401 if(!padr) padr = new TClonesArray("TVector",1000);
402 Int_t arg[6] = {0,0,0,0,0,0};
403 Double_t timeAmplitude;
406 // Fill maps from pList.
407 for(i=0;i<pList->GetMaxIndex();i++){
408 pList->GetMapIndex(i,arg[0],arg[1]);
409 for(j=0;j<pList->GetNEnteries();j++){
410 timeAmplitude = pList->GetTSignal(arg[0],arg[1],j);
411 if(timeAmplitude>0.0) continue;
412 arg[2] = pList->GetTrack(arg[0],arg[1],j);
413 arg[3] = pList->GetHit(arg[0],arg[1],j);
414 ListOfFiredCells(arg,timeAmplitude,alist,padr);
416 // Make sure map has full signal in it.
417 fHitMap2->SetHit(arg[0],arg[1],pList->GetSignal(arg[0],arg[1]));
426 fHitMap1->ClearMap();
427 fHitMap2->ClearMap();
429 //______________________________________________________________________
430 void AliITSsimulationSDD::FinishDigits(TObjArray *alist){
431 // introduce the electronics effects and do zero-suppression if required
432 Int_t nentries=alist->GetEntriesFast();
434 if(!nentries) return;
436 const char *kopt=fResponse->ZeroSuppOption();
437 ZeroSuppression(kopt);
439 //______________________________________________________________________
440 void AliITSsimulationSDD::HitsToAnalogDigits(AliITSmodule *mod,TObjArray *alst,
443 // create maps to build the lists of tracks for each digit
445 TObjArray *fHits = mod->GetHits();
446 Int_t nhits = fHits->GetEntriesFast();
447 Int_t arg[6] = {0,0,0,0,0,0};
449 Int_t nofAnodes = fNofMaps/2;
450 Float_t sddLength = fSegmentation->Dx();
451 Float_t sddWidth = fSegmentation->Dz();
452 Float_t anodePitch = fSegmentation->Dpz(dummy);
453 Float_t timeStep = fSegmentation->Dpx(dummy);
454 Float_t driftSpeed = fResponse->DriftSpeed();
455 Float_t maxadc = fResponse->MaxAdc();
456 Float_t topValue = fResponse->DynamicRange();
457 Float_t cHloss = fResponse->ChargeLoss();
458 Float_t norm = maxadc/topValue;
459 Float_t dfCoeff, s1; fResponse->DiffCoeff(dfCoeff,s1); // Signal 2d Shape
460 Double_t eVpairs = 3.6; // electron pair energy eV.
461 Float_t nsigma = fResponse->NSigmaIntegration(); //
462 Int_t nlookups = fResponse->GausNLookUp(); //
464 // Piergiorgio's part (apart for few variables which I made float
465 // when i thought that can be done
466 // Fill detector maps with GEANT hits
467 // loop over hits in the module
469 const Float_t kconv = 1.0e+6; // GeV->KeV
471 Int_t hitDetector; // detector number (lay,lad,hitDetector)
472 Int_t iWing; // which detector wing/side.
473 Int_t detector; // 2*(detector-1)+iWing
474 Int_t ii,kk,ka,kt; // loop indexs
475 Int_t ia,it,index; // sub-pixel integration indexies
476 Int_t iAnode; // anode number.
477 Int_t timeSample; // time buckett.
478 Int_t anodeWindow; // anode direction charge integration width
479 Int_t timeWindow; // time direction charge integration width
480 Int_t jamin,jamax; // anode charge integration window
481 Int_t jtmin,jtmax; // time charge integration window
482 Int_t ndiv; // Anode window division factor.
483 Int_t nsplit; // the number of splits in anode and time windows==1.
484 Int_t nOfSplits; // number of times track length is split into
485 Float_t nOfSplitsF; // Floating point version of nOfSplits.
486 Float_t kkF; // Floating point version of loop index kk.
487 Float_t pathInSDD; // Track length in SDD.
488 Float_t drPath; // average position of track in detector. in microns
489 Float_t drTime; // Drift time
490 Float_t nmul; // drift time window multiplication factor.
491 Float_t avDrft; // x position of path length segment in cm.
492 Float_t avAnode; // Anode for path length segment in Anode number (float)
493 Float_t xAnode; // Floating point anode number.
494 Float_t driftPath; // avDrft in microns.
495 Float_t width; // width of signal at anodes.
496 Double_t depEnergy; // Energy deposited in this GEANT step.
497 Double_t xL[3],dxL[3]; // local hit coordinates and diff.
498 Double_t sigA; // sigma of signal at anode.
499 Double_t sigT; // sigma in time/drift direction for track segment
500 Double_t aStep,aConst; // sub-pixel size and offset anode
501 Double_t tStep,tConst; // sub-pixel size and offset time
502 Double_t amplitude; // signal amplitude for track segment in nanoAmpere
503 Double_t chargeloss; // charge loss for track segment.
504 Double_t anodeAmplitude; // signal amplitude in anode direction
505 Double_t aExpo; // exponent of Gaussian anode direction
506 Double_t timeAmplitude; // signal amplitude in time direction
507 Double_t tExpo; // exponent of Gaussian time direction
508 // Double_t tof; // Time of flight in ns of this step.
510 for(ii=0; ii<nhits; ii++) {
511 if(!mod->LineSegmentL(ii,xL[0],dxL[0],xL[1],dxL[1],xL[2],dxL[2],
512 depEnergy,itrack)) continue;
514 hitDetector = mod->GetDet();
515 //tof = 1.E+09*(mod->GetHit(ii)->GetTOF()); // tof in ns.
516 //if(tof>sddLength/driftSpeed) continue; // hit happed too late.
518 // scale path to simulate a perpendicular track
519 // continue if the particle did not lose energy
520 // passing through detector
522 Warning("HitsToAnalogDigits",
523 "fTrack = %d hit=%d module=%d This particle has"
524 " passed without losing energy!",
525 itrack,ii,mod->GetIndex());
527 } // end if !depEnergy
529 pathInSDD = TMath::Sqrt(dxL[0]*dxL[0]+dxL[1]*dxL[1]+dxL[2]*dxL[2]);
531 if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
532 drPath = 10000.*(dxL[0]+2.*xL[0])*0.5;
533 if(drPath < 0) drPath = -drPath;
534 drPath = sddLength-drPath;
536 Warning("HitsToAnalogDigits",
537 "negative drift path drPath=%e sddLength=%e dxL[0]=%e "
539 drPath,sddLength,dxL[0],xL[0]);
541 } // end if drPath < 0
543 // Compute number of segments to brake step path into
544 drTime = drPath/driftSpeed; // Drift Time
545 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);// Sigma along the anodes
546 // calcuate the number of time the path length should be split into.
547 nOfSplits = (Int_t) (1. + 10000.*pathInSDD/sigA);
548 if(fFlag) nOfSplits = 1;
550 // loop over path segments, init. some variables.
551 depEnergy /= nOfSplits;
552 nOfSplitsF = (Float_t) nOfSplits;
553 for(kk=0;kk<nOfSplits;kk++) { // loop over path segments
554 kkF = (Float_t) kk + 0.5;
555 avDrft = xL[0]+dxL[0]*kkF/nOfSplitsF;
556 avAnode = xL[2]+dxL[2]*kkF/nOfSplitsF;
557 driftPath = 10000.*avDrft;
559 iWing = 2; // Assume wing is 2
560 if(driftPath < 0) { // if wing is not 2 it is 1.
562 driftPath = -driftPath;
563 } // end if driftPath < 0
564 driftPath = sddLength-driftPath;
565 detector = 2*(hitDetector-1) + iWing;
567 Warning("HitsToAnalogDigits","negative drift path "
568 "driftPath=%e sddLength=%e avDrft=%e dxL[0]=%e "
569 "xL[0]=%e",driftPath,sddLength,avDrft,dxL[0],xL[0]);
571 } // end if driftPath < 0
574 drTime = driftPath/driftSpeed; // drift time for segment.
575 timeSample = (Int_t) (fScaleSize*drTime/timeStep + 1);
576 // compute time Sample including tof information. The tof only
577 // effects the time of the signal is recoreded and not the
579 // timeSample = (Int_t) (fScaleSize*(drTime+tof)/timeStep + 1);
580 if(timeSample > fScaleSize*fMaxNofSamples) {
581 Warning("HItsToAnalogDigits","Wrong Time Sample: %e",
584 } // end if timeSample > fScaleSize*fMaxNoofSamples
587 xAnode = 10000.*(avAnode)/anodePitch + nofAnodes/2; // +1?
588 if(xAnode*anodePitch > sddWidth || xAnode*anodePitch < 0.)
589 Warning("HitsToAnalogDigits",
590 "Exceedubg sddWidth=%e Z = %e",
591 sddWidth,xAnode*anodePitch);
592 iAnode = (Int_t) (1.+xAnode); // xAnode?
593 if(iAnode < 1 || iAnode > nofAnodes) {
594 Warning("HitToAnalogDigits","Wrong iAnode: 1<%d>%d",
597 } // end if iAnode < 1 || iAnode > nofAnodes
599 // store straight away the particle position in the array
600 // of particles and take idhit=ii only when part is entering (this
601 // requires FillModules() in the macro for analysis) :
603 // Sigma along the anodes for track segment.
604 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);
605 sigT = sigA/driftSpeed;
606 // Peak amplitude in nanoAmpere
607 amplitude = fScaleSize*160.*depEnergy/
608 (timeStep*eVpairs*2.*acos(-1.)*sigT*sigA);
609 amplitude *= timeStep/25.; // WARNING!!!!! Amplitude scaling to
610 // account for clock variations
611 // (reference value: 40 MHz)
612 chargeloss = 1.-cHloss*driftPath/1000;
613 amplitude *= chargeloss;
614 width = 2.*nsigma/(nlookups-1);
622 } // end if drTime > 1200.
624 nsplit = 4; // hard-wired //nsplit=4;nsplit = (nsplit+1)/2*2;
625 // Sub-pixel size see computation of aExpo and tExpo.
626 aStep = anodePitch/(nsplit*fScaleSize*sigA);
627 aConst = xAnode*anodePitch/sigA;
628 tStep = timeStep/(nsplit*fScaleSize*sigT);
629 tConst = drTime/sigT;
630 // Define SDD window corresponding to the hit
631 anodeWindow = (Int_t)(fScaleSize*nsigma*sigA/anodePitch+1);
632 timeWindow = (Int_t) (fScaleSize*nsigma*sigT/timeStep+1.);
633 jamin = (iAnode - anodeWindow/ndiv - 1)*fScaleSize*nsplit +1;
634 jamax = (iAnode + anodeWindow/ndiv)*fScaleSize*nsplit;
635 if(jamin <= 0) jamin = 1;
636 if(jamax > fScaleSize*nofAnodes*nsplit)
637 jamax = fScaleSize*nofAnodes*nsplit;
638 // jtmin and jtmax are Hard-wired
639 jtmin = (Int_t)(timeSample-timeWindow*nmul-1)*nsplit+1;
640 jtmax = (Int_t)(timeSample+timeWindow*nmul)*nsplit;
641 if(jtmin <= 0) jtmin = 1;
642 if(jtmax > fScaleSize*fMaxNofSamples*nsplit)
643 jtmax = fScaleSize*fMaxNofSamples*nsplit;
644 // Spread the charge in the anode-time window
645 for(ka=jamin; ka <=jamax; ka++) {
646 ia = (ka-1)/(fScaleSize*nsplit) + 1;
648 Warning("HitsToAnalogDigits","ia < 1: ");
651 if(ia > nofAnodes) ia = nofAnodes;
652 aExpo = (aStep*(ka-0.5)-aConst);
653 if(TMath::Abs(aExpo) > nsigma) anodeAmplitude = 0.;
655 dummy = (Int_t) ((aExpo+nsigma)/width);
656 anodeAmplitude = amplitude*fResponse->GausLookUp(dummy);
657 } // end if TMath::Abs(aEspo) > nsigma
658 // index starts from 0
659 index = ((detector+1)%2)*nofAnodes+ia-1;
660 if(anodeAmplitude) for(kt=jtmin; kt<=jtmax; kt++) {
661 it = (kt-1)/nsplit+1; // it starts from 1
663 Warning("HitsToAnalogDigits","it < 1:");
666 if(it>fScaleSize*fMaxNofSamples)
667 it = fScaleSize*fMaxNofSamples;
668 tExpo = (tStep*(kt-0.5)-tConst);
669 if(TMath::Abs(tExpo) > nsigma) timeAmplitude = 0.;
671 dummy = (Int_t) ((tExpo+nsigma)/width);
672 timeAmplitude = anodeAmplitude*
673 fResponse->GausLookUp(dummy);
674 } // end if TMath::Abs(tExpo) > nsigma
675 // build the list of digits for this module
678 arg[2] = itrack; // track number
679 arg[3] = ii-1; // hit number.
680 timeAmplitude *= norm;
682 ListOfFiredCells(arg,timeAmplitude,alst,padr);
683 pList->AddSignal(index,it,itrack,ii-1,
684 mod->GetIndex(),timeAmplitude);
685 } // end if anodeAmplitude and loop over time in window
686 } // loop over anodes in window
687 } // end loop over "sub-hits"
688 } // end loop over hits
690 //______________________________________________________________________
691 void AliITSsimulationSDD::ListOfFiredCells(Int_t *arg,Double_t timeAmplitude,
692 TObjArray *alist,TClonesArray *padr){
693 // Returns the list of "fired" cells.
695 Int_t index = arg[0];
697 Int_t idtrack = arg[2];
698 Int_t idhit = arg[3];
699 Int_t counter = arg[4];
700 Int_t countadr = arg[5];
701 Double_t charge = timeAmplitude;
702 charge += fHitMap2->GetSignal(index,ik-1);
703 fHitMap2->SetHit(index, ik-1, charge);
706 Int_t it = (Int_t)((ik-1)/fScaleSize);
709 digits[2] = (Int_t)timeAmplitude;
711 if (idtrack >= 0) phys = (Float_t)timeAmplitude;
714 Double_t cellcharge = 0.;
715 AliITSTransientDigit* pdigit;
716 // build the list of fired cells and update the info
717 if (!fHitMap1->TestHit(index, it)) {
718 new((*padr)[countadr++]) TVector(3);
719 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
720 trinfo(0) = (Float_t)idtrack;
721 trinfo(1) = (Float_t)idhit;
722 trinfo(2) = (Float_t)timeAmplitude;
724 alist->AddAtAndExpand(new AliITSTransientDigit(phys,digits),counter);
725 fHitMap1->SetHit(index, it, counter);
727 pdigit=(AliITSTransientDigit*)alist->At(alist->GetLast());
729 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
730 trlist->Add(&trinfo);
732 pdigit = (AliITSTransientDigit*) fHitMap1->GetHit(index, it);
733 for(Int_t kk=0;kk<fScaleSize;kk++) {
734 cellcharge += fHitMap2->GetSignal(index,fScaleSize*it+kk);
737 (*pdigit).fSignal = (Int_t)cellcharge;
738 (*pdigit).fPhysics += phys;
739 // update list of tracks
740 TObjArray* trlist = (TObjArray*)pdigit->TrackList();
741 Int_t lastentry = trlist->GetLast();
742 TVector *ptrkp = (TVector*)trlist->At(lastentry);
743 TVector &trinfo = *ptrkp;
744 Int_t lasttrack = Int_t(trinfo(0));
745 Float_t lastcharge=(trinfo(2));
746 if (lasttrack==idtrack ) {
747 lastcharge += (Float_t)timeAmplitude;
748 trlist->RemoveAt(lastentry);
749 trinfo(0) = lasttrack;
751 trinfo(2) = lastcharge;
752 trlist->AddAt(&trinfo,lastentry);
754 new((*padr)[countadr++]) TVector(3);
755 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
756 trinfo(0) = (Float_t)idtrack;
757 trinfo(1) = (Float_t)idhit;
758 trinfo(2) = (Float_t)timeAmplitude;
759 trlist->Add(&trinfo);
760 } // end if lasttrack==idtrack
763 // check the track list - debugging
764 Int_t trk[20], htrk[20];
766 Int_t nptracks = trlist->GetEntriesFast();
769 for (tr=0;tr<nptracks;tr++) {
770 TVector *pptrkp = (TVector*)trlist->At(tr);
771 TVector &pptrk = *pptrkp;
772 trk[tr] = Int_t(pptrk(0));
773 htrk[tr] = Int_t(pptrk(1));
774 chtrk[tr] = (pptrk(2));
775 cout << "nptracks "<<nptracks << endl;
781 // update counter and countadr for next call.
785 //____________________________________________
787 void AliITSsimulationSDD::AddDigit(Int_t i, Int_t j, Int_t signal){
789 // tag with -1 signals coming from background tracks
790 // tag with -2 signals coming from pure electronic noise
792 Int_t digits[3], tracks[3], hits[3];
793 Float_t phys, charges[3];
795 Int_t trk[20], htrk[20];
798 Bool_t do10to8=fResponse->Do10to8();
800 if(do10to8) signal=Convert8to10(signal);
801 AliITSTransientDigit *obj = (AliITSTransientDigit*)fHitMap1->GetHit(i,j);
813 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
816 TObjArray* trlist=(TObjArray*)obj->TrackList();
817 Int_t nptracks=trlist->GetEntriesFast();
819 Warning("AddDigit","nptracks=%d > 20 nptracks set to 20",nptracks);
821 } // end if nptracks > 20
823 for (tr=0;tr<nptracks;tr++) {
824 TVector &pp =*((TVector*)trlist->At(tr));
825 trk[tr]=Int_t(pp(0));
826 htrk[tr]=Int_t(pp(1));
830 SortTracks(trk,chtrk,htrk,nptracks);
831 } // end if nptracks > 1
834 for (i=0; i<nptracks; i++) {
839 for (i=nptracks; i<3; i++) {
845 for (i=0; i<3; i++) {
850 } // end if/else nptracks < 3
852 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
854 } // end if/else !obj
856 //______________________________________________________________________
857 void AliITSsimulationSDD::SortTracks(Int_t *tracks,Float_t *charges,
858 Int_t *hits,Int_t ntr){
859 // Sort the list of tracks contributing to a given digit
860 // Only the 3 most significant tracks are acctually sorted
861 // Loop over signals, only 3 times
865 Int_t idx[3] = {-3,-3,-3};
866 Float_t jch[3] = {-3,-3,-3};
867 Int_t jtr[3] = {-3,-3,-3};
868 Int_t jhit[3] = {-3,-3,-3};
871 if (ntr<3) imax = ntr;
877 if((i == 1 && j == idx[i-1] )
878 ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
879 if(charges[j] > qmax) {
882 } // end if charges[j]>qmax
886 jch[i] = charges[jmax];
887 jtr[i] = tracks[jmax];
888 jhit[i] = hits[jmax];
901 } // end if jtr[i] == -3
904 //______________________________________________________________________
905 void AliITSsimulationSDD::ChargeToSignal() {
906 // add baseline, noise, electronics and ADC saturation effects
908 char opt1[20], opt2[20];
909 fResponse->ParamOptions(opt1,opt2);
910 char *read = strstr(opt1,"file");
911 Float_t baseline, noise;
914 static Bool_t readfile=kTRUE;
915 //read baseline and noise from file
916 if (readfile) ReadBaseline();
918 } else fResponse->GetNoiseParam(noise,baseline);
922 Float_t maxadc = fResponse->MaxAdc();
924 for (i=0;i<fNofMaps;i++) {
925 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
926 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
927 fInZR[k] = fHitMap2->GetSignal(i,k);
928 contrib = (baseline + noise*gRandom->Gaus());
931 for(k=0; k<fMaxNofSamples; k++) {
932 Double_t newcont = 0.;
933 Double_t maxcont = 0.;
934 for(kk=0;kk<fScaleSize;kk++) {
935 newcont = fInZR[fScaleSize*k+kk];
936 if(newcont > maxcont) maxcont = newcont;
939 if (newcont >= maxadc) newcont = maxadc -1;
940 if(newcont >= baseline){
941 Warning("","newcont=%d>=baseline=%d",newcont,baseline);
944 fHitMap2->SetHit(i,k,newcont);
946 } // end for i loop over anodes
950 for (i=0;i<fNofMaps;i++) {
951 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
952 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
953 fInZR[k] = fHitMap2->GetSignal(i,k);
954 contrib = (baseline + noise*gRandom->Gaus());
958 FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1);
959 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
960 Double_t rw = fElectronics->GetTraFunReal(k);
961 Double_t iw = fElectronics->GetTraFunImag(k);
962 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
963 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
965 FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1);
966 for(k=0; k<fMaxNofSamples; k++) {
967 Double_t newcont1 = 0.;
968 Double_t maxcont1 = 0.;
969 for(kk=0;kk<fScaleSize;kk++) {
970 newcont1 = fOutZR[fScaleSize*k+kk];
971 if(newcont1 > maxcont1) maxcont1 = newcont1;
974 if (newcont1 >= maxadc) newcont1 = maxadc -1;
975 fHitMap2->SetHit(i,k,newcont1);
977 } // end for i loop over anodes
980 //______________________________________________________________________
981 void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline,
983 // Returns the Baseline for a particular anode.
984 baseline = fBaseline[i];
987 //______________________________________________________________________
988 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl,
990 // Returns the compression alogirthm parameters
991 Int_t size = fD.GetSize();
993 db=fD[i]; tl=fT1[i]; th=fT2[i];
995 if (size <= 2 && i>=fNofMaps/2) {
996 db=fD[1]; tl=fT1[1]; th=fT2[1];
998 db=fD[0]; tl=fT1[0]; th=fT2[0];
999 } // end if size <=2 && i>=fNofMaps/2
1002 //______________________________________________________________________
1003 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){
1004 // returns the compression alogirthm parameters
1005 Int_t size = fD.GetSize();
1008 db=fD[i]; tl=fT1[i];
1010 if (size <= 2 && i>=fNofMaps/2) {
1011 db=fD[1]; tl=fT1[1];
1013 db=fD[0]; tl=fT1[0];
1014 } // end if size <=2 && i>=fNofMaps/2
1015 } // end if size > 2
1017 //______________________________________________________________________
1018 void AliITSsimulationSDD::SetCompressParam(){
1019 // Sets the compression alogirthm parameters
1022 fResponse->GiveCompressParam(cp);
1023 for (i=0; i<2; i++) {
1030 //______________________________________________________________________
1031 void AliITSsimulationSDD::ReadBaseline(){
1032 // read baseline and noise from file - either a .root file and in this
1033 // case data should be organised in a tree with one entry for each
1034 // module => reading should be done accordingly
1035 // or a classic file and do smth. like this:
1036 // Read baselines and noise for SDD
1040 char input[100], base[100], param[100];
1043 fResponse->Filenames(input,base,param);
1046 filtmp = gSystem->ExpandPathName(fFileName.Data());
1047 FILE *bline = fopen(filtmp,"r");
1051 while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) {
1053 Error("ReadBaseline","Anode number not in increasing order!",
1056 } // end if pos != na+1
1062 Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",filtmp);
1069 //______________________________________________________________________
1070 Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) const {
1071 // To the 10 to 8 bit lossive compression.
1072 // code from Davide C. and Albert W.
1074 if (signal < 128) return signal;
1075 if (signal < 256) return (128+((signal-128)>>1));
1076 if (signal < 512) return (192+((signal-256)>>3));
1077 if (signal < 1024) return (224+((signal-512)>>4));
1080 //______________________________________________________________________
1081 Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) const {
1082 // Undo the lossive 10 to 8 bit compression.
1083 // code from Davide C. and Albert W.
1084 if (signal < 0 || signal > 255) {
1085 Warning("Convert8to10","out of range signal=%d",signal);
1087 } // end if signal <0 || signal >255
1089 if (signal < 128) return signal;
1091 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
1092 else return (128+((signal-128)<<1)+1);
1093 } // end if signal < 192
1095 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
1096 else return (256+((signal-192)<<3)+4);
1097 } // end if signal < 224
1098 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
1099 return (512+((signal-224)<<4)+7);
1101 //______________________________________________________________________
1102 AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){
1103 //Return the correct map.
1105 return ((i==0)? fHitMap1 : fHitMap2);
1107 //______________________________________________________________________
1108 void AliITSsimulationSDD::ZeroSuppression(const char *option) {
1109 // perform the zero suppresion
1111 if (strstr(option,"2D")) {
1112 //Init2D(); // activate if param change module by module
1114 } else if (strstr(option,"1D")) {
1115 //Init1D(); // activate if param change module by module
1117 } else StoreAllDigits();
1119 //______________________________________________________________________
1120 void AliITSsimulationSDD::Init2D(){
1121 // read in and prepare arrays: fD, fT1, fT2
1122 // savemu[nanodes], savesigma[nanodes]
1123 // read baseline and noise from file - either a .root file and in this
1124 // case data should be organised in a tree with one entry for each
1125 // module => reading should be done accordingly
1126 // or a classic file and do smth. like this ( code from Davide C. and
1128 // Read 2D zero-suppression parameters for SDD
1130 if (!strstr(fParam,"file")) return;
1132 Int_t na,pos,tempTh;
1134 Float_t *savemu = new Float_t [fNofMaps];
1135 Float_t *savesigma = new Float_t [fNofMaps];
1136 char input[100],basel[100],par[100];
1138 Int_t minval = fResponse->MinVal();
1140 fResponse->Filenames(input,basel,par);
1143 filtmp = gSystem->ExpandPathName(fFileName.Data());
1144 FILE *param = fopen(filtmp,"r");
1148 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1150 Error("Init2D","Anode number not in increasing order!",filtmp);
1152 } // end if pos != na+1
1154 savesigma[na] = sigma;
1155 if ((2.*sigma) < mu) {
1156 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1159 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1160 if (tempTh < 0) tempTh=0;
1162 tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval;
1163 if (tempTh < 0) tempTh=0;
1168 Error("Init2D","THE FILE %s DOES NOT EXIST !",filtmp);
1175 delete [] savesigma;
1177 //______________________________________________________________________
1178 void AliITSsimulationSDD::Compress2D(){
1179 // simple ITS cluster finder -- online zero-suppression conditions
1182 Int_t minval = fResponse->MinVal();
1183 Bool_t write = fResponse->OutputOption();
1184 Bool_t do10to8 = fResponse->Do10to8();
1185 Int_t nz, nl, nh, low, i, j;
1187 for (i=0; i<fNofMaps; i++) {
1188 CompressionParam(i,db,tl,th);
1193 for (j=0; j<fMaxNofSamples; j++) {
1194 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1195 signal -= db; // if baseline eq. is done here
1196 if (signal <= 0) {nz++; continue;}
1197 if ((signal - tl) < minval) low++;
1198 if ((signal - th) >= minval) {
1201 FindCluster(i,j,signal,minval,cond);
1203 ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)){
1204 if(do10to8) signal = Convert10to8(signal);
1205 AddDigit(i,j,signal);
1206 } // end if cond&&j&&()
1207 } else if ((signal - tl) >= minval) nl++;
1208 } // end for j loop time samples
1209 if (write) TreeB()->Fill(nz,nl,nh,low,i+1);
1210 } //end for i loop anodes
1214 sprintf(hname,"TNtuple%d_%d",fModule,fEvent);
1215 TreeB()->Write(hname);
1220 //______________________________________________________________________
1221 void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal,
1222 Int_t minval,Bool_t &cond){
1223 // Find clusters according to the online 2D zero-suppression algorithm
1224 Bool_t do10to8 = fResponse->Do10to8();
1225 Bool_t high = kFALSE;
1227 fHitMap2->FlagHit(i,j);
1229 // check the online zero-suppression conditions
1231 const Int_t kMaxNeighbours = 4;
1234 Int_t xList[kMaxNeighbours], yList[kMaxNeighbours];
1235 fSegmentation->Neighbours(i,j,&nn,xList,yList);
1237 for (in=0; in<nn; in++) {
1240 if (fHitMap2->TestHit(ix,iy)==kUnused) {
1241 CompressionParam(ix,dbx,tlx,thx);
1242 Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy));
1243 qn -= dbx; // if baseline eq. is done here
1244 if ((qn-tlx) < minval) {
1245 fHitMap2->FlagHit(ix,iy);
1248 if ((qn - thx) >= minval) high=kTRUE;
1250 if(do10to8) signal = Convert10to8(signal);
1251 AddDigit(i,j,signal);
1253 if(do10to8) qns = Convert10to8(qn);
1255 if (!high) AddDigit(ix,iy,qns);
1257 if(!high) fHitMap2->FlagHit(ix,iy);
1258 } // end if qn-tlx < minval
1260 } // end for in loop over neighbours
1262 //______________________________________________________________________
1263 void AliITSsimulationSDD::Init1D(){
1264 // this is just a copy-paste of input taken from 2D algo
1265 // Torino people should give input
1266 // Read 1D zero-suppression parameters for SDD
1268 if (!strstr(fParam,"file")) return;
1270 Int_t na,pos,tempTh;
1272 Float_t *savemu = new Float_t [fNofMaps];
1273 Float_t *savesigma = new Float_t [fNofMaps];
1274 char input[100],basel[100],par[100];
1276 Int_t minval = fResponse->MinVal();
1278 fResponse->Filenames(input,basel,par);
1281 // set first the disable and tol param
1284 filtmp = gSystem->ExpandPathName(fFileName.Data());
1285 FILE *param = fopen(filtmp,"r");
1289 fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]);
1290 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1292 Error("Init1D","Anode number not in increasing order!",filtmp);
1294 } // end if pos != na+1
1296 savesigma[na]=sigma;
1297 if ((2.*sigma) < mu) {
1298 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1301 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1302 if (tempTh < 0) tempTh=0;
1307 Error("Init1D","THE FILE %s DOES NOT EXIST !",filtmp);
1314 delete [] savesigma;
1316 //______________________________________________________________________
1317 void AliITSsimulationSDD::Compress1D(){
1318 // 1D zero-suppression algorithm (from Gianluca A.)
1319 Int_t dis,tol,thres,decr,diff;
1320 UChar_t *str=fStream->Stream();
1322 Bool_t do10to8=fResponse->Do10to8();
1326 for (k=0; k<2; k++) {
1329 for (i=0; i<fNofMaps/2; i++) {
1330 Bool_t firstSignal=kTRUE;
1331 Int_t idx=i+k*fNofMaps/2;
1332 CompressionParam(idx,decr,thres);
1333 for (j=0; j<fMaxNofSamples; j++) {
1334 Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j));
1335 signal -= decr; // if baseline eq.
1336 if(do10to8) signal = Convert10to8(signal);
1337 if (signal <= thres) {
1341 // write diff in the buffer for HuffT
1342 str[counter]=(UChar_t)diff;
1345 } // end if signal <= thres
1347 if (diff > 127) diff=127;
1348 if (diff < -128) diff=-128;
1350 // tol has changed to 8 possible cases ? - one can write
1351 // this if(TMath::Abs(diff)<tol) ... else ...
1352 if(TMath::Abs(diff)<tol) diff=0;
1353 // or keep it as it was before
1355 if (tol==1 && (diff >= -2 && diff <= 1)) diff=0;
1356 if (tol==2 && (diff >= -4 && diff <= 3)) diff=0;
1357 if (tol==3 && (diff >= -16 && diff <= 15)) diff=0;
1359 AddDigit(idx,j,last+diff);
1361 AddDigit(idx,j,signal);
1362 } // end if singal < dis
1364 // write diff in the buffer used to compute Huffman tables
1365 if (firstSignal) str[counter]=(UChar_t)signal;
1366 else str[counter]=(UChar_t)diff;
1370 } // end for j loop time samples
1371 } // end for i loop anodes one half of detector
1375 fStream->CheckCount(counter);
1377 // open file and write out the stream of diff's
1378 static Bool_t open=kTRUE;
1379 static TFile *outFile;
1380 Bool_t write = fResponse->OutputOption();
1384 SetFileName("stream.root");
1385 cout<<"filename "<<fFileName<<endl;
1386 outFile=new TFile(fFileName,"recreate");
1387 cout<<"I have opened "<<fFileName<<" file "<<endl;
1394 fStream->ClearStream();
1396 // back to galice.root file
1398 TTree *fAli=gAlice->TreeK();
1401 if (fAli) file =fAli->GetCurrentFile();
1404 //______________________________________________________________________
1405 void AliITSsimulationSDD::StoreAllDigits(){
1406 // if non-zero-suppressed data
1407 Bool_t do10to8 = fResponse->Do10to8();
1408 Int_t i, j, digits[3];
1410 for (i=0; i<fNofMaps; i++) {
1411 for (j=0; j<fMaxNofSamples; j++) {
1412 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1413 if(do10to8) signal = Convert10to8(signal);
1414 if(do10to8) signal = Convert8to10(signal);
1418 fITS->AddRealDigit(1,digits);
1422 //______________________________________________________________________
1423 void AliITSsimulationSDD::CreateHistograms(Int_t scale){
1424 // Creates histograms of maps for debugging
1427 fHis=new TObjArray(fNofMaps);
1428 for (i=0;i<fNofMaps;i++) {
1429 TString sddName("sdd_");
1431 sprintf(candNum,"%d",i+1);
1432 sddName.Append(candNum);
1433 fHis->AddAt(new TH1F(sddName.Data(),"SDD maps",scale*fMaxNofSamples,
1434 0.,(Float_t) scale*fMaxNofSamples), i);
1437 //______________________________________________________________________
1438 void AliITSsimulationSDD::FillHistograms(){
1439 // fill 1D histograms from map
1443 for( Int_t i=0; i<fNofMaps; i++) {
1444 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
1445 Int_t nsamples = hist->GetNbinsX();
1446 for( Int_t j=0; j<nsamples; j++) {
1447 Double_t signal=fHitMap2->GetSignal(i,j);
1448 hist->Fill((Float_t)j,signal);
1452 //______________________________________________________________________
1453 void AliITSsimulationSDD::ResetHistograms(){
1454 // Reset histograms for this detector
1457 for (i=0;i<fNofMaps;i++ ) {
1458 if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset();
1461 //______________________________________________________________________
1462 TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
1463 // Fills a histogram from a give anode.
1465 if (!fHis) return 0;
1467 if(wing <=0 || wing > 2) {
1468 Warning("GetAnode","Wrong wing number: %d",wing);
1470 } // end if wing <=0 || wing >2
1471 if(anode <=0 || anode > fNofMaps/2) {
1472 Warning("GetAnode","Wrong anode number: %d",anode);
1474 } // end if ampde <=0 || andoe > fNofMaps/2
1476 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
1477 return (TH1F*)(fHis->At(index));
1479 //______________________________________________________________________
1480 void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
1481 // Writes the histograms to a file
1487 for(i=0; i<fNofMaps; i++) fHis->At(i)->Write(); //fAdcs[i]->Write();
1490 //______________________________________________________________________
1491 Float_t AliITSsimulationSDD::GetNoise() {
1492 // Returns the noise value
1493 //Bool_t do10to8=fResponse->Do10to8();
1494 //noise will always be in the liniar part of the signal
1496 Int_t threshold = fT1[0];
1497 char opt1[20], opt2[20];
1499 fResponse->ParamOptions(opt1,opt2);
1501 char *same = strstr(opt1,"same");
1502 Float_t noise,baseline;
1504 fResponse->GetNoiseParam(noise,baseline);
1506 static Bool_t readfile=kTRUE;
1507 //read baseline and noise from file
1508 if (readfile) ReadBaseline();
1512 TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2");
1513 if(c2) delete c2->GetPrimitive("noisehist");
1514 if(c2) delete c2->GetPrimitive("anode");
1515 else c2=new TCanvas("c2");
1517 c2->SetFillColor(0);
1519 TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold);
1520 TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0.,
1521 (float)fMaxNofSamples);
1523 for (i=0;i<fNofMaps;i++) {
1524 CompressionParam(i,decr,threshold);
1525 if (!same) GetAnodeBaseline(i,baseline,noise);
1527 for (k=0;k<fMaxNofSamples;k++) {
1528 Float_t signal=(Float_t)fHitMap2->GetSignal(i,k);
1529 //if (signal <= (float)threshold) noisehist->Fill(signal-baseline);
1530 if (signal <= (float)threshold) noisehist->Fill(signal);
1531 anode->Fill((float)k,signal);
1536 TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold);
1537 noisehist->Fit("gnoise","RQ");
1540 Float_t mnoise = gnoise->GetParameter(1);
1541 cout << "mnoise : " << mnoise << endl;
1542 Float_t rnoise = gnoise->GetParameter(2);
1543 cout << "rnoise : " << rnoise << endl;
1546 }//______________________________________________________________________
1547 void AliITSsimulationSDD::WriteSDigits(AliITSpList *pList){
1548 // Fills the Summable digits Tree
1550 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
1552 pList->GetMaxMapIndex(ni,nj);
1553 for(i=0;i<ni;i++)for(j=0;j<nj;j++){
1554 if(pList->GetSignalOnly(i,j)>0.5*fT1[0]){ // above small threshold.
1555 aliITS->AddSumDigit(*(pList->GetpListItem(i,j)));
1556 // cout << "pListSDD: " << *(pList->GetpListItem(i,j)) << endl;
1561 //______________________________________________________________________
1562 void AliITSsimulationSDD::Print() {
1563 // Print SDD simulation Parameters
1565 cout << "**************************************************" << endl;
1566 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1567 cout << "**************************************************" << endl;
1568 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1569 cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl;
1570 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1571 cout << "Number pf Anodes used: " << fNofMaps << endl;
1572 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1573 cout << "Scale size factor: " << fScaleSize << endl;
1574 cout << "**************************************************" << endl;