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
191 fHitMap1 = 0; // zero just in case of an error
192 fHitMap2 = 0; // zero just in case of an error
193 fElectronics = 0; // zero just in case of an error
194 fStream = 0; // zero just in case of an error
200 SetPerpendTracksFlag();
204 fHitMap2 = new AliITSMapA2(fSegmentation,fScaleSize,1);
205 fHitMap1 = new AliITSMapA1(fSegmentation);
208 fNofMaps=fSegmentation->Npz();
209 fMaxNofSamples=fSegmentation->Npx();
211 Float_t sddLength = fSegmentation->Dx();
212 Float_t sddWidth = fSegmentation->Dz();
215 Float_t anodePitch = fSegmentation->Dpz(dummy);
216 Double_t timeStep = (Double_t)fSegmentation->Dpx(dummy);
217 Float_t driftSpeed=fResponse->DriftSpeed();
219 if(anodePitch*(fNofMaps/2) > sddWidth) {
220 Warning("AliITSsimulationSDD",
221 "Too many anodes %d or too big pitch %f \n",
222 fNofMaps/2,anodePitch);
225 if(timeStep*fMaxNofSamples < sddLength/driftSpeed) {
226 Error("AliITSsimulationSDD",
227 "Time Interval > Allowed Time Interval: exit\n");
231 fElectronics = new AliITSetfSDD(timeStep/fScaleSize,
232 fResponse->Electronics());
234 char opt1[20], opt2[20];
235 fResponse->ParamOptions(opt1,opt2);
237 char *same = strstr(opt1,"same");
242 fNoise.Set(fNofMaps);
243 fBaseline.Set(fNofMaps);
247 const char *kopt=fResponse->ZeroSuppOption();
248 if (strstr(fParam,"file") ) {
251 if (strstr(kopt,"2D")) {
254 Init2D(); // desactivate if param change module by module
255 } else if(strstr(kopt,"1D")) {
258 Init1D(); // desactivate if param change module by module
266 } // end if else strstr
268 Bool_t write=fResponse->OutputOption();
269 if(write && strstr(kopt,"2D")) MakeTreeB();
271 // call here if baseline does not change by module
274 fITS = (AliITS*)gAlice->GetModule("ITS");
275 Int_t size=fNofMaps*fMaxNofSamples;
276 fStream = new AliITSInStream(size);
278 fInZR = new Double_t [fScaleSize*fMaxNofSamples];
279 fInZI = new Double_t [fScaleSize*fMaxNofSamples];
280 fOutZR = new Double_t [fScaleSize*fMaxNofSamples];
281 fOutZI = new Double_t [fScaleSize*fMaxNofSamples];
286 //_____________________________________________________________________________
288 AliITSsimulationSDD::~AliITSsimulationSDD() {
308 if(fTreeB) delete fTreeB;
309 if(fInZR) delete [] fInZR;
310 if(fInZI) delete [] fInZI;
311 if(fOutZR) delete [] fOutZR;
312 if(fOutZI) delete [] fOutZI;
314 //_____________________________________________________________________________
316 void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
317 // create maps to build the lists of tracks
319 //cout << "Module: " << md << endl;
323 TObjArray *fHits = mod->GetHits();
324 Int_t nhits = fHits->GetEntriesFast();
325 if (!nhits && fCheckNoise) {
328 fHitMap2->ClearMap();
330 } else if (!nhits) return;
332 //printf("simSDD: module nhits %d %d\n",md,nhits);
335 TObjArray *list=new TObjArray;
336 static TClonesArray *padr=0;
337 if(!padr) padr=new TClonesArray("TVector",1000);
338 Int_t arg[6] = {0,0,0,0,0,0};
339 fHitMap1->SetArray(list);
341 // cout << "set Parameters" << endl;
343 Int_t nofAnodes=fNofMaps/2;
345 Float_t sddLength = fSegmentation->Dx();
346 Float_t sddWidth = fSegmentation->Dz();
349 Float_t anodePitch = fSegmentation->Dpz(dummy);
350 Float_t timeStep = fSegmentation->Dpx(dummy);
352 Float_t driftSpeed=fResponse->DriftSpeed();
354 Float_t maxadc = fResponse->MaxAdc();
355 Float_t topValue = fResponse->DynamicRange();
356 Float_t CHloss = fResponse->ChargeLoss();
357 Float_t norm = maxadc/topValue;
359 // Piergiorgio's part (apart for few variables which I made float
360 // when i thought that can be done
362 // Fill detector maps with GEANT hits
363 // loop over hits in the module
368 const Float_t kconv=1.0e+6; // GeV->KeV
373 for(ii=0; ii<nhits; ii++) {
374 // cout << "hit: " << ii+1 << " of " << nhits << endl;
375 AliITShit *hit = (AliITShit*) fHits->At(ii);
378 // Take into account all hits when several GEANT steps are carried out
379 // inside the silicon
380 // Get and use the status of hit(track):
381 // 66 - for entering hit,
382 // 65 - for inside hit,
383 // 68 - for exiting hit,
384 // 33 - for stopping hit.
386 //Int_t status = hit->GetTrackStatus();
388 Int_t hitDetector = hit->GetDetector();
389 Float_t depEnergy = 0.;
390 if(hit->StatusEntering()) { // to be coupled to following hit
392 hit->GetPositionL(xL[0],xL[1],xL[2]);
394 hit1 = (AliITShit*) fHits->At(ii);
395 hit1->GetPositionL(xL1[0],xL1[1],xL1[2]);
396 status1 = hit1->GetTrackStatus();
397 depEnergy = kconv*hit1->GetIonization();
399 depEnergy = kconv*hit->GetIonization(); // Deposited energy in keV
400 hit->GetPositionL(xL1[0],xL1[1],xL1[2]);
402 // cout << "status: " << status << ", status1: " << status1 << ", dE: " << depEnergy << endl;
403 if(fFlag && status1 == 33) continue;
409 // Int_t status1 = -1;
411 //Take now the entering and inside hits only
412 // if(status == 66) {
414 // if(ii<nhits-1) ii++;
415 // hit1 = (AliITShit*) fHits->At(ii);
416 // hit1->GetPositionL(xL1[0],xL1[1],xL1[2]);
417 // status1 = hit1->GetTrackStatus();
418 // depEnergy += kconv*hit1->GetIonization();
419 // if(fFlag && status1 == 65) ctr++;
420 // } while(status1 != 68 && status1 != 33);
424 // scale path to simulate a perpendicular track
425 // continue if the particle did not lose energy
426 // passing through detector
428 printf("This particle has passed without losing energy!\n");
431 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]));
433 if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
434 Float_t Drft = (xL1[0]+xL[0])*0.5;
435 Float_t drPath = 10000.*Drft;
436 if(drPath < 0) drPath = -drPath;
437 drPath = sddLength-drPath;
439 cout << "Warning: negative drift path " << drPath << endl;
444 Float_t drTime = drPath/driftSpeed;
447 fResponse->DiffCoeff(dfCoeff,s1);
449 // Squared Sigma along the anodes
450 Double_t sig2A = 2.*dfCoeff*drTime+s1*s1;
451 Double_t sigA = TMath::Sqrt(sig2A);
453 nOfSplits = (Int_t) (1 + 10000.*pathInSDD/sigA);
454 //cout << "nOfSplits: " << nOfSplits << ", sigA: " << sigA << ", path: " << pathInSDD << endl;
456 if(fFlag) nOfSplits = 1;
457 depEnergy /= nOfSplits;
459 for(Int_t kk=0;kk<nOfSplits;kk++) {
461 xL[0]+(xL1[0]-xL[0])*((kk+0.5)/((Float_t) nOfSplits));
463 xL[2]+(xL1[2]-xL[2])*((kk+0.5)/((Float_t) nOfSplits));
464 Float_t driftPath = 10000.*avDrft;
469 driftPath = -driftPath;
471 driftPath = sddLength-driftPath;
472 Int_t detector = 2*(hitDetector-1) + iWing;
474 cout << "Warning: negative drift path " << driftPath << endl;
479 Float_t driftTime = driftPath/driftSpeed;
480 Int_t timeSample = (Int_t) (fScaleSize*driftTime/timeStep + 1);
481 if(timeSample > fScaleSize*fMaxNofSamples) {
482 cout << "Warning: Wrong Time Sample: " << timeSample << endl;
487 Float_t xAnode = 10000.*(avAnode)/anodePitch + nofAnodes/2; // +1?
488 if(xAnode*anodePitch > sddWidth || xAnode*anodePitch < 0.)
489 { cout << "Warning: Z = " << xAnode*anodePitch << endl; }
490 Int_t iAnode = (Int_t) (1.+xAnode); // xAnode?
491 if(iAnode < 1 || iAnode > nofAnodes) {
492 cout << "Warning: Wrong iAnode: " << iAnode << endl;
496 // work with the idtrack=entry number in the TreeH for the moment
497 //Int_t idhit,idtrack;
498 //mod->GetHitTrackAndHitIndex(ii,idtrack,idhit);
499 //Int_t idtrack=mod->GetHitTrackIndex(ii);
500 // or store straight away the particle position in the array
501 // of particles and take idhit=ii only when part is entering (this
502 // requires FillModules() in the macro for analysis) :
503 Int_t itrack = hit->GetTrack();
506 Float_t diffCoeff, s0;
507 fResponse->DiffCoeff(diffCoeff,s0);
509 // Squared Sigma along the anodes
510 Double_t sigma2A = 2.*diffCoeff*driftTime+s0*s0;
511 Double_t sigmaA = TMath::Sqrt(sigma2A);
512 Double_t sigmaT = sigmaA/driftSpeed;
513 // Peak amplitude in nanoAmpere
514 Double_t eVpairs = 3.6;
515 Double_t amplitude = fScaleSize*160.*depEnergy/(timeStep*eVpairs*2.*acos(-1.)*sigmaT*sigmaA);
516 amplitude *= timeStep/25.; // WARNING!!!!! Amplitude scaling to account for clock variations (reference value: 40 MHz)
517 Double_t chargeloss = 1.-CHloss*driftPath/1000;
518 amplitude *= chargeloss;
519 Float_t nsigma=fResponse->NSigmaIntegration();
520 Int_t nlookups = fResponse->GausNLookUp();
521 Float_t width = 2.*nsigma/(nlookups-1);
525 Int_t jt = timeSample;
528 if(driftTime > 1200.) {
533 Int_t nsplit = 4; // hard-wired
534 nsplit = (nsplit+1)/2*2;
536 Double_t aStep = anodePitch/(nsplit*fScaleSize);
537 Double_t tStep = timeStep/(nsplit*fScaleSize);
538 // Define SDD window corresponding to the hit
539 Int_t anodeWindow = (Int_t) (fScaleSize*nsigma*sigmaA/anodePitch + 1);
540 Int_t timeWindow = (Int_t) (fScaleSize*nsigma*sigmaT/timeStep + 1);
541 Int_t jamin = (ja - anodeWindow/ndiv - 1)*fScaleSize*nsplit + 1;
542 Int_t jamax = (ja + anodeWindow/ndiv)*fScaleSize*nsplit;
543 if(jamin <= 0) jamin = 1;
544 if(jamax > fScaleSize*nofAnodes*nsplit) jamax = fScaleSize*nofAnodes*nsplit;
545 Int_t jtmin = (Int_t) (jt - timeWindow*nmul - 1)*nsplit + 1; //hard-wired
546 Int_t jtmax = (Int_t) (jt + timeWindow*nmul)*nsplit; //hard-wired
547 if(jtmin <= 0) jtmin = 1;
548 if(jtmax > fScaleSize*fMaxNofSamples*nsplit) jtmax = fScaleSize*fMaxNofSamples*nsplit;
550 // Spread the charge in the anode-time window
552 //cout << "jamin: " << jamin << ", jamax: " << jamax << endl;
553 //cout << "jtmin: " << jtmin << ", jtmax: " << jtmax << endl;
554 for(ka=jamin; ka <=jamax; ka++) {
555 Int_t ia = (ka-1)/(fScaleSize*nsplit) + 1;
556 if(ia <= 0) { cout << "Warning: ia < 1: " << endl; continue; }
557 if(ia > nofAnodes) ia = nofAnodes;
558 Double_t aExpo = (aStep*(ka-0.5)-xAnode*anodePitch)/sigmaA;
559 Double_t anodeAmplitude = 0;
560 if(TMath::Abs(aExpo) > nsigma) {
562 //cout << "aExpo: " << aExpo << endl;
564 Int_t i = (Int_t) ((aExpo+nsigma)/width);
565 //cout << "eval ampl: " << i << ", " << amplitude << endl;
566 anodeAmplitude = amplitude*fResponse->GausLookUp(i);
567 //cout << "ampl: " << anodeAmplitude << endl;
569 Int_t index = ((detector+1)%2)*nofAnodes+ia-1; // index starts from 0
571 //Double_t rlTime = log(tStep*anodeAmplitude);
573 for(kt=jtmin; kt<=jtmax; kt++) {
574 Int_t it = (kt-1)/nsplit+1; // it starts from 1
575 if(it<=0) { cout << "Warning: it < 1: " << endl; continue; }
576 if(it>fScaleSize*fMaxNofSamples) it = fScaleSize*fMaxNofSamples;
577 Double_t tExpo = (tStep*(kt-0.5)-driftTime)/sigmaT;
578 Double_t timeAmplitude = 0.;
579 if(TMath::Abs(tExpo) > nsigma) {
581 //cout << "tExpo: " << tExpo << endl;
583 Int_t i = (Int_t) ((tExpo+nsigma)/width);
584 //cout << "eval ampl: " << i << ", " << anodeAmplitude << endl;
585 timeAmplitude = anodeAmplitude*fResponse->GausLookUp(i);
588 // build the list of digits for this module
593 timeAmplitude *= norm;
595 ListOfFiredCells(arg,timeAmplitude,list,padr);
596 //cout << "ampl: " << timeAmplitude << endl;
597 } // loop over time in window
598 } // end if anodeAmplitude
599 } // loop over anodes in window
600 } // end loop over "sub-hits"
601 for(Int_t ki=0; ki<3; ki++) xL[ki] = xL1[ki];
602 } // end loop over hits
604 // timer.Stop(); timer.Print();
606 // introduce the electronics effects and do zero-suppression if required
607 Int_t nentries=list->GetEntriesFast();
610 // TStopwatch timer1;
612 // timer1.Stop(); cout << "ele: "; timer1.Print();
614 const char *kopt=fResponse->ZeroSuppOption();
615 ZeroSuppression(kopt);
624 fHitMap1->ClearMap();
625 fHitMap2->ClearMap();
627 //gObjectTable->Print();
631 //____________________________________________
633 void AliITSsimulationSDD::ListOfFiredCells(Int_t *arg,Double_t timeAmplitude,
634 TObjArray *list,TClonesArray *padr){
635 // Returns the list of "fired" cells.
639 Int_t idtrack=arg[2];
641 Int_t counter=arg[4];
642 Int_t countadr=arg[5];
644 Double_t charge=timeAmplitude;
645 charge += fHitMap2->GetSignal(index,ik-1);
646 fHitMap2->SetHit(index, ik-1, charge);
649 Int_t it=(Int_t)((ik-1)/fScaleSize);
653 digits[2]=(Int_t)timeAmplitude;
655 if (idtrack >= 0) phys=(Float_t)timeAmplitude;
658 Double_t cellcharge=0.;
659 AliITSTransientDigit* pdigit;
660 // build the list of fired cells and update the info
661 if (!fHitMap1->TestHit(index, it)) {
663 new((*padr)[countadr++]) TVector(3);
664 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
665 trinfo(0)=(Float_t)idtrack;
666 trinfo(1)=(Float_t)idhit;
667 trinfo(2)=(Float_t)timeAmplitude;
669 list->AddAtAndExpand(
670 new AliITSTransientDigit(phys,digits),counter);
672 fHitMap1->SetHit(index, it, counter);
674 pdigit=(AliITSTransientDigit*)list->
677 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
678 trlist->Add(&trinfo);
682 (AliITSTransientDigit*) fHitMap1->GetHit(index, it);
683 for(Int_t kk=0;kk<fScaleSize;kk++) {
684 cellcharge += fHitMap2->GetSignal(index,fScaleSize*it+kk);
687 (*pdigit).fSignal=(Int_t)cellcharge;
688 (*pdigit).fPhysics+=phys;
689 // update list of tracks
690 TObjArray* trlist=(TObjArray*)pdigit->TrackList();
691 Int_t lastentry=trlist->GetLast();
692 TVector *ptrkp=(TVector*)trlist->At(lastentry);
693 TVector &trinfo=*ptrkp;
694 Int_t lasttrack=Int_t(trinfo(0));
695 //Int_t lasthit=Int_t(trinfo(1));
696 Float_t lastcharge=(trinfo(2));
698 if (lasttrack==idtrack ) {
699 lastcharge+=(Float_t)timeAmplitude;
700 trlist->RemoveAt(lastentry);
702 //trinfo(1)=lasthit; // or idhit
704 trinfo(2)=lastcharge;
705 trlist->AddAt(&trinfo,lastentry);
708 new((*padr)[countadr++]) TVector(3);
709 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
710 trinfo(0)=(Float_t)idtrack;
711 trinfo(1)=(Float_t)idhit;
712 trinfo(2)=(Float_t)timeAmplitude;
714 trlist->Add(&trinfo);
718 // check the track list - debugging
719 Int_t trk[20], htrk[20];
721 Int_t nptracks=trlist->GetEntriesFast();
724 for (tr=0;tr<nptracks;tr++) {
725 TVector *pptrkp=(TVector*)trlist->At(tr);
726 TVector &pptrk=*pptrkp;
727 trk[tr]=Int_t(pptrk(0));
728 htrk[tr]=Int_t(pptrk(1));
729 chtrk[tr]=(pptrk(2));
730 printf("nptracks %d \n",nptracks);
744 //____________________________________________
746 void AliITSsimulationSDD::AddDigit(Int_t i, Int_t j, Int_t signal){
748 // tag with -1 signals coming from background tracks
749 // tag with -2 signals coming from pure electronic noise
751 Int_t digits[3], tracks[3], hits[3];
752 Float_t phys, charges[3];
754 Int_t trk[20], htrk[20];
757 Bool_t do10to8=fResponse->Do10to8();
759 if(do10to8) signal=Convert8to10(signal);
760 AliITSTransientDigit *obj = (AliITSTransientDigit*)fHitMap1->GetHit(i,j);
772 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
775 TObjArray* trlist=(TObjArray*)obj->TrackList();
776 Int_t nptracks=trlist->GetEntriesFast();
779 cout<<"Attention - nptracks > 20 "<<nptracks<<endl;
783 for (tr=0;tr<nptracks;tr++) {
784 TVector &pp =*((TVector*)trlist->At(tr));
785 trk[tr]=Int_t(pp(0));
786 htrk[tr]=Int_t(pp(1));
790 //printf("nptracks > 2 -- %d\n",nptracks);
791 SortTracks(trk,chtrk,htrk,nptracks);
795 for (i=0; i<nptracks; i++) {
800 for (i=nptracks; i<3; i++) {
806 for (i=0; i<3; i++) {
813 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
819 //____________________________________________
821 void AliITSsimulationSDD::SortTracks(Int_t *tracks,Float_t *charges,Int_t *hits,Int_t ntr){
823 // Sort the list of tracks contributing to a given digit
824 // Only the 3 most significant tracks are acctually sorted
828 // Loop over signals, only 3 times
834 Int_t idx[3] = {-3,-3,-3};
835 Float_t jch[3] = {-3,-3,-3};
836 Int_t jtr[3] = {-3,-3,-3};
837 Int_t jhit[3] = {-3,-3,-3};
848 if((i == 1 && j == idx[i-1] )
849 ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
851 if(charges[j] > qmax) {
859 jch[i]=charges[jmax];
879 //____________________________________________
880 void AliITSsimulationSDD::ChargeToSignal() {
881 // add baseline, noise, electronics and ADC saturation effects
883 char opt1[20], opt2[20];
884 fResponse->ParamOptions(opt1,opt2);
885 char *read = strstr(opt1,"file");
887 Float_t baseline, noise;
890 static Bool_t readfile=kTRUE;
891 //read baseline and noise from file
892 if (readfile) ReadBaseline();
894 } else fResponse->GetNoiseParam(noise,baseline);
901 Float_t maxadc = fResponse->MaxAdc();
903 for (i=0;i<fNofMaps;i++) {
904 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
905 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
906 fInZR[k] = fHitMap2->GetSignal(i,k);
907 contrib = (baseline + noise*random.Gaus());
910 for(k=0; k<fMaxNofSamples; k++) {
911 Double_t newcont = 0.;
912 Double_t maxcont = 0.;
913 for(kk=0;kk<fScaleSize;kk++) {
915 newcont = fInZR[fScaleSize*k+kk];
916 if(newcont > maxcont) maxcont = newcont;
918 //newcont += (fInZR[fScaleSize*k+kk]/fScaleSize);
921 if (newcont >= maxadc) newcont = maxadc -1;
922 if(newcont >= baseline) cout << "newcont: " << newcont << endl;
924 fHitMap2->SetHit(i,k,newcont);
926 } // loop over anodes
930 for (i=0;i<fNofMaps;i++) {
931 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
932 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
933 fInZR[k] = fHitMap2->GetSignal(i,k);
934 contrib = (baseline + noise*random.Gaus());
938 FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1);
939 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
940 Double_t rw = fElectronics->GetTraFunReal(k);
941 Double_t iw = fElectronics->GetTraFunImag(k);
942 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
943 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
945 FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1);
946 for(k=0; k<fMaxNofSamples; k++) {
947 Double_t newcont1 = 0.;
948 Double_t maxcont1 = 0.;
949 for(kk=0;kk<fScaleSize;kk++) {
951 newcont1 = fOutZR[fScaleSize*k+kk];
952 if(newcont1 > maxcont1) maxcont1 = newcont1;
954 //newcont1 += (fInZR[fScaleSize*k+kk]/fScaleSize);
957 //cout << "newcont1: " << newcont1 << endl;
958 if (newcont1 >= maxadc) newcont1 = maxadc -1;
959 fHitMap2->SetHit(i,k,newcont1);
961 } // loop over anodes
966 //____________________________________________
967 void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline,
969 // Returns the Baseline for a particular anode.
970 baseline=fBaseline[i];
975 //____________________________________________
976 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl,
978 // Returns the compression alogirthm parameters
979 Int_t size = fD.GetSize();
981 db=fD[i]; tl=fT1[i]; th=fT2[i];
983 if (size <= 2 && i>=fNofMaps/2) {
984 db=fD[1]; tl=fT1[1]; th=fT2[1];
986 db=fD[0]; tl=fT1[0]; th=fT2[0];
990 //____________________________________________
991 void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){
992 // returns the compression alogirthm parameters
993 Int_t size = fD.GetSize();
997 if (size <= 2 && i>=fNofMaps/2) {
1000 db=fD[0]; tl=fT1[0];
1005 //____________________________________________
1006 void AliITSsimulationSDD::SetCompressParam(){
1007 // Sets the compression alogirthm parameters
1010 fResponse->GiveCompressParam(cp);
1011 for (i=0; i<2; i++) {
1016 printf("\n i, fD, fT1, fT2, fTol %d %d %d %d %d\n",
1017 i,fD[i],fT1[i],fT2[i],fTol[i]);
1021 //____________________________________________
1022 void AliITSsimulationSDD::ReadBaseline(){
1023 // read baseline and noise from file - either a .root file and in this
1024 // case data should be organised in a tree with one entry for each
1025 // module => reading should be done accordingly
1026 // or a classic file and do smth. like this:
1028 // 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");
1042 printf("filtmp %s\n",filtmp);
1046 while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) {
1048 Error("ReadBaseline","Anode number not in increasing order!",
1057 Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",
1066 //____________________________________________
1067 Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) {
1068 // To the 10 to 8 bit lossive compression.
1069 // code from Davide C. and Albert W.
1071 if (signal < 128) return signal;
1072 if (signal < 256) return (128+((signal-128)>>1));
1073 if (signal < 512) return (192+((signal-256)>>3));
1074 if (signal < 1024) return (224+((signal-512)>>4));
1079 //____________________________________________
1080 Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) {
1081 // Undo the lossive 10 to 8 bit compression.
1082 // code from Davide C. and Albert W.
1083 if (signal < 0 || signal > 255) {
1084 printf("<Convert8to10> out of range %d \n",signal);
1088 if (signal < 128) return signal;
1090 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
1091 else return (128+((signal-128)<<1)+1);
1094 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
1095 else return (256+((signal-192)<<3)+4);
1097 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
1098 else return (512+((signal-224)<<4)+7);
1103 //____________________________________________
1104 AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){
1105 //Return the correct map.
1106 return ((i==0)? fHitMap1 : fHitMap2);
1110 //____________________________________________
1111 void AliITSsimulationSDD::ZeroSuppression(const char *option) {
1112 // perform the zero suppresion
1113 if (strstr(option,"2D")) {
1114 //Init2D(); // activate if param change module by module
1116 } else if (strstr(option,"1D")) {
1117 //Init1D(); // activate if param change module by module
1119 } else StoreAllDigits();
1123 //____________________________________________
1124 void AliITSsimulationSDD::Init2D(){
1125 // read in and prepare arrays: fD, fT1, fT2
1126 // savemu[nanodes], savesigma[nanodes]
1127 // read baseline and noise from file - either a .root file and in this
1128 // case data should be organised in a tree with one entry for each
1129 // module => reading should be done accordingly
1130 // or a classic file and do smth. like this ( code from Davide C. and
1133 // Read 2D zero-suppression parameters for SDD
1136 if (!strstr(fParam,"file")) return;
1138 Int_t na,pos,tempTh;
1140 Float_t *savemu = new Float_t [fNofMaps];
1141 Float_t *savesigma = new Float_t [fNofMaps];
1142 char input[100],basel[100],par[100];
1146 Int_t minval = fResponse->MinVal();
1148 fResponse->Filenames(input,basel,par);
1152 filtmp = gSystem->ExpandPathName(fFileName.Data());
1153 FILE *param = fopen(filtmp,"r");
1157 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1159 Error("Init2D ","Anode number not in increasing order!",
1164 savesigma[na]=sigma;
1165 if ((2.*sigma) < mu) {
1166 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1169 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1170 if (tempTh < 0) tempTh=0;
1172 tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval;
1173 if (tempTh < 0) tempTh=0;
1179 Error("Init2D "," THE FILE %s DOES NOT EXIST !",
1187 delete [] savesigma;
1190 //____________________________________________
1191 void AliITSsimulationSDD::Compress2D(){
1193 // simple ITS cluster finder -- online zero-suppression conditions
1198 Int_t minval = fResponse->MinVal();
1199 Bool_t write=fResponse->OutputOption();
1200 Bool_t do10to8=fResponse->Do10to8();
1202 Int_t nz, nl, nh, low, i, j;
1204 for (i=0; i<fNofMaps; i++) {
1205 CompressionParam(i,db,tl,th);
1210 for (j=0; j<fMaxNofSamples; j++) {
1211 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1212 signal -= db; // if baseline eq. is done here
1213 if (signal <= 0) {nz++; continue;}
1214 if ((signal - tl) < minval) low++;
1215 if ((signal - th) >= minval) {
1218 FindCluster(i,j,signal,minval,cond);
1219 if (cond && j && ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)) {
1220 if(do10to8) signal = Convert10to8(signal);
1221 AddDigit(i,j,signal);
1223 } else if ((signal - tl) >= minval) nl++;
1224 } // loop time samples
1225 if (write) TreeB()->Fill(nz,nl,nh,low,i+1);
1230 sprintf(hname,"TNtuple%d_%d",fModule,fEvent);
1231 TreeB()->Write(hname);
1238 //_____________________________________________________________________________
1239 void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal,
1240 Int_t minval,Bool_t &cond){
1242 // Find clusters according to the online 2D zero-suppression algorithm
1245 Bool_t do10to8=fResponse->Do10to8();
1249 fHitMap2->FlagHit(i,j);
1251 // check the online zero-suppression conditions
1253 const Int_t maxNeighbours = 4;
1257 Int_t xList[maxNeighbours], yList[maxNeighbours];
1258 fSegmentation->Neighbours(i,j,&nn,xList,yList);
1260 for (in=0; in<nn; in++) {
1263 if (fHitMap2->TestHit(ix,iy)==kUnused) {
1264 CompressionParam(ix,dbx,tlx,thx);
1265 Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy));
1266 qn -= dbx; // if baseline eq. is done here
1267 if ((qn-tlx) < minval) {
1268 fHitMap2->FlagHit(ix,iy);
1271 if ((qn - thx) >= minval) high=kTRUE;
1273 if(do10to8) signal = Convert10to8(signal);
1274 AddDigit(i,j,signal);
1276 if(do10to8) qns = Convert10to8(qn);
1278 if (!high) AddDigit(ix,iy,qns);
1280 if(!high) fHitMap2->FlagHit(ix,iy);
1283 } // loop over neighbours
1287 //____________________________________________
1288 void AliITSsimulationSDD::Init1D(){
1289 // this is just a copy-paste of input taken from 2D algo
1290 // Torino people should give input
1292 // Read 1D zero-suppression parameters for SDD
1295 if (!strstr(fParam,"file")) return;
1297 Int_t na,pos,tempTh;
1299 Float_t *savemu = new Float_t [fNofMaps];
1300 Float_t *savesigma = new Float_t [fNofMaps];
1301 char input[100],basel[100],par[100];
1305 Int_t minval = fResponse->MinVal();
1306 fResponse->Filenames(input,basel,par);
1309 // set first the disable and tol param
1312 filtmp = gSystem->ExpandPathName(fFileName.Data());
1313 FILE *param = fopen(filtmp,"r");
1317 fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]);
1318 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1320 Error("Init1D ","Anode number not in increasing order!",
1325 savesigma[na]=sigma;
1326 if ((2.*sigma) < mu) {
1327 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1330 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1331 if (tempTh < 0) tempTh=0;
1336 Error("Init1D "," THE FILE %s DOES NOT EXIST !",
1344 delete [] savesigma;
1350 //____________________________________________
1351 void AliITSsimulationSDD::Compress1D(){
1352 // 1D zero-suppression algorithm (from Gianluca A.)
1354 Int_t dis,tol,thres,decr,diff;
1356 UChar_t *str=fStream->Stream();
1359 Bool_t do10to8=fResponse->Do10to8();
1363 for (k=0; k<2; k++) {
1366 for (i=0; i<fNofMaps/2; i++) {
1367 Bool_t firstSignal=kTRUE;
1368 Int_t idx=i+k*fNofMaps/2;
1369 CompressionParam(idx,decr,thres);
1370 for (j=0; j<fMaxNofSamples; j++) {
1371 Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j));
1372 signal -= decr; // if baseline eq.
1373 if(do10to8) signal = Convert10to8(signal);
1374 if (signal <= thres) {
1378 // write diff in the buffer for HuffT
1379 str[counter]=(UChar_t)diff;
1384 if (diff > 127) diff=127;
1385 if (diff < -128) diff=-128;
1388 // tol has changed to 8 possible cases ? - one can write
1389 // this if(TMath::Abs(diff)<tol) ... else ...
1390 if(TMath::Abs(diff)<tol) diff=0;
1391 // or keep it as it was before
1393 if (tol==1 && (diff >= -2 && diff <= 1)) diff=0;
1394 if (tol==2 && (diff >= -4 && diff <= 3)) diff=0;
1395 if (tol==3 && (diff >= -16 && diff <= 15)) diff=0;
1397 AddDigit(idx,j,last+diff);
1399 AddDigit(idx,j,signal);
1403 // write diff in the buffer used to compute Huffman tables
1404 if (firstSignal) str[counter]=(UChar_t)signal;
1405 else str[counter]=(UChar_t)diff;
1410 } // loop time samples
1411 } // loop anodes one half of detector
1415 fStream->CheckCount(counter);
1417 // open file and write out the stream of diff's
1419 static Bool_t open=kTRUE;
1420 static TFile *outFile;
1421 Bool_t write = fResponse->OutputOption();
1425 SetFileName("stream.root");
1426 cout<<"filename "<<fFileName<<endl;
1427 outFile=new TFile(fFileName,"recreate");
1428 cout<<"I have opened "<<fFileName<<" file "<<endl;
1435 fStream->ClearStream();
1437 // back to galice.root file
1439 TTree *fAli=gAlice->TreeK();
1442 if (fAli) file =fAli->GetCurrentFile();
1447 //____________________________________________
1448 void AliITSsimulationSDD::StoreAllDigits(){
1449 // if non-zero-suppressed data
1451 Bool_t do10to8=fResponse->Do10to8();
1453 Int_t i, j, digits[3];
1454 for (i=0; i<fNofMaps; i++) {
1455 for (j=0; j<fMaxNofSamples; j++) {
1456 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1457 if(do10to8) signal = Convert10to8(signal);
1458 if(do10to8) signal = Convert8to10(signal);
1462 fITS->AddRealDigit(1,digits);
1466 //____________________________________________
1468 void AliITSsimulationSDD::CreateHistograms(Int_t scale){
1469 // Creates histograms of maps for debugging
1473 fHis=new TObjArray(fNofMaps);
1474 for (i=0;i<fNofMaps;i++) {
1475 TString sddName("sdd_");
1477 sprintf(candNum,"%d",i+1);
1478 sddName.Append(candNum);
1479 (*fHis)[i] = new TH1F(sddName.Data(),"SDD maps",
1480 scale*fMaxNofSamples,0.,(Float_t) scale*fMaxNofSamples);
1484 //____________________________________________
1485 void AliITSsimulationSDD::FillHistograms(){
1486 // fill 1D histograms from map
1489 for( Int_t i=0; i<fNofMaps; i++) {
1490 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
1491 Int_t nsamples = hist->GetNbinsX();
1492 for( Int_t j=0; j<nsamples; j++) {
1493 Double_t signal=fHitMap2->GetSignal(i,j);
1494 hist->Fill((Float_t)j,signal);
1499 //____________________________________________
1501 void AliITSsimulationSDD::ResetHistograms(){
1503 // Reset histograms for this detector
1506 for (i=0;i<fNofMaps;i++ ) {
1507 if ((*fHis)[i]) ((TH1F*)(*fHis)[i])->Reset();
1513 //____________________________________________
1515 TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
1516 // Fills a histogram from a give anode.
1517 if (!fHis) return 0;
1519 if(wing <=0 || wing > 2) {
1520 cout << "Wrong wing number: " << wing << endl;
1523 if(anode <=0 || anode > fNofMaps/2) {
1524 cout << "Wrong anode number: " << anode << endl;
1528 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
1529 return (TH1F*)((*fHis)[index]);
1532 //____________________________________________
1534 void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
1535 // Writes the histograms to a file
1540 for(i=0; i<fNofMaps; i++) (*fHis)[i]->Write(); //fAdcs[i]->Write();
1543 //____________________________________________
1544 Float_t AliITSsimulationSDD::GetNoise() {
1545 // Returns the noise value
1547 //Bool_t do10to8=fResponse->Do10to8();
1548 //noise will always be in the liniar part of the signal
1551 Int_t threshold=fT1[0];
1553 char opt1[20], opt2[20];
1554 fResponse->ParamOptions(opt1,opt2);
1556 char *same = strstr(opt1,"same");
1557 Float_t noise,baseline;
1559 fResponse->GetNoiseParam(noise,baseline);
1561 static Bool_t readfile=kTRUE;
1562 //read baseline and noise from file
1563 if (readfile) ReadBaseline();
1567 TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2");
1568 if(c2) delete c2->GetPrimitive("noisehist");
1569 if(c2) delete c2->GetPrimitive("anode");
1570 else c2=new TCanvas("c2");
1572 c2->SetFillColor(0);
1574 TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold);
1575 TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0.,(float)fMaxNofSamples);
1577 for (i=0;i<fNofMaps;i++) {
1578 CompressionParam(i,decr,threshold);
1579 if (!same) GetAnodeBaseline(i,baseline,noise);
1581 for (k=0;k<fMaxNofSamples;k++) {
1582 Float_t signal=(Float_t)fHitMap2->GetSignal(i,k);
1583 //if (signal <= (float)threshold) noisehist->Fill(signal-baseline);
1584 if (signal <= (float)threshold) noisehist->Fill(signal);
1585 anode->Fill((float)k,signal);
1590 TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold);
1591 noisehist->Fit("gnoise","RQ");
1594 Float_t mnoise = gnoise->GetParameter(1);
1595 cout << "mnoise : " << mnoise << endl;
1596 Float_t rnoise = gnoise->GetParameter(2);
1597 cout << "rnoise : " << rnoise << endl;
1602 //____________________________________________
1604 void AliITSsimulationSDD::Print() {
1606 // Print SDD simulation Parameters
1608 cout << "**************************************************" << endl;
1609 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1610 cout << "**************************************************" << endl;
1611 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1612 cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl;
1613 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1614 cout << "Number pf Anodes used: " << fNofMaps << endl;
1615 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1616 cout << "Scale size factor: " << fScaleSize << endl;
1617 cout << "**************************************************" << endl;