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 **************************************************************************/
24 #include <TStopwatch.h>
36 #include "AliITShit.h"
37 #include "AliITSdigit.h"
38 #include "AliITSmodule.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 "AliITSsimulationSDD.h"
48 ClassImp(AliITSsimulationSDD)
49 ////////////////////////////////////////////////////////////////////////
51 // Written by Piergiorgio Cerello
54 // AliITSsimulationSDD is the simulation of SDDs.
58 <img src="picts/ITS/AliITShit_Class_Diagram.gif">
61 <font size=+2 color=red>
62 <p>This show the relasionships between the ITS hit class and the rest of Aliroot.
67 //_____________________________________________________________________________
69 Int_t power(Int_t b, Int_t e) {
70 // compute b to the e power, where both b and e are Int_ts.
72 for(i=0; i<e; i++) power *= b;
76 //_____________________________________________
78 void FastFourierTransform(AliITSetfSDD *alisddetf,Double_t *real,
79 Double_t *imag,Int_t direction) {
80 // Do a Fast Fourier Transform
81 //printf("FFT: direction %d\n",direction);
83 Int_t samples = alisddetf->GetSamples();
84 Int_t l = (Int_t) ((log((Float_t) samples)/log(2.))+0.5);
87 Int_t m2 = samples/m1;
90 for(j=0; j<samples; j += m1) {
92 for(k=j; k<= j+m-1; k++) {
93 Double_t wsr = alisddetf->GetWeightReal(p);
94 Double_t wsi = alisddetf->GetWeightImag(p);
95 if(direction == -1) wsi = -wsi;
96 Double_t xr = *(real+k+m);
97 Double_t xi = *(imag+k+m);
98 *(real+k+m) = wsr*(*(real+k)-xr) - wsi*(*(imag+k)-xi);
99 *(imag+k+m) = wsr*(*(imag+k)-xi) + wsi*(*(real+k)-xr);
110 for(j=0; j<samples; j++) {
114 for(i1=1; i1<=l; i1++) {
117 p = p + p + j2 - j1 - j1;
120 Double_t xr = *(real+j);
121 Double_t xi = *(imag+j);
122 *(real+j) = *(real+p);
123 *(imag+j) = *(imag+p);
128 if(direction == -1) {
129 for(i=0; i<samples; i++) {
130 *(real+i) /= samples;
131 *(imag+i) /= samples;
136 //_____________________________________________________________________________
138 AliITSsimulationSDD::AliITSsimulationSDD(){
139 // Default constructor
155 SetPerpendTracksFlag();
167 //_____________________________________________________________________________
168 AliITSsimulationSDD::AliITSsimulationSDD(AliITSsimulationSDD &source)
170 // Copy constructor to satify Coding roules only.
171 if(this==&source) return;
172 printf("Not allowed to make a copy of AliITSsimulationSDD "
173 "Using default creater instead\n");
174 AliITSsimulationSDD();
176 //_____________________________________________________________________________
177 AliITSsimulationSDD& AliITSsimulationSDD::operator=(AliITSsimulationSDD &source)
179 // Assignment operator to satify Coding roules only.
180 if(this==&source) return *this;
181 printf("Not allowed to make a = with AliITSsimulationSDD "
182 "Using default creater instead\n");
185 //_____________________________________________________________________________
187 AliITSsimulationSDD::AliITSsimulationSDD(AliITSsegmentation *seg,AliITSresponse *resp)
189 // Standard Constructor
196 SetPerpendTracksFlag();
200 fHitMap2 = new AliITSMapA2(fSegmentation,fScaleSize,1);
201 fHitMap1 = new AliITSMapA1(fSegmentation);
204 fNofMaps=fSegmentation->Npz();
205 fMaxNofSamples=fSegmentation->Npx();
207 Float_t sddLength = fSegmentation->Dx();
208 Float_t sddWidth = fSegmentation->Dz();
211 Float_t anodePitch = fSegmentation->Dpz(dummy);
212 Double_t timeStep = (Double_t)fSegmentation->Dpx(dummy);
213 Float_t driftSpeed=fResponse->DriftSpeed();
215 if(anodePitch*(fNofMaps/2) > sddWidth) {
216 Warning("AliITSsimulationSDD",
217 "Too many anodes %d or too big pitch %f \n",fNofMaps/2,anodePitch);
220 if(timeStep*fMaxNofSamples < sddLength/driftSpeed) {
221 Error("AliITSsimulationSDD",
222 "Time Interval > Allowed Time Interval: exit\n");
226 fElectronics = new AliITSetfSDD(timeStep/fScaleSize,fResponse->Electronics());
228 char opt1[20], opt2[20];
229 fResponse->ParamOptions(opt1,opt2);
231 char *same = strstr(opt1,"same");
236 fNoise.Set(fNofMaps);
237 fBaseline.Set(fNofMaps);
241 const char *kopt=fResponse->ZeroSuppOption();
242 if (strstr(fParam,"file") ) {
245 if (strstr(kopt,"2D")) {
248 Init2D(); // desactivate if param change module by module
249 } else if(strstr(kopt,"1D")) {
252 Init1D(); // desactivate if param change module by module
263 Bool_t write=fResponse->OutputOption();
264 if(write && strstr(kopt,"2D")) MakeTreeB();
266 // call here if baseline does not change by module
269 fITS = (AliITS*)gAlice->GetModule("ITS");
270 Int_t size=fNofMaps*fMaxNofSamples;
271 fStream = new AliITSInStream(size);
273 fInZR = new Double_t [fScaleSize*fMaxNofSamples];
274 fInZI = new Double_t [fScaleSize*fMaxNofSamples];
275 fOutZR = new Double_t [fScaleSize*fMaxNofSamples];
276 fOutZI = new Double_t [fScaleSize*fMaxNofSamples];
281 //_____________________________________________________________________________
283 AliITSsimulationSDD::~AliITSsimulationSDD() {
303 if(fTreeB) delete fTreeB;
304 if(fInZR) delete [] fInZR;
305 if(fInZI) delete [] fInZI;
306 if(fOutZR) delete [] fOutZR;
307 if(fOutZI) delete [] fOutZI;
309 //_____________________________________________________________________________
311 void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
312 // create maps to build the lists of tracks
314 cout << "Module: " << md << endl;
318 TObjArray *fHits = mod->GetHits();
319 Int_t nhits = fHits->GetEntriesFast();
320 if (!nhits && fCheckNoise) {
323 fHitMap2->ClearMap();
325 } else if (!nhits) return;
327 //printf("simSDD: module nhits %d %d\n",md,nhits);
330 TObjArray *list=new TObjArray;
331 static TClonesArray *padr=0;
332 if(!padr) padr=new TClonesArray("TVector",1000);
333 Int_t arg[6] = {0,0,0,0,0,0};
334 fHitMap1->SetArray(list);
336 // cout << "set Parameters" << endl;
338 Int_t nofAnodes=fNofMaps/2;
340 Float_t sddLength = fSegmentation->Dx();
341 Float_t sddWidth = fSegmentation->Dz();
344 Float_t anodePitch = fSegmentation->Dpz(dummy);
345 Float_t timeStep = fSegmentation->Dpx(dummy);
347 Float_t driftSpeed=fResponse->DriftSpeed();
349 Float_t maxadc = fResponse->MaxAdc();
350 Float_t topValue = fResponse->DynamicRange();
351 Float_t CHloss = fResponse->ChargeLoss();
352 Float_t norm = maxadc/topValue;
354 // Piergiorgio's part (apart for few variables which I made float
355 // when i thought that can be done
357 // Fill detector maps with GEANT hits
358 // loop over hits in the module
363 const Float_t kconv=1.0e+6; // GeV->KeV
368 for(ii=0; ii<nhits; ii++) {
369 // cout << "hit: " << ii+1 << " of " << nhits << endl;
370 AliITShit *hit = (AliITShit*) fHits->At(ii);
373 // Take into account all hits when several GEANT steps are carried out
374 // inside the silicon
375 // Get and use the status of hit(track):
376 // 66 - for entering hit,
377 // 65 - for inside hit,
378 // 68 - for exiting hit,
379 // 33 - for stopping hit.
381 Int_t status = hit->GetTrackStatus();
383 Int_t hitDetector = hit->GetDetector();
384 Float_t depEnergy = 0.;
385 if(hit->StatusEntering()) { // to be coupled to following hit
387 hit->GetPositionL(xL[0],xL[1],xL[2]);
389 hit1 = (AliITShit*) fHits->At(ii);
390 hit1->GetPositionL(xL1[0],xL1[1],xL1[2]);
391 status1 = hit1->GetTrackStatus();
392 depEnergy = kconv*hit1->GetIonization();
394 depEnergy = kconv*hit->GetIonization(); // Deposited energy in keV
395 hit->GetPositionL(xL1[0],xL1[1],xL1[2]);
397 // cout << "status: " << status << ", status1: " << status1 << ", dE: " << depEnergy << endl;
398 if(fFlag && status1 == 33) continue;
404 // Int_t status1 = -1;
406 //Take now the entering and inside hits only
407 // if(status == 66) {
409 // if(ii<nhits-1) ii++;
410 // hit1 = (AliITShit*) fHits->At(ii);
411 // hit1->GetPositionL(xL1[0],xL1[1],xL1[2]);
412 // status1 = hit1->GetTrackStatus();
413 // depEnergy += kconv*hit1->GetIonization();
414 // if(fFlag && status1 == 65) ctr++;
415 // } while(status1 != 68 && status1 != 33);
419 // scale path to simulate a perpendicular track
420 // continue if the particle did not lose energy
421 // passing through detector
423 printf("This particle has passed without losing energy!\n");
426 Float_t pathInSDD = TMath::Sqrt((xL[0]-xL1[0])*(xL[0]-xL1[0])+(xL[1]-xL1[1])*(xL[1]-xL1[1])+(xL[2]-xL1[2])*(xL[2]-xL1[2]));
428 if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
429 Float_t Drft = (xL1[0]+xL[0])*0.5;
430 Float_t drPath = 10000.*Drft;
431 if(drPath < 0) drPath = -drPath;
432 drPath = sddLength-drPath;
434 cout << "Warning: negative drift path " << drPath << endl;
439 Float_t drTime = drPath/driftSpeed;
442 fResponse->DiffCoeff(dfCoeff,s1);
444 // Squared Sigma along the anodes
445 Double_t sig2A = 2.*dfCoeff*drTime+s1*s1;
446 Double_t sigA = TMath::Sqrt(sig2A);
448 nOfSplits = (Int_t) (1 + 10000.*pathInSDD/sigA);
449 //cout << "nOfSplits: " << nOfSplits << ", sigA: " << sigA << ", path: " << pathInSDD << endl;
451 if(fFlag) nOfSplits = 1;
452 depEnergy /= nOfSplits;
454 for(Int_t kk=0;kk<nOfSplits;kk++) {
456 xL[0]+(xL1[0]-xL[0])*((kk+0.5)/((Float_t) nOfSplits));
458 xL[2]+(xL1[2]-xL[2])*((kk+0.5)/((Float_t) nOfSplits));
459 Float_t driftPath = 10000.*avDrft;
464 driftPath = -driftPath;
466 driftPath = sddLength-driftPath;
467 Int_t detector = 2*(hitDetector-1) + iWing;
469 cout << "Warning: negative drift path " << driftPath << endl;
474 Float_t driftTime = driftPath/driftSpeed;
475 Int_t timeSample = (Int_t) (fScaleSize*driftTime/timeStep + 1);
476 if(timeSample > fScaleSize*fMaxNofSamples) {
477 cout << "Warning: Wrong Time Sample: " << timeSample << endl;
482 Float_t xAnode = 10000.*(avAnode)/anodePitch + nofAnodes/2; // +1?
483 if(xAnode*anodePitch > sddWidth || xAnode*anodePitch < 0.)
484 { cout << "Warning: Z = " << xAnode*anodePitch << endl; }
485 Int_t iAnode = (Int_t) (1.+xAnode); // xAnode?
486 if(iAnode < 1 || iAnode > nofAnodes) {
487 cout << "Warning: Wrong iAnode: " << iAnode << endl;
491 // work with the idtrack=entry number in the TreeH for the moment
492 //Int_t idhit,idtrack;
493 //mod->GetHitTrackAndHitIndex(ii,idtrack,idhit);
494 //Int_t idtrack=mod->GetHitTrackIndex(ii);
495 // or store straight away the particle position in the array
496 // of particles and take idhit=ii only when part is entering (this
497 // requires FillModules() in the macro for analysis) :
498 Int_t itrack = hit->GetTrack();
501 Float_t diffCoeff, s0;
502 fResponse->DiffCoeff(diffCoeff,s0);
504 // Squared Sigma along the anodes
505 Double_t sigma2A = 2.*diffCoeff*driftTime+s0*s0;
506 Double_t sigmaA = TMath::Sqrt(sigma2A);
507 Double_t sigmaT = sigmaA/driftSpeed;
508 // Peak amplitude in nanoAmpere
509 Double_t eVpairs = 3.6;
510 Double_t amplitude = fScaleSize*160.*depEnergy/(timeStep*eVpairs*2.*acos(-1.)*sigmaT*sigmaA);
511 amplitude *= timeStep/25.; // WARNING!!!!! Amplitude scaling to account for clock variations (reference value: 40 MHz)
512 Double_t chargeloss = 1.-CHloss*driftPath/1000;
513 amplitude *= chargeloss;
514 Float_t nsigma=fResponse->NSigmaIntegration();
515 Int_t nlookups = fResponse->GausNLookUp();
516 Float_t width = 2.*nsigma/(nlookups-1);
520 Int_t jt = timeSample;
523 if(driftTime > 1200.) {
528 Int_t nsplit = 4; // hard-wired
529 nsplit = (nsplit+1)/2*2;
531 Double_t aStep = anodePitch/(nsplit*fScaleSize);
532 Double_t tStep = timeStep/(nsplit*fScaleSize);
533 // Define SDD window corresponding to the hit
534 Int_t anodeWindow = (Int_t) (fScaleSize*nsigma*sigmaA/anodePitch + 1);
535 Int_t timeWindow = (Int_t) (fScaleSize*nsigma*sigmaT/timeStep + 1);
536 Int_t jamin = (ja - anodeWindow/ndiv - 1)*fScaleSize*nsplit + 1;
537 Int_t jamax = (ja + anodeWindow/ndiv)*fScaleSize*nsplit;
538 if(jamin <= 0) jamin = 1;
539 if(jamax > fScaleSize*nofAnodes*nsplit) jamax = fScaleSize*nofAnodes*nsplit;
540 Int_t jtmin = (Int_t) (jt - timeWindow*nmul - 1)*nsplit + 1; //hard-wired
541 Int_t jtmax = (Int_t) (jt + timeWindow*nmul)*nsplit; //hard-wired
542 if(jtmin <= 0) jtmin = 1;
543 if(jtmax > fScaleSize*fMaxNofSamples*nsplit) jtmax = fScaleSize*fMaxNofSamples*nsplit;
545 // Spread the charge in the anode-time window
547 //cout << "jamin: " << jamin << ", jamax: " << jamax << endl;
548 //cout << "jtmin: " << jtmin << ", jtmax: " << jtmax << endl;
549 for(ka=jamin; ka <=jamax; ka++) {
550 Int_t ia = (ka-1)/(fScaleSize*nsplit) + 1;
551 if(ia <= 0) { cout << "Warning: ia < 1: " << endl; continue; }
552 if(ia > nofAnodes) ia = nofAnodes;
553 Double_t aExpo = (aStep*(ka-0.5)-xAnode*anodePitch)/sigmaA;
554 Double_t anodeAmplitude = 0;
555 if(TMath::Abs(aExpo) > nsigma) {
557 //cout << "aExpo: " << aExpo << endl;
559 Int_t i = (Int_t) ((aExpo+nsigma)/width);
560 //cout << "eval ampl: " << i << ", " << amplitude << endl;
561 anodeAmplitude = amplitude*fResponse->GausLookUp(i);
562 //cout << "ampl: " << anodeAmplitude << endl;
564 Int_t index = ((detector+1)%2)*nofAnodes+ia-1; // index starts from 0
566 //Double_t rlTime = log(tStep*anodeAmplitude);
568 for(kt=jtmin; kt<=jtmax; kt++) {
569 Int_t it = (kt-1)/nsplit+1; // it starts from 1
570 if(it<=0) { cout << "Warning: it < 1: " << endl; continue; }
571 if(it>fScaleSize*fMaxNofSamples) it = fScaleSize*fMaxNofSamples;
572 Double_t tExpo = (tStep*(kt-0.5)-driftTime)/sigmaT;
573 Double_t timeAmplitude = 0.;
574 if(TMath::Abs(tExpo) > nsigma) {
576 //cout << "tExpo: " << tExpo << endl;
578 Int_t i = (Int_t) ((tExpo+nsigma)/width);
579 //cout << "eval ampl: " << i << ", " << anodeAmplitude << endl;
580 timeAmplitude = anodeAmplitude*fResponse->GausLookUp(i);
583 // build the list of digits for this module
588 timeAmplitude *= norm;
590 ListOfFiredCells(arg,timeAmplitude,list,padr);
591 //cout << "ampl: " << timeAmplitude << endl;
592 } // loop over time in window
593 } // end if anodeAmplitude
594 } // loop over anodes in window
595 } // end loop over "sub-hits"
596 for(Int_t ki=0; ki<3; ki++) xL[ki] = xL1[ki];
597 } // end loop over hits
599 // timer.Stop(); timer.Print();
601 // introduce the electronics effects and do zero-suppression if required
602 Int_t nentries=list->GetEntriesFast();
605 // TStopwatch timer1;
607 // timer1.Stop(); cout << "ele: "; timer1.Print();
609 const char *kopt=fResponse->ZeroSuppOption();
610 ZeroSuppression(kopt);
619 fHitMap1->ClearMap();
620 fHitMap2->ClearMap();
622 //gObjectTable->Print();
626 //____________________________________________
628 void AliITSsimulationSDD::ListOfFiredCells(Int_t *arg,Double_t timeAmplitude,
629 TObjArray *list,TClonesArray *padr){
630 // Returns the list of "fired" cells.
634 Int_t idtrack=arg[2];
636 Int_t counter=arg[4];
637 Int_t countadr=arg[5];
639 Double_t charge=timeAmplitude;
640 charge += fHitMap2->GetSignal(index,ik-1);
641 fHitMap2->SetHit(index, ik-1, charge);
644 Int_t it=(Int_t)((ik-1)/fScaleSize);
648 digits[2]=(Int_t)timeAmplitude;
650 if (idtrack >= 0) phys=(Float_t)timeAmplitude;
653 Double_t cellcharge=0.;
654 AliITSTransientDigit* pdigit;
655 // build the list of fired cells and update the info
656 if (!fHitMap1->TestHit(index, it)) {
658 new((*padr)[countadr++]) TVector(3);
659 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
660 trinfo(0)=(Float_t)idtrack;
661 trinfo(1)=(Float_t)idhit;
662 trinfo(2)=(Float_t)timeAmplitude;
664 list->AddAtAndExpand(
665 new AliITSTransientDigit(phys,digits),counter);
667 fHitMap1->SetHit(index, it, counter);
669 pdigit=(AliITSTransientDigit*)list->
672 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
673 trlist->Add(&trinfo);
677 (AliITSTransientDigit*) fHitMap1->GetHit(index, it);
678 for(Int_t kk=0;kk<fScaleSize;kk++) {
679 cellcharge += fHitMap2->GetSignal(index,fScaleSize*it+kk);
682 (*pdigit).fSignal=(Int_t)cellcharge;
683 (*pdigit).fPhysics+=phys;
684 // update list of tracks
685 TObjArray* trlist=(TObjArray*)pdigit->TrackList();
686 Int_t lastentry=trlist->GetLast();
687 TVector *ptrkp=(TVector*)trlist->At(lastentry);
688 TVector &trinfo=*ptrkp;
689 Int_t lasttrack=Int_t(trinfo(0));
690 //Int_t lasthit=Int_t(trinfo(1));
691 Float_t lastcharge=(trinfo(2));
693 if (lasttrack==idtrack ) {
694 lastcharge+=(Float_t)timeAmplitude;
695 trlist->RemoveAt(lastentry);
697 //trinfo(1)=lasthit; // or idhit
699 trinfo(2)=lastcharge;
700 trlist->AddAt(&trinfo,lastentry);
703 new((*padr)[countadr++]) TVector(3);
704 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
705 trinfo(0)=(Float_t)idtrack;
706 trinfo(1)=(Float_t)idhit;
707 trinfo(2)=(Float_t)timeAmplitude;
709 trlist->Add(&trinfo);
713 // check the track list - debugging
714 Int_t trk[20], htrk[20];
716 Int_t nptracks=trlist->GetEntriesFast();
719 for (tr=0;tr<nptracks;tr++) {
720 TVector *pptrkp=(TVector*)trlist->At(tr);
721 TVector &pptrk=*pptrkp;
722 trk[tr]=Int_t(pptrk(0));
723 htrk[tr]=Int_t(pptrk(1));
724 chtrk[tr]=(pptrk(2));
725 printf("nptracks %d \n",nptracks);
739 //____________________________________________
741 void AliITSsimulationSDD::AddDigit(Int_t i, Int_t j, Int_t signal){
743 // tag with -1 signals coming from background tracks
744 // tag with -2 signals coming from pure electronic noise
746 Int_t digits[3], tracks[3], hits[3];
747 Float_t phys, charges[3];
749 Int_t trk[20], htrk[20];
752 Bool_t do10to8=fResponse->Do10to8();
754 if(do10to8) signal=Convert8to10(signal);
755 AliITSTransientDigit *obj = (AliITSTransientDigit*)fHitMap1->GetHit(i,j);
767 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
770 TObjArray* trlist=(TObjArray*)obj->TrackList();
771 Int_t nptracks=trlist->GetEntriesFast();
774 cout<<"Attention - nptracks > 20 "<<nptracks<<endl;
778 for (tr=0;tr<nptracks;tr++) {
779 TVector &pp =*((TVector*)trlist->At(tr));
780 trk[tr]=Int_t(pp(0));
781 htrk[tr]=Int_t(pp(1));
785 //printf("nptracks > 2 -- %d\n",nptracks);
786 SortTracks(trk,chtrk,htrk,nptracks);
790 for (i=0; i<nptracks; i++) {
795 for (i=nptracks; i<3; i++) {
801 for (i=0; i<3; i++) {
808 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
814 //____________________________________________
816 void AliITSsimulationSDD::SortTracks(Int_t *tracks,Float_t *charges,Int_t *hits,Int_t ntr){
818 // Sort the list of tracks contributing to a given digit
819 // Only the 3 most significant tracks are acctually sorted
823 // Loop over signals, only 3 times
829 Int_t idx[3] = {-3,-3,-3};
830 Float_t jch[3] = {-3,-3,-3};
831 Int_t jtr[3] = {-3,-3,-3};
832 Int_t jhit[3] = {-3,-3,-3};
843 if((i == 1 && j == idx[i-1] )
844 ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
846 if(charges[j] > qmax) {
854 jch[i]=charges[jmax];
874 //____________________________________________
875 void AliITSsimulationSDD::ChargeToSignal() {
876 // add baseline, noise, electronics and ADC saturation effects
878 char opt1[20], opt2[20];
879 fResponse->ParamOptions(opt1,opt2);
880 char *read = strstr(opt1,"file");
882 Float_t baseline, noise;
885 static Bool_t readfile=kTRUE;
886 //read baseline and noise from file
887 if (readfile) ReadBaseline();
889 } else fResponse->GetNoiseParam(noise,baseline);
896 Float_t maxadc = fResponse->MaxAdc();
898 for (i=0;i<fNofMaps;i++) {
899 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
900 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
901 fInZR[k] = fHitMap2->GetSignal(i,k);
902 contrib = (baseline + noise*random.Gaus());
905 for(k=0; k<fMaxNofSamples; k++) {
906 Double_t newcont = 0.;
907 Double_t maxcont = 0.;
908 for(kk=0;kk<fScaleSize;kk++) {
910 newcont = fInZR[fScaleSize*k+kk];
911 if(newcont > maxcont) maxcont = newcont;
913 //newcont += (fInZR[fScaleSize*k+kk]/fScaleSize);
916 if (newcont >= maxadc) newcont = maxadc -1;
917 if(newcont >= baseline) cout << "newcont: " << newcont << endl;
919 fHitMap2->SetHit(i,k,newcont);
921 } // loop over anodes
925 for (i=0;i<fNofMaps;i++) {
926 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
927 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
928 fInZR[k] = fHitMap2->GetSignal(i,k);
929 contrib = (baseline + noise*random.Gaus());
933 FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1);
934 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
935 Double_t rw = fElectronics->GetTraFunReal(k);
936 Double_t iw = fElectronics->GetTraFunImag(k);
937 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
938 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
940 FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1);
941 for(k=0; k<fMaxNofSamples; k++) {
942 Double_t newcont1 = 0.;
943 Double_t maxcont1 = 0.;
944 for(kk=0;kk<fScaleSize;kk++) {
946 newcont1 = fOutZR[fScaleSize*k+kk];
947 if(newcont1 > maxcont1) maxcont1 = newcont1;
949 //newcont1 += (fInZR[fScaleSize*k+kk]/fScaleSize);
952 //cout << "newcont1: " << newcont1 << endl;
953 if (newcont1 >= maxadc) newcont1 = maxadc -1;
954 fHitMap2->SetHit(i,k,newcont1);
956 } // loop over anodes
961 //____________________________________________
962 void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline,
964 // Returns the Baseline for a particular anode.
965 baseline=fBaseline[i];
970 //____________________________________________
971 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl,
973 // Returns the compression alogirthm parameters
974 Int_t size = fD.GetSize();
976 db=fD[i]; tl=fT1[i]; th=fT2[i];
978 if (size <= 2 && i>=fNofMaps/2) {
979 db=fD[1]; tl=fT1[1]; th=fT2[1];
981 db=fD[0]; tl=fT1[0]; th=fT2[0];
985 //____________________________________________
986 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){
987 // returns the compression alogirthm parameters
988 Int_t size = fD.GetSize();
992 if (size <= 2 && i>=fNofMaps/2) {
1000 //____________________________________________
1001 void AliITSsimulationSDD::SetCompressParam(){
1002 // Sets the compression alogirthm parameters
1005 fResponse->GiveCompressParam(cp);
1006 for (i=0; i<2; i++) {
1011 printf("\n i, fD, fT1, fT2, fTol %d %d %d %d %d\n",
1012 i,fD[i],fT1[i],fT2[i],fTol[i]);
1016 //____________________________________________
1017 void AliITSsimulationSDD::ReadBaseline(){
1018 // read baseline and noise from file - either a .root file and in this
1019 // case data should be organised in a tree with one entry for each
1020 // module => reading should be done accordingly
1021 // or a classic file and do smth. like this:
1023 // Read baselines and noise for SDD
1029 char input[100], base[100], param[100];
1032 fResponse->Filenames(input,base,param);
1035 filtmp = gSystem->ExpandPathName(fFileName.Data());
1036 FILE *bline = fopen(filtmp,"r");
1037 printf("filtmp %s\n",filtmp);
1041 while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) {
1043 Error("ReadBaseline","Anode number not in increasing order!",
1052 Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",
1061 //____________________________________________
1062 Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) {
1063 // To the 10 to 8 bit lossive compression.
1064 // code from Davide C. and Albert W.
1066 if (signal < 128) return signal;
1067 if (signal < 256) return (128+((signal-128)>>1));
1068 if (signal < 512) return (192+((signal-256)>>3));
1069 if (signal < 1024) return (224+((signal-512)>>4));
1074 //____________________________________________
1075 Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) {
1076 // Undo the lossive 10 to 8 bit compression.
1077 // code from Davide C. and Albert W.
1078 if (signal < 0 || signal > 255) {
1079 printf("<Convert8to10> out of range %d \n",signal);
1083 if (signal < 128) return signal;
1085 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
1086 else return (128+((signal-128)<<1)+1);
1089 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
1090 else return (256+((signal-192)<<3)+4);
1092 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
1093 else return (512+((signal-224)<<4)+7);
1098 //____________________________________________
1099 AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){
1100 //Return the correct map.
1101 return ((i==0)? fHitMap1 : fHitMap2);
1105 //____________________________________________
1106 void AliITSsimulationSDD::ZeroSuppression(const char *option) {
1107 // perform the zero suppresion
1108 if (strstr(option,"2D")) {
1109 //Init2D(); // activate if param change module by module
1111 } else if (strstr(option,"1D")) {
1112 //Init1D(); // activate if param change module by module
1114 } else StoreAllDigits();
1118 //____________________________________________
1119 void AliITSsimulationSDD::Init2D(){
1120 // read in and prepare arrays: fD, fT1, fT2
1121 // savemu[nanodes], savesigma[nanodes]
1122 // read baseline and noise from file - either a .root file and in this
1123 // case data should be organised in a tree with one entry for each
1124 // module => reading should be done accordingly
1125 // or a classic file and do smth. like this ( code from Davide C. and
1128 // Read 2D zero-suppression parameters for SDD
1131 if (!strstr(fParam,"file")) return;
1133 Int_t na,pos,tempTh;
1135 Float_t *savemu = new Float_t [fNofMaps];
1136 Float_t *savesigma = new Float_t [fNofMaps];
1137 char input[100],basel[100],par[100];
1141 Int_t minval = fResponse->MinVal();
1143 fResponse->Filenames(input,basel,par);
1147 filtmp = gSystem->ExpandPathName(fFileName.Data());
1148 FILE *param = fopen(filtmp,"r");
1152 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1154 Error("Init2D ","Anode number not in increasing order!",
1159 savesigma[na]=sigma;
1160 if ((2.*sigma) < mu) {
1161 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1164 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1165 if (tempTh < 0) tempTh=0;
1167 tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval;
1168 if (tempTh < 0) tempTh=0;
1174 Error("Init2D "," THE FILE %s DOES NOT EXIST !",
1182 delete [] savesigma;
1185 //____________________________________________
1186 void AliITSsimulationSDD::Compress2D(){
1188 // simple ITS cluster finder -- online zero-suppression conditions
1193 Int_t minval = fResponse->MinVal();
1194 Bool_t write=fResponse->OutputOption();
1195 Bool_t do10to8=fResponse->Do10to8();
1197 Int_t nz, nl, nh, low, i, j;
1199 for (i=0; i<fNofMaps; i++) {
1200 CompressionParam(i,db,tl,th);
1205 for (j=0; j<fMaxNofSamples; j++) {
1206 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1207 signal -= db; // if baseline eq. is done here
1208 if (signal <= 0) {nz++; continue;}
1209 if ((signal - tl) < minval) low++;
1210 if ((signal - th) >= minval) {
1213 FindCluster(i,j,signal,minval,cond);
1214 if (cond && j && ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)) {
1215 if(do10to8) signal = Convert10to8(signal);
1216 AddDigit(i,j,signal);
1218 } else if ((signal - tl) >= minval) nl++;
1219 } // loop time samples
1220 if (write) TreeB()->Fill(nz,nl,nh,low,i+1);
1225 sprintf(hname,"TNtuple%d_%d",fModule,fEvent);
1226 TreeB()->Write(hname);
1233 //_____________________________________________________________________________
1234 void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal,
1235 Int_t minval,Bool_t &cond){
1237 // Find clusters according to the online 2D zero-suppression algorithm
1240 Bool_t do10to8=fResponse->Do10to8();
1244 fHitMap2->FlagHit(i,j);
1246 // check the online zero-suppression conditions
1248 const Int_t maxNeighbours = 4;
1252 Int_t xList[maxNeighbours], yList[maxNeighbours];
1253 fSegmentation->Neighbours(i,j,&nn,xList,yList);
1255 for (in=0; in<nn; in++) {
1258 if (fHitMap2->TestHit(ix,iy)==kUnused) {
1259 CompressionParam(ix,dbx,tlx,thx);
1260 Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy));
1261 qn -= dbx; // if baseline eq. is done here
1262 if ((qn-tlx) < minval) {
1263 fHitMap2->FlagHit(ix,iy);
1266 if ((qn - thx) >= minval) high=kTRUE;
1268 if(do10to8) signal = Convert10to8(signal);
1269 AddDigit(i,j,signal);
1271 if(do10to8) qns = Convert10to8(qn);
1273 if (!high) AddDigit(ix,iy,qns);
1275 if(!high) fHitMap2->FlagHit(ix,iy);
1278 } // loop over neighbours
1282 //____________________________________________
1283 void AliITSsimulationSDD::Init1D(){
1284 // this is just a copy-paste of input taken from 2D algo
1285 // Torino people should give input
1287 // Read 1D zero-suppression parameters for SDD
1290 if (!strstr(fParam,"file")) return;
1292 Int_t na,pos,tempTh;
1294 Float_t *savemu = new Float_t [fNofMaps];
1295 Float_t *savesigma = new Float_t [fNofMaps];
1296 char input[100],basel[100],par[100];
1300 Int_t minval = fResponse->MinVal();
1301 fResponse->Filenames(input,basel,par);
1304 // set first the disable and tol param
1307 filtmp = gSystem->ExpandPathName(fFileName.Data());
1308 FILE *param = fopen(filtmp,"r");
1312 fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]);
1313 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1315 Error("Init1D ","Anode number not in increasing order!",
1320 savesigma[na]=sigma;
1321 if ((2.*sigma) < mu) {
1322 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1325 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1326 if (tempTh < 0) tempTh=0;
1331 Error("Init1D "," THE FILE %s DOES NOT EXIST !",
1339 delete [] savesigma;
1345 //____________________________________________
1346 void AliITSsimulationSDD::Compress1D(){
1347 // 1D zero-suppression algorithm (from Gianluca A.)
1349 Int_t dis,tol,thres,decr,diff;
1351 UChar_t *str=fStream->Stream();
1354 Bool_t do10to8=fResponse->Do10to8();
1358 for (k=0; k<2; k++) {
1361 for (i=0; i<fNofMaps/2; i++) {
1362 Bool_t firstSignal=kTRUE;
1363 Int_t idx=i+k*fNofMaps/2;
1364 CompressionParam(idx,decr,thres);
1365 for (j=0; j<fMaxNofSamples; j++) {
1366 Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j));
1367 signal -= decr; // if baseline eq.
1368 if(do10to8) signal = Convert10to8(signal);
1369 if (signal <= thres) {
1373 // write diff in the buffer for HuffT
1374 str[counter]=(UChar_t)diff;
1379 if (diff > 127) diff=127;
1380 if (diff < -128) diff=-128;
1383 // tol has changed to 8 possible cases ? - one can write
1384 // this if(TMath::Abs(diff)<tol) ... else ...
1385 if(TMath::Abs(diff)<tol) diff=0;
1386 // or keep it as it was before
1388 if (tol==1 && (diff >= -2 && diff <= 1)) diff=0;
1389 if (tol==2 && (diff >= -4 && diff <= 3)) diff=0;
1390 if (tol==3 && (diff >= -16 && diff <= 15)) diff=0;
1392 AddDigit(idx,j,last+diff);
1394 AddDigit(idx,j,signal);
1398 // write diff in the buffer used to compute Huffman tables
1399 if (firstSignal) str[counter]=(UChar_t)signal;
1400 else str[counter]=(UChar_t)diff;
1405 } // loop time samples
1406 } // loop anodes one half of detector
1410 fStream->CheckCount(counter);
1412 // open file and write out the stream of diff's
1414 static Bool_t open=kTRUE;
1415 static TFile *outFile;
1416 Bool_t write = fResponse->OutputOption();
1420 SetFileName("stream.root");
1421 cout<<"filename "<<fFileName<<endl;
1422 outFile=new TFile(fFileName,"recreate");
1423 cout<<"I have opened "<<fFileName<<" file "<<endl;
1430 fStream->ClearStream();
1432 // back to galice.root file
1434 TTree *fAli=gAlice->TreeK();
1437 if (fAli) file =fAli->GetCurrentFile();
1442 //____________________________________________
1443 void AliITSsimulationSDD::StoreAllDigits(){
1444 // if non-zero-suppressed data
1446 Bool_t do10to8=fResponse->Do10to8();
1448 Int_t i, j, digits[3];
1449 for (i=0; i<fNofMaps; i++) {
1450 for (j=0; j<fMaxNofSamples; j++) {
1451 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1452 if(do10to8) signal = Convert10to8(signal);
1453 if(do10to8) signal = Convert8to10(signal);
1457 fITS->AddRealDigit(1,digits);
1461 //____________________________________________
1463 void AliITSsimulationSDD::CreateHistograms(Int_t scale){
1464 // Creates histograms of maps for debugging
1468 fHis=new TObjArray(fNofMaps);
1469 for (i=0;i<fNofMaps;i++) {
1470 TString sddName("sdd_");
1472 sprintf(candNum,"%d",i+1);
1473 sddName.Append(candNum);
1474 (*fHis)[i] = new TH1F(sddName.Data(),"SDD maps",
1475 scale*fMaxNofSamples,0.,(Float_t) scale*fMaxNofSamples);
1479 //____________________________________________
1480 void AliITSsimulationSDD::FillHistograms(){
1481 // fill 1D histograms from map
1484 for( Int_t i=0; i<fNofMaps; i++) {
1485 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
1486 Int_t nsamples = hist->GetNbinsX();
1487 for( Int_t j=0; j<nsamples; j++) {
1488 Double_t signal=fHitMap2->GetSignal(i,j);
1489 hist->Fill((Float_t)j,signal);
1494 //____________________________________________
1496 void AliITSsimulationSDD::ResetHistograms(){
1498 // Reset histograms for this detector
1501 for (i=0;i<fNofMaps;i++ ) {
1502 if ((*fHis)[i]) ((TH1F*)(*fHis)[i])->Reset();
1508 //____________________________________________
1510 TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
1511 // Fills a histogram from a give anode.
1512 if (!fHis) return 0;
1514 if(wing <=0 || wing > 2) {
1515 cout << "Wrong wing number: " << wing << endl;
1518 if(anode <=0 || anode > fNofMaps/2) {
1519 cout << "Wrong anode number: " << anode << endl;
1523 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
1524 return (TH1F*)((*fHis)[index]);
1527 //____________________________________________
1529 void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
1530 // Writes the histograms to a file
1535 for(i=0; i<fNofMaps; i++) (*fHis)[i]->Write(); //fAdcs[i]->Write();
1538 //____________________________________________
1539 Float_t AliITSsimulationSDD::GetNoise() {
1540 // Returns the noise value
1542 //Bool_t do10to8=fResponse->Do10to8();
1543 //noise will always be in the liniar part of the signal
1546 Int_t threshold=fT1[0];
1548 char opt1[20], opt2[20];
1549 fResponse->ParamOptions(opt1,opt2);
1551 char *same = strstr(opt1,"same");
1552 Float_t noise,baseline;
1554 fResponse->GetNoiseParam(noise,baseline);
1556 static Bool_t readfile=kTRUE;
1557 //read baseline and noise from file
1558 if (readfile) ReadBaseline();
1562 TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2");
1563 if(c2) delete c2->GetPrimitive("noisehist");
1564 if(c2) delete c2->GetPrimitive("anode");
1565 else c2=new TCanvas("c2");
1567 c2->SetFillColor(0);
1569 TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold);
1570 TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0.,(float)fMaxNofSamples);
1572 for (i=0;i<fNofMaps;i++) {
1573 CompressionParam(i,decr,threshold);
1574 if (!same) GetAnodeBaseline(i,baseline,noise);
1576 for (k=0;k<fMaxNofSamples;k++) {
1577 Float_t signal=(Float_t)fHitMap2->GetSignal(i,k);
1578 //if (signal <= (float)threshold) noisehist->Fill(signal-baseline);
1579 if (signal <= (float)threshold) noisehist->Fill(signal);
1580 anode->Fill((float)k,signal);
1585 TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold);
1586 noisehist->Fit("gnoise","RQ");
1589 Float_t mnoise = gnoise->GetParameter(1);
1590 cout << "mnoise : " << mnoise << endl;
1591 Float_t rnoise = gnoise->GetParameter(2);
1592 cout << "rnoise : " << rnoise << endl;
1597 //____________________________________________
1599 void AliITSsimulationSDD::Print() {
1601 // Print SDD simulation Parameters
1603 cout << "**************************************************" << endl;
1604 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1605 cout << "**************************************************" << endl;
1606 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1607 cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl;
1608 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1609 cout << "Number pf Anodes used: " << fNofMaps << endl;
1610 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1611 cout << "Scale size factor: " << fScaleSize << endl;
1612 cout << "**************************************************" << endl;