]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ITS/AliITSsimulationFastPoints.cxx
bug fix
[u/mrichter/AliRoot.git] / ITS / AliITSsimulationFastPoints.cxx
index 9a51c3e12a9a606a0b2cf7314eded75067934e86..84afdd3d2b650c8a04cb8b58ca1d02e49026ba5d 100644 (file)
  * provided "as is" without express or implied warranty.                  *
  **************************************************************************/
 
-/*
-$Log$
-Revision 1.1.2.1  2000/06/11 20:16:05  barbera
-New: Fast simulation class for the ITS, class as part of new ITS code
-structure.
+/* $Id$ */
+//////////////////////////////////////////////////////////
+// implements fast simulation                           //
+//                                                      //
+//                                                      //
+//////////////////////////////////////////////////////////
 
-*/
 
-#include <TParticle.h>
+#include <TRandom.h>
+
 #include "AliITS.h"
+#include "AliITShit.h"
+#include "AliITSRecPoint.h"
+#include "AliITSmodule.h"
+#include "AliITSgeom.h"
+#include "AliRun.h"
 #include "AliITSsimulationFastPoints.h"
-#include "AliITSstatistics.h"
+
 
 ClassImp(AliITSsimulationFastPoints)
 
 AliITSsimulationFastPoints::AliITSsimulationFastPoints()
 {
   //constructor
-  fSx = new AliITSstatistics(2);
-  fSz = new AliITSstatistics(2);
+  fSigmaRPhi[0] = fSigmaRPhi[1] = 12e-4;
+  fSigmaRPhi[2] = fSigmaRPhi[3] = 38e-4;
+  fSigmaRPhi[4] = fSigmaRPhi[5] = 20e-4;
+  fSigmaZ[0] = fSigmaZ[1] = 120e-4;        // resolution for 425 micron pixels
+  fSigmaZ[2] = fSigmaZ[3] = 28e-4;
+  fSigmaZ[4] = fSigmaZ[5] = 830e-4;
+  fSigmaDe[0] = fSigmaDe[1] = 0.72e-6;
+  fSigmaDe[2] = fSigmaDe[3] = 0.90e-6;
+  fSigmaDe[4] = fSigmaDe[5] =  5e-6;
+  fThrDe[0] = fThrDe[1] = 7.2e-6;
+  fThrDe[2] = fThrDe[3] = 2.70e-6;
+  fThrDe[4] = fThrDe[5] = 10e-6;
 }
 
-//----------------------------------------------------------
-AliITSsimulationFastPoints::~AliITSsimulationFastPoints()
-{
-  //destructor
-  delete fSx;
-  delete fSz;
+//-------------------------------------------------------------
+void AliITSsimulationFastPoints::CreateFastRecPoints(Int_t module, TClonesArray* recp){
+    // Fast points simulator
+    AliITS *aliITS  = (AliITS*)gAlice->GetModule("ITS");
 
+    CreateFastRecPoints((AliITSmodule *)(aliITS->GetModule(module)),
+                       module,gRandom,recp);
 }
-
 //-------------------------------------------------------------
-void AliITSsimulationFastPoints::CreateFastRecPoints(AliITSmodule *mod){
-  // Fast points simulator for all of the ITS.
-  Int_t   nhit,h,trk,ifirst;
-  Float_t x,y,z,t,e;// local coordinate (cm) and time of flight, and dedx.
-  Float_t x1,y1,z1;
-  AliITShit *hit;
-
-  fSx->Reset(); // Start out with things clearly zeroed
-  fSz->Reset(); // Start out with things clearly zeroed
-  e = 0.; // Start out with things clearly zeroed
-  Double_t weight=1.;
-  nhit = mod->GetNhits();
-  ifirst = 1;
-  for(h=0;h<nhit;h++){
-    hit = mod->GetHit(h);
-    hit->GetPositionL(x,y,z,t);
-    if(ifirst) {x1=x;y1=y;z1=z;}
-    e += hit->GetIonization();
-    trk = hit->GetTrack();
-    fSx->AddValue((Double_t)x,weight);
-    fSz->AddValue((Double_t)z,weight);
-    ifirst = 0;
-    if(hit->StatusExiting()||  // leaving volume
-       hit->StatusDisappeared()|| // interacted/decayed...
-       hit->StatusStop() // dropped below E cuts.
-       ){ // exiting track, write out RecPoint.
-      //      if(fSz->GetRMS()>1.E-1) {
-      //       TParticle *part = hit->GetParticle();
-      //       printf("idpart %d energy %f \n",part->GetPdgCode(),part->Energy());
-      //       printf("diffx=%e diffy=%e diffz=%e\n",x-x1,y-y1,z-z1);
-      //      }
-      switch (mod->GetLayer()){
-      case 1: case 2:  // SPDs
-       AddSPD(e,mod,trk);
-       break;
-      case 3: case 4:  // SDDs
-       AddSDD(e,mod,trk);
-       break;
-      case 5: case 6:  // SSDs
-       AddSSD(e,mod,trk);
-       break;
-      } // end switch
-      fSx->Reset();
-      fSz->Reset();
-      e = 0.;
-      ifirst = 1;
-      continue;
-    }// end if
-  } // end for h
+void AliITSsimulationFastPoints::CreateFastRecPoints(AliITSmodule *mod,
+                                                    Int_t module,
+                                                    TRandom *random,
+                                                    TClonesArray* recp) {
+  // Fast points simulator 
+
+  TClonesArray &pt=*recp;
+  AliITS *aliITS  = (AliITS*)gAlice->GetModule("ITS");
+  AliITSgeom *gm = aliITS->GetITSgeom();
+  const Float_t kdEdXtoQ = 1.0e+6;  // GeV->KeV
+
+  Int_t lay,lad,det;
+  gm->GetModuleId(module,lay,lad,det);
+  Int_t ind=(lad-1)*gm->GetNdetectors(lay)+(det-1);
+  Int_t lyr=(lay-1);
+
+
+  Int_t ihit,flag,numofhits;
+  Float_t locals[3];
+  Float_t globals[3];
+  Double_t sigmarphi=0., sigmaz=0., sigmade=0., thrde=0.;
+  Float_t deltaXl,deltaZl,deltaDe;
+
+  Int_t hitlay, hitlad, hitdet, hitstatus;
+  Float_t hitpx, hitpy, hitpz, hitdestep;
+
+  Int_t   hitstatus1, hittrack1;
+  Float_t hitx1, hity1, hitz1;
+  Float_t hitdestep1;
+  Float_t xMg,yMg,zMg;
+  Int_t irecp=0;
+  numofhits = mod->GetNhits();
+  //printf("numofhits %d \n",numofhits);
+  for(ihit=0;ihit<numofhits;ihit++){
+    AliITShit *hit=mod->GetHit(ihit);
+    hit->GetPositionG(hitx1,hity1,hitz1);
+    hitstatus1 = hit->GetTrackStatus();
+    hitdestep1 = hit->GetIonization();
+    hittrack1 = hit->GetTrack();
+    
+    mod->MedianHit(module,hitx1,hity1,hitz1,hitstatus1,xMg,yMg,zMg,flag);
+    if (flag!=1) {
+      hitdestep = hit->GetIonization();
+      
+      if (hitdestep > 0) {
+       hit->GetDetectorID(hitlay,hitlad,hitdet);
+       hit->GetMomentumG(hitpx,hitpy,hitpz);            
+       hitstatus = hitstatus1;
+               // Transform to the module local frame
+       globals[0] = xMg; 
+       globals[1] = yMg;
+       globals[2] = zMg;
+       gm->GtoL(hitlay,hitlad,hitdet,globals,locals);
+       // Retrieve sigma values for position and energy, and energy
+       // threshold
+       sigmarphi = SigmaRPhi(hitlay);
+       sigmaz = SigmaZ(hitlay);
+       sigmade = SigmaDe(hitlay);
+       thrde = ThrDe(hitlay);
+       deltaXl = random->Gaus(0,sigmarphi);
+       deltaZl = random->Gaus(0,sigmaz);
+       deltaDe = random->Gaus(0,sigmade);
+       
+       // Apply energy threshold and trasform back to global reference
+       // system
+       
+       if ( (hitdestep+deltaDe) > thrde ){
+         locals[0] += deltaXl;
+         locals[2] += deltaZl;
+         Int_t lab[4] = {hit->GetTrack(),-3,-3,ind};
+         Float_t q=kdEdXtoQ*(hitdestep+deltaDe);
+         if(hitlay<3) q=1.; // SPD binary readout
+         Float_t hitv[6] = {locals[0],locals[2],sigmarphi*sigmarphi,sigmaz*sigmaz,q,q};
+         Int_t info[3] = {0,0,lyr};
+         AliITSRecPoint rp(lab,hitv,info,kTRUE);
+
+         new (pt[irecp]) AliITSRecPoint(rp);
+         irecp++;
+       } // end if ( (hitdestep+deltaDe)
+      } // end if (hitdestep > 0)
+    } // end if (flag!=1)
+  } // end for ihit
 }
 //_______________________________________________________________________
-void AliITSsimulationFastPoints::AddSPD(Float_t &e,
-                                        AliITSmodule *mod,Int_t trackNumber){
-  const Float_t kcmTomicron = 1.0e4;
-  //  const Float_t kdEdXtoQ = ;
-  const Float_t kRMSx = 12.0; // microns ITS TDR Table 1.3
-  const Float_t kRMSz = 70.0; // microns ITS TDR Table 1.3
-  Float_t a1,a2; // general float.
-  AliITSRecPoint rpSPD;
-  Int_t *trk = rpSPD.GetTracks();
-
-  trk[0] = trackNumber;
-  trk[1] = 0; trk[2] = 0;
-  rpSPD.SetX(kcmTomicron*fSx->GetMean());
-  rpSPD.SetZ(kcmTomicron*fSz->GetMean());
-  rpSPD.SetdEdX(0.0);
-  rpSPD.SetQ(1.0);
-  a1 = kcmTomicron*fSx->GetRMS(); a1 *= a1; a1 += kRMSx*kRMSx;
-  //  if(a1>1.E5) printf("addSPD: layer=%d track #%d dedx=%e sigmaX2= %e ",
-  //               mod->GetLayer(),trackNumber,e,a1);
-  rpSPD.SetSigmaX2(a1);
-  a2 = kcmTomicron*fSz->GetRMS(); a2 *= a2; a2 += kRMSz*kRMSz;
-  //  if(a1>1.E5) printf(" sigmaZ2= %e\n",a2);
-  rpSPD.SetSigmaZ2(a2);
-  rpSPD.SetProbability(1.0);
-
-  (mod->GetITS())->AddRecPoint(rpSPD);
+void AliITSsimulationFastPoints::SetSigmaRPhi(Double_t  srphi[6])
+{
+  // set sigmas in rphi
+
+    Int_t i;
+    for (i=0; i<6; i++) {
+       fSigmaRPhi[i]=srphi[i];
+    }
 }
 //_______________________________________________________________________
-void AliITSsimulationFastPoints::AddSDD(Float_t &e,
-                                        AliITSmodule *mod,Int_t trackNumber){
-
-  const Float_t kcmTomicron = 1.0e4;
-  const Float_t kdEdXtoQ = 2.778e+8; // Boris Batyuna June 10 2000.
-  const Float_t kRMSx = 38.0; // microns ITS TDR Table 1.3
-  const Float_t kRMSz = 28.0; // microns ITS TDR Table 1.3
-  Float_t a1,a2; // general float.
-  AliITSRecPoint rpSDD;
-  Int_t *trk = rpSDD.GetTracks();
-
-  trk[0] = trackNumber;
-  trk[1] = 0; trk[2] = 0;
-  rpSDD.SetX(kcmTomicron*fSx->GetMean());
-  rpSDD.SetZ(kcmTomicron*fSz->GetMean());
-  rpSDD.SetdEdX(e);
-  rpSDD.SetQ(kdEdXtoQ*e);
-  a1 = kcmTomicron*fSx->GetRMS(); a1 *= a1; a1 += kRMSx*kRMSx;
-  //  if(a1>1.E5) printf("addSDD: layer=%d track #%d dedx=%e sigmaX2= %e ",
-  //               mod->GetLayer(),trackNumber,e,a1);
-  rpSDD.SetSigmaX2(a1);
-  a2 = kcmTomicron*fSz->GetRMS(); a2 *= a2; a2 += kRMSz*kRMSz;
-  //  if(a1>1.E5) printf(" sigmaZ2= %e\n",a2);
-  rpSDD.SetSigmaZ2(a2);
-  rpSDD.SetProbability(1.0);
-
-  (mod->GetITS())->AddRecPoint(rpSDD);
+void AliITSsimulationFastPoints::SetSigmaZ(Double_t  sz[6])
+{
+  // set sigmas in z
+
+    Int_t i;
+    for (i=0; i<6; i++) {
+       fSigmaZ[i]=sz[i];
+    }
 }
 //_______________________________________________________________________
-void AliITSsimulationFastPoints::AddSSD(Float_t &e,
-                                        AliITSmodule *mod,Int_t trackNumber){
-
-  const Float_t kcmTomicron = 1.0e4;
-  const Float_t kdEdXtoQ = 2.778e+8; // Boris Batyuna June 10 2000.
-  const Float_t kRMSx = 20.0; // microns ITS TDR Table 1.3
-  const Float_t kRMSz = 830.0; // microns ITS TDR Table 1.3
-  Float_t a1,a2; // general float.
-  AliITSRecPoint rpSSD;
-  Int_t *trk = rpSSD.GetTracks();
-
-  trk[0] = trackNumber;
-  trk[1] = 0; trk[2] = 0;
-  rpSSD.SetX(kcmTomicron*fSx->GetMean());
-  rpSSD.SetZ(kcmTomicron*fSz->GetMean());
-  rpSSD.SetdEdX(e);
-  rpSSD.SetQ(kdEdXtoQ*e);
-  a1 = kcmTomicron*fSx->GetRMS(); a1 *= a1; a1 += kRMSx*kRMSx;
-  //  if(a1>1.E5) printf("addSSD: layer=%d track #%d dedx=%e sigmaX2= %e ",
-  //               mod->GetLayer(),trackNumber,e,a1);
-  rpSSD.SetSigmaX2(a1);
-  a2 = kcmTomicron*fSz->GetRMS(); a2 *= a2; a2 += kRMSz*kRMSz;
-  //  if(a1>1.E5) printf(" sigmaZ2= %e RMSx=%e RMSz=%e\n",a2,fSx->GetRMS(),fSz->GetRMS());
-  rpSSD.SetSigmaZ2(a2);
-  rpSSD.SetProbability(1.0);
-
-  (mod->GetITS())->AddRecPoint(rpSSD);
+void AliITSsimulationFastPoints::SetSigmaDe(Double_t  sde[6])
+{
+  // set sigmas in energy
+
+    Int_t i;
+    for (i=0; i<6; i++) {
+       fSigmaDe[i]=sde[i];
+    }
 }
 //_______________________________________________________________________
+void AliITSsimulationFastPoints::SetThrDe(Double_t  thrde[6])
+{
+  // set energy thersholds
+
+    Int_t i;
+    for (i=0; i<6; i++) {
+       fThrDe[i]=thrde[i];
+    }
+}
+