]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSsimulationSDD.cxx
New default value of noise for SDD simulations introduced.
[u/mrichter/AliRoot.git] / ITS / AliITSsimulationSDD.cxx
CommitLineData
b0f5e3fc 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
b0f5e3fc 16#include <iostream.h>
17#include <stdlib.h>
18#include <stdio.h>
1ca7869b 19#include <string.h>
20
94de3818 21#include <TSystem.h>
22#include <TROOT.h>
608f25d8 23#include <TStopwatch.h>
ece86d9a 24#include <TCanvas.h>
25#include <TF1.h>
26#include <TRandom.h>
1ca7869b 27#include <TH1.h>
28#include <TFile.h>
29#include <TVector.h>
30#include <TArrayI.h>
31#include <TArrayF.h>
ece86d9a 32
b0f5e3fc 33#include "AliRun.h"
e8189707 34#include "AliITS.h"
ece86d9a 35#include "AliITShit.h"
36#include "AliITSdigit.h"
37#include "AliITSmodule.h"
c7a4dac0 38#include "AliITSpList.h"
e8189707 39#include "AliITSMapA1.h"
40#include "AliITSMapA2.h"
e8189707 41#include "AliITSetfSDD.h"
42#include "AliITSRawData.h"
b0f5e3fc 43#include "AliITSHuffman.h"
50d05d7b 44#include "AliITSgeom.h"
ece86d9a 45#include "AliITSsegmentation.h"
46#include "AliITSresponse.h"
c7a4dac0 47#include "AliITSsegmentationSDD.h"
48#include "AliITSresponseSDD.h"
1ca7869b 49#include "AliITSsimulationSDD.h"
b0f5e3fc 50
b0f5e3fc 51ClassImp(AliITSsimulationSDD)
52////////////////////////////////////////////////////////////////////////
53// Version: 0
54// Written by Piergiorgio Cerello
55// November 23 1999
56//
57// AliITSsimulationSDD is the simulation of SDDs.
58 //
59//Begin_Html
60/*
61<img src="picts/ITS/AliITShit_Class_Diagram.gif">
62</pre>
63<br clear=left>
64<font size=+2 color=red>
65<p>This show the relasionships between the ITS hit class and the rest of Aliroot.
66</font>
67<pre>
68*/
69//End_Html
8a33ae9e 70//______________________________________________________________________
b0f5e3fc 71Int_t power(Int_t b, Int_t e) {
8a33ae9e 72 // compute b to the e power, where both b and e are Int_ts.
73 Int_t power = 1,i;
b0f5e3fc 74
8a33ae9e 75 for(i=0; i<e; i++) power *= b;
76 return power;
77}
78//______________________________________________________________________
b0f5e3fc 79void FastFourierTransform(AliITSetfSDD *alisddetf,Double_t *real,
80 Double_t *imag,Int_t direction) {
8a33ae9e 81 // Do a Fast Fourier Transform
8a33ae9e 82
83 Int_t samples = alisddetf->GetSamples();
84 Int_t l = (Int_t) ((log((Float_t) samples)/log(2.))+0.5);
85 Int_t m1 = samples;
86 Int_t m = samples/2;
87 Int_t m2 = samples/m1;
88 Int_t i,j,k;
89 for(i=1; i<=l; i++) {
50d05d7b 90 for(j=0; j<samples; j += m1) {
91 Int_t p = 0;
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);
100 *(real+k) += xr;
101 *(imag+k) += xi;
102 p += m2;
103 } // end for k
104 } // end for j
105 m1 = m;
106 m /= 2;
107 m2 += m2;
8a33ae9e 108 } // end for i
b0f5e3fc 109
8a33ae9e 110 for(j=0; j<samples; j++) {
50d05d7b 111 Int_t j1 = j;
112 Int_t p = 0;
113 Int_t i1;
114 for(i1=1; i1<=l; i1++) {
115 Int_t j2 = j1;
116 j1 /= 2;
117 p = p + p + j2 - j1 - j1;
118 } // end for i1
119 if(p >= j) {
120 Double_t xr = *(real+j);
121 Double_t xi = *(imag+j);
122 *(real+j) = *(real+p);
123 *(imag+j) = *(imag+p);
124 *(real+p) = xr;
125 *(imag+p) = xi;
126 } // end if p>=j
8a33ae9e 127 } // end for j
128 if(direction == -1) {
50d05d7b 129 for(i=0; i<samples; i++) {
130 *(real+i) /= samples;
131 *(imag+i) /= samples;
132 } // end for i
8a33ae9e 133 } // end if direction == -1
50d05d7b 134 return;
b0f5e3fc 135}
8a33ae9e 136//______________________________________________________________________
b0f5e3fc 137AliITSsimulationSDD::AliITSsimulationSDD(){
8a33ae9e 138 // Default constructor
139
140 fResponse = 0;
141 fSegmentation = 0;
142 fHis = 0;
50d05d7b 143// fpList = 0;
8a33ae9e 144 fHitMap2 = 0;
145 fElectronics = 0;
146 fStream = 0;
147 fInZR = 0;
148 fInZI = 0;
149 fOutZR = 0;
150 fOutZI = 0;
151 fNofMaps = 0;
152 fMaxNofSamples = 0;
153 fITS = 0;
8604b270 154 fTreeB = 0;
8a33ae9e 155 SetScaleFourier();
156 SetPerpendTracksFlag();
50d05d7b 157 SetCrosstalkFlag();
8a33ae9e 158 SetDoFFT();
159 SetCheckNoise();
b0f5e3fc 160}
8a33ae9e 161//______________________________________________________________________
162AliITSsimulationSDD::AliITSsimulationSDD(AliITSsimulationSDD &source){
163 // Copy constructor to satify Coding roules only.
164
165 if(this==&source) return;
c7a4dac0 166 Error("AliITSsimulationSSD","Not allowed to make a copy of "
50d05d7b 167 "AliITSsimulationSDD Using default creater instead");
8a33ae9e 168 AliITSsimulationSDD();
b0f5e3fc 169}
8a33ae9e 170//______________________________________________________________________
c7a4dac0 171AliITSsimulationSDD& AliITSsimulationSDD::operator=(AliITSsimulationSDD &src){
8a33ae9e 172 // Assignment operator to satify Coding roules only.
173
c7a4dac0 174 if(this==&src) return *this;
175 Error("AliITSsimulationSSD","Not allowed to make a = with "
50d05d7b 176 "AliITSsimulationSDD Using default creater instead");
8a33ae9e 177 return *this ;
b0f5e3fc 178}
8a33ae9e 179//______________________________________________________________________
c7a4dac0 180AliITSsimulationSDD::AliITSsimulationSDD(AliITSsegmentation *seg,
50d05d7b 181 AliITSresponse *resp){
c7a4dac0 182 // Standard Constructor
183
c7a4dac0 184 fResponse = 0;
185 fSegmentation = 0;
186 fHis = 0;
50d05d7b 187// fpList = 0;
c7a4dac0 188 fHitMap2 = 0;
189 fElectronics = 0;
190 fStream = 0;
191 fInZR = 0;
192 fInZI = 0;
193 fOutZR = 0;
194 fOutZI = 0;
195 fNofMaps = 0;
196 fMaxNofSamples = 0;
197 fITS = 0;
198 fTreeB = 0;
199
200 Init((AliITSsegmentationSDD*)seg,(AliITSresponseSDD*)resp);
201}
202//______________________________________________________________________
203void AliITSsimulationSDD::Init(AliITSsegmentationSDD *seg,
50d05d7b 204 AliITSresponseSDD *resp){
38867c90 205 // Standard Constructor
ece86d9a 206
38867c90 207 fResponse = resp;
208 fSegmentation = seg;
209 SetScaleFourier();
210 SetPerpendTracksFlag();
50d05d7b 211 SetCrosstalkFlag();
38867c90 212 SetDoFFT();
213 SetCheckNoise();
b0f5e3fc 214
50d05d7b 215 fpList = new AliITSpList( fSegmentation->Npz(),
216 fScaleSize*fSegmentation->Npx() );
38867c90 217 fHitMap2 = new AliITSMapA2(fSegmentation,fScaleSize,1);
b0f5e3fc 218
8a33ae9e 219 fNofMaps = fSegmentation->Npz();
220 fMaxNofSamples = fSegmentation->Npx();
221
38867c90 222 Float_t sddLength = fSegmentation->Dx();
8a33ae9e 223 Float_t sddWidth = fSegmentation->Dz();
b0f5e3fc 224
8a33ae9e 225 Int_t dummy = 0;
38867c90 226 Float_t anodePitch = fSegmentation->Dpz(dummy);
8a33ae9e 227 Double_t timeStep = (Double_t)fSegmentation->Dpx(dummy);
228 Float_t driftSpeed = fResponse->DriftSpeed();
b0f5e3fc 229
38867c90 230 if(anodePitch*(fNofMaps/2) > sddWidth) {
50d05d7b 231 Warning("AliITSsimulationSDD",
232 "Too many anodes %d or too big pitch %f \n",
233 fNofMaps/2,anodePitch);
38867c90 234 } // end if
b0f5e3fc 235
38867c90 236 if(timeStep*fMaxNofSamples < sddLength/driftSpeed) {
50d05d7b 237 Error("AliITSsimulationSDD",
238 "Time Interval > Allowed Time Interval: exit\n");
239 return;
38867c90 240 } // end if
b0f5e3fc 241
38867c90 242 fElectronics = new AliITSetfSDD(timeStep/fScaleSize,
50d05d7b 243 fResponse->Electronics());
b0f5e3fc 244
38867c90 245 char opt1[20], opt2[20];
246 fResponse->ParamOptions(opt1,opt2);
8a33ae9e 247 fParam = opt2;
38867c90 248 char *same = strstr(opt1,"same");
249 if (same) {
50d05d7b 250 fNoise.Set(0);
251 fBaseline.Set(0);
38867c90 252 } else {
50d05d7b 253 fNoise.Set(fNofMaps);
254 fBaseline.Set(fNofMaps);
38867c90 255 } // end if
b0f5e3fc 256
38867c90 257 const char *kopt=fResponse->ZeroSuppOption();
258 if (strstr(fParam,"file") ) {
50d05d7b 259 fD.Set(fNofMaps);
260 fT1.Set(fNofMaps);
261 if (strstr(kopt,"2D")) {
262 fT2.Set(fNofMaps);
b0f5e3fc 263 fTol.Set(0);
264 Init2D(); // desactivate if param change module by module
50d05d7b 265 } else if(strstr(kopt,"1D")) {
b0f5e3fc 266 fT2.Set(2);
267 fTol.Set(2);
268 Init1D(); // desactivate if param change module by module
50d05d7b 269 } // end if strstr
38867c90 270 } else {
50d05d7b 271 fD.Set(2);
272 fTol.Set(2);
273 fT1.Set(2);
274 fT2.Set(2);
275 SetCompressParam();
38867c90 276 } // end if else strstr
b0f5e3fc 277
8a33ae9e 278 Bool_t write = fResponse->OutputOption();
38867c90 279 if(write && strstr(kopt,"2D")) MakeTreeB();
b0f5e3fc 280
38867c90 281 // call here if baseline does not change by module
282 // ReadBaseline();
b0f5e3fc 283
8a33ae9e 284 fITS = (AliITS*)gAlice->GetModule("ITS");
285 Int_t size = fNofMaps*fMaxNofSamples;
286 fStream = new AliITSInStream(size);
b0f5e3fc 287
38867c90 288 fInZR = new Double_t [fScaleSize*fMaxNofSamples];
289 fInZI = new Double_t [fScaleSize*fMaxNofSamples];
290 fOutZR = new Double_t [fScaleSize*fMaxNofSamples];
291 fOutZI = new Double_t [fScaleSize*fMaxNofSamples];
5d18fa90 292
b0f5e3fc 293}
8a33ae9e 294//______________________________________________________________________
b0f5e3fc 295AliITSsimulationSDD::~AliITSsimulationSDD() {
8a33ae9e 296 // destructor
297
50d05d7b 298// delete fpList;
8a33ae9e 299 delete fHitMap2;
300 delete fStream;
301 delete fElectronics;
302
8a33ae9e 303 fITS = 0;
304
305 if (fHis) {
50d05d7b 306 fHis->Delete();
307 delete fHis;
8a33ae9e 308 } // end if fHis
309 if(fTreeB) delete fTreeB;
310 if(fInZR) delete [] fInZR;
50d05d7b 311 if(fInZI) delete [] fInZI;
8a33ae9e 312 if(fOutZR) delete [] fOutZR;
313 if(fOutZI) delete [] fOutZI;
b0f5e3fc 314}
8a33ae9e 315//______________________________________________________________________
50d05d7b 316void AliITSsimulationSDD::InitSimulationModule( Int_t module, Int_t event ) {
c7a4dac0 317 // create maps to build the lists of tracks for each summable digit
50d05d7b 318 fModule = module;
319 fEvent = event;
320 ClearMaps();
321}
322//______________________________________________________________________
323void AliITSsimulationSDD::ClearMaps() {
324 // clear maps
325 fpList->ClearMap();
326 fHitMap2->ClearMap();
327}
328//______________________________________________________________________
329void AliITSsimulationSDD::SDigitiseModule( AliITSmodule *mod, Int_t md, Int_t ev){
330 // digitize module using the "slow" detector simulator creating
331 // summable digits.
c7a4dac0 332
333 TObjArray *fHits = mod->GetHits();
334 Int_t nhits = fHits->GetEntriesFast();
50d05d7b 335 if( !nhits ) return;
c7a4dac0 336
50d05d7b 337 InitSimulationModule( md, ev );
338 HitsToAnalogDigits( mod );
339 WriteSDigits();
340 ClearMaps();
341}
342//______________________________________________________________________
343/*void AliITSsimulationSDD::AddSDigitsToModule( TClonesArray *pItemArray, Int_t mask ) {
344 // Add Summable digits to module maps.
345 Int_t nItems = pItemArray->GetEntries();
346
347 // cout << "Adding "<< nItems <<" SDigits to module " << fModule << endl;
348 for( Int_t i=0; i<nItems; i++ ) {
349 AliITSpListItem * pItem = (AliITSpListItem *)(pItemArray->At( i ));
350 if( pItem->GetModule() != fModule ) {
351 Error( "AliITSsimulationSDD",
352 "Error reading, SDigits module %d != current module %d: exit\n",
353 pItem->GetModule(), fModule );
354 return;
355 } // end if
356
357 fpList->AddItemTo( mask, pItem );
358 }
359}*/
360//______________________________________________________________________
361void AliITSsimulationSDD::FinishSDigitiseModule() {
362 // digitize module using the "slow" detector simulator from
363 // the sum of summable digits.
364 FinishDigits() ;
365 ClearMaps();
c7a4dac0 366}
367//______________________________________________________________________
b0f5e3fc 368void AliITSsimulationSDD::DigitiseModule(AliITSmodule *mod,Int_t md,Int_t ev){
8a33ae9e 369 // create maps to build the lists of tracks for each digit
b0f5e3fc 370
371 TObjArray *fHits = mod->GetHits();
8a33ae9e 372 Int_t nhits = fHits->GetEntriesFast();
8a33ae9e 373
50d05d7b 374 InitSimulationModule( md, ev );
375
376 if( !nhits && fCheckNoise ) {
ece86d9a 377 ChargeToSignal();
378 GetNoise();
50d05d7b 379 ClearMaps();
ece86d9a 380 return;
50d05d7b 381 } else
382 if( !nhits ) return;
383 HitsToAnalogDigits( mod );
384 FinishDigits();
385 ClearMaps();
c7a4dac0 386}
387//______________________________________________________________________
50d05d7b 388void AliITSsimulationSDD::FinishDigits() {
c7a4dac0 389 // introduce the electronics effects and do zero-suppression if required
8a33ae9e 390
50d05d7b 391 // Fill maps from fpList.
c5041e31 392 Int_t maxIndex = fpList->GetEntries();
393 for(Int_t i=0;i<maxIndex;i++) fHitMap2->SetHit(i,fpList->GetSignal(i));
50d05d7b 394
c7a4dac0 395 ChargeToSignal();
50d05d7b 396 ApplyDeadChannels();
397 if( fCrosstalkFlag ) ApplyCrosstalk();
398
399 const char *kopt = fResponse->ZeroSuppOption();
400 ZeroSuppression( kopt );
c7a4dac0 401}
402//______________________________________________________________________
50d05d7b 403void AliITSsimulationSDD::HitsToAnalogDigits( AliITSmodule *mod ) {
c7a4dac0 404 // create maps to build the lists of tracks for each digit
405
406 TObjArray *fHits = mod->GetHits();
407 Int_t nhits = fHits->GetEntriesFast();
50d05d7b 408// Int_t arg[6] = {0,0,0,0,0,0};
8a33ae9e 409 Int_t dummy = 0;
410 Int_t nofAnodes = fNofMaps/2;
411 Float_t sddLength = fSegmentation->Dx();
412 Float_t sddWidth = fSegmentation->Dz();
413 Float_t anodePitch = fSegmentation->Dpz(dummy);
414 Float_t timeStep = fSegmentation->Dpx(dummy);
415 Float_t driftSpeed = fResponse->DriftSpeed();
416 Float_t maxadc = fResponse->MaxAdc();
417 Float_t topValue = fResponse->DynamicRange();
418 Float_t cHloss = fResponse->ChargeLoss();
419 Float_t norm = maxadc/topValue;
420 Float_t dfCoeff, s1; fResponse->DiffCoeff(dfCoeff,s1); // Signal 2d Shape
421 Double_t eVpairs = 3.6; // electron pair energy eV.
422 Float_t nsigma = fResponse->NSigmaIntegration(); //
423 Int_t nlookups = fResponse->GausNLookUp(); //
44a312c3 424
b0f5e3fc 425 // Piergiorgio's part (apart for few variables which I made float
426 // when i thought that can be done
b0f5e3fc 427 // Fill detector maps with GEANT hits
428 // loop over hits in the module
429
8a33ae9e 430 const Float_t kconv = 1.0e+6; // GeV->KeV
431 Int_t itrack = 0;
432 Int_t hitDetector; // detector number (lay,lad,hitDetector)
433 Int_t iWing; // which detector wing/side.
434 Int_t detector; // 2*(detector-1)+iWing
435 Int_t ii,kk,ka,kt; // loop indexs
436 Int_t ia,it,index; // sub-pixel integration indexies
437 Int_t iAnode; // anode number.
438 Int_t timeSample; // time buckett.
439 Int_t anodeWindow; // anode direction charge integration width
440 Int_t timeWindow; // time direction charge integration width
441 Int_t jamin,jamax; // anode charge integration window
442 Int_t jtmin,jtmax; // time charge integration window
443 Int_t ndiv; // Anode window division factor.
444 Int_t nsplit; // the number of splits in anode and time windows==1.
445 Int_t nOfSplits; // number of times track length is split into
446 Float_t nOfSplitsF; // Floating point version of nOfSplits.
447 Float_t kkF; // Floating point version of loop index kk.
448 Float_t pathInSDD; // Track length in SDD.
449 Float_t drPath; // average position of track in detector. in microns
450 Float_t drTime; // Drift time
451 Float_t nmul; // drift time window multiplication factor.
452 Float_t avDrft; // x position of path length segment in cm.
453 Float_t avAnode; // Anode for path length segment in Anode number (float)
454 Float_t xAnode; // Floating point anode number.
455 Float_t driftPath; // avDrft in microns.
456 Float_t width; // width of signal at anodes.
457 Double_t depEnergy; // Energy deposited in this GEANT step.
458 Double_t xL[3],dxL[3]; // local hit coordinates and diff.
459 Double_t sigA; // sigma of signal at anode.
460 Double_t sigT; // sigma in time/drift direction for track segment
461 Double_t aStep,aConst; // sub-pixel size and offset anode
462 Double_t tStep,tConst; // sub-pixel size and offset time
463 Double_t amplitude; // signal amplitude for track segment in nanoAmpere
464 Double_t chargeloss; // charge loss for track segment.
465 Double_t anodeAmplitude; // signal amplitude in anode direction
466 Double_t aExpo; // exponent of Gaussian anode direction
467 Double_t timeAmplitude; // signal amplitude in time direction
468 Double_t tExpo; // exponent of Gaussian time direction
469// Double_t tof; // Time of flight in ns of this step.
44a312c3 470
b0f5e3fc 471 for(ii=0; ii<nhits; ii++) {
50d05d7b 472 if(!mod->LineSegmentL(ii,xL[0],dxL[0],xL[1],dxL[1],xL[2],dxL[2],
473 depEnergy,itrack)) continue;
474 depEnergy *= kconv;
475 hitDetector = mod->GetDet();
476 //tof = 1.E+09*(mod->GetHit(ii)->GetTOF()); // tof in ns.
477 //if(tof>sddLength/driftSpeed) continue; // hit happed too late.
478
479 // scale path to simulate a perpendicular track
480 // continue if the particle did not lose energy
481 // passing through detector
482 if (!depEnergy) {
483 Warning("HitsToAnalogDigits",
484 "fTrack = %d hit=%d module=%d This particle has"
485 " passed without losing energy!",
486 itrack,ii,mod->GetIndex());
487 continue;
488 } // end if !depEnergy
489
490 pathInSDD = TMath::Sqrt(dxL[0]*dxL[0]+dxL[1]*dxL[1]+dxL[2]*dxL[2]);
491
492 if (fFlag && pathInSDD) { depEnergy *= (0.03/pathInSDD); }
493 drPath = 10000.*(dxL[0]+2.*xL[0])*0.5;
494 if(drPath < 0) drPath = -drPath;
495 drPath = sddLength-drPath;
496 if(drPath < 0) {
497 Warning("HitsToAnalogDigits",
498 "negative drift path drPath=%e sddLength=%e dxL[0]=%e "
499 "xL[0]=%e",
500 drPath,sddLength,dxL[0],xL[0]);
501 continue;
502 } // end if drPath < 0
503
504 // Compute number of segments to brake step path into
505 drTime = drPath/driftSpeed; // Drift Time
506 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);// Sigma along the anodes
507 // calcuate the number of time the path length should be split into.
508 nOfSplits = (Int_t) (1. + 10000.*pathInSDD/sigA);
509 if(fFlag) nOfSplits = 1;
510
511 // loop over path segments, init. some variables.
512 depEnergy /= nOfSplits;
513 nOfSplitsF = (Float_t) nOfSplits;
514 for(kk=0;kk<nOfSplits;kk++) { // loop over path segments
515 kkF = (Float_t) kk + 0.5;
516 avDrft = xL[0]+dxL[0]*kkF/nOfSplitsF;
517 avAnode = xL[2]+dxL[2]*kkF/nOfSplitsF;
518 driftPath = 10000.*avDrft;
519
520 iWing = 2; // Assume wing is 2
521 if(driftPath < 0) { // if wing is not 2 it is 1.
522 iWing = 1;
523 driftPath = -driftPath;
524 } // end if driftPath < 0
525 driftPath = sddLength-driftPath;
526 detector = 2*(hitDetector-1) + iWing;
527 if(driftPath < 0) {
528 Warning("HitsToAnalogDigits","negative drift path "
529 "driftPath=%e sddLength=%e avDrft=%e dxL[0]=%e "
530 "xL[0]=%e",driftPath,sddLength,avDrft,dxL[0],xL[0]);
531 continue;
532 } // end if driftPath < 0
533
534 // Drift Time
535 drTime = driftPath/driftSpeed; // drift time for segment.
536 timeSample = (Int_t) (fScaleSize*drTime/timeStep + 1);
537 // compute time Sample including tof information. The tof only
538 // effects the time of the signal is recoreded and not the
539 // the defusion.
540 // timeSample = (Int_t) (fScaleSize*(drTime+tof)/timeStep + 1);
541 if(timeSample > fScaleSize*fMaxNofSamples) {
542 Warning("HitsToAnalogDigits","Wrong Time Sample: %e",
543 timeSample);
544 continue;
545 } // end if timeSample > fScaleSize*fMaxNoofSamples
546
547 // Anode
548 xAnode = 10000.*(avAnode)/anodePitch + nofAnodes/2; // +1?
549 if(xAnode*anodePitch > sddWidth || xAnode*anodePitch < 0.)
550 Warning("HitsToAnalogDigits",
551 "Exceedubg sddWidth=%e Z = %e",
552 sddWidth,xAnode*anodePitch);
553 iAnode = (Int_t) (1.+xAnode); // xAnode?
554 if(iAnode < 1 || iAnode > nofAnodes) {
555 Warning("HitToAnalogDigits","Wrong iAnode: 1<%d>%d",
556 iAnode,nofAnodes);
557 continue;
558 } // end if iAnode < 1 || iAnode > nofAnodes
559
560 // store straight away the particle position in the array
561 // of particles and take idhit=ii only when part is entering (this
562 // requires FillModules() in the macro for analysis) :
b0f5e3fc 563
50d05d7b 564 // Sigma along the anodes for track segment.
565 sigA = TMath::Sqrt(2.*dfCoeff*drTime+s1*s1);
566 sigT = sigA/driftSpeed;
567 // Peak amplitude in nanoAmpere
568 amplitude = fScaleSize*160.*depEnergy/
569 (timeStep*eVpairs*2.*acos(-1.)*sigT*sigA);
570 amplitude *= timeStep/25.; // WARNING!!!!! Amplitude scaling to
571 // account for clock variations
8a33ae9e 572 // (reference value: 40 MHz)
50d05d7b 573 chargeloss = 1.-cHloss*driftPath/1000;
574 amplitude *= chargeloss;
575 width = 2.*nsigma/(nlookups-1);
576 // Spread the charge
577 // Pixel index
578 ndiv = 2;
579 nmul = 3.;
580 if(drTime > 1200.) {
581 ndiv = 4;
582 nmul = 1.5;
583 } // end if drTime > 1200.
584 // Sub-pixel index
585 nsplit = 4; // hard-wired //nsplit=4;nsplit = (nsplit+1)/2*2;
586 // Sub-pixel size see computation of aExpo and tExpo.
587 aStep = anodePitch/(nsplit*fScaleSize*sigA);
588 aConst = xAnode*anodePitch/sigA;
589 tStep = timeStep/(nsplit*fScaleSize*sigT);
590 tConst = drTime/sigT;
591 // Define SDD window corresponding to the hit
592 anodeWindow = (Int_t)(fScaleSize*nsigma*sigA/anodePitch+1);
593 timeWindow = (Int_t) (fScaleSize*nsigma*sigT/timeStep+1.);
594 jamin = (iAnode - anodeWindow/ndiv - 1)*fScaleSize*nsplit +1;
595 jamax = (iAnode + anodeWindow/ndiv)*fScaleSize*nsplit;
596 if(jamin <= 0) jamin = 1;
597 if(jamax > fScaleSize*nofAnodes*nsplit)
598 jamax = fScaleSize*nofAnodes*nsplit;
599 // jtmin and jtmax are Hard-wired
600 jtmin = (Int_t)(timeSample-timeWindow*nmul-1)*nsplit+1;
601 jtmax = (Int_t)(timeSample+timeWindow*nmul)*nsplit;
602 if(jtmin <= 0) jtmin = 1;
603 if(jtmax > fScaleSize*fMaxNofSamples*nsplit)
604 jtmax = fScaleSize*fMaxNofSamples*nsplit;
605 // Spread the charge in the anode-time window
606 for(ka=jamin; ka <=jamax; ka++) {
607 ia = (ka-1)/(fScaleSize*nsplit) + 1;
608 if(ia <= 0) {
609 Warning("HitsToAnalogDigits","ia < 1: ");
610 continue;
611 } // end if
612 if(ia > nofAnodes) ia = nofAnodes;
613 aExpo = (aStep*(ka-0.5)-aConst);
614 if(TMath::Abs(aExpo) > nsigma) anodeAmplitude = 0.;
615 else {
616 dummy = (Int_t) ((aExpo+nsigma)/width);
617 anodeAmplitude = amplitude*fResponse->GausLookUp(dummy);
618 } // end if TMath::Abs(aEspo) > nsigma
619 // index starts from 0
620 index = ((detector+1)%2)*nofAnodes+ia-1;
621 if(anodeAmplitude) for(kt=jtmin; kt<=jtmax; kt++) {
622 it = (kt-1)/nsplit+1; // it starts from 1
623 if(it<=0){
624 Warning("HitsToAnalogDigits","it < 1:");
625 continue;
626 } // end if
627 if(it>fScaleSize*fMaxNofSamples)
628 it = fScaleSize*fMaxNofSamples;
629 tExpo = (tStep*(kt-0.5)-tConst);
630 if(TMath::Abs(tExpo) > nsigma) timeAmplitude = 0.;
631 else {
632 dummy = (Int_t) ((tExpo+nsigma)/width);
633 timeAmplitude = anodeAmplitude*
634 fResponse->GausLookUp(dummy);
635 } // end if TMath::Abs(tExpo) > nsigma
636 // build the list of Sdigits for this module
637// arg[0] = index;
638// arg[1] = it;
639// arg[2] = itrack; // track number
640// arg[3] = ii-1; // hit number.
641 timeAmplitude *= norm;
642 timeAmplitude *= 10;
643// ListOfFiredCells(arg,timeAmplitude,alst,padr);
644 fpList->AddSignal(index,it-1,itrack,ii-1,
645 mod->GetIndex(),timeAmplitude);
646 } // end if anodeAmplitude and loop over time in window
647 } // loop over anodes in window
648 } // end loop over "sub-hits"
44a312c3 649 } // end loop over hits
b0f5e3fc 650}
50d05d7b 651
652/*
8a33ae9e 653//______________________________________________________________________
b0f5e3fc 654void AliITSsimulationSDD::ListOfFiredCells(Int_t *arg,Double_t timeAmplitude,
8a33ae9e 655 TObjArray *alist,TClonesArray *padr){
656 // Returns the list of "fired" cells.
657
658 Int_t index = arg[0];
659 Int_t ik = arg[1];
660 Int_t idtrack = arg[2];
661 Int_t idhit = arg[3];
662 Int_t counter = arg[4];
663 Int_t countadr = arg[5];
664 Double_t charge = timeAmplitude;
665 charge += fHitMap2->GetSignal(index,ik-1);
666 fHitMap2->SetHit(index, ik-1, charge);
667
668 Int_t digits[3];
669 Int_t it = (Int_t)((ik-1)/fScaleSize);
670 digits[0] = index;
671 digits[1] = it;
672 digits[2] = (Int_t)timeAmplitude;
673 Float_t phys;
674 if (idtrack >= 0) phys = (Float_t)timeAmplitude;
675 else phys = 0;
676
677 Double_t cellcharge = 0.;
678 AliITSTransientDigit* pdigit;
679 // build the list of fired cells and update the info
680 if (!fHitMap1->TestHit(index, it)) {
50d05d7b 681 new((*padr)[countadr++]) TVector(3);
682 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
683 trinfo(0) = (Float_t)idtrack;
684 trinfo(1) = (Float_t)idhit;
685 trinfo(2) = (Float_t)timeAmplitude;
686
687 alist->AddAtAndExpand(new AliITSTransientDigit(phys,digits),counter);
688 fHitMap1->SetHit(index, it, counter);
689 counter++;
690 pdigit=(AliITSTransientDigit*)alist->At(alist->GetLast());
691 // list of tracks
692 TObjArray *trlist=(TObjArray*)pdigit->TrackList();
693 trlist->Add(&trinfo);
8a33ae9e 694 } else {
50d05d7b 695 pdigit = (AliITSTransientDigit*) fHitMap1->GetHit(index, it);
696 for(Int_t kk=0;kk<fScaleSize;kk++) {
697 cellcharge += fHitMap2->GetSignal(index,fScaleSize*it+kk);
698 } // end for kk
699 // update charge
700 (*pdigit).fSignal = (Int_t)cellcharge;
701 (*pdigit).fPhysics += phys;
702 // update list of tracks
703 TObjArray* trlist = (TObjArray*)pdigit->TrackList();
704 Int_t lastentry = trlist->GetLast();
705 TVector *ptrkp = (TVector*)trlist->At(lastentry);
706 TVector &trinfo = *ptrkp;
707 Int_t lasttrack = Int_t(trinfo(0));
708 Float_t lastcharge=(trinfo(2));
709 if (lasttrack==idtrack ) {
710 lastcharge += (Float_t)timeAmplitude;
711 trlist->RemoveAt(lastentry);
712 trinfo(0) = lasttrack;
713 trinfo(1) = idhit;
714 trinfo(2) = lastcharge;
715 trlist->AddAt(&trinfo,lastentry);
716 } else {
717 new((*padr)[countadr++]) TVector(3);
718 TVector &trinfo=*((TVector*) (*padr)[countadr-1]);
719 trinfo(0) = (Float_t)idtrack;
720 trinfo(1) = (Float_t)idhit;
721 trinfo(2) = (Float_t)timeAmplitude;
722 trlist->Add(&trinfo);
723 } // end if lasttrack==idtrack
b0f5e3fc 724
725#ifdef print
50d05d7b 726 // check the track list - debugging
727 Int_t trk[20], htrk[20];
728 Float_t chtrk[20];
729 Int_t nptracks = trlist->GetEntriesFast();
730 if (nptracks > 2) {
731 Int_t tr;
732 for (tr=0;tr<nptracks;tr++) {
733 TVector *pptrkp = (TVector*)trlist->At(tr);
734 TVector &pptrk = *pptrkp;
735 trk[tr] = Int_t(pptrk(0));
736 htrk[tr] = Int_t(pptrk(1));
737 chtrk[tr] = (pptrk(2));
738 cout << "nptracks "<<nptracks << endl;
739 } // end for tr
740 } // end if nptracks
b0f5e3fc 741#endif
8a33ae9e 742 } // end if pdigit
b0f5e3fc 743
8a33ae9e 744 // update counter and countadr for next call.
745 arg[4] = counter;
746 arg[5] = countadr;
b0f5e3fc 747}
50d05d7b 748*/
749
b0f5e3fc 750//____________________________________________
50d05d7b 751void AliITSsimulationSDD::AddDigit( Int_t i, Int_t j, Int_t signal ) {
752 // Adds a Digit.
753 Int_t digits[3], tracks[3], hits[3];
754 Float_t phys, charges[3];
755
756 if( fResponse->Do10to8() ) signal = Convert8to10( signal );
757 digits[0] = i;
758 digits[1] = j;
759 digits[2] = signal;
760
761 Int_t jdx = j*fScaleSize;
762 AliITSpListItem pItem;
763
764 // put the fScaleSize analog digits in only one
765 for( Int_t ik=0; ik<fScaleSize; ik++ ) {
766 AliITSpListItem* pItemTmp = fpList->GetpListItem( i, jdx+ik );
767 if( pItemTmp == 0 ) continue;
768 Double_t sig = 0.0;
769 for( Int_t l=0; l<pItemTmp->GetNsignals(); l++ ) {
770 Double_t signalT = pItemTmp->GetSignal( l );
771 // if( signalT <= 0.0 ) break; // no more signals
772 sig += signalT;
773 Int_t track = pItemTmp->GetTrack( l );
774 Int_t hit = pItemTmp->GetHit( l );
775 pItem.AddSignal( track, hit, -1, -1, signalT );
776 }
777 // add to noise : total signal - sum of signal tracks
778 pItem.AddNoise( -1, -1, pItemTmp->GetSignal() - sig );
779 }
780
781 Int_t idtrack = pItem.GetTrack( 0 );
782 if( idtrack >= 0 ) phys = pItem.GetSumSignal();
783 else phys = 0;
784
785 for( Int_t l=0; l<3; l++ ) {
786 tracks[l] = pItem.GetTrack( l );
787 hits[l] = pItem.GetHit( l );
788 charges[l] = pItem.GetSignal( l );
789 }
b0f5e3fc 790
50d05d7b 791 fITS->AddSimDigit( 1, phys, digits, tracks, hits, charges );
792}
793
794/*
795//____________________________________________
b0f5e3fc 796void AliITSsimulationSDD::AddDigit(Int_t i, Int_t j, Int_t signal){
8a33ae9e 797 // Adds a Digit.
b0f5e3fc 798 // tag with -1 signals coming from background tracks
799 // tag with -2 signals coming from pure electronic noise
800
e8189707 801 Int_t digits[3], tracks[3], hits[3];
b0f5e3fc 802 Float_t phys, charges[3];
803
e8189707 804 Int_t trk[20], htrk[20];
b0f5e3fc 805 Float_t chtrk[20];
806
ece86d9a 807 Bool_t do10to8=fResponse->Do10to8();
808
809 if(do10to8) signal=Convert8to10(signal);
b0f5e3fc 810 AliITSTransientDigit *obj = (AliITSTransientDigit*)fHitMap1->GetHit(i,j);
8a33ae9e 811 digits[0] = i;
812 digits[1] = j;
813 digits[2] = signal;
b0f5e3fc 814 if (!obj) {
815 phys=0;
816 Int_t k;
e8189707 817 for (k=0;k<3;k++) {
50d05d7b 818 tracks[k]=-2;
819 charges[k]=0;
820 hits[k]=-1;
821 } // end for k
e8189707 822 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
b0f5e3fc 823 } else {
50d05d7b 824 phys=obj->fPhysics;
825 TObjArray* trlist=(TObjArray*)obj->TrackList();
826 Int_t nptracks=trlist->GetEntriesFast();
827 if (nptracks > 20) {
828 Warning("AddDigit","nptracks=%d > 20 nptracks set to 20",nptracks);
829 nptracks=20;
830 } // end if nptracks > 20
831 Int_t tr;
832 for (tr=0;tr<nptracks;tr++) {
833 TVector &pp =*((TVector*)trlist->At(tr));
834 trk[tr]=Int_t(pp(0));
835 htrk[tr]=Int_t(pp(1));
836 chtrk[tr]=(pp(2));
837 } // end for tr
838 if (nptracks > 1) {
839 SortTracks(trk,chtrk,htrk,nptracks);
840 } // end if nptracks > 1
841 Int_t i;
842 if (nptracks < 3 ) {
843 for (i=0; i<nptracks; i++) {
844 tracks[i]=trk[i];
845 charges[i]=chtrk[i];
846 hits[i]=htrk[i];
847 } // end for i
848 for (i=nptracks; i<3; i++) {
849 tracks[i]=-3;
850 hits[i]=-1;
851 charges[i]=0;
852 } // end for i
853 } else {
854 for (i=0; i<3; i++) {
855 tracks[i]=trk[i];
856 charges[i]=chtrk[i];
857 hits[i]=htrk[i];
858 } // end for i
859 } // end if/else nptracks < 3
860
861 fITS->AddSimDigit(1,phys,digits,tracks,hits,charges);
b0f5e3fc 862
8a33ae9e 863 } // end if/else !obj
b0f5e3fc 864}
50d05d7b 865
866
8a33ae9e 867//______________________________________________________________________
868void AliITSsimulationSDD::SortTracks(Int_t *tracks,Float_t *charges,
50d05d7b 869 Int_t *hits,Int_t ntr){
8a33ae9e 870 // Sort the list of tracks contributing to a given digit
871 // Only the 3 most significant tracks are acctually sorted
872 // Loop over signals, only 3 times
873
874 Float_t qmax;
875 Int_t jmax;
876 Int_t idx[3] = {-3,-3,-3};
877 Float_t jch[3] = {-3,-3,-3};
878 Int_t jtr[3] = {-3,-3,-3};
879 Int_t jhit[3] = {-3,-3,-3};
880 Int_t i,j,imax;
881
882 if (ntr<3) imax = ntr;
883 else imax = 3;
884 for(i=0;i<imax;i++){
50d05d7b 885 qmax = 0;
886 jmax = 0;
887 for(j=0;j<ntr;j++){
888 if((i == 1 && j == idx[i-1] )
889 ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
890 if(charges[j] > qmax) {
891 qmax = charges[j];
892 jmax=j;
893 } // end if charges[j]>qmax
894 } // end for j
895 if(qmax > 0) {
896 idx[i] = jmax;
897 jch[i] = charges[jmax];
898 jtr[i] = tracks[jmax];
899 jhit[i] = hits[jmax];
900 } // end if qmax > 0
8a33ae9e 901 } // end for i
902
903 for(i=0;i<3;i++){
50d05d7b 904 if (jtr[i] == -3) {
905 charges[i] = 0;
906 tracks[i] = -3;
907 hits[i] = -1;
908 } else {
909 charges[i] = jch[i];
910 tracks[i] = jtr[i];
911 hits[i] = jhit[i];
912 } // end if jtr[i] == -3
8a33ae9e 913 } // end for i
b0f5e3fc 914}
50d05d7b 915*/
8a33ae9e 916//______________________________________________________________________
b0f5e3fc 917void AliITSsimulationSDD::ChargeToSignal() {
8a33ae9e 918 // add baseline, noise, electronics and ADC saturation effects
b0f5e3fc 919
8a33ae9e 920 char opt1[20], opt2[20];
921 fResponse->ParamOptions(opt1,opt2);
922 char *read = strstr(opt1,"file");
923 Float_t baseline, noise;
924
925 if (read) {
50d05d7b 926 static Bool_t readfile=kTRUE;
927 //read baseline and noise from file
928 if (readfile) ReadBaseline();
929 readfile=kFALSE;
8a33ae9e 930 } else fResponse->GetNoiseParam(noise,baseline);
931
932 Float_t contrib=0;
933 Int_t i,k,kk;
934 Float_t maxadc = fResponse->MaxAdc();
935 if(!fDoFFT) {
50d05d7b 936 for (i=0;i<fNofMaps;i++) {
937 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
938 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
939 fInZR[k] = fHitMap2->GetSignal(i,k);
940 contrib = (baseline + noise*gRandom->Gaus());
941 fInZR[k] += contrib;
942 } // end for k
943 for(k=0; k<fMaxNofSamples; k++) {
944 Double_t newcont = 0.;
945 Double_t maxcont = 0.;
946 for(kk=0;kk<fScaleSize;kk++) {
947 newcont = fInZR[fScaleSize*k+kk];
948 if(newcont > maxcont) maxcont = newcont;
949 } // end for kk
950 newcont = maxcont;
951 if (newcont >= maxadc) newcont = maxadc -1;
952 if(newcont >= baseline){
953 Warning("","newcont=%d>=baseline=%d",newcont,baseline);
954 } // end if
955 // back to analog: ?
956 fHitMap2->SetHit(i,k,newcont);
957 } // end for k
958 } // end for i loop over anodes
959 return;
8a33ae9e 960 } // end if DoFFT
ece86d9a 961
ece86d9a 962 for (i=0;i<fNofMaps;i++) {
50d05d7b 963 if (read && i<fNofMaps) GetAnodeBaseline(i,baseline,noise);
964 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
965 fInZR[k] = fHitMap2->GetSignal(i,k);
966 contrib = (baseline + noise*gRandom->Gaus());
967 fInZR[k] += contrib;
968 fInZI[k] = 0.;
969 } // end for k
970 FastFourierTransform(fElectronics,&fInZR[0],&fInZI[0],1);
971 for(k=0; k<fScaleSize*fMaxNofSamples; k++) {
972 Double_t rw = fElectronics->GetTraFunReal(k);
973 Double_t iw = fElectronics->GetTraFunImag(k);
974 fOutZR[k] = fInZR[k]*rw - fInZI[k]*iw;
975 fOutZI[k] = fInZR[k]*iw + fInZI[k]*rw;
976 } // end for k
977 FastFourierTransform(fElectronics,&fOutZR[0],&fOutZI[0],-1);
978 for(k=0; k<fMaxNofSamples; k++) {
979 Double_t newcont1 = 0.;
980 Double_t maxcont1 = 0.;
981 for(kk=0;kk<fScaleSize;kk++) {
982 newcont1 = fOutZR[fScaleSize*k+kk];
983 if(newcont1 > maxcont1) maxcont1 = newcont1;
984 } // end for kk
985 newcont1 = maxcont1;
986 if (newcont1 >= maxadc) newcont1 = maxadc -1;
987 fHitMap2->SetHit(i,k,newcont1);
988 } // end for k
8a33ae9e 989 } // end for i loop over anodes
ece86d9a 990 return;
b0f5e3fc 991}
50d05d7b 992//____________________________________________________________________
993void AliITSsimulationSDD::ApplyDeadChannels() {
994 // Set dead channel signal to zero
995 AliITSresponseSDD * response = (AliITSresponseSDD *)fResponse;
996
997 // nothing to do
998 if( response->GetDeadModules() == 0 &&
999 response->GetDeadChips() == 0 &&
1000 response->GetDeadChannels() == 0 )
1001 return;
1002
1003 static AliITS *iTS = (AliITS*)gAlice->GetModule( "ITS" );
1004
1005 Int_t fMaxNofSamples = fSegmentation->Npx();
1006 AliITSgeom *geom = iTS->GetITSgeom();
1007 Int_t firstSDDMod = geom->GetStartDet( 1 );
1008 // loop over wings
1009 for( Int_t j=0; j<2; j++ ) {
1010 Int_t mod = (fModule-firstSDDMod)*2 + j;
1011 for( Int_t u=0; u<response->Chips(); u++ )
1012 for( Int_t v=0; v<response->Channels(); v++ ) {
1013 Float_t Gain = response->Gain( mod, u, v );
1014 for( Int_t k=0; k<fMaxNofSamples; k++ ) {
1015 Int_t i = j*response->Chips()*response->Channels() +
1016 u*response->Channels() +
1017 v;
1018 Double_t signal = Gain * fHitMap2->GetSignal( i, k );
1019 fHitMap2->SetHit( i, k, signal ); ///
1020 }
1021 }
1022 }
1023}
1024//______________________________________________________________________
1025void AliITSsimulationSDD::ApplyCrosstalk() {
1026 // function add the crosstalk effect to signal
1027 // temporal function, should be checked...!!!
1028
1029 Int_t fNofMaps = fSegmentation->Npz();
1030 Int_t fMaxNofSamples = fSegmentation->Npx();
1031
1032 // create and inizialice crosstalk map
1033 Float_t* ctk = new Float_t[fNofMaps*fMaxNofSamples+1];
1034 if( ctk == NULL ) {
1035 Error( "ApplyCrosstalk", "no memory for temporal map: exit \n" );
1036 return;
1037 }
1038 memset( ctk, 0, sizeof(Float_t)*(fNofMaps*fMaxNofSamples+1) );
1039
1040 Float_t noise, baseline;
1041 fResponse->GetNoiseParam( noise, baseline );
1042
1043 for( Int_t z=0; z<fNofMaps; z++ ) {
1044 Bool_t on = kFALSE;
1045 Int_t tstart = 0;
1046 Int_t tstop = 0;
1047 Int_t nTsteps = 0;
1048
1049 for( Int_t l=0; l<fMaxNofSamples; l++ ) {
1050 Float_t fadc = (Float_t)fHitMap2->GetSignal( z, l );
1051 if( fadc > baseline ) {
1052 if( on == kFALSE && l<fMaxNofSamples-4 ) {
1053 Float_t fadc1 = (Float_t)fHitMap2->GetSignal( z, l+1 );
1054 if( fadc1 < fadc ) continue;
1055 on = kTRUE;
1056 nTsteps = 0;
1057 tstart = l;
1058 }
1059 nTsteps++;
1060 }
1061 else { // end fadc > baseline
1062 if( on == kTRUE ) {
1063 if( nTsteps > 2 ) {
1064 tstop = l;
1065 // make smooth derivative
1066 Float_t* dev = new Float_t[fMaxNofSamples+1];
1067 memset( dev, 0, sizeof(Float_t)*(fMaxNofSamples+1) );
1068 if( ctk == NULL ) {
1069 Error( "ApplyCrosstalk",
1070 "no memory for temporal array: exit \n" );
1071 return;
1072 }
1073 for( Int_t i=tstart; i<tstop; i++ ) {
1074 if( i > 2 && i < fMaxNofSamples-2 )
1075 dev[i] = -0.2*fHitMap2->GetSignal( z,i-2 )
1076 -0.1*fHitMap2->GetSignal( z,i-1 )
1077 +0.1*fHitMap2->GetSignal( z,i+1 )
1078 +0.2*fHitMap2->GetSignal( z,i+2 );
1079 }
1080
1081 // add crosstalk contribution to neibourg anodes
1082 for( Int_t i=tstart; i<tstop; i++ ) {
1083 Int_t anode = z - 1;
1084 Int_t i1 = (Int_t)((i-tstart)*.61+tstart+0.5); //
1085 Float_t ctktmp = -dev[i1] * 0.25;
1086 if( anode > 0 ) {
1087 ctk[anode*fMaxNofSamples+i] += ctktmp;
1088 }
1089 anode = z + 1;
1090 if( anode < fNofMaps ) {
1091 ctk[anode*fMaxNofSamples+i] += ctktmp;
1092 }
1093 }
1094 delete [] dev;
1095
1096 } // if( nTsteps > 2 )
1097 on = kFALSE;
1098 } // if( on == kTRUE )
1099 } // else
1100 }
1101 }
1102
1103 for( Int_t a=0; a<fNofMaps; a++ )
1104 for( Int_t t=0; t<fMaxNofSamples; t++ ) {
1105 Float_t signal = fHitMap2->GetSignal( a, t ) + ctk[a*fMaxNofSamples+t];
1106 fHitMap2->SetHit( a, t, signal );
1107 }
1108
1109 delete [] ctk;
1110}
8a33ae9e 1111//______________________________________________________________________
b0f5e3fc 1112void AliITSsimulationSDD::GetAnodeBaseline(Int_t i,Float_t &baseline,
1113 Float_t &noise){
8a33ae9e 1114 // Returns the Baseline for a particular anode.
1115 baseline = fBaseline[i];
1116 noise = fNoise[i];
b0f5e3fc 1117}
8a33ae9e 1118//______________________________________________________________________
b0f5e3fc 1119void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl,
1120 Int_t &th){
8a33ae9e 1121 // Returns the compression alogirthm parameters
1122 Int_t size = fD.GetSize();
1123 if (size > 2 ) {
50d05d7b 1124 db=fD[i]; tl=fT1[i]; th=fT2[i];
8a33ae9e 1125 } else {
50d05d7b 1126 if (size <= 2 && i>=fNofMaps/2) {
1127 db=fD[1]; tl=fT1[1]; th=fT2[1];
1128 } else {
1129 db=fD[0]; tl=fT1[0]; th=fT2[0];
1130 } // end if size <=2 && i>=fNofMaps/2
8a33ae9e 1131 } // end if size >2
b0f5e3fc 1132}
8a33ae9e 1133//______________________________________________________________________
b0f5e3fc 1134void AliITSsimulationSDD::CompressionParam(Int_t i,Int_t &db,Int_t &tl){
8a33ae9e 1135 // returns the compression alogirthm parameters
1136 Int_t size = fD.GetSize();
b0f5e3fc 1137
8a33ae9e 1138 if (size > 2 ) {
50d05d7b 1139 db=fD[i]; tl=fT1[i];
8a33ae9e 1140 } else {
50d05d7b 1141 if (size <= 2 && i>=fNofMaps/2) {
1142 db=fD[1]; tl=fT1[1];
1143 } else {
1144 db=fD[0]; tl=fT1[0];
1145 } // end if size <=2 && i>=fNofMaps/2
8a33ae9e 1146 } // end if size > 2
b0f5e3fc 1147}
8a33ae9e 1148//______________________________________________________________________
b0f5e3fc 1149void AliITSsimulationSDD::SetCompressParam(){
8a33ae9e 1150 // Sets the compression alogirthm parameters
1151 Int_t cp[8],i;
1152
1153 fResponse->GiveCompressParam(cp);
1154 for (i=0; i<2; i++) {
50d05d7b 1155 fD[i] = cp[i];
1156 fT1[i] = cp[i+2];
1157 fT2[i] = cp[i+4];
1158 fTol[i] = cp[i+6];
8a33ae9e 1159 } // end for i
b0f5e3fc 1160}
8a33ae9e 1161//______________________________________________________________________
b0f5e3fc 1162void AliITSsimulationSDD::ReadBaseline(){
8a33ae9e 1163 // read baseline and noise from file - either a .root file and in this
1164 // case data should be organised in a tree with one entry for each
1165 // module => reading should be done accordingly
1166 // or a classic file and do smth. like this:
1167 // Read baselines and noise for SDD
b0f5e3fc 1168
1169 Int_t na,pos;
1170 Float_t bl,n;
e8189707 1171 char input[100], base[100], param[100];
b0f5e3fc 1172 char *filtmp;
1173
e8189707 1174 fResponse->Filenames(input,base,param);
b0f5e3fc 1175 fFileName=base;
1176//
1177 filtmp = gSystem->ExpandPathName(fFileName.Data());
1178 FILE *bline = fopen(filtmp,"r");
b0f5e3fc 1179 na = 0;
1180
1181 if(bline) {
50d05d7b 1182 while(fscanf(bline,"%d %f %f",&pos, &bl, &n) != EOF) {
1183 if (pos != na+1) {
1184 Error("ReadBaseline","Anode number not in increasing order!",
1185 filtmp);
1186 exit(1);
1187 } // end if pos != na+1
1188 fBaseline[na]=bl;
1189 fNoise[na]=n;
1190 na++;
1191 } // end while
b0f5e3fc 1192 } else {
50d05d7b 1193 Error("ReadBaseline"," THE BASELINE FILE %s DOES NOT EXIST !",filtmp);
1194 exit(1);
b0f5e3fc 1195 } // end if(bline)
8a33ae9e 1196
b0f5e3fc 1197 fclose(bline);
1198 delete [] filtmp;
b0f5e3fc 1199}
8a33ae9e 1200//______________________________________________________________________
1201Int_t AliITSsimulationSDD::Convert10to8(Int_t signal) const {
1202 // To the 10 to 8 bit lossive compression.
1203 // code from Davide C. and Albert W.
1204
1205 if (signal < 128) return signal;
1206 if (signal < 256) return (128+((signal-128)>>1));
1207 if (signal < 512) return (192+((signal-256)>>3));
1208 if (signal < 1024) return (224+((signal-512)>>4));
1209 return 0;
b0f5e3fc 1210}
8a33ae9e 1211//______________________________________________________________________
1212Int_t AliITSsimulationSDD::Convert8to10(Int_t signal) const {
1213 // Undo the lossive 10 to 8 bit compression.
1214 // code from Davide C. and Albert W.
1215 if (signal < 0 || signal > 255) {
50d05d7b 1216 Warning("Convert8to10","out of range signal=%d",signal);
1217 return 0;
8a33ae9e 1218 } // end if signal <0 || signal >255
1219
1220 if (signal < 128) return signal;
1221 if (signal < 192) {
50d05d7b 1222 if (TMath::Odd(signal)) return (128+((signal-128)<<1));
1223 else return (128+((signal-128)<<1)+1);
8a33ae9e 1224 } // end if signal < 192
1225 if (signal < 224) {
50d05d7b 1226 if (TMath::Odd(signal)) return (256+((signal-192)<<3)+3);
1227 else return (256+((signal-192)<<3)+4);
8a33ae9e 1228 } // end if signal < 224
1229 if (TMath::Odd(signal)) return (512+((signal-224)<<4)+7);
9226cc6a 1230 return (512+((signal-224)<<4)+7);
8a33ae9e 1231}
50d05d7b 1232
1233/*
8a33ae9e 1234//______________________________________________________________________
b0f5e3fc 1235AliITSMap* AliITSsimulationSDD::HitMap(Int_t i){
8a33ae9e 1236 //Return the correct map.
1237
b0f5e3fc 1238 return ((i==0)? fHitMap1 : fHitMap2);
50d05d7b 1239}*/
1240
8a33ae9e 1241//______________________________________________________________________
e8189707 1242void AliITSsimulationSDD::ZeroSuppression(const char *option) {
8a33ae9e 1243 // perform the zero suppresion
1244
1245 if (strstr(option,"2D")) {
50d05d7b 1246 //Init2D(); // activate if param change module by module
1247 Compress2D();
8a33ae9e 1248 } else if (strstr(option,"1D")) {
50d05d7b 1249 //Init1D(); // activate if param change module by module
1250 Compress1D();
8a33ae9e 1251 } else StoreAllDigits();
b0f5e3fc 1252}
8a33ae9e 1253//______________________________________________________________________
b0f5e3fc 1254void AliITSsimulationSDD::Init2D(){
8a33ae9e 1255 // read in and prepare arrays: fD, fT1, fT2
1256 // savemu[nanodes], savesigma[nanodes]
1257 // read baseline and noise from file - either a .root file and in this
1258 // case data should be organised in a tree with one entry for each
1259 // module => reading should be done accordingly
1260 // or a classic file and do smth. like this ( code from Davide C. and
1261 // Albert W.) :
1262 // Read 2D zero-suppression parameters for SDD
b0f5e3fc 1263
1264 if (!strstr(fParam,"file")) return;
1265
1266 Int_t na,pos,tempTh;
1267 Float_t mu,sigma;
8a33ae9e 1268 Float_t *savemu = new Float_t [fNofMaps];
e8189707 1269 Float_t *savesigma = new Float_t [fNofMaps];
1270 char input[100],basel[100],par[100];
b0f5e3fc 1271 char *filtmp;
b0f5e3fc 1272 Int_t minval = fResponse->MinVal();
1273
e8189707 1274 fResponse->Filenames(input,basel,par);
8a33ae9e 1275 fFileName = par;
b0f5e3fc 1276//
1277 filtmp = gSystem->ExpandPathName(fFileName.Data());
1278 FILE *param = fopen(filtmp,"r");
1279 na = 0;
1280
1281 if(param) {
50d05d7b 1282 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1283 if (pos != na+1) {
1284 Error("Init2D","Anode number not in increasing order!",filtmp);
1285 exit(1);
1286 } // end if pos != na+1
1287 savemu[na] = mu;
8a33ae9e 1288 savesigma[na] = sigma;
b0f5e3fc 1289 if ((2.*sigma) < mu) {
1290 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1291 mu = 2.0 * sigma;
50d05d7b 1292 } else fD[na] = 0;
b0f5e3fc 1293 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1294 if (tempTh < 0) tempTh=0;
1295 fT1[na] = tempTh;
1296 tempTh = (Int_t)floor(mu+3.0*sigma+0.5) - minval;
1297 if (tempTh < 0) tempTh=0;
1298 fT2[na] = tempTh;
1299 na++;
50d05d7b 1300 } // end while
b0f5e3fc 1301 } else {
50d05d7b 1302 Error("Init2D","THE FILE %s DOES NOT EXIST !",filtmp);
1303 exit(1);
b0f5e3fc 1304 } // end if(param)
8a33ae9e 1305
b0f5e3fc 1306 fclose(param);
1307 delete [] filtmp;
5d18fa90 1308 delete [] savemu;
e8189707 1309 delete [] savesigma;
8a33ae9e 1310}
1311//______________________________________________________________________
b0f5e3fc 1312void AliITSsimulationSDD::Compress2D(){
8a33ae9e 1313 // simple ITS cluster finder -- online zero-suppression conditions
b0f5e3fc 1314
b0f5e3fc 1315 Int_t db,tl,th;
8a33ae9e 1316 Int_t minval = fResponse->MinVal();
1317 Bool_t write = fResponse->OutputOption();
1318 Bool_t do10to8 = fResponse->Do10to8();
b0f5e3fc 1319 Int_t nz, nl, nh, low, i, j;
1320
e8189707 1321 for (i=0; i<fNofMaps; i++) {
b0f5e3fc 1322 CompressionParam(i,db,tl,th);
8a33ae9e 1323 nz = 0;
1324 nl = 0;
1325 nh = 0;
1326 low = 0;
50d05d7b 1327 for (j=0; j<fMaxNofSamples; j++) {
1328 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1329 signal -= db; // if baseline eq. is done here
b0f5e3fc 1330 if (signal <= 0) {nz++; continue;}
50d05d7b 1331 if ((signal - tl) < minval) low++;
b0f5e3fc 1332 if ((signal - th) >= minval) {
50d05d7b 1333 nh++;
1334 Bool_t cond=kTRUE;
1335 FindCluster(i,j,signal,minval,cond);
1336 if(cond && j &&
1337 ((TMath::Abs(fHitMap2->GetSignal(i,j-1))-th)>=minval)){
1338 if(do10to8) signal = Convert10to8(signal);
1339 AddDigit(i,j,signal);
1340 } // end if cond&&j&&()
1341 } else if ((signal - tl) >= minval) nl++;
1342 } // end for j loop time samples
1343 if (write) TreeB()->Fill(nz,nl,nh,low,i+1);
8a33ae9e 1344 } //end for i loop anodes
b0f5e3fc 1345
8a33ae9e 1346 char hname[30];
1347 if (write) {
50d05d7b 1348 sprintf(hname,"TNtuple%d_%d",fModule,fEvent);
1349 TreeB()->Write(hname);
1350 // reset tree
b0f5e3fc 1351 TreeB()->Reset();
8a33ae9e 1352 } // end if write
1353}
1354//______________________________________________________________________
b0f5e3fc 1355void AliITSsimulationSDD::FindCluster(Int_t i,Int_t j,Int_t signal,
ece86d9a 1356 Int_t minval,Bool_t &cond){
8a33ae9e 1357 // Find clusters according to the online 2D zero-suppression algorithm
1358 Bool_t do10to8 = fResponse->Do10to8();
1359 Bool_t high = kFALSE;
b0f5e3fc 1360
1361 fHitMap2->FlagHit(i,j);
1362//
1363// check the online zero-suppression conditions
1364//
8a33ae9e 1365 const Int_t kMaxNeighbours = 4;
b0f5e3fc 1366 Int_t nn;
1367 Int_t dbx,tlx,thx;
8a33ae9e 1368 Int_t xList[kMaxNeighbours], yList[kMaxNeighbours];
e8189707 1369 fSegmentation->Neighbours(i,j,&nn,xList,yList);
ece86d9a 1370 Int_t in,ix,iy,qns;
e8189707 1371 for (in=0; in<nn; in++) {
50d05d7b 1372 ix=xList[in];
e8189707 1373 iy=yList[in];
b0f5e3fc 1374 if (fHitMap2->TestHit(ix,iy)==kUnused) {
50d05d7b 1375 CompressionParam(ix,dbx,tlx,thx);
1376 Int_t qn = (Int_t)(fHitMap2->GetSignal(ix,iy));
1377 qn -= dbx; // if baseline eq. is done here
1378 if ((qn-tlx) < minval) {
1379 fHitMap2->FlagHit(ix,iy);
1380 continue;
1381 } else {
1382 if ((qn - thx) >= minval) high=kTRUE;
1383 if (cond) {
1384 if(do10to8) signal = Convert10to8(signal);
1385 AddDigit(i,j,signal);
1386 } // end if cond
1387 if(do10to8) qns = Convert10to8(qn);
1388 else qns=qn;
1389 if (!high) AddDigit(ix,iy,qns);
1390 cond=kFALSE;
1391 if(!high) fHitMap2->FlagHit(ix,iy);
1392 } // end if qn-tlx < minval
1393 } // end if TestHit
8a33ae9e 1394 } // end for in loop over neighbours
b0f5e3fc 1395}
8a33ae9e 1396//______________________________________________________________________
b0f5e3fc 1397void AliITSsimulationSDD::Init1D(){
8a33ae9e 1398 // this is just a copy-paste of input taken from 2D algo
1399 // Torino people should give input
1400 // Read 1D zero-suppression parameters for SDD
b0f5e3fc 1401
1402 if (!strstr(fParam,"file")) return;
1403
1404 Int_t na,pos,tempTh;
1405 Float_t mu,sigma;
8a33ae9e 1406 Float_t *savemu = new Float_t [fNofMaps];
e8189707 1407 Float_t *savesigma = new Float_t [fNofMaps];
1408 char input[100],basel[100],par[100];
b0f5e3fc 1409 char *filtmp;
b0f5e3fc 1410 Int_t minval = fResponse->MinVal();
8a33ae9e 1411
e8189707 1412 fResponse->Filenames(input,basel,par);
1413 fFileName=par;
b0f5e3fc 1414
1415// set first the disable and tol param
1416 SetCompressParam();
1417//
1418 filtmp = gSystem->ExpandPathName(fFileName.Data());
1419 FILE *param = fopen(filtmp,"r");
1420 na = 0;
1421
1422 if (param) {
50d05d7b 1423 fscanf(param,"%d %d %d %d ", &fT2[0], &fT2[1], &fTol[0], &fTol[1]);
1424 while(fscanf(param,"%d %f %f",&pos, &mu, &sigma) != EOF) {
1425 if (pos != na+1) {
1426 Error("Init1D","Anode number not in increasing order!",filtmp);
1427 exit(1);
1428 } // end if pos != na+1
1429 savemu[na]=mu;
1430 savesigma[na]=sigma;
1431 if ((2.*sigma) < mu) {
1432 fD[na] = (Int_t)floor(mu - 2.0*sigma + 0.5);
1433 mu = 2.0 * sigma;
1434 } else fD[na] = 0;
1435 tempTh = (Int_t)floor(mu+2.25*sigma+0.5) - minval;
1436 if (tempTh < 0) tempTh=0;
1437 fT1[na] = tempTh;
1438 na++;
1439 } // end while
b0f5e3fc 1440 } else {
50d05d7b 1441 Error("Init1D","THE FILE %s DOES NOT EXIST !",filtmp);
1442 exit(1);
b0f5e3fc 1443 } // end if(param)
8a33ae9e 1444
b0f5e3fc 1445 fclose(param);
1446 delete [] filtmp;
749bd21a 1447 delete [] savemu;
e8189707 1448 delete [] savesigma;
8a33ae9e 1449}
1450//______________________________________________________________________
b0f5e3fc 1451void AliITSsimulationSDD::Compress1D(){
1452 // 1D zero-suppression algorithm (from Gianluca A.)
8a33ae9e 1453 Int_t dis,tol,thres,decr,diff;
b0f5e3fc 1454 UChar_t *str=fStream->Stream();
8a33ae9e 1455 Int_t counter=0;
1456 Bool_t do10to8=fResponse->Do10to8();
1457 Int_t last=0;
1458 Int_t k,i,j;
ece86d9a 1459
ece86d9a 1460 for (k=0; k<2; k++) {
50d05d7b 1461 tol = Tolerance(k);
1462 dis = Disable(k);
1463 for (i=0; i<fNofMaps/2; i++) {
1464 Bool_t firstSignal=kTRUE;
1465 Int_t idx=i+k*fNofMaps/2;
1466 CompressionParam(idx,decr,thres);
1467 for (j=0; j<fMaxNofSamples; j++) {
1468 Int_t signal=(Int_t)(fHitMap2->GetSignal(idx,j));
1469 signal -= decr; // if baseline eq.
1470 if(do10to8) signal = Convert10to8(signal);
1471 if (signal <= thres) {
1472 signal=0;
1473 diff=128;
1474 last=0;
1475 // write diff in the buffer for HuffT
1476 str[counter]=(UChar_t)diff;
1477 counter++;
1478 continue;
1479 } // end if signal <= thres
1480 diff=signal-last;
1481 if (diff > 127) diff=127;
1482 if (diff < -128) diff=-128;
1483 if (signal < dis) {
1484 // tol has changed to 8 possible cases ? - one can write
1485 // this if(TMath::Abs(diff)<tol) ... else ...
1486 if(TMath::Abs(diff)<tol) diff=0;
1487 // or keep it as it was before
1488 /*
b0f5e3fc 1489 if (tol==1 && (diff >= -2 && diff <= 1)) diff=0;
1490 if (tol==2 && (diff >= -4 && diff <= 3)) diff=0;
1491 if (tol==3 && (diff >= -16 && diff <= 15)) diff=0;
50d05d7b 1492 */
ece86d9a 1493 AddDigit(idx,j,last+diff);
50d05d7b 1494 } else {
1495 AddDigit(idx,j,signal);
1496 } // end if singal < dis
1497 diff += 128;
1498 // write diff in the buffer used to compute Huffman tables
1499 if (firstSignal) str[counter]=(UChar_t)signal;
1500 else str[counter]=(UChar_t)diff;
1501 counter++;
1502 last=signal;
1503 firstSignal=kFALSE;
1504 } // end for j loop time samples
1505 } // end for i loop anodes one half of detector
8a33ae9e 1506 } // end for k
b0f5e3fc 1507
1508 // check
1509 fStream->CheckCount(counter);
1510
1511 // open file and write out the stream of diff's
b0f5e3fc 1512 static Bool_t open=kTRUE;
e8189707 1513 static TFile *outFile;
b0f5e3fc 1514 Bool_t write = fResponse->OutputOption();
9ad8b5dd 1515 TDirectory *savedir = gDirectory;
b0f5e3fc 1516
1517 if (write ) {
50d05d7b 1518 if(open) {
1519 SetFileName("stream.root");
1520 cout<<"filename "<<fFileName<<endl;
1521 outFile=new TFile(fFileName,"recreate");
1522 cout<<"I have opened "<<fFileName<<" file "<<endl;
1523 } // end if open
1524 open = kFALSE;
1525 outFile->cd();
b0f5e3fc 1526 fStream->Write();
50d05d7b 1527 } // endif write
b0f5e3fc 1528
8a33ae9e 1529 fStream->ClearStream();
b0f5e3fc 1530
8a33ae9e 1531 // back to galice.root file
9ad8b5dd 1532 if(savedir) savedir->cd();
8a33ae9e 1533}
1534//______________________________________________________________________
b0f5e3fc 1535void AliITSsimulationSDD::StoreAllDigits(){
8a33ae9e 1536 // if non-zero-suppressed data
1537 Bool_t do10to8 = fResponse->Do10to8();
ece86d9a 1538 Int_t i, j, digits[3];
8a33ae9e 1539
e8189707 1540 for (i=0; i<fNofMaps; i++) {
1541 for (j=0; j<fMaxNofSamples; j++) {
50d05d7b 1542 Int_t signal=(Int_t)(fHitMap2->GetSignal(i,j));
1543 if(do10to8) signal = Convert10to8(signal);
1544 if(do10to8) signal = Convert8to10(signal);
1545 digits[0] = i;
1546 digits[1] = j;
1547 digits[2] = signal;
1548 fITS->AddRealDigit(1,digits);
1549 } // end for j
8a33ae9e 1550 } // end for i
b0f5e3fc 1551}
8a33ae9e 1552//______________________________________________________________________
ece86d9a 1553void AliITSsimulationSDD::CreateHistograms(Int_t scale){
8a33ae9e 1554 // Creates histograms of maps for debugging
1555 Int_t i;
ece86d9a 1556
1557 fHis=new TObjArray(fNofMaps);
e8189707 1558 for (i=0;i<fNofMaps;i++) {
50d05d7b 1559 TString sddName("sdd_");
1560 Char_t candNum[4];
1561 sprintf(candNum,"%d",i+1);
1562 sddName.Append(candNum);
1563 fHis->AddAt(new TH1F(sddName.Data(),"SDD maps",scale*fMaxNofSamples,
1564 0.,(Float_t) scale*fMaxNofSamples), i);
8a33ae9e 1565 } // end for i
b0f5e3fc 1566}
8a33ae9e 1567//______________________________________________________________________
ece86d9a 1568void AliITSsimulationSDD::FillHistograms(){
8a33ae9e 1569 // fill 1D histograms from map
1570
1571 if (!fHis) return;
1572
1573 for( Int_t i=0; i<fNofMaps; i++) {
50d05d7b 1574 TH1F *hist =(TH1F *)fHis->UncheckedAt(i);
1575 Int_t nsamples = hist->GetNbinsX();
1576 for( Int_t j=0; j<nsamples; j++) {
1577 Double_t signal=fHitMap2->GetSignal(i,j);
1578 hist->Fill((Float_t)j,signal);
1579 } // end for j
8a33ae9e 1580 } // end for i
ece86d9a 1581}
8a33ae9e 1582//______________________________________________________________________
b0f5e3fc 1583void AliITSsimulationSDD::ResetHistograms(){
b0f5e3fc 1584 // Reset histograms for this detector
b0f5e3fc 1585 Int_t i;
8a33ae9e 1586
e8189707 1587 for (i=0;i<fNofMaps;i++ ) {
50d05d7b 1588 if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset();
8a33ae9e 1589 } // end for i
b0f5e3fc 1590}
8a33ae9e 1591//______________________________________________________________________
b0f5e3fc 1592TH1F *AliITSsimulationSDD::GetAnode(Int_t wing, Int_t anode) {
8a33ae9e 1593 // Fills a histogram from a give anode.
1594
1595 if (!fHis) return 0;
1596
1597 if(wing <=0 || wing > 2) {
50d05d7b 1598 Warning("GetAnode","Wrong wing number: %d",wing);
1599 return NULL;
8a33ae9e 1600 } // end if wing <=0 || wing >2
1601 if(anode <=0 || anode > fNofMaps/2) {
50d05d7b 1602 Warning("GetAnode","Wrong anode number: %d",anode);
1603 return NULL;
8a33ae9e 1604 } // end if ampde <=0 || andoe > fNofMaps/2
1605
1606 Int_t index = (wing-1)*fNofMaps/2 + anode-1;
8a33ae9e 1607 return (TH1F*)(fHis->At(index));
b0f5e3fc 1608}
8a33ae9e 1609//______________________________________________________________________
b0f5e3fc 1610void AliITSsimulationSDD::WriteToFile(TFile *hfile) {
8a33ae9e 1611 // Writes the histograms to a file
b0f5e3fc 1612
8a33ae9e 1613 if (!fHis) return;
1614
1615 hfile->cd();
1616 Int_t i;
8a33ae9e 1617 for(i=0; i<fNofMaps; i++) fHis->At(i)->Write(); //fAdcs[i]->Write();
1618 return;
b0f5e3fc 1619}
8a33ae9e 1620//______________________________________________________________________
ece86d9a 1621Float_t AliITSsimulationSDD::GetNoise() {
8a33ae9e 1622 // Returns the noise value
1623 //Bool_t do10to8=fResponse->Do10to8();
1624 //noise will always be in the liniar part of the signal
1625 Int_t decr;
1626 Int_t threshold = fT1[0];
1627 char opt1[20], opt2[20];
1628
1629 fResponse->ParamOptions(opt1,opt2);
1630 fParam=opt2;
1631 char *same = strstr(opt1,"same");
1632 Float_t noise,baseline;
1633 if (same) {
50d05d7b 1634 fResponse->GetNoiseParam(noise,baseline);
8a33ae9e 1635 } else {
50d05d7b 1636 static Bool_t readfile=kTRUE;
1637 //read baseline and noise from file
1638 if (readfile) ReadBaseline();
1639 readfile=kFALSE;
8a33ae9e 1640 } // end if same
1641
1642 TCanvas *c2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c2");
1643 if(c2) delete c2->GetPrimitive("noisehist");
1644 if(c2) delete c2->GetPrimitive("anode");
1645 else c2=new TCanvas("c2");
1646 c2->cd();
1647 c2->SetFillColor(0);
1648
1649 TH1F *noisehist = new TH1F("noisehist","noise",100,0.,(float)2*threshold);
1650 TH1F *anode = new TH1F("anode","Anode Projection",fMaxNofSamples,0.,
50d05d7b 1651 (float)fMaxNofSamples);
8a33ae9e 1652 Int_t i,k;
1653 for (i=0;i<fNofMaps;i++) {
50d05d7b 1654 CompressionParam(i,decr,threshold);
1655 if (!same) GetAnodeBaseline(i,baseline,noise);
1656 anode->Reset();
1657 for (k=0;k<fMaxNofSamples;k++) {
1658 Float_t signal=(Float_t)fHitMap2->GetSignal(i,k);
1659 //if (signal <= (float)threshold) noisehist->Fill(signal-baseline);
1660 if (signal <= (float)threshold) noisehist->Fill(signal);
1661 anode->Fill((float)k,signal);
1662 } // end for k
1663 anode->Draw();
1664 c2->Update();
8a33ae9e 1665 } // end for i
1666 TF1 *gnoise = new TF1("gnoise","gaus",0.,threshold);
1667 noisehist->Fit("gnoise","RQ");
1668 noisehist->Draw();
ece86d9a 1669 c2->Update();
8a33ae9e 1670 Float_t mnoise = gnoise->GetParameter(1);
1671 cout << "mnoise : " << mnoise << endl;
1672 Float_t rnoise = gnoise->GetParameter(2);
1673 cout << "rnoise : " << rnoise << endl;
1674 delete noisehist;
1675 return rnoise;
50d05d7b 1676}
1677//______________________________________________________________________
1678void AliITSsimulationSDD::WriteSDigits(){
c7a4dac0 1679 // Fills the Summable digits Tree
1680 Int_t i,ni,j,nj;
1681 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
1682
50d05d7b 1683 fpList->GetMaxMapIndex( ni, nj );
1684 for( i=0; i<ni; i++ )
1685 for( j=0; j<nj; j++ ) {
1686// if( fpList->GetSignalOnly( i, j ) > 0.5*fT1[0] ) {
1687 if( fpList->GetpListItem( i, j ) != 0 )
1688 // above small threshold.
1689 aliITS->AddSumDigit( *(fpList->GetpListItem( i, j ) ) );
1690// cout << "pListSDD: " << *(pList->GetpListItem(i,j)) << endl;
1691// } // end if
1692 } // end for i,j
c7a4dac0 1693 return;
b0f5e3fc 1694}
8a33ae9e 1695//______________________________________________________________________
44a312c3 1696void AliITSsimulationSDD::Print() {
8a33ae9e 1697 // Print SDD simulation Parameters
1698
1699 cout << "**************************************************" << endl;
1700 cout << " Silicon Drift Detector Simulation Parameters " << endl;
1701 cout << "**************************************************" << endl;
1702 cout << "Flag for Perpendicular tracks: " << (Int_t) fFlag << endl;
1703 cout << "Flag for noise checking: " << (Int_t) fCheckNoise << endl;
1704 cout << "Flag to switch off electronics: " << (Int_t) fDoFFT << endl;
1705 cout << "Number pf Anodes used: " << fNofMaps << endl;
1706 cout << "Number of Time Samples: " << fMaxNofSamples << endl;
1707 cout << "Scale size factor: " << fScaleSize << endl;
1708 cout << "**************************************************" << endl;
44a312c3 1709}