1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
19 #include <Riostream.h>
22 #include <TParticle.h>
27 #include "AliITSMapA2.h"
28 #include "AliITSdigitSPD.h"
29 #include "AliITShit.h"
30 #include "AliITSmodule.h"
31 #include "AliITSpList.h"
32 #include "AliITSresponseSPD.h"
33 #include "AliITSsegmentationSPD.h"
34 #include "AliITSsimulationSPDdubna.h"
40 ClassImp(AliITSsimulationSPDdubna)
41 ////////////////////////////////////////////////////////////////////////
43 // Modified by Bjorn S. Nilsen
45 // Written by Boris Batyunya
48 // AliITSsimulationSPDdubna is to do the simulation of SPDs.
49 //______________________________________________________________________
50 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna():
55 // Default constructor.
61 // A default constructed AliITSsimulationSPDdubna class.
63 AliDebug(1,Form("Calling degault constructor"));
65 //______________________________________________________________________
66 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(AliITSsegmentation *seg,
67 AliITSresponse *resp,Int_t cup):
68 AliITSsimulation(seg,resp),
72 // standard constructor
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
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.
86 // A default constructed AliITSsimulationSPDdubna class.
89 Form("Calling degault constructor seg=%p resp=%p cup=%d",seg,resp,cup));
90 if(cup==1||cup==2){ // For the moment, remove defusion if Coupling is
92 resp->SetTemperature(0.0);
93 resp->SetDistanceOverVoltage(0.0);
97 //______________________________________________________________________
98 void AliITSsimulationSPDdubna::Init(){
106 const Double_t kmictocm = 1.0e-4; // convert microns to cm.
110 SetMap(new AliITSpList(GetNPixelsZ(),GetNPixelsX()));
111 GetResp(0,0)->SetDistanceOverVoltage(kmictocm*GetSeg()->Dy(),50.0);
113 //______________________________________________________________________
114 AliITSsimulationSPDdubna::~AliITSsimulationSPDdubna(){
128 //______________________________________________________________________
129 AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(const
130 AliITSsimulationSPDdubna
131 &s) : AliITSsimulation(s){
134 // AliITSsimulationSPDdubna &s The original class for which
135 // this class is a copy of
143 //______________________________________________________________________
144 AliITSsimulationSPDdubna& AliITSsimulationSPDdubna::operator=(const
145 AliITSsimulationSPDdubna &s){
146 // Assignment operator
148 // AliITSsimulationSPDdubna &s The original class for which
149 // this class is a copy of
154 if(&s == this) return *this;
156 fCoupling = s.fCoupling;
157 fSPDname = s.fSPDname;
160 //______________________________________________________________________
161 AliITSsimulation& AliITSsimulationSPDdubna::operator=(const
162 AliITSsimulation &s){
163 // Assignment operator
165 // AliITSsimulationSPDdubna &s The original class for which
166 // this class is a copy of
171 if(&s == this) return *this;
172 Error("AliITSsimulationSPDdubna","Not allowed to make a = with "
173 "AliITSsimulationSPDdubna","Using default creater instead");
177 //______________________________________________________________________
178 void AliITSsimulationSPDdubna::InitSimulationModule(Int_t module, Int_t event){
179 // This function creates maps to build the list of tracks for each
180 // summable digit. Inputs defined by base class.
182 // Int_t module // Module number to be simulated
183 // Int_t event // Event number to be simulated
189 AliDebug(1,Form("(module=%d,event=%d)",module,event));
190 SetModuleNumber(module);
191 SetEventNumber(event);
194 //_____________________________________________________________________
195 void AliITSsimulationSPDdubna::SDigitiseModule(AliITSmodule *mod,Int_t,
197 // This function begins the work of creating S-Digits. Inputs defined
200 // AliITSmodule *mod // module
202 // Int_t event // Event number
206 // test // test returns kTRUE if the module contained hits
207 // // test returns kFALSE if it did not contain hits
209 AliDebug(1,Form("(mod=%p, ,event=%d)",mod,event));
210 if(!(mod->GetNhits())){
211 AliDebug(1,Form("In event %d module %d there are %d hits returning.",
212 event, mod->GetIndex(),mod->GetNhits()));
213 return;// if module has no hits don't create Sdigits
215 SetModuleNumber(mod->GetIndex());
216 SetEventNumber(event);
221 //______________________________________________________________________
222 void AliITSsimulationSPDdubna::WriteSDigits(){
223 // This function adds each S-Digit to pList
230 Int_t ix, nix, iz, niz;
231 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
233 AliDebug(1,Form("Writing SDigits for module %d",GetModuleNumber()));
234 GetMap()->GetMaxMapIndex(niz, nix);
235 for(iz=0; iz<niz; iz++)for(ix=0; ix<nix; ix++){
236 if(GetMap()->GetSignalOnly(iz,ix)>0.0){
237 aliITS->AddSumDigit(*(GetMap()->GetpListItem(iz,ix)));
238 if(AliDebugLevel()>0) {
239 AliDebug(1,Form("%d, %d",iz,ix));
240 cout << *(GetMap()->GetpListItem(iz,ix)) << endl;
242 } // end if GetMap()->GetSignalOnly(iz,ix)>0.0
246 //______________________________________________________________________
247 void AliITSsimulationSPDdubna::FinishSDigitiseModule(){
248 // This function calls SDigitsToDigits which creates Digits from SDigits
257 pListToDigits(); // Charge To Signal both adds noise and
261 //______________________________________________________________________
262 void AliITSsimulationSPDdubna::DigitiseModule(AliITSmodule *mod,Int_t,
264 // This function creates Digits straight from the hits and then adds
265 // electronic noise to the digits before adding them to pList
266 // Each of the input variables is passed along to HitToSDigit
268 // AliITSmodule *mod module
276 AliDebug(1,Form("(mod=%p,,)",mod));
281 //______________________________________________________________________
282 void AliITSsimulationSPDdubna::HitToSDigit(AliITSmodule *mod){
283 // Does the charge distributions using Gaussian diffusion charge charing.
285 // AliITSmodule *mod Pointer to this module
290 const Double_t kmictocm = 1.0e-4; // convert microns to cm.
291 TObjArray *hits = mod->GetHits();
292 Int_t nhits = hits->GetEntriesFast();
295 Double_t x0=0.0,x1=0.0,y0=0.0,y1=0.0,z0=0.0,z1=0.0,de=0.0;
296 Double_t x,y,z,t,tp,st,dt=0.2,el,sig;
297 Double_t thick = kmictocm*GetSeg()->Dy();
299 AliDebug(1,Form("(mod=%p) fCoupling=%d",mod,fCoupling));
301 for(h=0;h<nhits;h++){
302 if(AliDebugLevel()>0) {
303 AliDebug(1,Form("Hits, %d", h));
304 cout << *(mod->GetHit(h)) << endl;
306 if(!mod->LineSegmentL(h,x0,x1,y0,y1,z0,z1,de,idtrack)) continue;
307 st = TMath::Sqrt(x1*x1+y1*y1+z1*z1);
309 st = (Double_t)((Int_t)(st/kmictocm)); // number of microns
310 if(st<=1.0) st = 1.0;
312 for(t=0.0;t<1.0;t+=dt){ // Integrate over t
317 if(!(GetSeg()->LocalToDet(x,z,ix,iz))) continue; // outside
318 el = GetResp(ix,iz)->GeVToCharge((Double_t)(dt*de));
320 if(el<=0.0) cout<<"el="<<el<<" dt="<<dt
323 sig = GetResp(ix,iz)->SigmaDiffusion1D(thick + y);
324 SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
326 } else { // st == 0.0 deposit it at this point
330 if(!(GetSeg()->LocalToDet(x,z,ix,iz))) continue; // outside
331 el = GetResp(ix,iz)->GeVToCharge((Double_t)de);
332 sig = GetResp(ix,iz)->SigmaDiffusion1D(thick + y);
333 SpreadCharge(x,z,ix,iz,el,sig,idtrack,h);
340 // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
341 for(i=0;i<GetMap()->GetEntries();i++)
342 if(GetMap()->GetpListItem(i)==0) continue;
344 GetMap()->GetMapIndex(
345 GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
346 SetCoupling(iz,ix,idtrack,h);
350 // x is column and z is row (see AliITSsegmentationSPD::GetPadIxz)
351 for(i=0;i<GetMap()->GetEntries();i++)
352 if(GetMap()->GetpListItem(i)==0) continue;
354 GetMap()->GetMapIndex(
355 GetMap()->GetpListItem(i)->GetIndex(),iz,ix);
356 SetCouplingOld(iz,ix,idtrack,h);
360 } // Loop over all hits h
361 if(GetDebug(2))Info("HitToSDigit","Finished fCoupling=%d",fCoupling);
363 //______________________________________________________________________
364 void AliITSsimulationSPDdubna::SpreadCharge(Double_t x0,Double_t z0,
366 Double_t el,Double_t sig,Int_t t,
368 // Spreads the charge over neighboring cells. Assume charge is distributed
369 // as charge(x,z) = (el/2*pi*sig*sig)*exp(-arg)
370 // arg=((x-x0)*(x-x0)/2*sig*sig)+((z-z0*z-z0)/2*sig*sig)
371 // Defined this way, the integral over all x and z is el.
373 // Double_t x0 x position of point where charge is liberated
374 // Double_t y0 y position of point where charge is liberated
375 // Double_t z0 z position of point where charge is liberated
376 // Int_t ix0 row of cell corresponding to point x0
377 // Int_t iz0 columb of cell corresponding to point z0
378 // Double_t el number of electrons liberated in this step
379 // Double_t sig Sigma difusion for this step (y0 dependent)
380 // Int_t t track number
381 // Int_t ti hit track index number
382 // Int_t hi hit "hit" index number
387 const Int_t knx = 3,knz = 2;
388 const Double_t kRoot2 = 1.414213562; // Sqrt(2).
389 const Double_t kmictocm = 1.0e-4; // convert microns to cm.
390 Int_t ix,iz,ixs,ixe,izs,ize;
392 Double_t x1,x2,z1,z2,s,sp;
394 if(GetDebug(4)) Info("SpreadCharge","(x0=%e,z0=%e,ix0=%d,iz0=%d,el=%e,"
395 "sig=%e,t=%d,i=%d)",x0,z0,ix0,iz0,el,sig,t,hi);
396 if(sig<=0.0) { // if sig<=0 No diffusion to simulate.
397 GetMap()->AddSignal(iz0,ix0,t,hi,GetModuleNumber(),el);
399 cout << "sig<=0.0=" << sig << endl;
403 sp = 1.0/(sig*kRoot2);
405 cout << "sig=" << sig << " sp=" << sp << endl;
407 ixs = TMath::Max(-knx+ix0,0);
408 ixe = TMath::Min(knx+ix0,GetSeg()->Npx()-1);
409 izs = TMath::Max(-knz+iz0,0);
410 ize = TMath::Min(knz+iz0,GetSeg()->Npz()-1);
411 for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){
412 GetSeg()->DetToLocal(ix,iz,x,z); // pixel center
415 x2 = x1 + 0.5*kmictocm*GetSeg()->Dpx(ix); // Upper
416 x1 -= 0.5*kmictocm*GetSeg()->Dpx(ix); // Lower
417 z2 = z1 + 0.5*kmictocm*GetSeg()->Dpz(iz); // Upper
418 z1 -= 0.5*kmictocm*GetSeg()->Dpz(iz); // Lower
419 x1 -= x0; // Distance from where track traveled
420 x2 -= x0; // Distance from where track traveled
421 z1 -= z0; // Distance from where track traveled
422 z2 -= z0; // Distance from where track traveled
423 s = 0.25; // Correction based on definision of Erfc
424 s *= TMath::Erfc(sp*x1) - TMath::Erfc(sp*x2);
426 cout <<"el="<<el<<" ix0="<<ix0<<" ix="<<ix<<" x0="<<x<<
427 " iz0="<<iz0<<" iz="<<iz<<" z0="<<z<<
428 " sp*x1="<<sp*x1<<" sp*x2="<<sp*x2<<" s="<<s;
430 s *= TMath::Erfc(sp*z1) - TMath::Erfc(sp*z2);
432 cout<<" sp*z1="<<sp*z1<<" sp*z2="<<sp*z2<<" s="<<s<< endl;
434 GetMap()->AddSignal(iz,ix,t,hi,GetModuleNumber(),s*el);
437 //______________________________________________________________________
438 void AliITSsimulationSPDdubna::pListToDigits(){
439 // add noise and electronics, perform the zero suppression and add the
447 static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS");
449 Double_t electronics;
451 const Int_t nmaxtrk=AliITSdigitSPD::GetNTracks();
452 static AliITSdigitSPD dig;
454 if(GetDebug(1)) Info("pListToDigits","()");
455 for(iz=0; iz<GetNPixelsZ(); iz++) for(ix=0; ix<GetNPixelsX(); ix++){
456 // Apply Noise/Dead channals and the like
457 if(GetResp(ix,iz)->IsPixelDead(GetModuleNumber(),ix,iz)) continue;
458 electronics = GetResp(ix,iz)->ApplyBaselineAndNoise();
459 UpdateMapNoise(ix,iz,electronics);
461 // Apply Threshold and write Digits.
462 sig = GetMap()->GetSignalOnly(iz,ix);
463 FillHistograms(ix,iz,sig+electronics);
465 cout<<sig<<"+"<<electronics<<">threshold("<<ix<<","<<iz
466 <<")="<<GetThreshold(ix,iz) <<endl;
468 if (sig+electronics <= GetThreshold(ix,iz)) continue;
472 dig.SetSignalSPD((Int_t) GetMap()->GetSignal(iz,ix));
473 for(j=0;j<nmaxtrk;j++){
474 if (j<GetMap()->GetNEnteries()) {
475 dig.SetTrack(j,GetMap()->GetTrack(iz,ix,j));
476 dig.SetHit(j,GetMap()->GetHit(iz,ix,j));
477 }else { // Default values
483 cout<<iz<<","<<ix<<","<<*(GetMap()->GetpListItem(iz,ix))<<endl;
485 aliITS->AddSimDigit(0,&dig);
488 //______________________________________________________________________
489 void AliITSsimulationSPDdubna::CreateHistograms(){
490 // create 1D histograms for tests
498 if(GetDebug(1)) Info("CreateHistograms","create histograms");
500 fHis = new TObjArray(GetNPixelsZ());
501 TString fSPDname("spd_");
502 for(Int_t i=0;i<GetNPixelsZ();i++) {
504 sprintf(pixelz,"%d",i);
505 fSPDname.Append(pixelz);
506 fHis->AddAt(new TH1F(fSPDname.Data(),"SPD maps",
507 GetNPixelsX(),0.,(Double_t)GetNPixelsX()),i);
510 //______________________________________________________________________
511 void AliITSsimulationSPDdubna::FillHistograms(Int_t ix,Int_t iz,Double_t v){
512 // Fill the histogram
520 if(!GetHistArray()) return; // Only fill if setup.
521 if(GetDebug(2)) Info("FillHistograms","fill histograms");
522 GetHistogram(iz)->Fill(ix,v);
524 //______________________________________________________________________
525 void AliITSsimulationSPDdubna::ResetHistograms(){
526 // Reset histograms for this detector
534 if(!GetHistArray()) return; // Only fill if setup.
535 if(GetDebug(2)) Info("FillHistograms","fill histograms");
536 for ( int i=0;i<GetNPixelsZ();i++ ) {
537 if (fHis->At(i)) ((TH1F*)fHis->At(i))->Reset();
541 //______________________________________________________________________
542 void AliITSsimulationSPDdubna::SetCoupling(Int_t row, Int_t col, Int_t ntrack,
544 // Take into account the coupling between adiacent pixels.
545 // The parameters probcol and probrow are the probability of the
546 // signal in one pixel shared in the two adjacent pixels along
547 // the column and row direction, respectively.
548 // Note pList is goten via GetMap() and module is not need any more.
549 // Otherwise it is identical to that coded by Tiziano Virgili (BSN).
552 <img src="picts/ITS/barimodel_3.gif">
555 <font size=+2 color=red>
556 <a href="mailto:tiziano.virgili@cern.ch"></a>.
562 // Int_t row z cell index
563 // Int_t col x cell index
564 // Int_t ntrack track incex number
565 // Int_t idhit hit index number
571 Double_t pulse1,pulse2;
572 Double_t couplR=0.0,couplC=0.0;
575 GetCouplings(couplR,couplC);
576 if(GetDebug(3)) Info("SetCoupling","(row=%d,col=%d,ntrack=%d,idhit=%d) "
577 "Calling SetCoupling couplR=%e couplC=%e",
578 row,col,ntrack,idhit,couplR,couplC);
581 pulse1 = GetMap()->GetSignalOnly(row,col);
583 for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
587 xr = gRandom->Rndm();
588 //if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold(j1,col))){
589 if ((j1<0) || (j1>GetNPixelsZ()-1) || (xr>couplR)){
593 UpdateMapSignal(col,j1,ntrack,idhit,pulse1);
595 flag = 1; // only first next!!
598 // loop in column direction
602 xr = gRandom->Rndm();
603 //if((j2<0)||j2>(GetNPixelsX()-1)||pulse2<GetThreshold(row,j2)){
604 if ((j2<0) || (j2>GetNPixelsX()-1) || (xr>couplC)){
608 UpdateMapSignal(j2,row,ntrack,idhit,pulse2);
610 flag = 1; // only first next!!
615 //______________________________________________________________________
616 void AliITSsimulationSPDdubna::SetCouplingOld(Int_t row, Int_t col,
617 Int_t ntrack,Int_t idhit) {
618 // Take into account the coupling between adiacent pixels.
619 // The parameters probcol and probrow are the fractions of the
620 // signal in one pixel shared in the two adjacent pixels along
621 // the column and row direction, respectively.
624 <img src="picts/ITS/barimodel_3.gif">
627 <font size=+2 color=red>
628 <a href="mailto:Rocco.Caliandro@ba.infn.it"></a>.
634 // Int_t row z cell index
635 // Int_t col x cell index
636 // Int_t ntrack track incex number
637 // Int_t idhit hit index number
638 // Int_t module module number
644 Double_t pulse1,pulse2;
645 Double_t couplR=0.0,couplC=0.0;
647 GetCouplings(couplR,couplC);
648 if(GetDebug(3)) Info("SetCouplingOld","(row=%d,col=%d,ntrack=%d,idhit=%d) "
649 "Calling SetCoupling couplR=%e couplC=%e",
650 row,col,ntrack,idhit,couplR,couplC);
653 pulse1 = GetMap()->GetSignalOnly(row,col);
655 for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
659 if ((j1<0)||(j1>GetNPixelsZ()-1)||(pulse1<GetThreshold(j1,col))){
660 pulse1 = GetMap()->GetSignalOnly(row,col);
664 UpdateMapSignal(col,j1,ntrack,idhit,pulse1);
668 // loop in column direction
672 if((j2<0)||(j2>(GetNPixelsX()-1))||(pulse2<GetThreshold(row,j2))){
673 pulse2 = GetMap()->GetSignalOnly(row,col);
677 UpdateMapSignal(j2,row,ntrack,idhit,pulse2);