]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSsimulationSPDdubna.cxx
Removing GetDebug and SetDebug from AliRun and AliModule. Using AliLog for the messages
[u/mrichter/AliRoot.git] / ITS / AliITSsimulationSPDdubna.cxx
CommitLineData
f74211b0 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
12e7c97c 16/*
17$Id$
18*/
4ae5bbc4 19#include <Riostream.h>
409f8c84 20#include <TRandom.h>
21#include <TH1.h>
22#include <TMath.h>
23#include <TString.h>
24#include <TParticle.h>
25
409f8c84 26#include "AliRun.h"
27#include "AliITS.h"
28#include "AliITShit.h"
e869281d 29#include "AliITSdigitSPD.h"
409f8c84 30#include "AliITSmodule.h"
31#include "AliITSMapA2.h"
2cc6b29a 32#include "AliITSpList.h"
409f8c84 33#include "AliITSsimulationSPDdubna.h"
f74211b0 34#include "AliITSsegmentationSPD.h"
aacedc3e 35#include "AliITSresponseSPD.h"
409f8c84 36
f74211b0 37//#define DEBUG
409f8c84 38
39ClassImp(AliITSsimulationSPDdubna)
40////////////////////////////////////////////////////////////////////////
12e7c97c 41// Version: 1
42// Modified by Bjorn S. Nilsen
409f8c84 43// Version: 0
44// Written by Boris Batyunya
45// December 20 1999
46//
12e7c97c 47// AliITSsimulationSPDdubna is to do the simulation of SPDs.
2cc6b29a 48//______________________________________________________________________
aacedc3e 49AliITSsimulationSPDdubna::AliITSsimulationSPDdubna():
50AliITSsimulation(),
51fHis(0),
52fSPDname(),
53fCoupling(0){
12e7c97c 54 // Default constructor.
55 // Inputs:
56 // none.
57 // Outputs:
58 // none.
59 // Return:
60 // A default constructed AliITSsimulationSPDdubna class.
2cc6b29a 61
aacedc3e 62 if(GetDebug(1)) Info("AliITSsimulationSPDdubda()",
63 "Calling degault constructor");
409f8c84 64}
2cc6b29a 65//______________________________________________________________________
66AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(AliITSsegmentation *seg,
aacedc3e 67 AliITSresponse *resp,Int_t cup):
68AliITSsimulation(seg,resp),
69fHis(0),
70fSPDname(),
71fCoupling(cup){
2cc6b29a 72 // standard constructor
12e7c97c 73 // Inputs:
74 // AliITSsegmentation *seg A pointer to the segmentation class
75 // to be used for this simulation
76 // AliITSresponse *resp A pointer to the responce class to
77 // be used for this simulation
aacedc3e 78 // Int_t cup The type of coupling to be used
79 // =1 uses SetCoupling, =2 uses SetCouplingOld
80 // With diffusion tured off
81 // =3 uses SetCoupling, =4 uses SetCouplingOld
82 // with diffusion on other, no coupling.
12e7c97c 83 // Outputs:
84 // none.
85 // Return:
86 // A default constructed AliITSsimulationSPDdubna class.
409f8c84 87
aacedc3e 88 if(GetDebug(1)) Info("AliITSsimulationSPDdubda",
89 "Calling degault constructor seg=%p resp=%p cup=%d",
90 seg,resp,cup);
91 if(cup==1||cup==2){ // For the moment, remove defusion if Coupling is
92 // set.
93 resp->SetTemperature(0.0);
94 resp->SetDistanceOverVoltage(0.0);
95 } // end if
96 Init();
12e7c97c 97}
98//______________________________________________________________________
aacedc3e 99void AliITSsimulationSPDdubna::Init(){
12e7c97c 100 // Initilization
101 // Inputs:
aacedc3e 102 // none.
12e7c97c 103 // Outputs:
104 // none.
105 // Return:
106 // none.
107 const Double_t kmictocm = 1.0e-4; // convert microns to cm.
409f8c84 108
12e7c97c 109 SetModuleNumber(0);
110 SetEventNumber(0);
111 SetMap(new AliITSpList(GetNPixelsZ(),GetNPixelsX()));
aacedc3e 112 GetResp(0,0)->SetDistanceOverVoltage(kmictocm*GetSeg()->Dy(),50.0);
2cc6b29a 113}
114//______________________________________________________________________
115AliITSsimulationSPDdubna::~AliITSsimulationSPDdubna(){
116 // destructor
12e7c97c 117 // Inputs:
118 // none.
119 // Outputs:
120 // none.
121 // Return:
122 // none.
409f8c84 123
2cc6b29a 124 if (fHis) {
aacedc3e 125 fHis->Delete();
126 delete fHis;
2cc6b29a 127 } // end if fHis
409f8c84 128}
2cc6b29a 129//______________________________________________________________________
130AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(const
131 AliITSsimulationSPDdubna
12e7c97c 132 &s) : AliITSsimulation(s){
133 // Copy Constructor
134 // Inputs:
135 // AliITSsimulationSPDdubna &s The original class for which
136 // this class is a copy of
137 // Outputs:
138 // none.
139 // Return:
140
141 *this = s;
2cc6b29a 142 return;
143}
144//______________________________________________________________________
145AliITSsimulationSPDdubna& AliITSsimulationSPDdubna::operator=(const
12e7c97c 146 AliITSsimulationSPDdubna &s){
2cc6b29a 147 // Assignment operator
12e7c97c 148 // Inputs:
149 // AliITSsimulationSPDdubna &s The original class for which
150 // this class is a copy of
151 // Outputs:
152 // none.
153 // Return:
154
155 if(&s == this) return *this;
aacedc3e 156 this->fHis = s.fHis;
157 fCoupling = s.fCoupling;
158 fSPDname = s.fSPDname;
2cc6b29a 159 return *this;
160}
161//______________________________________________________________________
d2f55a22 162AliITSsimulationSPDdubna& AliITSsimulationSPDdubna::operator=(const
163 AliITSsimulation &s){
164 // Assignment operator
165 // Inputs:
166 // AliITSsimulationSPDdubna &s The original class for which
167 // this class is a copy of
168 // Outputs:
169 // none.
170 // Return:
171
172 if(&s == this) return *this;
173 Error("AliITSsimulationSPDdubna","Not allowed to make a = with "
174 "AliITSsimulationSPDdubna","Using default creater instead");
175
176 return *this;
177}
178//______________________________________________________________________
2cc6b29a 179void AliITSsimulationSPDdubna::InitSimulationModule(Int_t module, Int_t event){
180 // This function creates maps to build the list of tracks for each
12e7c97c 181 // summable digit. Inputs defined by base class.
2cc6b29a 182 // Inputs:
183 // Int_t module // Module number to be simulated
184 // Int_t event // Event number to be simulated
2cc6b29a 185 // Outputs:
186 // none
2cc6b29a 187 // Returns:
188 // none
409f8c84 189
aacedc3e 190 if(GetDebug(1)) Info("InitSimulationModule","(module=%d,event=%d)",
191 module,event);
12e7c97c 192 SetModuleNumber(module);
193 SetEventNumber(event);
194 ClearMap();
2cc6b29a 195}
196//_____________________________________________________________________
aacedc3e 197void AliITSsimulationSPDdubna::SDigitiseModule(AliITSmodule *mod,Int_t,
198 Int_t event){
12e7c97c 199 // This function begins the work of creating S-Digits. Inputs defined
200 // by base class.
2cc6b29a 201 // Inputs:
202 // AliITSmodule *mod // module
aacedc3e 203 // Int_t // not used
204 // Int_t event // Event number
2cc6b29a 205 // Outputs:
206 // none
2cc6b29a 207 // Return:
208 // test // test returns kTRUE if the module contained hits
209 // // test returns kFALSE if it did not contain hits
409f8c84 210
aacedc3e 211 if(GetDebug(1)) Info("SDigitiseModule","(mod=%p, ,event=%d)",mod,event);
212 if(!(mod->GetNhits())){
213 if(GetDebug(1)) Info("SDigitiseModule","In event %d module %d there "
214 "are %d hits returning.",event,
215 mod->GetIndex(),mod->GetNhits());
216 return;// if module has no hits don't create Sdigits
217 } // end if
12e7c97c 218 SetModuleNumber(mod->GetIndex());
219 SetEventNumber(event);
220 HitToSDigit(mod);
221 WriteSDigits();
222 ClearMap();
409f8c84 223}
2cc6b29a 224//______________________________________________________________________
12e7c97c 225void AliITSsimulationSPDdubna::WriteSDigits(){
2cc6b29a 226 // This function adds each S-Digit to pList
2cc6b29a 227 // Inputs:
12e7c97c 228 // none.
2cc6b29a 229 // Outputs:
12e7c97c 230 // none.
2cc6b29a 231 // Return:
232 // none
7d50ea21 233 Int_t ix, nix, iz, niz;
2cc6b29a 234 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
235
aacedc3e 236 if(GetDebug(1))Info("WriteSDigits","Writing SDigits for module %d",
237 GetModuleNumber());
12e7c97c 238 GetMap()->GetMaxMapIndex(niz, nix);
239 for(iz=0; iz<niz; iz++)for(ix=0; ix<nix; ix++){
aacedc3e 240 if(GetMap()->GetSignalOnly(iz,ix)>0.0){
241 aliITS->AddSumDigit(*(GetMap()->GetpListItem(iz,ix)));
242 if(GetDebug(1)){
243 cout <<"AliITSsimulationSPDdubna:WriteSDigits " << iz << ","
244 << ix << "," << *(GetMap()->GetpListItem(iz,ix)) << endl;
245 } // end if GetDebug
246 } // end if GetMap()->GetSignalOnly(iz,ix)>0.0
7d50ea21 247 } // end for iz,ix
2cc6b29a 248 return;
249}
250//______________________________________________________________________
251void AliITSsimulationSPDdubna::FinishSDigitiseModule(){
252 // This function calls SDigitsToDigits which creates Digits from SDigits
2cc6b29a 253 // Inputs:
254 // none
2cc6b29a 255 // Outputs:
256 // none
257 // Return
258 // none
409f8c84 259
aacedc3e 260 if(GetDebug(1)) Info("SDigitiseModule","()");
261 pListToDigits(); // Charge To Signal both adds noise and
12e7c97c 262 ClearMap();
aacedc3e 263 return;
409f8c84 264}
2cc6b29a 265//______________________________________________________________________
aacedc3e 266void AliITSsimulationSPDdubna::DigitiseModule(AliITSmodule *mod,Int_t,
267 Int_t){
2cc6b29a 268 // This function creates Digits straight from the hits and then adds
269 // electronic noise to the digits before adding them to pList
12e7c97c 270 // Each of the input variables is passed along to HitToSDigit
2cc6b29a 271 // Inputs:
12e7c97c 272 // AliITSmodule *mod module
aacedc3e 273 // Int_t Dummy.
274 // Int_t Dummy
2cc6b29a 275 // Outputs:
12e7c97c 276 // none.
2cc6b29a 277 // Return:
12e7c97c 278 // none.
409f8c84 279
aacedc3e 280 if(GetDebug(1)) Info("DigitiseModule","(mod=%p,,)",mod);
12e7c97c 281 HitToSDigit(mod);
aacedc3e 282 pListToDigits();
12e7c97c 283 ClearMap();
2cc6b29a 284}
285//______________________________________________________________________
12e7c97c 286void AliITSsimulationSPDdubna::HitToSDigit(AliITSmodule *mod){
f74211b0 287 // Does the charge distributions using Gaussian diffusion charge charing.
12e7c97c 288 // Inputs:
289 // AliITSmodule *mod Pointer to this module
290 // Output:
291 // none.
292 // Return:
293 // none.
f74211b0 294 const Double_t kmictocm = 1.0e-4; // convert microns to cm.
295 TObjArray *hits = mod->GetHits();
296 Int_t nhits = hits->GetEntriesFast();
aacedc3e 297 Int_t h,ix,iz,i;
f74211b0 298 Int_t idtrack;
299 Double_t x0=0.0,x1=0.0,y0=0.0,y1=0.0,z0=0.0,z1=0.0,de=0.0;
300 Double_t x,y,z,t,tp,st,dt=0.2,el,sig;
301 Double_t thick = kmictocm*GetSeg()->Dy();
302
aacedc3e 303 if(GetDebug(1)) Info("HitsToSDigits","(mod=%p) fCoupling=%d",
304 mod,fCoupling);
f74211b0 305 if(nhits<=0) return;
306 for(h=0;h<nhits;h++){
aacedc3e 307 if(GetDebug(1)){
308 cout << "Hits=" << h << "," << *(mod->GetHit(h)) << endl;
309 } // end if GetDebug
310 if(!mod->LineSegmentL(h,x0,x1,y0,y1,z0,z1,de,idtrack)) continue;
311 st = TMath::Sqrt(x1*x1+y1*y1+z1*z1);
312 if(st>0.0){
c789ee28 313 st = (Double_t)((Int_t)(st/kmictocm)); // number of microns
314 if(st<=1.0) st = 1.0;
aacedc3e 315 dt = 1.0/st;
316 for(t=0.0;t<1.0;t+=dt){ // Integrate over t
317 tp = t+0.5*dt;
318 x = x0+x1*tp;
319 y = y0+y1*tp;
320 z = z0+z1*tp;
321 if(!(GetSeg()->LocalToDet(x,z,ix,iz))) continue; // outside
322 el = GetResp(ix,iz)->GeVToCharge((Double_t)(dt*de));
323 if(GetDebug(1)){
324 if(el<=0.0) cout<<"el="<<el<<" dt="<<dt
325 <<" de="<<de<<endl;
326 } // end if GetDebug
327 sig = GetResp(ix,iz)->SigmaDiffusion1D(thick + y);
328 SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
329 } // end for t
330 } else { // st == 0.0 deposit it at this point
331 x = x0;
332 y = y0;
333 z = z0;
334 if(!(GetSeg()->LocalToDet(x,z,ix,iz))) continue; // outside
335 el = GetResp(ix,iz)->GeVToCharge((Double_t)de);
336 sig = GetResp(ix,iz)->SigmaDiffusion1D(thick + y);
337 SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
338 } // end if st>0.0
339 // Coupling
340 switch (fCoupling) {
341 default:
342 break;
343 case 1: case 3:
344 // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
345 for(i=0;i<GetMap()->GetEntries();i++)
346 if(GetMap()->GetpListItem(i)==0) continue;
347 else{
348 GetMap()->GetMapIndex(
349 GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
350 SetCoupling(iz,ix,idtrack,h);
351 } // end for i
352 break;
353 case 2: case 4:
354 // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
355 for(i=0;i<GetMap()->GetEntries();i++)
356 if(GetMap()->GetpListItem(i)==0) continue;
357 else{
358 GetMap()->GetMapIndex(
359 GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
360 SetCouplingOld(iz,ix,idtrack,h);
361 } // end for i
362 break;
363 } // end switch
12e7c97c 364 } // Loop over all hits h
aacedc3e 365 if(GetDebug(2))Info("HitToSDigit","Finished fCoupling=%d",fCoupling);
12e7c97c 366}
f74211b0 367//______________________________________________________________________
12e7c97c 368void AliITSsimulationSPDdubna::SpreadCharge(Double_t x0,Double_t z0,
369 Int_t ix0,Int_t iz0,
f74211b0 370 Double_t el,Double_t sig,Int_t t,
12e7c97c 371 Int_t hi){
f74211b0 372 // Spreads the charge over neighboring cells. Assume charge is distributed
373 // as charge(x,z) = (el/2*pi*sig*sig)*exp(-arg)
374 // arg=((x-x0)*(x-x0)/2*sig*sig)+((z-z0*z-z0)/2*sig*sig)
375 // Defined this way, the integral over all x and z is el.
12e7c97c 376 // Inputs:
377 // Double_t x0 x position of point where charge is liberated
378 // Double_t y0 y position of point where charge is liberated
379 // Double_t z0 z position of point where charge is liberated
380 // Int_t ix0 row of cell corresponding to point x0
381 // Int_t iz0 columb of cell corresponding to point z0
382 // Double_t el number of electrons liberated in this step
383 // Double_t sig Sigma difusion for this step (y0 dependent)
384 // Int_t t track number
385 // Int_t ti hit track index number
386 // Int_t hi hit "hit" index number
387 // Outputs:
388 // none.
389 // Return:
390 // none.
f74211b0 391 const Int_t knx = 3,knz = 2;
392 const Double_t kRoot2 = 1.414213562; // Sqrt(2).
393 const Double_t kmictocm = 1.0e-4; // convert microns to cm.
394 Int_t ix,iz,ixs,ixe,izs,ize;
395 Float_t x,z;
396 Double_t x1,x2,z1,z2,s,sp;
397
aacedc3e 398 if(GetDebug(4)) Info("SpreadCharge","(x0=%e,z0=%e,ix0=%d,iz0=%d,el=%e,"
399 "sig=%e,t=%d,i=%d)",x0,z0,ix0,iz0,el,sig,t,hi);
400 if(sig<=0.0) { // if sig<=0 No diffusion to simulate.
401 GetMap()->AddSignal(iz0,ix0,t,hi,GetModuleNumber(),el);
402 if(GetDebug(2)){
403 cout << "sig<=0.0=" << sig << endl;
404 } // end if GetDebug
405 return;
f74211b0 406 } // end if
407 sp = 1.0/(sig*kRoot2);
aacedc3e 408 if(GetDebug(2)){
409 cout << "sig=" << sig << " sp=" << sp << endl;
410 } // end if GetDebug
f74211b0 411 ixs = TMath::Max(-knx+ix0,0);
412 ixe = TMath::Min(knx+ix0,GetSeg()->Npx()-1);
413 izs = TMath::Max(-knz+iz0,0);
414 ize = TMath::Min(knz+iz0,GetSeg()->Npz()-1);
415 for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){
aacedc3e 416 GetSeg()->DetToLocal(ix,iz,x,z); // pixel center
417 x1 = x;
418 z1 = z;
419 x2 = x1 + 0.5*kmictocm*GetSeg()->Dpx(ix); // Upper
420 x1 -= 0.5*kmictocm*GetSeg()->Dpx(ix); // Lower
421 z2 = z1 + 0.5*kmictocm*GetSeg()->Dpz(iz); // Upper
422 z1 -= 0.5*kmictocm*GetSeg()->Dpz(iz); // Lower
423 x1 -= x0; // Distance from where track traveled
424 x2 -= x0; // Distance from where track traveled
425 z1 -= z0; // Distance from where track traveled
426 z2 -= z0; // Distance from where track traveled
427 s = 0.25; // Correction based on definision of Erfc
428 s *= TMath::Erfc(sp*x1) - TMath::Erfc(sp*x2);
429 if(GetDebug(3)){
430 cout <<"el="<<el<<" ix0="<<ix0<<" ix="<<ix<<" x0="<<x<<
431 " iz0="<<iz0<<" iz="<<iz<<" z0="<<z<<
432 " sp*x1="<<sp*x1<<" sp*x2="<<sp*x2<<" s="<<s;
433 } // end if GetDebug
434 s *= TMath::Erfc(sp*z1) - TMath::Erfc(sp*z2);
435 if(GetDebug(3)){
436 cout<<" sp*z1="<<sp*z1<<" sp*z2="<<sp*z2<<" s="<<s<< endl;
437 } // end if GetDebug
438 GetMap()->AddSignal(iz,ix,t,hi,GetModuleNumber(),s*el);
f74211b0 439 } // end for ix, iz
440}
441//______________________________________________________________________
aacedc3e 442void AliITSsimulationSPDdubna::pListToDigits(){
2cc6b29a 443 // add noise and electronics, perform the zero suppression and add the
444 // digit to the list
12e7c97c 445 // Inputs:
446 // none.
447 // Outputs:
448 // none.
449 // Return:
450 // none.
f74211b0 451 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
12e7c97c 452 Int_t j,ix,iz;
aacedc3e 453 Double_t electronics;
2cc6b29a 454 Double_t sig;
ee86d557 455 const Int_t nmaxtrk=AliITSdigitSPD::GetNTracks();
f74211b0 456 static AliITSdigitSPD dig;
457
aacedc3e 458 if(GetDebug(1)) Info("pListToDigits","()");
12e7c97c 459 for(iz=0; iz<GetNPixelsZ(); iz++) for(ix=0; ix<GetNPixelsX(); ix++){
aacedc3e 460 // Apply Noise/Dead channals and the like
461 if(GetResp(ix,iz)->IsPixelDead(GetModuleNumber(),ix,iz)) continue;
462 electronics = GetResp(ix,iz)->ApplyBaselineAndNoise();
463 UpdateMapNoise(ix,iz,electronics);
464 //
465 // Apply Threshold and write Digits.
466 sig = GetMap()->GetSignalOnly(iz,ix);
467 FillHistograms(ix,iz,sig+electronics);
468 if(GetDebug(3)){
469 cout<<sig<<"+"<<electronics<<">threshold("<<ix<<","<<iz
470 <<")="<<GetThreshold(ix,iz) <<endl;
471 } // end if GetDebug
472 if (sig+electronics <= GetThreshold(ix,iz)) continue;
473 dig.SetCoord1(iz);
474 dig.SetCoord2(ix);
475 dig.SetSignal(1);
476 dig.SetSignalSPD((Int_t) GetMap()->GetSignal(iz,ix));
477 for(j=0;j<nmaxtrk;j++){
478 if (j<GetMap()->GetNEnteries()) {
479 dig.SetTrack(j,GetMap()->GetTrack(iz,ix,j));
480 dig.SetHit(j,GetMap()->GetHit(iz,ix,j));
481 }else { // Default values
482 dig.SetTrack(j,-3);
483 dig.SetHit(j,-1);
484 } // end if GetMap()
485 } // end for j
486 if(GetDebug(3)){
487 cout<<iz<<","<<ix<<","<<*(GetMap()->GetpListItem(iz,ix))<<endl;
488 } // end if GetDebug
489 aliITS->AddSimDigit(0,&dig);
12e7c97c 490 } // for ix/iz
409f8c84 491}
2cc6b29a 492//______________________________________________________________________
493void AliITSsimulationSPDdubna::CreateHistograms(){
494 // create 1D histograms for tests
12e7c97c 495 // Inputs:
496 // none.
497 // Outputs:
498 // none.
499 // Return:
500 // none.
2cc6b29a 501
aacedc3e 502 if(GetDebug(1)) Info("CreateHistograms","create histograms");
2cc6b29a 503
aacedc3e 504 fHis = new TObjArray(GetNPixelsZ());
505 TString fSPDname("spd_");
506 for(Int_t i=0;i<GetNPixelsZ();i++) {
507 Char_t pixelz[4];
508 sprintf(pixelz,"%d",i);
509 fSPDname.Append(pixelz);
510 fHis->AddAt(new TH1F(fSPDname.Data(),"SPD maps",
511 GetNPixelsX(),0.,(Double_t)GetNPixelsX()),i);
2cc6b29a 512 } // end for i
409f8c84 513}
2cc6b29a 514//______________________________________________________________________
aacedc3e 515void AliITSsimulationSPDdubna::FillHistograms(Int_t ix,Int_t iz,Double_t v){
516 // Fill the histogram
517 // Inputs:
518 // none.
519 // Outputs:
520 // none.
521 // Return:
522 // none.
523
524 if(!GetHistArray()) return; // Only fill if setup.
525 if(GetDebug(2)) Info("FillHistograms","fill histograms");
526 GetHistogram(iz)->Fill(ix,v);
527}
528//______________________________________________________________________
2cc6b29a 529void AliITSsimulationSPDdubna::ResetHistograms(){
409f8c84 530 // Reset histograms for this detector
12e7c97c 531 // Inputs:
532 // none.
533 // Outputs:
534 // none.
535 // Return:
536 // none.
409f8c84 537
aacedc3e 538 if(!GetHistArray()) return; // Only fill if setup.
539 if(GetDebug(2)) Info("FillHistograms","fill histograms");
12e7c97c 540 for ( int i=0;i<GetNPixelsZ();i++ ) {
aacedc3e 541 if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset();
2cc6b29a 542 } // end for i
409f8c84 543}
12e7c97c 544
545//______________________________________________________________________
546void AliITSsimulationSPDdubna::SetCoupling(Int_t row, Int_t col, Int_t ntrack,
547 Int_t idhit) {
548 // Take into account the coupling between adiacent pixels.
549 // The parameters probcol and probrow are the probability of the
550 // signal in one pixel shared in the two adjacent pixels along
551 // the column and row direction, respectively.
aacedc3e 552 // Note pList is goten via GetMap() and module is not need any more.
553 // Otherwise it is identical to that coded by Tiziano Virgili (BSN).
12e7c97c 554 //Begin_Html
555 /*
556 <img src="picts/ITS/barimodel_3.gif">
557 </pre>
558 <br clear=left>
559 <font size=+2 color=red>
560 <a href="mailto:tiziano.virgili@cern.ch"></a>.
561 </font>
562 <pre>
563 */
564 //End_Html
565 // Inputs:
566 // Int_t row z cell index
567 // Int_t col x cell index
568 // Int_t ntrack track incex number
569 // Int_t idhit hit index number
12e7c97c 570 // Outputs:
571 // none.
572 // Return:
573 // none.
574 Int_t j1,j2,flag=0;
575 Double_t pulse1,pulse2;
aacedc3e 576 Double_t couplR=0.0,couplC=0.0;
12e7c97c 577 Double_t xr=0.;
12e7c97c 578
579 GetCouplings(couplR,couplC);
aacedc3e 580 if(GetDebug(3)) Info("SetCoupling","(row=%d,col=%d,ntrack=%d,idhit=%d) "
581 "Calling SetCoupling couplR=%e couplC=%e",
582 row,col,ntrack,idhit,couplR,couplC);
12e7c97c 583 j1 = row;
584 j2 = col;
aacedc3e 585 pulse1 = GetMap()->GetSignalOnly(row,col);
12e7c97c 586 pulse2 = pulse1;
587 for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
aacedc3e 588 do{
589 j1 += isign;
590 // pulse1 *= couplR;
591 xr = gRandom->Rndm();
592 //if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold(j1,col))){
593 if ((j1<0) || (j1>GetNPixelsZ()-1) || (xr>couplR)){
594 j1 = row;
595 flag = 1;
596 }else{
597 UpdateMapSignal(col,j1,ntrack,idhit,pulse1);
598 // flag = 0;
599 flag = 1; // only first next!!
600 } // end if
601 } while(flag == 0);
602 // loop in column direction
603 do{
604 j2 += isign;
605 // pulse2 *= couplC;
606 xr = gRandom->Rndm();
607 //if((j2<0)||j2>(GetNPixelsX()-1)||pulse2<GetThreshold(row,j2)){
608 if ((j2<0) || (j2>GetNPixelsX()-1) || (xr>couplC)){
609 j2 = col;
610 flag = 1;
611 }else{
612 UpdateMapSignal(j2,row,ntrack,idhit,pulse2);
613 // flag = 0;
614 flag = 1; // only first next!!
615 } // end if
616 } while(flag == 0);
12e7c97c 617 } // for isign
618}
619//______________________________________________________________________
620void AliITSsimulationSPDdubna::SetCouplingOld(Int_t row, Int_t col,
621 Int_t ntrack,Int_t idhit) {
622 // Take into account the coupling between adiacent pixels.
623 // The parameters probcol and probrow are the fractions of the
624 // signal in one pixel shared in the two adjacent pixels along
625 // the column and row direction, respectively.
626 //Begin_Html
627 /*
628 <img src="picts/ITS/barimodel_3.gif">
629 </pre>
630 <br clear=left>
631 <font size=+2 color=red>
632 <a href="mailto:Rocco.Caliandro@ba.infn.it"></a>.
633 </font>
634 <pre>
635 */
636 //End_Html
637 // Inputs:
638 // Int_t row z cell index
639 // Int_t col x cell index
640 // Int_t ntrack track incex number
641 // Int_t idhit hit index number
642 // Int_t module module number
643 // Outputs:
644 // none.
645 // Return:
646 // none.
647 Int_t j1,j2,flag=0;
648 Double_t pulse1,pulse2;
aacedc3e 649 Double_t couplR=0.0,couplC=0.0;
12e7c97c 650
651 GetCouplings(couplR,couplC);
aacedc3e 652 if(GetDebug(3)) Info("SetCouplingOld","(row=%d,col=%d,ntrack=%d,idhit=%d) "
653 "Calling SetCoupling couplR=%e couplC=%e",
654 row,col,ntrack,idhit,couplR,couplC);
12e7c97c 655 j1 = row;
656 j2 = col;
aacedc3e 657 pulse1 = GetMap()->GetSignalOnly(row,col);
12e7c97c 658 pulse2 = pulse1;
659 for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
aacedc3e 660 do{
661 j1 += isign;
662 pulse1 *= couplR;
663 if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold(j1,col))){
664 pulse1 = GetMap()->GetSignalOnly(row,col);
665 j1 = row;
666 flag = 1;
667 }else{
668 UpdateMapSignal(col,j1,ntrack,idhit,pulse1);
669 flag = 0;
670 } // end if
671 } while(flag == 0);
672 // loop in column direction
673 do{
674 j2 += isign;
675 pulse2 *= couplC;
676 if((j2<0)||(j2>(GetNPixelsX()-1))||(pulse2<GetThreshold(row,j2))){
677 pulse2 = GetMap()->GetSignalOnly(row,col);
678 j2 = col;
679 flag = 1;
680 }else{
681 UpdateMapSignal(j2,row,ntrack,idhit,pulse2);
682 flag = 0;
683 } // end if
684 } while(flag == 0);
12e7c97c 685 } // for isign
686}