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