/[PAMELA software]/trieste/pamVMC/trk/src/PamVMCTrkDig.cxx_good
ViewVC logotype

Contents of /trieste/pamVMC/trk/src/PamVMCTrkDig.cxx_good

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1.1.1 - (show annotations) (download) (vendor branch)
Wed Mar 4 12:51:30 2009 UTC (15 years, 11 months ago) by pamelats
Branch: MAIN, pamVMC
CVS Tags: start, v0r00, HEAD
Changes since 1.1: +0 -0 lines
Test pamVMC

1 #include "PamVMCTrkDig.h"
2 #include <TMath.h>
3
4 ClassImp(PamVMCTrkDig)
5
6 using TMath::Abs;
7
8 void PamVMCTrkDig::LoadFile(){
9
10 cout<<"Loading Tracker Calibrations..."<<endl;
11
12 Int_t time=0;
13 Int_t type=8;
14
15 fdberr = fsql->Query_GL_PARAM(time,type);
16 //we hardcoded values of time and type of calibrations
17 //later we'll introduce an options which allow us to load
18 //from different sources
19
20 if(fdberr<0){
21 TString crfname = fpath+"/CalibTrk_00110_000_000.root";
22 cout<<"No such record in DB for TRK: time="<<time
23 <<" type="<<type<<endl
24 <<"We will use file:"<<crfname<<endl;
25
26 fcrfile = new TFile(crfname);
27
28 } else {
29
30 fquery.str("");
31 fquery << fsql->GetPAR()->PATH.Data() << "/";
32 fquery << fsql->GetPAR()->NAME.Data();
33
34 ThrowCalFileUsage("TRK",fquery.str().c_str());
35
36 fcrfile = new TFile(fquery.str().c_str());
37 }
38
39 if(!fcrfile) ThrowCalFileWarning("TRK"); else
40 if(fcrfile->IsZombie()){
41 ThrowCalFileWarning("TRK"); //critical error
42 return;
43 }
44
45 //Load calibrations. I don't know, some empiric numbers of
46 //broken VA chips introduced. Should be load as options..
47
48 //------X-plane------//
49 TTree *tr1 = (TTree*)fcrfile->Get("CalibTrk1");
50 CalibTrk1Event *caldata1 = 0;
51 tr1->SetBranchAddress("CalibTrk1",&caldata1);
52 tr1->GetEntry(1);
53
54 Int_t Kview,jj;
55
56 for (Int_t i=0; i<6;i++) {
57 for (Int_t j=0; j<3072;j++) {
58 jj=j;
59 // broken va1 replaced
60 BrokenStrip(0,4,i,j,jj);
61 BrokenStrip(5,12,i,j,jj);
62
63 Kview=caldata1->DSPnumber[i]-1;
64 fPedeTrack[Kview][j]=caldata1->DSPped_par[i][jj];
65 fSigmaTrack[Kview][j]=caldata1->DSPsig_par[i][jj];
66 if(caldata1->DSPbad_par[i][jj]==1) fSigmaTrack[Kview][j]=-fSigmaTrack[Kview][j];
67 };
68 };
69
70 //------Y-plane------//
71 TTree *tr2 = (TTree*)fcrfile->Get("CalibTrk2");
72 CalibTrk2Event *caldata2 = 0;
73 tr2->SetBranchAddress("CalibTrk2",&caldata2);
74 tr2->GetEntry(1);
75
76 for (Int_t i=0; i<6;i++) {
77 for (Int_t j=0; j<3072;j++) {
78 jj=j;
79 // broken va1 replaced
80 BrokenStrip(2,3,i,j,jj);
81 BrokenStrip(2,5,i,j,jj);
82 BrokenStrip(2,6,i,j,jj);
83 BrokenStrip(2,7,i,j,jj);
84 BrokenStrip(2,11,i,j,jj);
85
86 Kview=caldata2->DSPnumber[i]-1;
87 fPedeTrack[Kview][j]=caldata2->DSPped_par[i][jj];
88 fSigmaTrack[Kview][j]=caldata2->DSPsig_par[i][jj];
89 if(caldata2->DSPbad_par[i][jj]==1) fSigmaTrack[Kview][j]=-fSigmaTrack[Kview][j];
90 };
91 };
92
93 fcrfile->Close();
94 }
95
96
97 void PamVMCTrkDig::DigitizeTrackCalib(Int_t n){
98
99 cout<<"Starting Tracker Calibrations..."<<endl;
100 if( (n!=1)&&(n!=2) ) {
101 cout << "!!!ERROR: Wrong DigitizeTrackCalib argument!!!" << endl;
102 return;
103 };
104
105 fData.clear();
106
107 UShort_t Dato;
108
109 USBuffer tempdat; //temporary buffer to keep data before writing DSP data
110 USBuffer::const_iterator p; //iterator;
111
112 Float_t dato1,dato2,dato3,dato4;
113 UShort_t DatoDec,DatoDec1,DatoDec2,DatoDec3,DatoDec4;
114 UShort_t EVENT_CAL, PED_L1, ReLength,OveCheckCode;
115 UShort_t CkSum;
116
117 for (Int_t j=n-1; j<fNviews;j+=2) { //0(1)...12
118 CkSum=0;
119 // here skip the dsp header and his trailer , to be written later
120 tempdat.clear(); //clearing temporary buffer
121 for (Int_t i=0; i<fNladder;i++) { //0...2
122 for (Int_t k=0; k<fNstrips_ladder;k++) { //0...1023
123 // write in buffer the current LADDER
124 Dato=(UShort_t)fPedeTrack[j][i*fNstrips_ladder+k];
125 dato1=fPedeTrack[j][i*fNstrips_ladder+k]-Dato;
126
127 DatoDec1=(UShort_t)(dato1*2);
128 dato2=dato1*2-DatoDec1;
129
130 DatoDec2=(UShort_t)(dato2*2);
131 dato3=dato2*2-DatoDec2;
132
133 DatoDec3=(UShort_t)(dato3*2);
134 dato4=dato3*2-DatoDec3;
135
136 DatoDec4=(UShort_t)(dato4*2);
137
138 DatoDec=DatoDec1*0x0008+DatoDec2*0x0004+DatoDec3*0x0002+DatoDec4*0x0001;
139
140 tempdat.push_back( ((Dato << 4) | (DatoDec & 0x000F)) );
141 CkSum^=tempdat.back();
142 };
143
144 for (Int_t k=0; k<fNstrips_ladder;k++) { //0...1023
145 // write in buffer the current LADDER
146 Dato=(UShort_t)fabs(fSigmaTrack[j][i*fNstrips_ladder+k]);
147 dato1=fabs(fSigmaTrack[j][i*fNstrips_ladder+k])-Dato;
148
149 DatoDec1=(UShort_t)(dato1*2);
150 dato2=dato1*2-DatoDec1;
151
152 DatoDec2=(UShort_t)(dato2*2);
153 dato3=dato2*2-DatoDec2;
154
155 DatoDec3=(UShort_t)(dato3*2);
156 dato4=dato3*2-DatoDec3;
157
158 DatoDec4=(UShort_t)(dato4*2);
159
160 DatoDec=DatoDec1*0x0008+DatoDec2*0x0004+DatoDec3*0x0002+DatoDec4*0x0001;
161
162 tempdat.push_back( ((Dato << 4) | (DatoDec & 0x000F)) );
163 CkSum^=tempdat.back();
164 };
165
166 for (Int_t k=0; k<64;k++) { //0...63
167 UShort_t DatoBad=0x0000;
168 for (Int_t nb=0; nb<16;nb++) {
169 if( fSigmaTrack[j][i*fNstrips_ladder+k*16+nb]<0. ) DatoBad=( DatoBad | (0x8000 >> nb) );
170 };
171
172 tempdat.push_back(DatoBad);
173 CkSum^=tempdat.back();
174 };
175 // end ladder
176
177 // write in buffer the end ladder word
178
179 switch(i){
180 case 0:
181 tempdat.push_back(0x1807);
182 break;
183 case 1:
184 tempdat.push_back(0x1808);
185 break;
186 case 2:
187 tempdat.push_back(0x1808);
188 break;
189 default:
190 break;
191 }
192 CkSum^=tempdat.back();
193
194 // write in buffer the TRAILER
195 ReLength=(UShort_t)((fNstrips_ladder*2+64+1)*2+3);
196 OveCheckCode=0x0000;
197
198
199 tempdat.push_back(0x0000);
200 tempdat.push_back((ReLength >> 8));
201 tempdat.push_back(( (ReLength << 8) | (OveCheckCode & 0x00FF) ));
202
203 }; // end TRAILER
204
205 cout<<"LENGTH OF TEMPDAT"<<tempdat.size()<<endl;
206 // write in buffer the DSP header
207
208 fData.push_back((0xE800 | ( ((j+1) << 3) | 0x0005) ));
209 fData.push_back(0x01A9);
210 fData.push_back(0x8740);
211 EVENT_CAL=0;
212 fData.push_back((0x1A00 | ( (0x03FF & EVENT_CAL)>> 1) ));
213 PED_L1=0;
214 fData.push_back(( ((EVENT_CAL << 15) | 0x5002 ) | ((0x03FF & PED_L1) << 2) ));
215 fData.push_back(0x8014);
216 fData.push_back(0x00A0);
217 fData.push_back(0x0500);
218 fData.push_back(0x2801);
219 fData.push_back(0x400A);
220 fData.push_back(0x0050);
221 CkSum=(CkSum >> 8)^(CkSum&0x00FF);
222 fData.push_back((0x0280 | (CkSum >> 3)));
223 fData.push_back((0x1FFF | (CkSum << 13) ));
224
225 ReLength=(UShort_t)((13*2)+3);
226 OveCheckCode=0x0000;
227 fData.push_back(0x0000);
228 fData.push_back((ReLength >> 8));
229 fData.push_back(( (ReLength << 8) | (OveCheckCode & 0x00FF) ));
230 cout<<"fDATA length before adding:"<<fData.size()<<endl;
231 // Now we will copy to fData vector data from tempdat:
232 p = tempdat.begin();
233 while( p != tempdat.end() ){
234 fData.push_back(*p);
235 p++;
236 }
237 cout<<"fDATA length after adding:"<<fData.size()<<endl;
238 };
239 }
240
241 void PamVMCTrkDig::WriteCalib(){
242 cout<<"Writing Tracker Calibrations..."<<endl;
243 fraw->WritePSCU(&fDataPSCU);
244 fraw->CopyUShortToBuff(&fData);
245 if(fPadding) fraw->WritePadding(&fDataPadding);
246 }
247
248 void PamVMCTrkDig::LoadMipCor(){
249
250 std:: cout << "Entering LoadMipCor" << endl;
251
252 Float_t xfactor=1./151.6*1.04;
253 Float_t yfactor=1./152.1;
254
255 fMipCor[0][0]=140.02*yfactor;
256 fMipCor[0][1]=140.99*xfactor;
257 fMipCor[0][2]=134.48*yfactor;
258 fMipCor[0][3]=144.41*xfactor;
259 fMipCor[0][4]=140.74*yfactor;
260 fMipCor[0][5]=142.28*xfactor;
261 fMipCor[0][6]=134.53*yfactor;
262 fMipCor[0][7]=140.63*xfactor;
263 fMipCor[0][8]=135.55*yfactor;
264 fMipCor[0][9]=138.00*xfactor;
265 fMipCor[0][10]=154.95*yfactor;
266 fMipCor[0][11]=158.44*xfactor;
267
268
269 fMipCor[1][0]=136.07*yfactor;
270 fMipCor[1][1]=135.59*xfactor;
271 fMipCor[1][2]=142.69*yfactor;
272 fMipCor[1][3]=138.19*xfactor;
273 fMipCor[1][4]=137.35*yfactor;
274 fMipCor[1][5]=140.23*xfactor;
275 fMipCor[1][6]=153.15*yfactor;
276 fMipCor[1][7]=151.42*xfactor;
277 fMipCor[1][8]=129.76*yfactor;
278 fMipCor[1][9]=140.63*xfactor;
279 fMipCor[1][10]=157.87*yfactor;
280 fMipCor[1][11]=153.64*xfactor;
281
282 fMipCor[2][0]=134.98*yfactor;
283 fMipCor[2][1]=143.95*xfactor;
284 fMipCor[2][2]=140.23*yfactor;
285 fMipCor[2][3]=138.88*xfactor;
286 fMipCor[2][4]=137.95*yfactor;
287 fMipCor[2][5]=134.87*xfactor;
288 fMipCor[2][6]=157.56*yfactor;
289 fMipCor[2][7]=157.31*xfactor;
290 fMipCor[2][8]=141.37*yfactor;
291 fMipCor[2][9]=143.39*xfactor;
292 fMipCor[2][10]=156.15*yfactor;
293 fMipCor[2][11]=158.79*xfactor;
294
295 }
296
297
298 void PamVMCTrkDig::Digitize(){
299
300 cout<<"For now I Digitized DUMMY for Tracker"<<endl;
301
302 fData.clear();
303 //fData.push_back(0xCAAA);
304
305
306 //fhits->Print();
307 //for (Int_t i =1; i<64; i++) fData.push_back(0xFFFF);
308
309
310
311 Int_t Iview;
312 Int_t Nstrip;
313
314 for (Int_t j=0; j<fNviews;j++) {
315
316 for (Int_t i=0; i<fNladder;i++) {
317
318 Float_t commonN1=gRandom->Gaus(0.,fSigmaCommon);
319 Float_t commonN2=gRandom->Gaus(0.,fSigmaCommon);
320 for (Int_t k=0; k<fNstrips_ladder;k++) {
321 Nstrip=i*fNstrips_ladder+k;
322 Float_t Sigma=Abs(fSigmaTrack[j][Nstrip]);
323 AdcTrack[j][Nstrip]=gRandom->Gaus(fPedeTrack[j][Nstrip],Sigma );
324 if(k<4*128) {AdcTrack[j][Nstrip] += commonN1;} // full correlation of 4 VA1 Com. Noise
325 else {AdcTrack[j][Nstrip] += commonN2;} // full correlation of 4 VA1 Com. Noise
326 if(AdcTrack[j][Nstrip] < 0. ) AdcTrack[j][Nstrip]=0.;
327 if(AdcTrack[j][Nstrip] > 4095.) AdcTrack[j][Nstrip]=4095.;
328 };
329 };
330 };
331
332 if (fhits){
333 Float_t ADCfull;
334 Int_t iladd=0;
335 for (Int_t ix=0; ix<fhits->GetNXHit();ix++) {
336 Iview=((fhits->GetXHit(ix))->fnpstrip)*2-1;
337 Nstrip=(Int_t)((fhits->GetXHit(ix))->fistrip)-1;
338 if(Nstrip<fNstrips_ladder) iladd=0;
339 if((Nstrip>=fNstrips_ladder)&&(Nstrip<2*fNstrips_ladder)) iladd=1;
340 if((Nstrip>=2*fNstrips_ladder)&&(Nstrip<3*fNstrips_ladder)) iladd=2;
341 ADCfull=AdcTrack[Iview][Nstrip] += ((fhits->GetXHit(ix))->fqstrip)*fMipCor[iladd][Iview];
342 AdcTrack[Iview][Nstrip] *= SaturationTrackx(ADCfull);
343 };
344
345
346 for (Int_t iy=0; iy<fhits->GetNYHit();iy++) {
347 Iview=((fhits->GetYHit(iy))->fnpstrip)*2-2;
348 Nstrip=(Int_t)((fhits->GetYHit(iy))->fistrip)-1;
349 if(Nstrip<fNstrips_ladder) iladd=0;
350 if((Nstrip>=fNstrips_ladder)&&(Nstrip<2*fNstrips_ladder)) iladd=1;
351 if((Nstrip>=2*fNstrips_ladder)&&(Nstrip<3*fNstrips_ladder)) iladd=2;
352 ADCfull=AdcTrack[Iview][Nstrip] -= ((fhits->GetYHit(iy))->fqstrip)*fMipCor[iladd][Iview];
353 AdcTrack[Iview][Nstrip] *= SaturationTracky(ADCfull);
354 };
355
356 CompressTrackData();
357 };
358 }
359
360
361 void PamVMCTrkDig::CompressTrackData(){
362
363 // copy of the corresponding compression fortran routine + new digitization
364
365 Int_t oldval=0;
366 Int_t newval=0;
367 Int_t trasmesso=0;
368 Int_t ntrastot=0;
369 Float_t real, inte;
370 Int_t cercacluster=0;
371 Int_t kt=0;
372 static const int DSPbufferSize = 4000; // 13 bit buffer to be rearranged in 16 bit Track buffer
373 UShort_t DataDSP[DSPbufferSize]; // 13 bit buffer to be rearranged in 16 bit Track buffer
374 UShort_t DSPlength; // 13 bit buffer to be rearranged in 16 bit Track buffer
375 UShort_t CheckSum, Nword, Dato, DATA, ReLength, OveCheckCode;
376 Int_t k, diff, clval, clvalp, klp, kl1, kl2, Bit16free, Bit13ToWrite;
377
378 //To be removed:
379 const Int_t TRACKbuffer = 50000;
380 UShort_t DataTrack[TRACKbuffer];
381 memset(DataTrack,0,sizeof(UShort_t)*TRACKbuffer);
382 Int_t Tracklength=0;
383 UShort_t Ch2;
384
385 for (Int_t iv=0; iv<fNviews;iv++) {
386 memset(DataDSP,0,sizeof(UShort_t)*DSPbufferSize);
387 DSPlength=16; // skip the header, to be written later
388 CheckSum=0;
389
390 Ch2=0;
391
392 for (Int_t ladder=0; ladder<fNladder;ladder++) {
393 k=0;
394 while (k<fNstrips_ladder) {
395 // compress write in buffer the current LADDER
396 if ( k == 0) {
397 real=modff(AdcTrack[iv][ladder*fNstrips_ladder+k],&inte);
398 if (real > 0.5) inte=inte+1;
399 newval=(Int_t)inte -(Int_t)fPedeTrack[iv][ladder*fNstrips_ladder+k];
400 // first strip of ladder is transmitted
401 DataDSP[DSPlength]=( ((UShort_t)inte) & 0x0FFF);
402 DSPlength++;
403 ntrastot++;
404 trasmesso=1;
405 oldval=newval;
406 kt=k;
407 k++;
408 continue;
409 };
410 real=modff(AdcTrack[iv][ladder*fNstrips_ladder+k],&inte);
411 if (real > 0.5) inte=inte+1;
412 newval=(Int_t)inte -(Int_t)(fPedeTrack[iv][ladder*fNstrips_ladder+k]);
413
414 cercacluster=1; //PAR?
415
416 if (cercacluster==1) {
417
418 diff=0;
419 switch ((iv+1)%2) {
420 case 0: diff=newval-oldval;
421 break;
422 case 1: diff=oldval-newval;
423 break;
424 };
425
426
427 if (diff>fCutclu*(Int_t)fabs(fSigmaTrack[iv][ladder*fNstrips_ladder+k]) ) {
428 clval=newval;
429 klp=k; // go on to search for maximum
430 klp++;
431
432 while(klp<fNstrips_ladder) {
433 real=modff(AdcTrack[iv][ladder*fNstrips_ladder+klp],&inte);
434 if (real > 0.5) inte=inte+1;
435 clvalp=(Int_t)inte -(Int_t)fPedeTrack[iv][ladder*fNstrips_ladder+klp];
436 if((iv+1)%2==0) {
437
438 if(clvalp>clval) {
439 clval=clvalp;
440 k=klp;}
441 else break; // max of cluster found
442
443 } else {
444
445 if(clvalp<clval) {
446 clval=clvalp;
447 k=klp;}
448 else break; // max of cluster found
449 };
450 klp++;
451 };
452
453 kl1=k-fNclst; // max of cluster (or end of ladder ?)
454 trasmesso=0;
455 if(kl1<0) kl1=0;
456 if(kt>=kl1) kl1=kt+1;
457 if( (kt+1)==kl1 ) trasmesso=1;
458
459 kl2=k+fNclst;
460 if(kl2>=fNstrips_ladder) kl2=fNstrips_ladder-1;
461
462 for(Int_t klt=kl1 ; klt<=kl2 ; klt++) {
463 if(trasmesso==0) {
464
465 DataDSP[DSPlength]=( ((UShort_t)klt) | 0x1000);
466 DSPlength++;
467 ntrastot++;
468
469 real=modff(AdcTrack[iv][ladder*fNstrips_ladder+klt],&inte);
470 if (real > 0.5) inte=inte+1;
471 DataDSP[DSPlength]=( ((UShort_t)inte) & 0x0FFF);
472 DSPlength++;
473 ntrastot++;
474
475 }
476 else {
477
478 real=modff(AdcTrack[iv][ladder*fNstrips_ladder+klt],&inte);
479 if (real > 0.5) inte=inte+1;
480 DataDSP[DSPlength]=( ((UShort_t)inte) & 0x0FFF);
481 DSPlength++;
482 ntrastot++;
483 };
484 trasmesso=1;
485 }; // end trasmission
486 kt=kl2;
487 k=kl2;
488 real=modff(AdcTrack[iv][ladder*fNstrips_ladder+kt],&inte);
489 if (real > 0.5) inte=inte+1;
490 oldval=(Int_t)inte -(Int_t)fPedeTrack[iv][ladder*fNstrips_ladder+kt];
491 k++;
492 continue;
493 }
494 }//END cercacluster
495
496 // start ZOP check for strips no
497
498 if(abs(newval-oldval)>=fCutzop*(Int_t)fabs(fSigmaTrack[iv][ladder*fNstrips_ladder+k]) ) { if(trasmesso==0) {
499
500 DataDSP[DSPlength]=( ((UShort_t)k) | 0x1000);
501 DSPlength++;
502 ntrastot++;
503
504 real=modff(AdcTrack[iv][ladder*fNstrips_ladder+k],&inte);
505 if (real > 0.5) inte=inte+1;
506 DataDSP[DSPlength]=( ((UShort_t)inte) & 0x0FFF);
507 DSPlength++;
508 ntrastot++;
509
510 } else {
511 real=modff(AdcTrack[iv][ladder*fNstrips_ladder+k],&inte);
512 if (real > 0.5) inte=inte+1;
513 DataDSP[DSPlength]=( ((UShort_t)inte) & 0x0FFF);
514 DSPlength++;
515 ntrastot++;
516 };
517 trasmesso=1;
518 oldval=newval;
519 kt=k;
520
521 }
522 else trasmesso=0;
523 // end zop
524 k++;
525 };
526
527 DataDSP[DSPlength]=( ((UShort_t)(ladder+1)) | 0x1800);
528 DSPlength++;
529 ntrastot++;
530 trasmesso=0;
531
532 }; //end cycle inside dsp
533
534 // here put DSP header
535 DataDSP[0]=(0x1CA0 | ((UShort_t)(iv+1)) );
536 Nword=(DSPlength*13)/16;
537 if( ((DSPlength*13)%16)!=0) Nword++;
538 DataDSP[1]=(0x1400 | ( Nword >> 10));
539 DataDSP[2]=(0x1400 | ( Nword & 0x03FF) );
540 DataDSP[3]=(0x1400 | (( (UShort_t)(fraw->GetCounter() >> 10) ) & 0x03FF) );
541 DataDSP[4]=(0x1400 | (( (UShort_t)(fraw->GetCounter()) ) & 0x03FF) );
542 DataDSP[5]=(0x1400 | ( (UShort_t)(fNclst << 7) ) | ( (UShort_t)(fCutzop << 4) )
543 | ( (UShort_t)fCutzop ) );
544 DataDSP[6]=0x1400;
545 DataDSP[7]=0x1400;
546 DataDSP[8]=0x1400;
547 DataDSP[9]=0x1400;
548 DataDSP[10]=0x1400;
549 DataDSP[11]=0x1400;
550 DataDSP[12]=0x1400;
551 DataDSP[13]=0x1400;
552 DataDSP[14]=(0x1400 | (CheckSum & 0x00FF) );
553 DataDSP[15]=0x1C00;
554 // end DSP header
555
556 // write 13 bit DataDSP bufer inside 16 bit fDataTrack buffer
557 cout<<"++++Beginning cycle for DSP "<<iv<<endl;
558 cout<<" Tracklength="<<fData.size()<<endl;
559 Bit16free=16;
560 DATA = 0;
561 for (Int_t NDSP=0; NDSP<DSPlength;NDSP++) {
562 Bit13ToWrite=13;
563 while(Bit13ToWrite>0) {
564 if(Bit13ToWrite<=Bit16free) {
565 Dato=((DataDSP[NDSP]&(0xFFFF >> (16-Bit13ToWrite)))<<(Bit16free-Bit13ToWrite));
566 //cout<<" before DATA:"<<DATA<<" DataTrack: "<<DataTrack[Tracklength]<<endl;
567 DATA = DATA | Dato ;
568 DataTrack[Tracklength]=DataTrack[Tracklength] | Dato ;
569 //cout<<" after DATA:"<<DATA<<" DataTrack: "<<DataTrack[Tracklength]<<endl;
570 Bit16free=Bit16free-Bit13ToWrite;
571 Bit13ToWrite=0;
572 if(Bit16free==0) {
573 if(NDSP>15) Ch2=Ch2^DataTrack[Tracklength];
574 //cout<<"Datatrack"<<Int_t(DataTrack[Tracklength])<<" DATA"<<Int_t(DATA)<<" Tracklength:"<< Tracklength<<endl;
575 Tracklength++;
576 if(NDSP>15) CheckSum=CheckSum^DATA;
577 fData.push_back(DATA);
578 cout<<"i="<<Tracklength<<" Inside first DATA:"<<hex<<DATA<<" DataTrack: "<<hex<<DataTrack[Tracklength]<<dec<<endl;
579 DATA = 0;
580 //cout<<" Writing Data["<<fData.size()<<"]"<<endl;
581 Bit16free=16;
582 };
583 }
584 else if(Bit13ToWrite>Bit16free) {
585 Dato=( (DataDSP[NDSP]&(0xFFFF >> (16-Bit13ToWrite) ) ) >> (Bit13ToWrite-Bit16free) );
586 DataTrack[Tracklength]=DataTrack[Tracklength] | Dato ;
587 if(NDSP>15) Ch2=Ch2^DataTrack[Tracklength];
588 DATA = DATA | Dato ;
589 cout<<"i="<<Tracklength<<" Inside second DATA:"<<hex<<DATA<<" DataTrack: "<<hex<<DataTrack[Tracklength]<<dec<<endl;
590 // cout<<"Datatrack"<<Int_t(DataTrack[Tracklength])<<" DATA"<<Int_t(DATA)<<" Tracklength:"<< Tracklength<<endl;
591 Tracklength++;
592
593 fData.push_back(DATA);
594 if(NDSP>15) CheckSum=CheckSum^DATA;
595 DATA = 0;
596 //cout<<" Writing Data["<<fData.size()<<"]"<<endl;
597 Bit13ToWrite=Bit13ToWrite-Bit16free;
598 Bit16free=16;
599 };
600
601 }; // end cycle while(Bit13ToWrite>0)
602
603 }; // end cycle DataDSP
604
605 cout<<"!!!i=71"<<" DataTrack: "<<hex<<DataTrack[71]<<dec<<endl;
606 cout<<"Checksum:"<<Ch2<<" "<<CheckSum<<endl;
607 if(Bit16free!=16) {
608 fData.push_back(DATA);
609 DATA= 0;
610 CheckSum=CheckSum^DATA;
611 };
612 if(Bit16free!=16) { Tracklength++; CheckSum=CheckSum^DataTrack[Tracklength]; };
613 CheckSum=(CheckSum >> 8)^(CheckSum&0x00FF);
614 Ch2=(Ch2 >> 8)^(Ch2&0x00FF);
615 cout<<"Checksum:"<<Ch2<<" "<<CheckSum<<endl;
616 cout<<"+++Out of cycle++"<<endl;
617 cout<<" Writing Data["<<fData.size()-Nword+12<<"]"<<endl; //at(
618 //DataTrack[Tracklength-Nword+11]=(0x0280 | (CheckSum >> 3));
619 //DataTrack[Tracklength-Nword+12]=(0x1C00 | (CheckSum << 13) );
620 cout<<"Checksum1:"<<UShort_t((0x0280 | (CheckSum >> 3)))<<" Checksum2"<<UShort_t((0x1C00 | (CheckSum << 13)) )<<endl;
621 fData.at(fData.size()-Nword+11)=(0x0280 | (CheckSum >> 3));
622 //cout<<" Writing Data["<<fData.size()-Nword+13<<"]"<<endl;
623
624 fData.at(fData.size()-Nword+12)=(0x1C00 | (CheckSum << 13) );
625
626 // end write 13 bit DataDSP bufer inside 16 bit fDataTrack buffer
627
628 cout<<"+++Writing Trailer+++"<<endl;
629 //write trailer on buffer
630 ReLength=(UShort_t)((Nword+13)*2+3);
631 OveCheckCode=0x0000;
632
633 cout<<"i="<<Tracklength<<" DataTrack: "<<hex<<DataTrack[Tracklength]<<dec<<endl;
634 DataTrack[Tracklength]=0x0000;
635 Tracklength++;
636 cout<<"i="<<Tracklength<<" DataTrack: "<<hex<<DataTrack[Tracklength]<<dec<<endl;
637 //DataTrack[Tracklength]=(ReLength >> 8);
638 //Tracklength++;
639
640 //DataTrack[Tracklength]=( (ReLength << 8) | (OveCheckCode & 0x00FF) );
641 //Tracklength++;
642 //cout<<" Writing Data["<<fData.size()<<"]"<<endl;
643 //fData.push_back(0x0000);
644 fData.push_back(0x0000);
645 cout<<" Writing Data["<<fData.size()<<"]"<<endl;
646 fData.push_back((ReLength >> 8));
647
648 //cout<<" Writing Data["<<fData.size()<<"]"<<endl;
649 fData.push_back(( (ReLength << 8) | (OveCheckCode & 0x00FF) ));
650
651 // end trailer
652
653
654 }//END VIEWS CYCLE
655
656 // for(Int_t i=0; i<Tracklength; i++) fData.push_back(DataTrack[i]);
657 }
658
659
660 Float_t PamVMCTrkDig::SaturationTrackx(Float_t ADC) {
661 Float_t SatFact=1.;
662 if(ADC<1.) { SatFact=1./ADC; };
663 if(ADC>3000.) { SatFact=3000./ADC; };
664 return SatFact;
665 };
666
667
668 Float_t PamVMCTrkDig::SaturationTracky(Float_t ADC) {
669 Float_t SatFact=1.;
670 if(ADC<70.) { SatFact=70./ADC; };
671 if(ADC>4095.) { SatFact=4095./ADC; };
672 return SatFact;
673 };

  ViewVC Help
Powered by ViewVC 1.1.23