]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Adding digitize function (Marian)
authormarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 24 Feb 2009 11:39:23 +0000 (11:39 +0000)
committermarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 24 Feb 2009 11:39:23 +0000 (11:39 +0000)
TPC/fastSimul/AliTPCclusterFast.cxx

index bd5871fedf9260ebca3146eedb955aa94ce36f1f..b3051856819483786344690c85fd4a6fbb798cbd 100644 (file)
@@ -1,17 +1,29 @@
 /*
   gSystem->Load("libSTAT.so");
-.x ~/UliStyle.C
+  .x ~/NimStyle.C
  
  .L $ALICE_ROOT/TPC/fastSimul/AliTPCclusterFast.cxx+
- AliTPCclusterFast::Simul("aaa.root",100000);
-  gSystem->Load("libSTAT.so");
+ //
+ AliTPCclusterFast::fPRF = new TF1("fprf","gausn");
+ AliTPCclusterFast::fTRF = new TF1("ftrf","gausn");
+ AliTPCclusterFast::fPRF->SetParameters(1,0,0.5);
+ AliTPCclusterFast::fTRF->SetParameters(1,0,0.5);
+ //
+
+ AliTPCclusterFast::Simul("aaa.root",50000); 
+ gSystem->Load("libSTAT.so");
+
+ TFile f("aaa.root");
+ TTree * tree = (TTree*)f.Get("simul");
 
 */
 
 #include "TObject.h"
+#include "TF1.h"
 #include "TMath.h"
 #include "TRandom.h"
 #include "TVectorD.h"
+#include "TMatrixD.h"
 #include "TTreeStream.h"
 
 class AliTPCclusterFast: public TObject {
@@ -20,33 +32,49 @@ public:
   virtual ~AliTPCclusterFast();
   void SetParam(Float_t mnprim, Float_t diff, Float_t y, Float_t z, Float_t ky, Float_t kz);
   void GenerElectrons();
+  void Digitize();
+  Double_t GetNsec();
   static void Simul(const char* simul, Int_t npoints);
 public:
   Float_t fMNprim;     // mean number of primary electrons
   //                   //electrons part input
-  Float_t fNprim;      // mean number of primary electrons
+  Int_t   fNprim;      // mean number of primary electrons
   Int_t   fNtot;       // total number of primary electrons 
   Float_t fDiff;       // diffusion sigma
   Float_t fY;          // y position 
   Float_t fZ;          // z postion 
   Float_t fAngleY;     // y angle - tan(y)
   Float_t fAngleZ;     // z angle - tan z
+  //
+  //
   //                   // electron part simul
+  TVectorD fSec;       // number of secondary electrons
   TVectorD fPosY;      //! position y for each electron
   TVectorD fPosZ;      //! position z for each electron
   TVectorD fGain;      //! gg for each electron
-  TVectorD fStatY;     // stat Y
+  //
+  TVectorD fStatY;     // stat Y  
   TVectorD fStatZ;     // stat Y
   //
-  //  
+  // digitization part
+  //
+  TMatrixD fDigits;    // response matrix
+  static TF1* fPRF;    // Pad response
+  static TF1* fTRF;    // Time response function 
   ClassDef(AliTPCclusterFast,1)  // container for
 };
 ClassImp(AliTPCclusterFast)
 
+TF1 *AliTPCclusterFast::fPRF=0;
+TF1 *AliTPCclusterFast::fTRF=0;
 
 
 AliTPCclusterFast::AliTPCclusterFast(){
+  //
+  //
+  fDigits.ResizeTo(5,7);
 }
+
 AliTPCclusterFast::~AliTPCclusterFast(){
 }
 
@@ -58,24 +86,33 @@ void AliTPCclusterFast::SetParam(Float_t mnprim, Float_t diff, Float_t y, Float_
   fY=y; fZ=z; 
   fAngleY=ky; fAngleZ=kz;
 }
-
-void AliTPCclusterFast::GenerElectrons(){
-  //
+Double_t AliTPCclusterFast::GetNsec(){
   //
+  // Generate number of secondary electrons
+  // copy of procedure implemented in geant
   //
-  const Int_t knMax=1000;
   const Double_t FPOT=20.77E-9, EEND=10E-6, EEXPO=2.2, EEND1=1E-6;
   const Double_t XEXPO=-EEXPO+1, YEXPO=1/XEXPO;
   const Double_t W=20.77E-9;
+  Float_t RAN = gRandom->Rndm();
+  return TMath::Nint(TMath::Power((TMath::Power(FPOT,XEXPO)*(1-RAN)+TMath::Power(EEND,XEXPO)*RAN),YEXPO)/W);
+}
+
+void AliTPCclusterFast::GenerElectrons(){
   //
+  //
+  //
+  //
+  const Int_t knMax=1000;
   if (fPosY.GetNrows()<knMax){
     fPosY.ResizeTo(knMax);
     fPosZ.ResizeTo(knMax);
     fGain.ResizeTo(knMax);
+    fSec.ResizeTo(knMax);
     fStatY.ResizeTo(3);
     fStatZ.ResizeTo(3);
   }
-  fNprim =-fMNprim*TMath::Log(gRandom->Rndm());  //number of primary electrons
+  fNprim = gRandom->Poisson(fMNprim);  //number of primary electrons
   fNtot=0;
   //
   Double_t sumQ=0;
@@ -83,15 +120,15 @@ void AliTPCclusterFast::GenerElectrons(){
   Double_t sumZQ=0;
   Double_t sumY2Q=0;
   Double_t sumZ2Q=0;
-
-
+  for (Int_t i=0;i<knMax;i++){ 
+    fSec[i]=0;
+  }
   for (Int_t iprim=0; iprim<fNprim;iprim++){
-    Float_t RAN = gRandom->Rndm();
-    //      Float_t dE   =  (FPOT**XEXPO*(1-RAN)+EEND**XEXPO*RAN)**YEXPO;
-    Float_t dN   =  TMath::Power((TMath::Power(FPOT,XEXPO)*(1-RAN)+TMath::Power(EEND,XEXPO)*RAN),YEXPO)/W;
+    Float_t dN   =  GetNsec();
+    fSec[iprim]=dN;
     Double_t yc = fY+(gRandom->Rndm()-0.5)*fAngleY;
     Double_t zc = fZ+(gRandom->Rndm()-0.5)*fAngleZ;
-    for (Int_t isec=0;isec<dN;isec++){
+    for (Int_t isec=0;isec<=dN;isec++){
       //
       //
       Double_t y = gRandom->Gaus(0,fDiff)+yc;
@@ -120,12 +157,34 @@ void AliTPCclusterFast::GenerElectrons(){
   }
 }
 
+void AliTPCclusterFast::Digitize(){
+  //
+  //
+  //
+  // 1. Clear digits
+  for (Int_t i=0; i<5;i++)
+    for (Int_t j=0; j<7;j++){
+      fDigits(i,j)=0;
+    }
+  //
+  // Fill digits
+  for (Int_t iel = 0; iel<fNtot; iel++){
+    for (Int_t di=-2; di<=2;di++)
+      for (Int_t dj=-2; dj<=2;dj++){
+       Float_t fac = fPRF->Eval(di-fPosY[iel])*fTRF->Eval(dj-fPosZ[iel]);
+       fac*=fGain[iel];
+       fDigits(2+di,3+dj)+=fac;
+      }
+  }
+}
+
+
 
 void AliTPCclusterFast::Simul(const char* fname, Int_t npoints){
   AliTPCclusterFast fast;
   TTreeSRedirector cstream(fname);
   for (Int_t icl=0; icl<npoints; icl++){
-    Float_t nprim=(10+20*gRandom->Rndm());
+    Float_t nprim=(10+40*gRandom->Rndm());
     Float_t diff =0.01 +0.3*gRandom->Rndm();
     Float_t posY = gRandom->Rndm()-0.5;
     Float_t posZ = gRandom->Rndm()-0.5;
@@ -133,6 +192,7 @@ void AliTPCclusterFast::Simul(const char* fname, Int_t npoints){
     Float_t kz   = 1.*(gRandom->Rndm()-0.5);
     fast.SetParam(nprim,diff,posY,posZ,ky,kz);
     fast.GenerElectrons();
+    fast.Digitize();
     cstream<<"simul"<<
       "s.="<<&fast<<
       "\n";
@@ -141,6 +201,6 @@ void AliTPCclusterFast::Simul(const char* fname, Int_t npoints){
 
 
 /*
-
+  TH2F *hisL = new TH2F("hisL","hisL",10,10,50,100,0,10)
 
 */