/[PAMELA software]/PamVMC_update/s4/src/PamVMCS4Dig.cxx
ViewVC logotype

Annotation of /PamVMC_update/s4/src/PamVMCS4Dig.cxx

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1.1.1 - (hide annotations) (download) (vendor branch)
Tue Oct 15 15:52:32 2013 UTC (11 years, 1 month ago) by formato
Branch: MAIN, rel
CVS Tags: reltag, HEAD
Changes since 1.1: +0 -0 lines
PamVMC update

1 formato 1.1 #include "PamVMCS4Dig.h"
2     #include <TRandom.h>
3    
4     ClassImp(PamVMCS4Dig)
5    
6    
7     void PamVMCS4Dig::Digitize(){
8    
9     #ifdef DIG_DEBUG
10     cout<<"Starting S4 Digitization...."<<endl;
11     #endif
12     //This is simple approximation, reflection + atennuation..
13     //I hope in future we will use for generation and propagation of
14     //optical photons abilitiees of montecarlo...
15    
16    
17    
18    
19     Int_t DEBUG=0;
20     // creato: S. Borisov, INFN Roma2 e MEPHI, Sett 2007
21     Int_t i,j,t,NdF,pmt,S4;
22     Float_t E0,X,Y,Z,x,y,z,V[3],q,w,l;
23    
24     const Float_t Xs[2] = {-24.1, 24.1};
25     const Float_t Ys[2] = {-24.1, 24.1};
26     const Float_t Zs[2] = {-0.5, 0.5};
27     const Float_t Yp[6] = {-20., -17., -1., 2., 17., 20.}; //PMT positions
28    
29     const Float_t E1 = 1e-6; //Energy which we spend to produce 1 photon
30     const Float_t p = 0.1;
31     const Int_t l0 = 500;
32    
33     const Int_t S4p = 32; //pedestal
34    
35     UShort_t S4v = 0; // ADC count
36    
37     Int_t NdFT = 0;
38    
39     Float_t Ert = 0.;
40    
41     TClonesArray * hc = (TClonesArray*)fhitscolmap.GetValue("S4");
42     PamVMCDetectorHit * hit = 0;
43    
44     if(hc){
45     for(i=0;i<hc->GetEntriesFast();i++){
46    
47     hit = (PamVMCDetectorHit*)hc->At(i);
48    
49     Ert += hit->GetEREL();
50     NdF=Int_t(hit->GetEREL()/E1);
51     NdFT=0;
52     X=hit->GetXIN();
53     Y=hit->GetYIN();
54     Z=(Float_t)(random())/(Float_t)(0x7fffffff)-0.5;
55     for(j=0;j<NdF;j++){
56     q=(Float_t)random()/(Float_t)0x7fffffff;
57     w=(Float_t)random()/(Float_t)0x7fffffff;
58     V[0]=p*cos(6.28318*q);
59     V[1]=p*sin(6.28318*q);
60     V[2]=p*(2.*w-1.);
61     pmt=0;
62     x=X;
63     y=Y;
64     z=Z;
65     while(pmt==0 && (x>Xs[0] && x<Xs[1])&&(y>Ys[0] && y<Ys[1])&&(z>Zs[0] && z<Zs[1])){
66     l=0;
67     while(pmt==0 && (x>Xs[0] && x<Xs[1])&&(y>Ys[0] && y<Ys[1])&&(z>Zs[0] && z<Zs[1])){
68     x+=V[0];
69     y+=V[1];
70     z+=V[2];
71     l+=p;
72     }
73     if((x<Xs[0]+p || x>Xs[1]-p)&&(y>Ys[0]+p && y<Ys[1]-p)&&(z>Zs[0]+p && z<Zs[1]-p)){
74     for(t=0;t<3;t++){
75     if(y>=Yp[2*t] && y<Yp[2*t+1]){
76     if(pmt==0)NdFT++;
77     pmt=1;
78     break;
79     }
80     }
81     if(pmt==1)break;
82     V[0]=-V[0];
83     }
84     q=(Float_t)random()/(Float_t)0x7fffffff;
85     w=1-exp(-l/l0);
86     if(q<w)break;
87     q=(Float_t)random()/(Float_t)0x7fffffff;
88     w=0.5;
89     if(q<w)break;
90     if((x>Xs[0]+p && x<Xs[1]-p)&&(y<Ys[0]+p || y>Ys[1]-p)&&(z>Zs[0]+p && z<Zs[1]-p))V[1]=-V[1];
91     if((x>Xs[0]+p && x<Xs[1]-p)&&(y>Ys[0]+p && y<Ys[1]-p)&&(z<Zs[0]+p || z>Zs[1]-p))V[2]=-V[2];
92     x+=V[0];
93     y+=V[1];
94     z+=V[2];
95     l=0;
96     }
97     }
98     }
99     Ert=Ert/0.002;
100     q=(Float_t)(random())/(Float_t)0x7fffffff;
101     w=0.7;
102     E0=4064./7.;
103     if(Ert<1) S4=0;
104     else S4=(Int_t)(4064.*(1.-exp(-(Ert-1.)/E0)));
105     i=S4/4;
106     if(S4%4==0)
107     S4v=S4+S4p;
108     else if(S4%4==1){
109     if(q<w) S4v=S4-1+S4p;
110     else S4v=S4+1+S4p;
111     } else if(S4%4==2) S4v=S4+S4p;
112     else if(S4%4==3){
113     if(q<w) S4v=S4+1+S4p;
114     else S4v=S4-1+S4p;
115     }
116     if (DEBUG)
117     cout<<"Ert_S4 = " << Ert << " --- S4v = " << S4v << endl;
118    
119     fData.clear();
120    
121     fData.push_back(S4v);
122     fData.push_back(0xd800);
123     fData.push_back(0x0300);
124     }
125     }

  ViewVC Help
Powered by ViewVC 1.1.23