/[PAMELA software]/PamelaDigitizer/DigitizeTOF.cxx
ViewVC logotype

Contents of /PamelaDigitizer/DigitizeTOF.cxx

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.8 - (show annotations) (download)
Fri Jan 17 13:08:40 2014 UTC (10 years, 10 months ago) by mocchiut
Branch: MAIN
CVS Tags: HEAD
Changes since 1.7: +2 -2 lines
Compilation warnings using GCC4.7 fixed

1 #include "Digitizer.h"
2
3 void Digitizer::DigitizeTOF(int np,float *atte1,float *atte2,float *lambda1,float *lambda2){
4 //fDataTof: 12 x 23 bytes (=276 bytes)
5 UChar_t *pTof=fDataTof;
6 Bool_t DEBUG=false;
7
8 Int_t cdp[75] = {0,1,1,0,1,1,0,1,1,0,1,1,0,1,1, //0-14
9 0,0,0,1,0,1,0,1,1,0,0,1,0,1,0, //15-29
10 1,1,1,1,2,2,2,3,3,3,3,4,4,4,1, //30-44
11 1,2,0,2,0,0,5,5,5,5,6,6,6,6,7, //45-59
12 3,3,4,4,5,5,6,7,8,9,10,11,12,13,14 }; //60-74
13
14 int Z = cdp[Ipa-1];
15
16 float time_res[8] = {425.,210.,170.,130.,120.,120.,120.,120.};
17 for(Int_t i=0;i<8;i++)time_res[i]/=1.4;//1.17;1.5;1.3*/
18 Float_t dt1 = 0.;// = 1.e-12*time_res[0]; // single PMT resolution for Z=1 (WM, Nov'07)
19
20 if ((Z > 1) && (Z < 9)) dt1=1.e-12*time_res[(Z-1)];
21 if (Z > 8) dt1=120.e-12;
22
23
24 // ------ evaluate energy in each pmt: ------
25 // strip geometry (lenght/width)
26 Float_t dimel[6] = {33.0, 40.8 ,18.0, 15.0, 15.0, 18.0};
27 // S11 8 paddles 33.0 x 5.1 cm
28 // S12 6 paddles 40.8 x 5.5 cm
29 // S21 2 paddles 18.0 x 7.5 cm
30 // S22 2 paddles 15.0 x 9.0 cm
31 // S31 3 paddles 15.0 x 6.0 cm
32 // S32 3 paddles 18.0 x 5.0 cm
33 Float_t FGeo[2]={0., 0.}; /* geometrical factor */
34 const Float_t Pho_keV = 10.; // photons per keV in scintillator
35 const Float_t echarge = 1.6e-19; // electron charge
36 Float_t Npho=0.;
37 Float_t QevePmt_pC[48];
38 Float_t QhitPad_pC[2]={0.,0.};
39 Float_t QhitPmt_pC[2]={0.,0.};
40 Float_t pmGain = 3.5e6; /* PMT Gain: the same for all PMTs */
41 Float_t effi=0.21; /* Efficienza di fotocatodo */
42 // pC < 800
43 Float_t ADC_pC0A =-4.437616e+01;
44 Float_t ADC_pC1A = 1.573329e+00;
45 Float_t ADC_pC2A = 2.780518e-04;
46 Float_t ADC_pC3A =-2.302160e-07;
47 // pC > 800:
48 Float_t ADC_pC0B =-2.245756e+02;
49 Float_t ADC_pC1B = 2.184156e+00;
50 Float_t ADC_pC2B =-4.171825e-04;
51 Float_t ADC_pC3B = 3.789715e-08;
52
53 Float_t pCthres=40.; // threshold in charge
54 Int_t ADClast=4095; // no signal --> ADC ch=4095
55 Int_t ADCsat=3100; // saturation value for the ADCs
56 Int_t ADCtof[48];
57 Float_t ScaleFact[48]={0.39, 0.49, 0.38, 0.40, 0.65, 0.51, 0.43, 0.49,
58 0.58, 0.38, 0.53, 0.57, 0.53, 0.45, 0.49, 0.22,
59 0.21, 0.44, 0.28, 0.57, 0.26, 0.72, 0.37, 0.29,
60 0.30, 0.89, 0.37, 0.12, 0.27, 0.23, 0.15, 0.22,
61 0.19, 0.20, 0.21, 0.19, 0.41, 0.32, 0.39, 0.38,
62 0.28, 0.66, 0.28, 0.40, 0.39, 0.40, 0.37, 0.35};//15:0.7--0.95, 16:0.9--1.25, 27:0.9--1.3, 30:0.9--1.15, 32:0.85--1.05, 33:0.85--1.05
63 for(Int_t i=0; i<48; i++){
64 QevePmt_pC[i] = 0;
65 ADCtof[i]=0;
66 }
67 Int_t ip,ipad,pmtleft=0,pmtright=0;
68 // TDC variables:
69 Int_t TDClast=4095,TDCint[48];
70 Float_t tdc[48],tdc1[48],tdcpmt[48];
71 for(Int_t i=0; i<48; i++) {
72 tdcpmt[i] = 1000.;
73 tdc[i] = 0.; // 18-oct WM
74 tdc1[i] = 0.; // 18-oct WM
75 }
76 Float_t thresh=20.; // to be defined better... (Wolfgang)
77 // === TDC: simulate timing for each paddle
78 Float_t tdcres[50],c1_S[50],c2_S[50],c3_S[50];
79 for(Int_t j=0;j<48;j++){
80 tdcres[j] = 50.E-12; // TDC resolution 50 picosec
81 c1_S[j] = 500.; // cable length in channels
82 c2_S[j] = 0.;
83 c3_S[j] = 1000.;
84 c1_S[j] = c1_S[j]*tdcres[j]; // cable length in sec
85 c2_S[j] = c2_S[j]*tdcres[j];
86 }
87 /* ********************************** start loop over hits */
88 if(Nthtof>*ntof)cout<<"NTHTOF > "<<*ntof<<" , event rejected ! "<<Nthtof<<endl;
89 else{
90 for(Int_t nh=0; nh<Nthtof; nh++){
91 Float_t s_l_g[6] = {8.0, 8.0, 20.9, 22.0, 9.8, 8.3 }; // length of the lightguide
92 Float_t t1,t2,veff,veff1,veff0 ;
93 veff0 = 100.*1.0e8 ; // light velocity in the scintillator in m/sec
94 veff1 = 100.*1.5e8; // light velocity in the lightguide in m/sec
95 veff=veff0; // signal velocity in the paddle
96 t1 = Timetof[nh] ; // Start
97 t2 = Timetof[nh] ;
98 // Donatella: redefinition plane and pad for vectors in C
99 ip = Ipltof[nh]-1;
100 ipad = Ipaddle[nh]-1;
101 pmtleft=0;
102 pmtright=0;
103 // WM: S12 paddles are "reversed" (Nov'07)
104 if (ip==2){
105 if (ipad==0)
106 ipad=1;
107 else
108 ipad=0;
109 }
110 if ((ip>-1)&&(ip<6)) { //ToF paddles only, not S4
111 Paddle2Pmt(ip, ipad, &pmtleft, &pmtright);
112 // DC: evaluates mean position and path inside the paddle
113 Float_t tpos=0.;
114 Float_t path[2] = {0., 0.};
115 //--- Strip in Y = S11,S22,S31 ------
116 if(ip==0 || ip==3 || ip==4)
117 tpos = (Yintof[nh]+Youttof[nh])/2.;
118 else
119 if(ip==1 || ip==2 || ip==5) //--- Strip in X for S12,S21,S32
120 tpos = (Xintof[nh]+Xouttof[nh])/2.;
121 else //if (ip!=6)
122 printf("*** WARNING TOF: this option should never occur! (ip=%2i, nh=%2i)\n",ip,nh);
123 path[0]= tpos + dimel[ip]/2.; // path to left PMT
124 path[1]= dimel[ip]/2.- tpos; // path to right PMT
125 if (DEBUG) {
126 cout <<" plane "<<ip<<" strip # ="<< ipad <<" tpos "<< tpos <<"\n";
127 cout <<"pmtleft, pmtright "<<pmtleft<<" "<<pmtright<<endl;
128 }
129 // constant geometric factor, the rest is handled by the scaling factor
130 FGeo[0] =0.5;
131 FGeo[1] =0.5;
132 Npho = Ereltof[nh]*Pho_keV*1.0e6; // Eloss in GeV
133
134 Float_t knorm[2]={0., 0.}; // Donatella
135 Float_t Atten[2]={0., 0.}; // Donatella
136 for(Int_t j=0; j<2; j++){
137 QhitPad_pC[j]= Npho*FGeo[j]*effi*pmGain*echarge*1.E12*ScaleFact[pmtleft+j];
138 // WM
139 knorm[j]=atte1[pmtleft+j]*exp(lambda1[pmtleft+j]*dimel[ip]/2.*pow(-1,j+1)) + atte2[pmtleft+j]*exp(lambda2[pmtleft+j]*dimel[ip]/2.*pow(-1,j+1));
140 Atten[j]=atte1[pmtleft+j]*exp(tpos*lambda1[pmtleft+j]) + atte2[pmtleft+j]*exp(tpos*lambda2[pmtleft+j]) ;
141 QhitPmt_pC[j]= QhitPad_pC[j]*Atten[j]/knorm[j];
142 if (DEBUG) {
143 cout<<"pmtleft "<<pmtleft<<" j "<<j<<endl;
144 cout<<" atte1 "<<atte1[pmtleft+j]<<"lambda1 "<<lambda1[pmtleft+j]<<" atte2 "<<atte2[pmtleft+j]<<"lambda2 "<<lambda2[pmtleft+j] <<endl;
145 cout<<j<<" tpos "<<tpos<<" knorm "<<knorm[j]<<" "<<Atten[j]<<" "<<"QhitPmt_pC "<<QhitPmt_pC[j]<<endl;
146 }
147 }
148 if(DEBUG)cout<<"Npho "<<Npho<<" QhitPmt_pC "<<QhitPmt_pC[0]<<" "<<QhitPmt_pC[1]<<endl;
149 QevePmt_pC[pmtleft] += QhitPmt_pC[0];
150 QevePmt_pC[pmtright] += QhitPmt_pC[1];
151 //TDC
152 // WM right and left <->
153 t1 = t1 + fabs(path[0]/veff) + s_l_g[ip]/veff1;
154 t2 = t2 + fabs(path[1]/veff) + s_l_g[ip]/veff1 ; // Signal reaches PMT
155 t1 = gRandom->Gaus(t1,dt1); //apply gaussian error dt
156 t2 = gRandom->Gaus(t2,dt1); //apply gaussian error dt
157 t1 = t1 + c1_S[pmtleft] ; // Signal reaches Discriminator ,TDC starts to run
158 t2 = t2 + c1_S[pmtright] ;
159 // check if signal is above threshold
160 // then check if tdcpmt is already filled by another hit...
161 // only re-fill if time is smaller
162 if (QhitPmt_pC[0] > thresh) {
163 if (tdcpmt[pmtleft] == 1000.) { // fill for the first time
164 tdcpmt[pmtleft] = t1;
165 tdc[pmtleft] = t1 + c2_S[pmtleft] ; // Signal reaches Coincidence
166 }
167 if (tdcpmt[pmtleft] < 1000.) // is already filled!
168 if (t1 < tdcpmt[pmtleft]) {
169 tdcpmt[pmtleft] = t1;
170 t1 = t1 + c2_S[pmtleft] ; // Signal reaches Coincidence
171 tdc[pmtleft] = t1;
172 }
173 }
174 if (QhitPmt_pC[1] > thresh) {
175 if (tdcpmt[pmtright] == 1000.) { // fill for the first time
176 tdcpmt[pmtright] = t2;
177 tdc[pmtright] = t2 + c2_S[pmtright] ; // Signal reaches Coincidence
178 }
179 if (tdcpmt[pmtright] < 1000.) // is already filled!
180 if (t2 < tdcpmt[pmtright]) {
181 tdcpmt[pmtright] = t2;
182 t2 = t2 + c2_S[pmtright] ;
183 tdc[pmtright] = t2;
184 }
185 }
186 if(DEBUG)cout<<nh<<" "<<Timetof[nh]<<" "<<t1<<" "<<t2<<endl;
187 } // ip > -1 && ip < 6
188 } // **************************************** end loop over hits
189 } // NTHTOF < 200
190 // ====== ADC ======
191 for(Int_t i=0; i<48; i++){
192 if (QevePmt_pC[i] <= 800.) ADCtof[i]= (Int_t)(ADC_pC0A + ADC_pC1A*QevePmt_pC[i] + ADC_pC2A*pow(QevePmt_pC[i],2) + ADC_pC3A*pow(QevePmt_pC[i],3));
193 if (QevePmt_pC[i] > 800.) ADCtof[i]= (Int_t)(ADC_pC0B + ADC_pC1B*QevePmt_pC[i] + ADC_pC2B*pow(QevePmt_pC[i],2) + ADC_pC3B*pow(QevePmt_pC[i],3));
194 if (QevePmt_pC[i] > 2485.) ADCtof[i]= (Int_t)(1758. + 0.54*QevePmt_pC[i]); //assuming a fictional 0.54 ch/pC above ADCsat
195 if (ADCtof[i]>ADCsat) ADCtof[i]=ADCsat;
196 if (QevePmt_pC[i] < pCthres) ADCtof[i]= ADClast;
197 if (ADCtof[i] < 0) ADCtof[i]=ADClast;
198 if (ADCtof[i] > ADClast) ADCtof[i]=ADClast;
199 //if(ADCtof[i]!=4095)cout<<ADCtof[i]<<" ";
200 //if((i+1)%4==0)cout<<endl;
201 }
202 // cin>>ciao;
203
204 // ====== build TDC coincidence ======
205
206 //
207 for(Int_t i=0; i<48; i++) {
208 if((tdcpmt[i] - c1_S[i]) > 1e-7) {
209 tdcpmt[i] = 0.;
210 tdc[i] = 0.;
211 }
212 }// cycle to introduce a window for tdc
213
214
215 Float_t t_coinc = 0;
216 // Int_t ilast = 100;
217 for (Int_t ii=0; ii<48;ii++)
218 if (tdc[ii] > t_coinc) {
219 t_coinc = tdc[ii];
220 // ilast = ii;
221 }
222
223 // cout<<ilast<<" "<<t_coinc<<endl;
224 // At t_coinc trigger condition is fulfilled
225
226 for (Int_t ii=0; ii<48;ii++){
227 // if (tdc[ii] != 0) tdc1[ii] = t_coinc - tdc[ii]; // test 1
228 if (tdc[ii] != 0) tdc1[ii] = t_coinc - tdcpmt[ii]; // test 2
229 tdc1[ii] = tdc1[ii]/tdcres[ii]; // divide by TDC resolution
230 if (tdc[ii] != 0) tdc1[ii] = tdc1[ii] + c3_S[ii]; // add cable length c3
231 } // missing parenthesis inserted! (Silvio)
232
233 for(Int_t i=0; i<48; i++){
234 if(tdc1[i] != 0.){
235 TDCint[i]=(Int_t)tdc1[i];
236 if (TDCint[i]>4093) TDCint[i]=TDClast; // 18-oct WM
237 if (DEBUG)cout<<i<<" "<<TDCint[i]<<endl;
238 } else
239 TDCint[i]= TDClast;
240 }
241 if (DEBUG)cout<<"-----------"<<endl;
242 //------ use channelmap for ToF: 18-oct WM
243 Int_t channelmap[] = {3,21,11,29,19,45,27,37,36,28,44,20,5,12,13,4,
244 6,47,14,39,22,31,30,23,38,15,46,7,0,33,16,24,
245 8,41,32,40,25,17,34,9,42,1,2,10,18,26,35,43};
246 Int_t ADChelp[48],TDChelp[48];
247 for(Int_t i=0; i<48; i++){
248 Int_t ii=channelmap[i];
249 ADChelp[ii]= ADCtof[i];
250 TDChelp[ii]= TDCint[i];
251 }
252 for(Int_t i=0; i<48; i++){
253 ADCtof[i]= ADChelp[i];
254 TDCint[i]= TDChelp[i];
255 }
256 // ====== write fDataTof =======
257 UChar_t Ctrl3bit[8]={32,0,96,64,160,128,224,192}; // DC (msb in 8 bit word )
258 UChar_t tofBin;
259 for (Int_t j=0; j < 12; j++){ // loop on TDC #12
260 Int_t j12=j*23; // for each TDC 23 bytes (8 bits)
261 fDataTof[j12+0]=0x00; // TDC_ID
262 fDataTof[j12+1]=0x00; // EV_COUNT
263 fDataTof[j12+2]=0x00; // TDC_MASK (1)
264 fDataTof[j12+3]=0x00; // TDC_MASK (2)
265 for (Int_t k=0; k < 4; k++){ // for each TDC 4 channels (ADC+TDC)
266 Int_t jk12=j12+4*k; // ADC,TDC channel (0-47)
267 tofBin =(UChar_t)(ADCtof[k+4*j]/256); // ADC# (msb)
268 fDataTof[jk12+4] = Bin2GrayTof(tofBin,fDataTof[jk12+4]);
269 /* control bits inserted here, after the bin to gray conv - DC*/
270 fDataTof[jk12+4] = Ctrl3bit[2*k] | fDataTof[jk12+4];
271 tofBin=(UChar_t)(ADCtof[k+4*j]%256); // ADC# (lsb)
272 fDataTof[jk12+5] = Bin2GrayTof(tofBin,fDataTof[jk12+5]);
273 tofBin=(UChar_t)(TDCint[k+4*j]/256); // TDC# (msb)
274 fDataTof[jk12+6]=Bin2GrayTof(tofBin,fDataTof[jk12+6]);
275 /* control bits inserted here, after the bin to gray conv - DC*/
276 fDataTof[jk12+6] = Ctrl3bit[2*k+1] | fDataTof[jk12+6];
277 tofBin=(UChar_t)(TDCint[k+4*j]%256); // TDC# (lsb)
278 fDataTof[jk12+7]=Bin2GrayTof(tofBin,fDataTof[jk12+7]);
279 }
280 fDataTof[j12+20]=0x00; // TEMP1
281 fDataTof[j12+21]=0x00; // TEMP2
282 fDataTof[j12+22]= EvaluateCrcTof(pTof); // CRC
283 pTof+=23;
284 }
285 // ====== evaluate trigger variables =======
286 //fDataTrigger: 152 bytes (corrected 30/11/'07 SO - it was 153)
287 // initialization:
288 for (Int_t j=0; j < 152; j++)fDataTrigger[j]=0x00;
289 UChar_t *pTrg=fDataTrigger;
290 // Only the variables with a (*) are modified; the others are set to 0
291 // info given in #bites data + #bites crc
292 // TB_READ_PMT_PLANE : 6 + 1
293 // TB_READ_EVENT_COUNT : 3 + 1 (*)
294 // TB_READ_TRIGGER_RATE : 12 + 1
295 // TB_READ_D_L_TIME : 4 + 1
296 // TB_READ_S4_CAL_COUNT : 4 + 1
297 // TB_READ_PMT_COUNT1 : 48 + 1
298 // TB_READ_PMT_COUNT2 : 48 + 1
299 // TB_READ_PATTERN_BUSY : 8 + 1
300 // TB_READ_PATTERN_TRIGGER: 7 + 1 (*)
301 // TB_READ_TRIGGER_CONF : 2 + 1 (*)
302
303 // TB_READ_EVENT_COUNT
304 fhBookTree->SetBranchStatus("Ievnt",&Ievnt);
305 UInt_t cTrg = (UInt_t)Ievnt; //counter
306 UInt_t cTrg2 = 0; //counter with bits inverted, according to document
307 //"formato dati provenienti dalla trigger board"
308 for (Int_t i=0; i < 24; i++){ // Use the first 24 bits
309 if (cTrg & (0x1 << i) )
310 cTrg2 = cTrg2 | (0x1 << (24-i));
311 }
312 fDataTrigger[7] = (UChar_t)(cTrg2 >> 16); // 8 MSbits (out of 24)
313 fDataTrigger[8] = (UChar_t)(cTrg2 >> 8); // 8 "middle" bits
314 fDataTrigger[9] = (UChar_t)(cTrg2); // 8 LSbits
315 pTrg=fDataTrigger+7;
316 fDataTrigger[10]=EvaluateCrcTrigger(pTrg, 3);
317
318 // TB_READ_PATTERN_TRIGGER: bytes 141-148:
319 // PatternTrigMap[i] corresponds to bit i in TB_READ_PATTERN_TRIGGER:
320 // mapping according to documents:
321 // 1. "formato dati provenienti dalla trigger board"
322 // 2. "The ToF quicklook software", Appendix A (Campana, Nagni)
323 Int_t PatternTrigMap[]={29,42,43,1,16,7,17,28,33,41,46,2,15,8,18,27,
324 30,40,44,3,14,9,19,26,32,37,47,4,13,10,20,25,
325 34,31,38,45,5,12,21,24,36,35,39,48,6,11,22,23};
326 for (Int_t i=0; i < 48; i++)
327 //if (ADCtof[i]>thrTrg)
328 if (tdc1[channelmap[i]]!=0)
329 fDataTrigger[147-(Int_t)((PatternTrigMap[i]+1)/8)]=fDataTrigger[147-(Int_t)((PatternTrigMap[i]+1)/8)] | (0x1 << (PatternTrigMap[i]%8));
330 pTrg=fDataTrigger+141;
331 fDataTrigger[148]=EvaluateCrcTrigger(pTrg, 7);
332
333 // TB_READ_TRIGGER_CONF : set always acq.mode TOF4
334 //
335 // TOF1: S1-S2-S3 (&,|)
336 // TOF4: S2-S3 (&,&)
337 fDataTrigger[149]=0x02;
338 fDataTrigger[150]=0x0;
339 pTrg=fDataTrigger+149;
340 fDataTrigger[151]=EvaluateCrcTrigger(pTrg, 2);
341 }
342
343
344 UChar_t Digitizer::Bin2GrayTof(UChar_t binaTOF,UChar_t grayTOF){
345 union graytof_data {
346 UChar_t word;
347 struct bit_field {
348 unsigned b0:1;
349 unsigned b1:1;
350 unsigned b2:1;
351 unsigned b3:1;
352 unsigned b4:1;
353 unsigned b5:1;
354 unsigned b6:1;
355 unsigned b7:1;
356 } bit;
357 } bi,gr;
358 //
359 bi.word = binaTOF;
360 gr.word = grayTOF;
361 //
362 gr.bit.b0 = bi.bit.b1 ^ bi.bit.b0;
363 gr.bit.b1 = bi.bit.b2 ^ bi.bit.b1;
364 gr.bit.b2 = bi.bit.b3 ^ bi.bit.b2;
365 gr.bit.b3 = bi.bit.b3;
366 //
367 /* bin to gray conversion 4 bit per time*/
368 //
369 gr.bit.b4 = bi.bit.b5 ^ bi.bit.b4;
370 gr.bit.b5 = bi.bit.b6 ^ bi.bit.b5;
371 gr.bit.b6 = bi.bit.b7 ^ bi.bit.b6;
372 gr.bit.b7 = bi.bit.b7;
373 //
374 return(gr.word);
375 }
376
377 UChar_t Digitizer::EvaluateCrcTof(UChar_t *pTof) {
378 Bool_t DEBUG=false;
379 if (DEBUG)
380 return(0x00);
381
382 UChar_t crcTof=0x00;
383 UChar_t *pc=&crcTof, *pc2;
384 pc2=pTof;
385 for (Int_t jp=0; jp < 22; jp++){ // cecilia: fixed 23->22
386 //crcTof = crc8(...)
387 Crc8Tof(pc2++,pc);
388 // printf("%2i --- %x\n",jp,crcTof);
389 }
390 return(crcTof);
391 }
392
393 UChar_t Digitizer::EvaluateCrcTrigger(UChar_t *pTrg, Int_t nb) {
394 Bool_t DEBUG=false;
395 if (DEBUG)
396 return(0x00);
397
398 UChar_t crcTrg=0x00;
399 UChar_t *pc=&crcTrg, *pc2;
400 pc2=pTrg;
401 for (Int_t jp=0; jp < nb; jp++)
402 Crc8Tof(pc2++,pc);
403 return(crcTrg);
404 }
405
406 void Digitizer::Crc8Tof(UChar_t *oldCRC, UChar_t *crcTof){
407 union crctof_data {
408 UChar_t word;
409 struct bit_field {
410 unsigned b0:1;
411 unsigned b1:1;
412 unsigned b2:1;
413 unsigned b3:1;
414 unsigned b4:1;
415 unsigned b5:1;
416 unsigned b6:1;
417 unsigned b7:1;
418 } bit;
419 } c,d,r;
420
421 c.word = *oldCRC;
422 //d.word = *newCRC;
423 d.word = *crcTof;
424 r.word = 0;
425
426 r.bit.b0 = c.bit.b7 ^ c.bit.b6 ^ c.bit.b0 ^
427 d.bit.b0 ^ d.bit.b6 ^ d.bit.b7;
428
429 r.bit.b1 = c.bit.b6 ^ c.bit.b1 ^ c.bit.b0 ^
430 d.bit.b0 ^ d.bit.b1 ^ d.bit.b6;
431
432 r.bit.b2 = c.bit.b6 ^ c.bit.b2 ^ c.bit.b1 ^ c.bit.b0 ^
433 d.bit.b0 ^ d.bit.b1 ^ d.bit.b2 ^ d.bit.b6;
434
435 r.bit.b3 = c.bit.b7 ^ c.bit.b3 ^ c.bit.b2 ^ c.bit.b1 ^
436 d.bit.b1 ^ d.bit.b2 ^ d.bit.b3 ^ d.bit.b7;
437
438 r.bit.b4 = c.bit.b4 ^ c.bit.b3 ^ c.bit.b2 ^
439 d.bit.b2 ^ d.bit.b3 ^ d.bit.b4;
440
441 r.bit.b5 = c.bit.b5 ^ c.bit.b4 ^ c.bit.b3 ^
442 d.bit.b3 ^ d.bit.b4 ^ d.bit.b5;
443
444 r.bit.b6 = c.bit.b6 ^ c.bit.b5 ^ c.bit.b4 ^
445 d.bit.b4 ^ d.bit.b5 ^ d.bit.b6;
446
447 r.bit.b7 = c.bit.b7 ^ c.bit.b6 ^ c.bit.b5 ^
448 d.bit.b5 ^ d.bit.b6 ^ d.bit.b7 ;
449
450 *crcTof=r.word;
451 //return r.word;
452 };
453
454 void Digitizer::Paddle2Pmt(Int_t plane, Int_t paddle, Int_t *pl, Int_t *pr){
455 //* @param plane (0 - 5)
456 //* @param paddle (plane=0, paddle = 0,...5)
457 //* @param padid (0 - 23)
458 //
459 Int_t padid=-1;
460 Int_t pads[6]={8,6,2,2,3,3};
461 //
462 Int_t somma=0;
463 Int_t np=plane;
464 for(Int_t j=0; j<np; j++)somma+=pads[j];
465 padid=paddle+somma;
466 *pl = padid*2;
467 // *pr = *pr + 1;
468 *pr = *pl + 1; // WM
469 };
470
471 void Digitizer::LoadTOFCalib(int np,float *atte1,float *atte2,float *lambda1,float *lambda2){
472 stringstream calfile;
473 Int_t error = 0,temp=0;
474 GL_PARAM *glparam = new GL_PARAM();
475 error = glparam->Query_GL_PARAM(3,202,fDbc);
476 if(!error){
477 calfile.str("");
478 calfile << glparam->PATH.Data() << "/";
479 calfile << glparam->NAME.Data();
480 printf("\n Using TOF calibration file: \n %s\n",calfile.str().c_str());
481 ifstream fileTriggerCalib;
482 fileTriggerCalib.open(calfile.str().c_str());
483 for(Int_t i=0; i<np; i++){
484 fileTriggerCalib >> temp;
485 fileTriggerCalib >> atte1[i];
486 fileTriggerCalib >> lambda1[i];
487 fileTriggerCalib >> atte2[i];
488 fileTriggerCalib >> lambda2[i];
489 fileTriggerCalib >> temp;
490 }
491 fileTriggerCalib.close();
492 }
493 else{
494 cout<<endl<<" *********** ATTENTION ***********"<<endl;
495 cout<<endl<<" TOF: NO trigger calib file!"<<endl<<endl;
496 cout<<endl<<" TOF digitized data will be wrong!"<<endl<<endl;
497 for(Int_t i=0; i<np; i++){
498 atte1[i]=0.;
499 lambda1[i]=0.;
500 atte2[i]=0.;
501 lambda2[i]=0.;
502 }
503 }
504 //end tof calib
505 }

  ViewVC Help
Powered by ViewVC 1.1.23