/************************************************************************** * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ /* $Id$ */ #include #include #include #include #include #include #include "AliRun.h" #include "AliITS.h" #include "AliITShit.h" #include "AliITSdigit.h" #include "AliITSmodule.h" #include "AliITSMapA2.h" #include "AliITSpList.h" #include "AliITSsimulationSPDdubna.h" #include "AliITSsegmentationSPD.h" #include "AliITSresponseSPDdubna.h" //#define DEBUG ClassImp(AliITSsimulationSPDdubna) //////////////////////////////////////////////////////////////////////// // Version: 1 // Modified by Bjorn S. Nilsen // Version: 0 // Written by Boris Batyunya // December 20 1999 // // AliITSsimulationSPDdubna is to do the simulation of SPDs. //______________________________________________________________________ AliITSsimulationSPDdubna::AliITSsimulationSPDdubna() : AliITSsimulation(){ // Default constructor. // Inputs: // none. // Outputs: // none. // Return: // A default constructed AliITSsimulationSPDdubna class. fHis = 0; } //______________________________________________________________________ AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(AliITSsegmentation *seg, AliITSresponse *resp){ // standard constructor // Inputs: // AliITSsegmentation *seg A pointer to the segmentation class // to be used for this simulation // AliITSresponse *resp A pointer to the responce class to // be used for this simulation // Outputs: // none. // Return: // A default constructed AliITSsimulationSPDdubna class. fHis = 0; Init(seg,resp); } //______________________________________________________________________ void AliITSsimulationSPDdubna::Init(AliITSsegmentation *seg, AliITSresponse *resp){ // Initilization // Inputs: // AliITSsegmentation *seg A pointer to the segmentation class // to be used for this simulation // AliITSresponse *resp A pointer to the responce class to // be used for this simulation // Outputs: // none. // Return: // none. const Double_t kmictocm = 1.0e-4; // convert microns to cm. SetResponseModel(resp); SetSegmentationModel(seg); SetModuleNumber(0); SetEventNumber(0); SetMap(new AliITSpList(GetNPixelsZ(),GetNPixelsX())); GetResp()->SetDistanceOverVoltage(kmictocm*GetSeg()->Dy(),50.0); } //______________________________________________________________________ AliITSsimulationSPDdubna::~AliITSsimulationSPDdubna(){ // destructor // Inputs: // none. // Outputs: // none. // Return: // none. if (fHis) { fHis->Delete(); delete fHis; } // end if fHis } //______________________________________________________________________ AliITSsimulationSPDdubna::AliITSsimulationSPDdubna(const AliITSsimulationSPDdubna &s) : AliITSsimulation(s){ // Copy Constructor // Inputs: // AliITSsimulationSPDdubna &s The original class for which // this class is a copy of // Outputs: // none. // Return: *this = s; return; } //______________________________________________________________________ AliITSsimulationSPDdubna& AliITSsimulationSPDdubna::operator=(const AliITSsimulationSPDdubna &s){ // Assignment operator // Inputs: // AliITSsimulationSPDdubna &s The original class for which // this class is a copy of // Outputs: // none. // Return: if(&s == this) return *this; this->fHis = s.fHis; return *this; } //______________________________________________________________________ void AliITSsimulationSPDdubna::InitSimulationModule(Int_t module, Int_t event){ // This function creates maps to build the list of tracks for each // summable digit. Inputs defined by base class. // Inputs: // Int_t module // Module number to be simulated // Int_t event // Event number to be simulated // Outputs: // none // Returns: // none SetModuleNumber(module); SetEventNumber(event); ClearMap(); } //_____________________________________________________________________ void AliITSsimulationSPDdubna::SDigitiseModule(AliITSmodule *mod, Int_t mask, Int_t event){ // This function begins the work of creating S-Digits. Inputs defined // by base class. // Inputs: // AliITSmodule *mod // module // Int_t mask // mask to be applied to the module // Int_t event // Event number // Outputs: // none // Return: // test // test returns kTRUE if the module contained hits // // test returns kFALSE if it did not contain hits mask = mod->GetNhits(); if(!mask) return;// if module has no hits don't create Sdigits SetModuleNumber(mod->GetIndex()); SetEventNumber(event); HitToSDigit(mod); WriteSDigits(); ClearMap(); } //______________________________________________________________________ void AliITSsimulationSPDdubna::WriteSDigits(){ // This function adds each S-Digit to pList // Inputs: // none. // Outputs: // none. // Return: // none Int_t ix, nix, iz, niz; static AliITS *aliITS = (AliITS*)gAlice->GetModule("ITS"); GetMap()->GetMaxMapIndex(niz, nix); for(iz=0; izGetSignalOnly(iz,ix)>0.0){ aliITS->AddSumDigit(*(GetMap()->GetpListItem(iz,ix))); #ifdef DEBUG cout <<"SDigits " << iz << "," << ix << "," << *(GetMap()->GetpListItem(iz,ix)) << endl; #endif } // end if GetMap()->GetSignalOnly(iz,ix)>0.0 } // end for iz,ix return; } //______________________________________________________________________ void AliITSsimulationSPDdubna::FinishSDigitiseModule(){ // This function calls SDigitsToDigits which creates Digits from SDigits // Inputs: // none // Outputs: // none // Return // none SDigitsToDigits(); return; } //______________________________________________________________________ void AliITSsimulationSPDdubna::SDigitsToDigits(){ // This function adds electronic noise to the S-Digits and then adds them // to a new pList // Inputs: // none. // Outputs: // none. // Return: // none ChargeToSignal(); // Charge To Signal both adds noise and ClearMap(); } //______________________________________________________________________ void AliITSsimulationSPDdubna::DigitiseModule(AliITSmodule *mod, Int_t module, Int_t dummy){ // This function creates Digits straight from the hits and then adds // electronic noise to the digits before adding them to pList // Each of the input variables is passed along to HitToSDigit // Inputs: // AliITSmodule *mod module // Int_t module module number Dummy. // Int_t dummy // Outputs: // none. // Return: // none. //This calls the module for HitToSDigit fModule = dummy = module = mod->GetIndex(); HitToSDigit(mod); ChargeToSignal(); ClearMap(); } //______________________________________________________________________ void AliITSsimulationSPDdubna::UpdateMapSignal(Int_t iz,Int_t ix,Int_t trk, Int_t ht,Double_t signal){ // This function adds a signal to the pList from the pList class // Inputs: // Int_t iz // row number // Int_t ix // column number // Int_t trk // track number // Int_t ht // hit number // Double_t signal // signal strength // Outputs: // none. // Return: // none GetMap()->AddSignal(iz,ix,trk,ht,GetModuleNumber(),signal); } //______________________________________________________________________ void AliITSsimulationSPDdubna::UpdateMapNoise(Int_t iz,Int_t ix,Float_t noise){ // This function adds noise to data in the MapA2 as well as the pList // Inputs: // Int_t iz // row number // Int_t ix // column number // Int_t mod // module number // Double_t sig // signal strength // Double_t noise // electronic noise generated by ChargeToSignal // Outputs: // All of the inputs are passed to AliITSMapA2::AddSignal or // AliITSpList::AddNoise // Return: // none GetMap()->AddNoise(iz,ix,GetModuleNumber(),noise); } //______________________________________________________________________ void AliITSsimulationSPDdubna::HitToDigit(AliITSmodule *mod){ // Standard interface to DigitiseModule // Inputs: // AliITSmodule *mod Pointer to this module // Outputs: // none. // Return: // none. DigitiseModule(mod,GetModuleNumber(),0); } //______________________________________________________________________ void AliITSsimulationSPDdubna::HitToSDigit(AliITSmodule *mod){ // Does the charge distributions using Gaussian diffusion charge charing. // Inputs: // AliITSmodule *mod Pointer to this module // Output: // none. // Return: // none. const Double_t kmictocm = 1.0e-4; // convert microns to cm. TObjArray *hits = mod->GetHits(); Int_t nhits = hits->GetEntriesFast(); Int_t h,ix,iz; Int_t idtrack; Double_t x0=0.0,x1=0.0,y0=0.0,y1=0.0,z0=0.0,z1=0.0,de=0.0; Double_t x,y,z,t,tp,st,dt=0.2,el,sig; Double_t thick = kmictocm*GetSeg()->Dy(); if(nhits<=0) return; for(h=0;hGetHit(h)) << endl; #endif if(mod->LineSegmentL(h,x0,x1,y0,y1,z0,z1,de,idtrack)){ st = TMath::Sqrt(x1*x1+y1*y1+z1*z1); if(st>0.0){ st = (Double_t)((Int_t)(1.0E+04*st)); // number of microns if(st<=0.0) st = 1.0; dt = 1.0/st; for(t=0;t<1.0;t+=dt){ // Integrate over t tp = t+0.5*dt; el = GetResp()->GeVToCharge((Float_t)(dt*de)); #ifdef DEBUG if(el<=0.0) cout << "el="<AddSignal(iz0,ix0,t,hi,GetModuleNumber(),el); return; } // end if sp = 1.0/(sig*kRoot2); #ifdef DEBUG cout << "sig=" << sig << " sp=" << sp << endl; #endif ixs = TMath::Max(-knx+ix0,0); ixe = TMath::Min(knx+ix0,GetSeg()->Npx()-1); izs = TMath::Max(-knz+iz0,0); ize = TMath::Min(knz+iz0,GetSeg()->Npz()-1); for(ix=ixs;ix<=ixe;ix++) for(iz=izs;iz<=ize;iz++){ GetSeg()->DetToLocal(ix,iz,x,z); // pixel center x1 = x; z1 = z; x2 = x1 + 0.5*kmictocm*GetSeg()->Dpx(ix); // Upper x1 -= 0.5*kmictocm*GetSeg()->Dpx(ix); // Lower z2 = z1 + 0.5*kmictocm*GetSeg()->Dpz(iz); // Upper z1 -= 0.5*kmictocm*GetSeg()->Dpz(iz); // Lower x1 -= x0; // Distance from where track traveled x2 -= x0; // Distance from where track traveled z1 -= z0; // Distance from where track traveled z2 -= z0; // Distance from where track traveled s = 0.25; // Correction based on definision of Erfc s *= TMath::Erfc(sp*x1) - TMath::Erfc(sp*x2); #ifdef DEBUG cout << "el=" << el << " ix0=" << ix0 << " ix=" << ix << " x0="<< x << " iz0=" << iz0 << " iz=" << iz << " z0=" << z << " sp*x1=" << sp*x1 <<" sp*x2=" << sp*x2 << " s=" << s; #endif s *= TMath::Erfc(sp*z1) - TMath::Erfc(sp*z2); #ifdef DEBUG cout << " sp*z1=" << sp*z1 <<" sp*z2=" << sp*z2 << " s=" << s << endl; #endif GetMap()->AddSignal(iz,ix,t,hi,GetModuleNumber(),s*el); } // end for ix, iz } //______________________________________________________________________ /* Double_t *AliITSsimulationSPDdubna::CreateFindCellEdges(Double_t x0, Double_t x1,Double_t z0,Double_t z1,Int_t &n){ // Note: This function is a potensial source for a memory leak. The memory // pointed to in its return, must be deleted. // Inputs: // Double_t x0 The starting location of the track step in x // Double_t x1 The distance allong x for the track step // Double_t z0 The starting location of the track step in z // Double_t z1 The distance allong z for the track step // Output: // Int_t &n The size of the array returned. Minimal n=2. // Return: // The pointer to the array of track steps. Int_t ix0,ix1,ix,iz0,iz1,iz,i; Double_t x,z,lx,ux,lz,uz,a,b,c,d; Double_t *t; GetSeg()->LocalToDet(x0,z0,ix0,iz0); GetSeg()->LocalToDet(x1,z1,ix1,iz1); n = 2 + TMath::Abs(ix1-ix0) + TMath::Abs(iz1-iz0); t = new Double_t[n]; t[0] = 0.0; t[n-1] = 1.0; x = x0; z = z0; for(i=1;iLocalToDet(x,z,ix,iz); GetSeg()->CellBoundries(ix,iz,lx,ux,lz,uz); a = (lx-x0)/x1; if(a<=t[i-1]) a = 1.0; b = (ux-x0)/x1; if(b<=t[i-1]) b = 1.0; c = (lz-z0)/z1; if(c<=t[i-1]) c = 1.0; d = (uz-z0)/z1; if(d<=t[i-1]) d = 1.0; t[i] = TMath::Min(TMath::Min(TMath::Min(a,b),c),d); x = x0+x1*(t[i]*1.00000001); z = z0+z1*(t[i]*1.00000001); i++; } // end for i return t; } */ //______________________________________________________________________ void AliITSsimulationSPDdubna::HitToSDigitOld(AliITSmodule *mod){ // digitize module Old method // Inputs: // AliITSmodule *mod Pointer to this module // Output: // none. // Return: // none. const Float_t kEnToEl = 2.778e+8; // GeV->charge in electrons // for 3.6 eV/pair const Float_t kconv = 10000.; // cm -> microns Float_t spdLength = GetSeg()->Dz(); Float_t spdWidth = GetSeg()->Dx(); Float_t spdThickness = GetSeg()->Dy(); Float_t difCoef=GetResp()->SigmaDiffusion1D(1.0); Float_t zPix0 = 1e+6; Float_t xPix0 = 1e+6; Float_t yPrev = 1e+6; Float_t zPitch = GetSeg()->Dpz(0); Float_t xPitch = GetSeg()->Dpx(0); // Array of pointers to the label-signal list Int_t indexRange[4] = {0,0,0,0}; // Fill detector maps with GEANT hits // loop over hits in the module static Bool_t first; Int_t lasttrack=-2; Int_t hit, iZi, jz, jx; Int_t idhit=-1; //! TObjArray *fHits = mod->GetHits(); Int_t nhits = fHits->GetEntriesFast(); Float_t yPix0 = -spdThickness/2; if (!nhits) return; #ifdef DEBUG cout<<"len,wid,thickness,nx,nz,pitchx,pitchz,difcoef ="<At(hit); #ifdef DEBUG cout << "Hits=" << hit << "," << *iHit << endl; #endif //Int_t layer = iHit->GetLayer(); // work with the idtrack=entry number in the TreeH //Int_t idhit,idtrack; //! //mod->GetHitTrackAndHitIndex(hit,idtrack,idhit); //! //Int_t idtrack=mod->GetHitTrackIndex(hit); // or store straight away the particle position in the array // of particles : if(iHit->StatusEntering()) idhit=hit; Int_t itrack = iHit->GetTrack(); Int_t dray = 0; if (lasttrack != itrack || hit==(nhits-1)) first = kTRUE; //Int_t parent = iHit->GetParticle()->GetFirstMother(); Int_t partcode = iHit->GetParticle()->GetPdgCode(); // partcode (pdgCode): 11 - e-, 13 - mu-, 22 - gamma, 111 - pi0, // 211 - pi+, 310 - K0s, 321 - K+, 2112 - n, 2212 - p, 3122 - lambda Float_t pmod = iHit->GetParticle()->P(); // total momentum at the // vertex pmod *= 1000; if(partcode == 11 && pmod < 6) dray = 1; // delta ray is e- // at p < 6 MeV/c // Get hit z and x(r*phi) cordinates for each module (detector) // in local system. Float_t zPix = kconv*iHit->GetZL(); Float_t xPix = kconv*iHit->GetXL(); Float_t yPix = kconv*iHit->GetYL(); // Get track status Int_t status = iHit->GetTrackStatus(); // Check boundaries if(zPix > spdLength/2) { #ifdef DEBUG cout<<"!!! SPD: z outside ="<GetThreshold(); Int_t j,ix,iz; Float_t electronics; Double_t sig; const Int_t nmaxtrk=AliITSdigitSPD::GetNTracks(); static AliITSdigitSPD dig; for(iz=0; izIsPixelDead(GetModuleNumber(),ix,iz)) continue; electronics = GetResp()->ApplyBaselineAndNoise(); sig = GetMap()->GetSignalOnly(iz,ix); UpdateMapNoise(iz,ix,electronics); #ifdef DEBUG cout << sig << "+" << electronics <<">threshold=" << threshold << endl; #endif if (sig+electronics <= threshold) continue; dig.SetCoord1(iz); dig.SetCoord2(ix); dig.SetSignal(1); dig.SetSignalSPD((Int_t) GetMap()->GetSignal(iz,ix)); for(j=0;jGetNEnteries()) { dig.SetTrack(j,GetMap()->GetTrack(iz,ix,j)); dig.SetHit(j,GetMap()->GetHit(iz,ix,j)); }else { // Default values dig.SetTrack(j,-3); dig.SetHit(j,-1); } // end if GetMap() } // end for j #ifdef DEBUG cout<GetpListItem(iz,ix)) << endl; #endif aliITS->AddSimDigit(0,&dig); } // for ix/iz } //______________________________________________________________________ void AliITSsimulationSPDdubna::CreateHistograms(){ // create 1D histograms for tests // Inputs: // none. // Outputs: // none. // Return: // none. printf("SPD - create histograms\n"); fHis=new TObjArray(GetNPixelsZ()); TString spdName("spd_"); for (Int_t i=0;iAddAt(new TH1F(spdName.Data(),"SPD maps", GetNPixelsX(),0.,(Float_t) GetNPixelsX()), i); } // end for i } //______________________________________________________________________ void AliITSsimulationSPDdubna::ResetHistograms(){ // Reset histograms for this detector // Inputs: // none. // Outputs: // none. // Return: // none. for ( int i=0;iAt(i)) ((TH1F*)fHis->At(i))->Reset(); } // end for i } //______________________________________________________________________ void AliITSsimulationSPDdubna::SetCoupling(Int_t row, Int_t col, Int_t ntrack, Int_t idhit) { // Take into account the coupling between adiacent pixels. // The parameters probcol and probrow are the probability of the // signal in one pixel shared in the two adjacent pixels along // the column and row direction, respectively. //Begin_Html /*
.
    */
    //End_Html
    // Inputs:
    //    Int_t row            z cell index
    //    Int_t col            x cell index
    //    Int_t ntrack         track incex number
    //    Int_t idhit          hit index number
    //    Int_t module         module number
    // Outputs:
    //    none.
    // Return:
    //     none.
    Int_t j1,j2,flag=0;
    Double_t pulse1,pulse2;
    Float_t couplR=0.0,couplC=0.0;
    Double_t xr=0.;
    AliITSpList *pList = GetMap();

    GetCouplings(couplR,couplC);
    j1 = row;
    j2 = col;
    pulse1 = pList->GetSignalOnly(row,col);
    pulse2 = pulse1;
    for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
      do{
	j1 += isign;
	//   pulse1 *= couplR; 
	xr = gRandom->Rndm();
	//   if ((j1<0) || (j1>GetNPixelsZ()-1) || (pulse1GetNPixelsZ()-1) || (xr>couplR)){
	  j1 = row;
	  flag = 1;
	}else{
	  UpdateMapSignal(j1,col,ntrack,idhit,pulse1);
	  //  flag = 0;
	  flag = 1; // only first next!!
	} // end if
      } while(flag == 0);
      // loop in column direction
      do{
	j2 += isign;
	// pulse2 *= couplC; 
	xr = gRandom->Rndm();
	//  if ((j2<0) || (j2>(GetNPixelsX()-1)) || (pulse2GetNPixelsX()-1) || (xr>couplC)){
	  j2 = col;
	  flag = 1;
	}else{
	  UpdateMapSignal(row,j2,ntrack,idhit,pulse2);
	  //  flag = 0;
	  flag = 1; // only first next!!
	} // end if
      } while(flag == 0);
    } // for isign
}
//______________________________________________________________________
void AliITSsimulationSPDdubna::SetCouplingOld(Int_t row, Int_t col,
                Int_t ntrack,Int_t idhit) {
    //  Take into account the coupling between adiacent pixels.
    //  The parameters probcol and probrow are the fractions of the
    //  signal in one pixel shared in the two adjacent pixels along
    //  the column and row direction, respectively.
    //Begin_Html
    /*
      
      

.
    */
    //End_Html
    // Inputs:
    //    Int_t row            z cell index
    //    Int_t col            x cell index
    //    Int_t ntrack         track incex number
    //    Int_t idhit          hit index number
    //    Int_t module         module number
    // Outputs:
    //    none.
    // Return:
    //     none.
    Int_t j1,j2,flag=0;
    Double_t pulse1,pulse2;
    Float_t couplR=0.0,couplC=0.0;
    AliITSpList *pList = GetMap();

    GetCouplings(couplR,couplC);
    j1 = row;
    j2 = col;
    pulse1 = pList->GetSignalOnly(row,col);
    pulse2 = pulse1;
    for (Int_t isign=-1;isign<=1;isign+=2){// loop in row direction
      do{
	j1 += isign;
	pulse1 *= couplR;
	if ((j1<0) || (j1>GetNPixelsZ()-1) || (pulse1GetSignalOnly(row,col);
	  j1 = row;
	  flag = 1;
	}else{
	  UpdateMapSignal(j1,col,ntrack,idhit,pulse1);
	  flag = 0;
	} // end if
      } while(flag == 0);
      // loop in column direction
      do{
	j2 += isign;
	pulse2 *= couplC;
	if ((j2<0) || (j2>(GetNPixelsX()-1)) || (pulse2GetSignalOnly(row,col);
	  j2 = col;
	  flag = 1;
	}else{
	  UpdateMapSignal(row,j2,ntrack,idhit,pulse2);
	  flag = 0;
	} // end if
      } while(flag == 0);
    } // for isign
}