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 **************************************************************************/
18 #include <Riostream.h>
30 #include "AliITSMapA2.h"
31 #include "AliITSRawData.h"
32 #include "AliITSdigitSPD.h"
33 #include "AliITSetfSDD.h"
34 #include "AliITSmodule.h"
35 #include "AliITShit.h"
36 #include "AliITSpList.h"
37 #include "AliITSCalibrationSDD.h"
38 #include "AliITSsimulationSDD.h"
42 ClassImp(AliITSsimulationSDD)
43 ////////////////////////////////////////////////////////////////////////
45 // Written by Piergiorgio Cerello //
46 // November 23 1999 //
48 // AliITSsimulationSDD is the simulation of SDDs. //
49 ////////////////////////////////////////////////////////////////////////
51 //______________________________________________________________________
52 AliITSsimulationSDD::AliITSsimulationSDD():
66 fCrosstalkFlag(kFALSE),
71 // Default constructor
72 SetPerpendTracksFlag();
76 //______________________________________________________________________
77 AliITSsimulationSDD::AliITSsimulationSDD(const AliITSsimulationSDD &source) :
78 AliITSsimulation(source),
80 fHitMap2(source.fHitMap2),
81 fHitSigMap2(source.fHitSigMap2),
82 fHitNoiMap2(source.fHitNoiMap2),
83 fElectronics(source.fElectronics),
86 fOutZR(source.fOutZR),
87 fOutZI(source.fOutZI),
88 fAnodeFire(source.fAnodeFire),
91 fCrosstalkFlag(source.fCrosstalkFlag),
92 fDoFFT(source.fDoFFT),
93 fNofMaps(source.fNofMaps),
94 fMaxNofSamples(source.fMaxNofSamples),
95 fScaleSize(source.fScaleSize){
96 // Copy constructor to satify Coding roules only.
99 //______________________________________________________________________
100 AliITSsimulationSDD& AliITSsimulationSDD::operator=(const AliITSsimulationSDD &src){
101 // Assignment operator to satify Coding roules only.
103 if(this==&src) return *this;
104 Error("AliITSsimulationSDD","Not allowed to make a = with "
105 "AliITSsimulationSDD Using default creater instead");
109 //______________________________________________________________________
110 AliITSsimulation& AliITSsimulationSDD::operator=(const AliITSsimulation &src){
111 // Assignment operator to satify Coding roules only.
113 if(this==&src) return *this;
114 Error("AliITSsimulationSSD","Not allowed to make a = with "
115 "AliITSsimulationSDD Using default creater instead");
119 //______________________________________________________________________
120 AliITSsimulationSDD::AliITSsimulationSDD(AliITSDetTypeSim* dettyp):
121 AliITSsimulation(dettyp),
134 fCrosstalkFlag(kFALSE),
139 // Default Constructor
142 //______________________________________________________________________
143 void AliITSsimulationSDD::Init(){
144 // Standard Constructor
146 AliITSsegmentationSDD* seg = (AliITSsegmentationSDD*)GetSegmentationModel(1);
147 fScaleSize = ScaleFourier(seg);
148 SetPerpendTracksFlag();
152 AliITSSimuParam* simpar = fDetType->GetSimuParam();
153 fpList = new AliITSpList( seg->Npz(),
154 fScaleSize*seg->Npx() );
155 fHitSigMap2 = new AliITSMapA2(seg,fScaleSize,1);
156 fHitNoiMap2 = new AliITSMapA2(seg,fScaleSize,1);
157 fHitMap2 = fHitSigMap2;
159 fNofMaps = seg->Npz();
160 fMaxNofSamples = seg->Npx();
161 fAnodeFire = new Bool_t [fNofMaps];
163 Float_t sddWidth = seg->Dz();
164 Float_t anodePitch = seg->Dpz(0);
165 Double_t timeStep = (Double_t)seg->Dpx(0);
167 if(anodePitch*(fNofMaps/2) > sddWidth) {
168 AliWarning(Form("Too many anodes %d or too big pitch %f ",
169 fNofMaps/2,anodePitch));
173 fElectronics = new AliITSetfSDD(timeStep/fScaleSize,
174 simpar->GetSDDElectronics());
177 fITS = (AliITS*)gAlice->GetModule("ITS");
179 fInZR = new Double_t [fScaleSize*fMaxNofSamples];
180 fInZI = new Double_t [fScaleSize*fMaxNofSamples];
181 fOutZR = new Double_t [fScaleSize*fMaxNofSamples];
182 fOutZI = new Double_t [fScaleSize*fMaxNofSamples];
184 //______________________________________________________________________
185 AliITSsimulationSDD::~AliITSsimulationSDD() {
199 if(fInZR) delete [] fInZR;
200 if(fInZI) delete [] fInZI;
201 if(fOutZR) delete [] fOutZR;
202 if(fOutZI) delete [] fOutZI;
203 if(fAnodeFire) delete [] fAnodeFire;
205 //______________________________________________________________________
206 void AliITSsimulationSDD::InitSimulationModule( Int_t module, Int_t event ) {
207 // create maps to build the lists of tracks for each summable digit
211 memset(fAnodeFire,0,sizeof(Bool_t)*fNofMaps);
213 //______________________________________________________________________
214 void AliITSsimulationSDD::ClearMaps() {
217 fHitSigMap2->ClearMap();
218 fHitNoiMap2->ClearMap();
220 //______________________________________________________________________
221 void AliITSsimulationSDD::FastFourierTransform(Double_t *real,
222 Double_t *imag,Int_t direction) {
223 // Do a Fast Fourier Transform
225 Int_t samples = fElectronics->GetSamples();
226 Int_t l = (Int_t) ((log((Float_t) samples)/log(2.))+0.5);
229 Int_t m2 = samples/m1;
231 for(i=1; i<=l; i++) {
232 for(j=0; j<samples; j += m1) {
234 for(k=j; k<= j+m-1; k++) {
235 Double_t wsr = fElectronics->GetWeightReal(p);
236 Double_t wsi = fElectronics->GetWeightImag(p);
237 if(direction == -1) wsi = -wsi;
238 Double_t xr = *(real+k+m);
239 Double_t xi = *(imag+k+m);
240 *(real+k+m) = wsr*(*(real+k)-xr) - wsi*(*(imag+k)-xi);
241 *(imag+k+m) = wsr*(*(imag+k)-xi) + wsi*(*(real+k)-xr);
251 for(j=0; j<samples; j++) {
255 for(i1=1; i1<=l; i1++) {
258 p = p + p + j2 - j1 - j1;
261 Double_t xr = *(real+j);
262 Double_t xi = *(imag+j);
263 *(real+j) = *(real+p);
264 *(imag+j) = *(imag+p);
269 if(direction == -1) {
270 for(i=0; i<samples; i++) {
271 *(real+i) /= samples;
272 *(imag+i) /= samples;
274 } // end if direction == -1
278 //______________________________________________________________________
279 void AliITSsimulationSDD::SDigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
280 // digitize module using the "slow" detector simulator creating
283 TObjArray *fHits = mod->GetHits();
284 Int_t nhits = fHits->GetEntriesFast();
287 InitSimulationModule( md, ev );
288 HitsToAnalogDigits( mod ); // fills fHitMap2 which is = fHitSigmap2
289 ChargeToSignal( fModule,kFALSE,kTRUE ); // - Process signal adding gain without adding noise
290 fHitMap2 = fHitNoiMap2; // - Swap to noise map
291 ChargeToSignal( fModule,kTRUE,kFALSE ); // - Process only noise
292 fHitMap2 = fHitSigMap2; // - Return to signal map
296 //______________________________________________________________________
297 Bool_t AliITSsimulationSDD::AddSDigitsToModule(TClonesArray *pItemArray,
299 // Add Summable digits to module maps.
300 AliITSSimuParam* simpar = fDetType->GetSimuParam();
301 Int_t nItems = pItemArray->GetEntries();
302 Double_t maxadc = simpar->GetSDDMaxAdc();
305 // cout << "Adding "<< nItems <<" SDigits to module " << fModule << endl;
306 for( Int_t i=0; i<nItems; i++ ) {
307 AliITSpListItem * pItem = (AliITSpListItem *)(pItemArray->At( i ));
308 if( pItem->GetModule() != fModule ) {
309 Error( "AliITSsimulationSDD","Error reading, SDigits module "
310 "%d != current module %d: exit",
311 pItem->GetModule(), fModule );
315 if(pItem->GetSignal()>0.0 ) sig = kTRUE;
317 fpList->AddItemTo( mask, pItem ); // Add SignalAfterElect + noise
318 AliITSpListItem * pItem2 = fpList->GetpListItem( pItem->GetIndex() );
319 Double_t sigAE = pItem2->GetSignalAfterElect();
320 if( sigAE >= maxadc ) sigAE = maxadc-1; // avoid overflow signal
323 fpList->GetMapIndex( pItem->GetIndex(), ia, it );
324 fHitMap2->SetHit( ia, it, sigAE );
325 fAnodeFire[ia] = kTRUE;
329 //______________________________________________________________________
330 void AliITSsimulationSDD::FinishSDigitiseModule() {
331 // digitize module using the "slow" detector simulator from
332 // the sum of summable digits.
336 //______________________________________________________________________
337 void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
338 // create maps to build the lists of tracks for each digit
340 TObjArray *fHits = mod->GetHits();
341 Int_t nhits = fHits->GetEntriesFast();
343 InitSimulationModule( md, ev );
346 HitsToAnalogDigits( mod );
347 ChargeToSignal( fModule,kTRUE,kTRUE ); // process signal + noise
349 for( Int_t i=0; i<fNofMaps; i++ ) {
350 for( Int_t j=0; j<fMaxNofSamples; j++ ) {
351 Int_t jdx = j*fScaleSize;
352 Int_t index = fpList->GetHitIndex( i, j );
353 AliITSpListItem pItemTmp2( fModule, index, 0. );
354 // put the fScaleSize analog digits in only one
355 for( Int_t ik=0; ik<fScaleSize; ik++ ) {
356 AliITSpListItem *pItemTmp = fpList->GetpListItem( i, jdx+ik );
357 if( pItemTmp == 0 ) continue;
358 pItemTmp2.Add( pItemTmp );
360 fpList->DeleteHit( i, j );
361 fpList->AddItemTo( 0, &pItemTmp2 );
367 //______________________________________________________________________
368 void AliITSsimulationSDD::FinishDigits() {
369 // introduce the electronics effects and do zero-suppression if required
371 if( fCrosstalkFlag ) ApplyCrosstalk(fModule);
373 AliITSCalibrationSDD* res = (AliITSCalibrationSDD*)GetCalibrationModel(fModule);
374 Bool_t isZeroSupp = res->GetZeroSupp();
375 if (isZeroSupp) Compress2D();
376 else StoreAllDigits();
378 //______________________________________________________________________
379 void AliITSsimulationSDD::HitsToAnalogDigits( AliITSmodule *mod ) {
380 // create maps to build the lists of tracks for each digit
381 AliITSsegmentationSDD* seg = (AliITSsegmentationSDD*)GetSegmentationModel(1);
382 AliITSCalibrationSDD* res = (AliITSCalibrationSDD*)GetCalibrationModel(fModule);
383 AliITSSimuParam* simpar = fDetType->GetSimuParam();
384 TObjArray *hits = mod->GetHits();
385 Int_t nhits = hits->GetEntriesFast();
387 // Int_t arg[6] = {0,0,0,0,0,0};
388 Int_t nofAnodes = fNofMaps/2;
389 Double_t sddLength = seg->Dx();
390 Double_t anodePitch = seg->Dpz(0);
391 Double_t timeStep = seg->Dpx(0);
392 Double_t driftSpeed ; // drift velocity (anode dependent)
393 Double_t nanoampToADC = simpar->GetSDDMaxAdc()/simpar->GetSDDDynamicRange(); // maxadc/topValue;
394 Double_t cHloss = simpar->GetSDDChargeLoss();
396 simpar->GetSDDDiffCoeff(dfCoeff,s1); // Signal 2d Shape
397 Double_t eVpairs = simpar->GetGeVToCharge()*1.0E9; // 3.6 eV by def.
398 Double_t nsigma = simpar->GetNSigmaIntegration(); //
399 Int_t nlookups = simpar->GetGausNLookUp(); //
400 Float_t jitter = simpar->GetSDDJitterError(); //
402 // Piergiorgio's part (apart for few variables which I made float
403 // when i thought that can be done
404 // Fill detector maps with GEANT hits
405 // loop over hits in the module
407 const Float_t kconv = 1.0e+6; // GeV->KeV
409 Int_t iWing; // which detector wing/side.
410 Int_t ii,kk,ka,kt; // loop indexs
411 Int_t ia,it,index; // sub-pixel integration indexies
412 Int_t iAnode; // anode number.
413 Int_t timeSample; // time buckett.
414 Int_t anodeWindow; // anode direction charge integration width
415 Int_t timeWindow; // time direction charge integration width
416 Int_t jamin,jamax; // anode charge integration window
417 Int_t jtmin,jtmax; // time charge integration window
418 Int_t nsplitAn; // the number of splits in anode and time windows
419 Int_t nsplitTb; // the number of splits in anode and time windows
420 Int_t nOfSplits; // number of times track length is split into
421 Float_t nOfSplitsF; // Floating point version of nOfSplits.
422 Float_t kkF; // Floating point version of loop index kk.
423 Double_t pathInSDD; // Track length in SDD.
424 Double_t drPath; // average position of track in detector. in microns
425 Double_t drTime; // Drift time
426 Double_t avDrft; // x position of path length segment in cm.
427 Double_t avAnode; // Anode for path length segment in Anode number (float)
428 Double_t zAnode; // Floating point anode number.
429 Double_t driftPath; // avDrft in microns.
430 Double_t width; // width of signal at anodes.
431 Double_t depEnergy; // Energy deposited in this GEANT step.
432 Double_t xL[3],dxL[3]; // local hit coordinates and diff.
433 Double_t sigA; // sigma of signal at anode.
434 Double_t sigT; // sigma in time/drift direction for track segment
435 Double_t aStep,aConst; // sub-pixel size and offset anode
436 Double_t tStep,tConst; // sub-pixel size and offset time
437 Double_t amplitude; // signal amplitude for track segment in nanoAmpere
438 Double_t chargeloss; // charge loss for track segment.
439 Double_t anodeAmplitude; // signal amplitude in anode direction
440 Double_t aExpo; // exponent of Gaussian anode direction
441 Double_t timeAmplitude; // signal amplitude in time direction
442 Double_t tExpo; // exponent of Gaussian time direction
443 Double_t tof; // Time of flight in ns of this step.
445 for(ii=0; ii<nhits; ii++) {
446 if(!mod->LineSegmentL(ii,xL[0],dxL[0],xL[1],dxL[1],xL[2],dxL[2],
447 depEnergy,itrack)) continue;
449 if(xloc>0) iWing=0; // left side, carlos channel 0
450 else iWing=1; // right side
452 Float_t zloc=xL[2]+0.5*dxL[2];
453 zAnode=seg->GetAnodeFromLocal(xloc,zloc); // anode number in the range 0.-511.
454 driftSpeed = res->GetDriftSpeedAtAnode(zAnode);
455 if(timeStep*fMaxNofSamples < sddLength/driftSpeed) {
456 AliWarning("Time Interval > Allowed Time Interval");
461 Form("fTrack = %d hit=%d module=%d This particle has passed without losing energy!",
462 itrack,ii,mod->GetIndex()));
464 // continue if the particle did not lose energy
465 // passing through detector
466 } // end if !depEnergy
469 AliITShit* h=(AliITShit*)hits->At(ii);
470 if(h) tof=h->GetTOF()*1E9;
471 AliDebug(1,Form("TOF for hit %d on mod %d (particle %d)=%g",ii,fModule,h->Track(),tof));
473 xL[0] += 0.0001*gRandom->Gaus( 0, jitter ); //
474 pathInSDD = TMath::Sqrt(dxL[0]*dxL[0]+dxL[1]*dxL[1]+dxL[2]*dxL[2]);
476 if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
477 drPath = TMath::Abs(10000.*(dxL[0]+2.*xL[0])*0.5);
478 drPath = sddLength-drPath;
480 AliDebug(1, // this should be fixed at geometry level
481 Form("negative drift path drPath=%e sddLength=%e dxL[0]=%e xL[0]=%e",
482 drPath,sddLength,dxL[0],xL[0]));
484 } // end if drPath < 0
486 // Compute number of segments to brake step path into
487 drTime = drPath/driftSpeed; // Drift Time
488 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);// Sigma along the anodes
489 // calcuate the number of time the path length should be split into.
490 nOfSplits = (Int_t) (1. + 10000.*pathInSDD/sigA);
491 if(fFlag) nOfSplits = 1;
493 // loop over path segments, init. some variables.
494 depEnergy /= nOfSplits;
495 nOfSplitsF = (Float_t) nOfSplits;
496 Float_t theAverage=0.,theSteps=0.;
497 for(kk=0;kk<nOfSplits;kk++) { // loop over path segments
498 kkF = (Float_t) kk + 0.5;
499 avDrft = xL[0]+dxL[0]*kkF/nOfSplitsF;
500 avAnode = xL[2]+dxL[2]*kkF/nOfSplitsF;
503 zAnode = seg->GetAnodeFromLocal(avDrft,avAnode);
504 driftSpeed = res->GetDriftSpeedAtAnode(zAnode);
505 driftPath = TMath::Abs(10000.*avDrft);
506 driftPath = sddLength-driftPath;
508 AliDebug(1, // this should be fixed at geometry level
509 Form("negative drift path driftPath=%e sddLength=%e avDrft=%e dxL[0]=%e xL[0]=%e",
510 driftPath,sddLength,avDrft,dxL[0],xL[0]));
512 } // end if driftPath < 0
513 drTime = driftPath/driftSpeed; // drift time for segment.
514 // Sigma along the anodes for track segment.
515 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);
516 sigT = sigA/driftSpeed;
518 drTime+=tof; // take into account Time Of Flight from production point
519 timeSample = (Int_t) (fScaleSize*drTime/timeStep + 1.001); // time bin in range 1-256 !!!
520 if(zAnode>nofAnodes) zAnode-=nofAnodes; // to have the anode number between 0. and 256.
521 iAnode = (Int_t) (1.001+zAnode); // iAnode in range 1-256 !!!!
523 // Peak amplitude in nanoAmpere
524 amplitude = fScaleSize*160.*depEnergy/
525 (timeStep*eVpairs*2.*acos(-1.));
526 chargeloss = 1.-cHloss*driftPath/1000.;
527 amplitude *= chargeloss;
528 width = 2.*nsigma/(nlookups-1);
532 aStep = anodePitch/(nsplitAn*sigA);
533 aConst = zAnode*anodePitch/sigA;
534 tStep = timeStep/(nsplitTb*fScaleSize*sigT);
535 tConst = drTime/sigT;
536 // Define SDD window corresponding to the hit
537 anodeWindow = (Int_t)(nsigma*sigA/anodePitch+1);
538 timeWindow = (Int_t) (fScaleSize*nsigma*sigT/timeStep+1.);
539 jamin = (iAnode - anodeWindow - 2)*nsplitAn+1;
540 if(jamin <= 0) jamin = 1;
541 if(jamin > nofAnodes*nsplitAn){
542 AliDebug(1,Form("Energy deposition completely outside anode acceptance: anode min=%d",jamin));
545 jamax = (iAnode + anodeWindow + 2)*nsplitAn;
546 if(jamax > nofAnodes*nsplitAn) jamax = nofAnodes*nsplitAn;
548 AliDebug(1,Form("Energy deposition completely outside anode acceptance: anode max=%d",jamax));
551 jtmin = (Int_t)(timeSample-timeWindow-2)*nsplitTb+1;
552 if(jtmin <= 0) jtmin = 1;
553 if(jtmin > fScaleSize*fMaxNofSamples*nsplitTb){
554 AliDebug(1,Form("Energy deposition completely outside time acceptance: time sample min=%d tof=%f",jtmin,tof));
557 jtmax = (Int_t)(timeSample+timeWindow+2)*nsplitTb;
558 if(jtmax > fScaleSize*fMaxNofSamples*nsplitTb) jtmax = fScaleSize*fMaxNofSamples*nsplitTb;
560 AliDebug(1,Form("Energy deposition completely outside time acceptance: time sample max=%d tof=%f",jtmax,tof));
564 // Spread the charge in the anode-time window
565 for(ka=jamin; ka <=jamax; ka++) {
566 ia = (ka-1)/nsplitAn + 1;
568 if(ia > nofAnodes) ia = nofAnodes;
569 aExpo = (aStep*(ka-0.5)-aConst);
570 if(TMath::Abs(aExpo) > nsigma) anodeAmplitude = 0.;
572 Int_t theBin = (Int_t) ((aExpo+nsigma)/width+0.5);
573 anodeAmplitude = amplitude*simpar->GetGausLookUp(theBin);
575 // index starts from 0
576 index = iWing*nofAnodes+ia-1;
578 for(kt=jtmin; kt<=jtmax; kt++) {
579 it = (kt-1)/nsplitTb+1; // it starts from 1
581 if(it>fScaleSize*fMaxNofSamples)
582 it = fScaleSize*fMaxNofSamples;
583 tExpo = (tStep*(kt-0.5)-tConst);
584 if(TMath::Abs(tExpo) > nsigma) timeAmplitude = 0.;
586 Int_t theBin = (Int_t) ((tExpo+nsigma)/width+0.5);
587 timeAmplitude = anodeAmplitude*simpar->GetGausLookUp(theBin)*aStep*tStep;
589 timeAmplitude *= nanoampToADC;
590 // ListOfFiredCells(arg,timeAmplitude,alst,padr);
591 Double_t charge = timeAmplitude;
592 charge += fHitMap2->GetSignal(index,it-1);
593 fHitMap2->SetHit(index, it-1, charge);
594 fpList->AddSignal(index,it-1,itrack,ii-1,
595 mod->GetIndex(),timeAmplitude);
596 fAnodeFire[index] = kTRUE;
597 } // end loop over time in window
598 } // end if anodeAmplitude
599 } // loop over anodes in window
600 } // end loop over "sub-hits"
601 } // end loop over hits
604 //____________________________________________
605 void AliITSsimulationSDD::AddDigit( Int_t i, Int_t j, Int_t signalc, Int_t signale) {
607 Int_t size = AliITSdigit::GetNTracks();
610 Int_t * tracks = new Int_t[size];
611 Int_t * hits = new Int_t[size];
613 Float_t * charges = new Float_t[size];
619 AliITSpListItem *pItem = fpList->GetpListItem( i, j );
622 for( Int_t l=0; l<size; l++ ) {
628 Int_t idtrack = pItem->GetTrack( 0 );
629 if( idtrack >= 0 ) phys = pItem->GetSignal();
632 for( Int_t l=0; l<size; l++ ) if(l<pItem->GetMaxKept()) {
633 tracks[l] = pItem->GetTrack( l );
634 hits[l] = pItem->GetHit( l );
635 charges[l] = pItem->GetSignal( l );
643 fITS->AddSimDigit( 1, phys, digits, tracks, hits, charges, signale );
648 //______________________________________________________________________
649 void AliITSsimulationSDD::ChargeToSignal(Int_t mod,Bool_t bAddNoise, Bool_t bAddGain) {
650 // add baseline, noise, gain, electronics and ADC saturation effects
651 // apply dead channels
653 AliITSCalibrationSDD* res = (AliITSCalibrationSDD*)GetCalibrationModel(mod);
659 AliITSSimuParam* simpar = fDetType->GetSimuParam();
660 Float_t maxadc = simpar->GetSDDMaxAdc();
661 Int_t nGroup=fScaleSize;
662 if(res->IsAMAt20MHz()){
666 for (i=0;i<fNofMaps;i++) {
667 if( !fAnodeFire[i] ) continue;
668 baseline = res->GetBaseline(i);
669 noise = res->GetNoise(i);
670 gain = res->GetChannelGain(i)/fDetType->GetAverageGainSDD();
671 if(res->IsBad()) gain=0.;
672 if( res->IsChipBad(res->GetChip(i)) )gain=0.;
673 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
674 fInZR[k] = fHitMap2->GetSignal(i,k);
675 if(bAddGain) fInZR[k]*=gain;
677 contrib = (baseline + noise*gRandom->Gaus());
683 for(k=0; k<fMaxNofSamples; k++) {
684 Double_t newcont = 0.;
685 Double_t maxcont = 0.;
686 for(kk=0;kk<fScaleSize;kk++) {
687 newcont = fInZR[fScaleSize*k+kk];
688 if(newcont > maxcont) maxcont = newcont;
691 if (newcont >= maxadc) newcont = maxadc -1;
692 if(newcont >= baseline){
693 Warning("","newcont=%f>=baseline=%f",newcont,baseline);
696 fHitMap2->SetHit(i,k,newcont);
699 FastFourierTransform(&fInZR[0],&fInZI[0],1);
700 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
701 Double_t rw = fElectronics->GetTraFunReal(k);
702 Double_t iw = fElectronics->GetTraFunImag(k);
703 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
704 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
706 FastFourierTransform(&fOutZR[0],&fOutZI[0],-1);
707 for(k=0; k<fMaxNofSamples; k++) {
708 Double_t newcont1 = 0.;
709 Double_t maxcont1 = 0.;
710 for(kk=0;kk<nGroup;kk++) {
711 newcont1 = fOutZR[fScaleSize*k+kk];
712 if(newcont1 > maxcont1) maxcont1 = newcont1;
715 if (newcont1 >= maxadc) newcont1 = maxadc -1;
716 fHitMap2->SetHit(i,k,newcont1);
719 } // end for i loop over anodes
723 //______________________________________________________________________
724 void AliITSsimulationSDD::ApplyCrosstalk(Int_t mod) {
725 // function add the crosstalk effect to signal
726 // temporal function, should be checked...!!!
728 // create and inizialice crosstalk map
729 Float_t* ctk = new Float_t[fNofMaps*fMaxNofSamples+1];
731 Error( "ApplyCrosstalk", "no memory for temporal map: exit " );
734 memset( ctk, 0, sizeof(Float_t)*(fNofMaps*fMaxNofSamples+1) );
735 AliITSCalibrationSDD* calibr = (AliITSCalibrationSDD*)GetCalibrationModel(mod);
736 for( Int_t z=0; z<fNofMaps; z++ ) {
737 Double_t baseline = calibr->GetBaseline(z);
743 for( Int_t l=0; l<fMaxNofSamples; l++ ) {
744 Float_t fadc = (Float_t)fHitMap2->GetSignal( z, l );
745 if( fadc > baseline ) {
746 if( on == kFALSE && l<fMaxNofSamples-4 ) {
747 Float_t fadc1 = (Float_t)fHitMap2->GetSignal( z, l+1 );
748 if( fadc1 < fadc ) continue;
755 else { // end fadc > baseline
759 // make smooth derivative
760 Float_t* dev = new Float_t[fMaxNofSamples+1];
761 memset( dev, 0, sizeof(Float_t)*(fMaxNofSamples+1) );
763 Error( "ApplyCrosstalk",
764 "no memory for temporal array: exit " );
767 for( Int_t i=tstart; i<tstop; i++ ) {
768 if( i > 2 && i < fMaxNofSamples-2 )
769 dev[i] = -0.2*fHitMap2->GetSignal( z,i-2 )
770 -0.1*fHitMap2->GetSignal( z,i-1 )
771 +0.1*fHitMap2->GetSignal( z,i+1 )
772 +0.2*fHitMap2->GetSignal( z,i+2 );
775 // add crosstalk contribution to neibourg anodes
776 for( Int_t i=tstart; i<tstop; i++ ) {
778 Int_t i1 = (Int_t)((i-tstart)*.61+tstart+0.5); //
779 Float_t ctktmp = -dev[i1] * 0.25;
781 ctk[anode*fMaxNofSamples+i] += ctktmp;
784 if( anode < fNofMaps ) {
785 ctk[anode*fMaxNofSamples+i] += ctktmp;
790 } // if( nTsteps > 2 )
792 } // if( on == kTRUE )
797 for( Int_t a=0; a<fNofMaps; a++ )
798 for( Int_t t=0; t<fMaxNofSamples; t++ ) {
799 Float_t signal = fHitMap2->GetSignal(a,t)+ctk[a*fMaxNofSamples+t];
800 fHitMap2->SetHit( a, t, signal );
806 //______________________________________________________________________
807 Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) const {
808 // To the 10 to 8 bit lossive compression.
809 // code from Davide C. and Albert W.
811 if (signal < 128) return signal;
812 if (signal < 256) return (128+((signal-128)>>1));
813 if (signal < 512) return (192+((signal-256)>>3));
814 if (signal < 1024) return (224+((signal-512)>>4));
817 //______________________________________________________________________
818 Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) const {
819 // Decompression from 8 to 10 bit
821 if (signal < 0 || signal > 255) {
822 AliWarning(Form("Signal value %d out of range",signal));
824 } // end if signal <0 || signal >255
826 if (signal < 128) return signal;
828 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
829 else return (128+((signal-128)<<1)+1);
830 } // end if signal < 192
832 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
833 else return (256+((signal-192)<<3)+4);
834 } // end if signal < 224
835 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
836 return (512+((signal-224)<<4)+8);
838 //______________________________________________________________________
839 void AliITSsimulationSDD::Compress2D(){
840 // 2D zero-suppression algorithm as described in ALICE-INT-1999-28 V10
841 AliITSCalibrationSDD* res = (AliITSCalibrationSDD*)GetCalibrationModel(fModule);
842 for (Int_t iWing=0; iWing<2; iWing++) {
843 Int_t tL=res->GetZSLowThreshold(iWing);
844 Int_t tH=res->GetZSHighThreshold(iWing);
845 for (Int_t i=0; i<fNofMaps/2; i++) {
846 Int_t ian=i+iWing*fNofMaps/2;
847 if( !fAnodeFire[ian] ) continue;
848 for (Int_t itb=0; itb<fMaxNofSamples; itb++) {
849 Int_t nLow=0, nHigh=0;
850 Float_t cC=fHitMap2->GetSignal(ian,itb);
852 nLow++; // cC is greater than tL
855 // Get "quintuple": WCE
858 if(itb>0) wW=fHitMap2->GetSignal(ian,itb-1);
862 if(itb<fMaxNofSamples-1) eE=fHitMap2->GetSignal(ian,itb+1);
866 if(i<(fNofMaps/2-1)) nN=fHitMap2->GetSignal(ian+1,itb);
870 if(i>0) sS=fHitMap2->GetSignal(ian-1,itb);
874 if(nLow>=2 && nHigh>=1){
875 Int_t signal=(Int_t)cC;
876 Int_t signalc = Convert10to8(signal);
877 Int_t signale = Convert8to10(signalc);
878 signalc-=tL; // subtract low threshold after 10 to 8 bit compression
879 if(signalc>=4) AddDigit(ian,itb,signalc,signale); // store C
887 //______________________________________________________________________
888 void AliITSsimulationSDD::StoreAllDigits(){
889 // store digits for non-zero-suppressed data
890 for (Int_t ian=0; ian<fNofMaps; ian++) {
891 for (Int_t itb=0; itb<fMaxNofSamples; itb++){
892 Int_t signal=(Int_t)(fHitMap2->GetSignal(ian,itb));
893 Int_t signalc = Convert10to8(signal);
894 Int_t signale = Convert8to10(signalc);
895 AddDigit(ian,itb,signalc,signale);
899 //______________________________________________________________________
900 void AliITSsimulationSDD::CreateHistograms(Int_t scale){
901 // Creates histograms of maps for debugging
904 fHis=new TObjArray(fNofMaps);
905 for (i=0;i<fNofMaps;i++) {
906 TString sddName("sdd_");
908 sprintf(candNum,"%d",i+1);
909 sddName.Append(candNum);
910 fHis->AddAt(new TH1F(sddName.Data(),"SDD maps",scale*fMaxNofSamples,
911 0.,(Float_t) scale*fMaxNofSamples), i);
914 //______________________________________________________________________
915 void AliITSsimulationSDD::FillHistograms(){
916 // fill 1D histograms from map
920 for( Int_t i=0; i<fNofMaps; i++) {
921 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
922 Int_t nsamples = hist->GetNbinsX();
923 for( Int_t j=0; j<nsamples; j++) {
924 Double_t signal=fHitMap2->GetSignal(i,j);
925 hist->Fill((Float_t)j,signal);
929 //______________________________________________________________________
930 void AliITSsimulationSDD::ResetHistograms(){
931 // Reset histograms for this detector
934 for (i=0;i<fNofMaps;i++ ) {
935 if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset();
938 //______________________________________________________________________
939 TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
940 // Fills a histogram from a give anode.
944 if(wing <=0 || wing > 2) {
945 Warning("GetAnode","Wrong wing number: %d",wing);
947 } // end if wing <=0 || wing >2
948 if(anode <=0 || anode > fNofMaps/2) {
949 Warning("GetAnode","Wrong anode number: %d",anode);
951 } // end if ampde <=0 || andoe > fNofMaps/2
953 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
954 return (TH1F*)(fHis->At(index));
956 //______________________________________________________________________
957 void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
958 // Writes the histograms to a file
964 for(i=0; i<fNofMaps; i++) fHis->At(i)->Write(); //fAdcs[i]->Write();
967 //______________________________________________________________________
968 void AliITSsimulationSDD::WriteSDigits(){
969 // Fills the Summable digits Tree
970 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
972 for( Int_t i=0; i<fNofMaps; i++ ) {
973 if( !fAnodeFire[i] ) continue;
974 for( Int_t j=0; j<fMaxNofSamples; j++ ) {
975 Double_t sig = fHitMap2->GetSignal( i, j );
977 Int_t jdx = j*fScaleSize;
978 Int_t index = fpList->GetHitIndex( i, j );
979 AliITSpListItem pItemTmp2( fModule, index, 0. );
980 // put the fScaleSize analog digits in only one
981 for( Int_t ik=0; ik<fScaleSize; ik++ ) {
982 AliITSpListItem *pItemTmp = fpList->GetpListItem(i,jdx+ik);
983 if( pItemTmp == 0 ) continue;
984 pItemTmp2.Add( pItemTmp );
986 pItemTmp2.AddSignalAfterElect( fModule, index, sig );
987 pItemTmp2.AddNoise(fModule,index,fHitNoiMap2->GetSignal(i,j));
988 aliITS->AddSumDigit( pItemTmp2 );
989 } // end if (sig > 0.2)
994 //______________________________________________________________________
995 void AliITSsimulationSDD::PrintStatus() const {
996 // Print SDD simulation Parameters
998 cout << "**************************************************" << endl;
999 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1000 cout << "**************************************************" << endl;
1001 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1002 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1003 cout << "Number of Anodes used: " << fNofMaps << endl;
1004 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1005 cout << "Scale size factor: " << fScaleSize << endl;
1006 cout << "**************************************************" << endl;