/[PAMELA software]/PamVMC_update/cal/src/PamVMCCaloDig.cxx
ViewVC logotype

Annotation of /PamVMC_update/cal/src/PamVMCCaloDig.cxx

Parent Directory Parent Directory | Revision Log Revision Log


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

1 formato 1.1 #include "PamVMCCaloDig.h"
2     #include <TTree.h>
3     #include "CalibCalPedEvent.h"
4    
5     using namespace pamela;
6    
7     ClassImp(PamVMCCaloDig)
8    
9     extern "C"{
10     short crc(short, short);
11     };
12    
13    
14     void PamVMCCaloDig::LoadCalib(){
15    
16     cout<<"Loading Calorimeter Calibrations..."<<endl;
17    
18    
19    
20     fhc =(TClonesArray*)fhitscolmap.GetValue("CAST");
21     if (!fhc) cout <<"!!!WARNING Hit Collection of CAL not found by PamVMCCaloDig!!!!"<<endl;
22    
23     Int_t GivenCaloCalib = 0; //@@@ should be given as input par @@@
24    
25     Int_t time=3;
26     Int_t type=101;
27    
28     fdberr = fsql->Query_GL_PARAM(time,type);
29    
30     if(fdberr<0){
31     cout<<"No such record in DB for CAL: time="<<time<<" type="<<type<<endl;
32     ThrowCalFileWarning("CAL");
33    
34     } else {
35    
36     fquery.str("");
37     fquery << fsql->GetPAR()->PATH.Data() << "/";
38     fquery << fsql->GetPAR()->NAME.Data();
39    
40     ThrowCalFileUsage("CAL",fquery.str().c_str());
41    
42    
43     FILE *f;
44     f = fopen(fquery.str().c_str(),"rb");
45    
46     if(f==NULL) ThrowCalFileWarning("CAL");
47     else{
48     memset(fCalomip,0,4224*sizeof(fCalomip[0][0][0]));
49    
50     for (Int_t m = 0; m < 2 ; m++ ){
51     for (Int_t k = 0; k < 22; k++ ){
52     for (Int_t l = 0; l < 96; l++ ){
53     fread(&fCalomip[m][k][l],sizeof(fCalomip[m][k][l]),1,f);
54     };
55     };
56     };
57     fclose(f);
58    
59    
60     // determine which calibration has to be used
61     // and load it for each section
62    
63     UInt_t idcalib, calibno;
64     UInt_t utime = 0;
65     TString calname;
66     for (UInt_t s=0; s<4; s++){
67    
68     // clear calo calib variables for section s
69    
70     ClearCaloCalib(s);
71    
72     if ( GivenCaloCalib ){
73    
74     // a time has been given as input on the command line
75     // so retrieve the calibration that preceed that time
76    
77     fsql->Query_GL_CALO_CALIB(GivenCaloCalib,utime,s);
78    
79     calibno = fsql->GetCaloCalib()->EV_ROOT;
80     idcalib = fsql->GetCaloCalib()->ID_ROOT_L0;
81    
82     // determine path and name and entry of the calibration file
83    
84     printf("\n");
85     printf(" ** SECTION %i **\n",s);
86    
87     fsql->Query_GL_ROOT(idcalib);
88    
89     fquery.str("");
90     fquery << fsql->GetROOT()->PATH.Data() << "/";
91     fquery << fsql->GetROOT()->NAME.Data();
92    
93     calname = (TString)fquery.str().c_str();
94    
95     printf("\n Section %i : using file %s calibration at entry %i: \n",s,calname.Data(),calibno);
96    
97     } else {
98    
99     fdberr = 0;
100     fdberr = fsql->Query_GL_PARAM(1,104);
101    
102     fquery.str("");
103     fquery << fsql->GetPAR()->PATH.Data() << "/";
104     fquery << fsql->GetPAR()->NAME.Data();
105    
106     printf("\n Section %i : using default calorimeter calibration: \n %s \n",s,fquery.str().c_str());
107    
108     calname = (TString)fquery.str().c_str();
109     calibno = s;
110    
111     };
112    
113     // load calibration variables in memory
114    
115     CaloLoadCalib(s,calname,calibno);
116    
117     };
118     //
119     // at this point we have in memory the calorimeter calibration
120     // and we can save it to disk in the correct format and use it to digitize the data
121     //
122    
123     DigitizeCALOCALIB();
124     }
125     }
126     }
127    
128    
129     void PamVMCCaloDig::ClearCaloCalib(Int_t s){
130    
131     fcstwerr[s] = 0;
132     fcperror[s] = 0.;
133     for ( Int_t d=0 ; d<11 ;d++ ){
134     Int_t pre = -1;
135     for ( Int_t j=0; j<96 ;j++){
136     if ( j%16 == 0 ) pre++;
137     fcalped[s][d][j] = 0.;
138     fcstwerr[s] = 0.;
139     fcperror[s] = 0.;
140     fcalgood[s][d][j] = 0.;
141     fcalthr[s][d][pre] = 0.;
142     fcalrms[s][d][j] = 0.;
143     fcalbase[s][d][pre] = 0.;
144     fcalvar[s][d][pre] = 0.;
145     };
146     };
147     return;
148     }
149    
150     Int_t PamVMCCaloDig::CaloLoadCalib(Int_t s,TString calname, UInt_t calibno){
151    
152     UInt_t e = 0;
153     switch(s){
154     case 0: e = 0;
155     break;
156     case 1: e = 2;
157     break;
158     case 2: e = 3;
159     break;
160     case 3: e = 1;
161     break;
162     default:
163     break;
164     }
165     fcfile.open(calname.Data());
166     ThrowCalFileUsage("CAL",calname.Data());
167     if ( !fcfile ){
168     ThrowCalFileWarning("CAL");
169     return(-107);
170     };
171     fcfile.close();
172    
173     fcrfile = new TFile(calname.Data());
174     if ( !fcrfile ){
175     ThrowCalFileWarning("CAL");
176     return(-108);
177     }
178     TTree *tr = (TTree*)fcrfile->Get("CalibCalPed");
179     if ( !tr ){
180     cout<<"!!!WARNING Tree CalibCalPed in file "<<calname.Data()
181     <<" was NOT found!!!"<<endl;
182     return(-109);
183     }
184     TBranch *calo = tr->GetBranch("CalibCalPed");
185     CalibCalPedEvent *ce = 0;
186     tr->SetBranchAddress("CalibCalPed", &ce);
187     Long64_t ncalibs = calo->GetEntries();
188     if ( !ncalibs ){
189     cout<<"!!!WARNING No entries in calibration file!!!"<<endl;
190     return(-110);
191     }
192    
193     calo->GetEntry(calibno);
194    
195     if (ce->cstwerr[s] != 0 && ce->cperror[s] == 0 ) {
196     fcstwerr[s] = ce->cstwerr[s];
197     fcperror[s] = ce->cperror[s];
198     for ( Int_t d=0 ; d<11 ;d++ ){
199     Int_t pre = -1;
200     for ( Int_t j=0; j<96 ;j++){
201     if ( j%16 == 0 ) pre++;
202     fcalped[s][d][j] = ce->calped[e][d][j];
203     fcalgood[s][d][j] = ce->calgood[e][d][j];
204     fcalthr[s][d][pre] = ce->calthr[e][d][pre];
205     fcalrms[s][d][j] = ce->calrms[e][d][j];
206     fcalbase[s][d][pre] = ce->calbase[e][d][pre];
207     fcalvar[s][d][pre] = ce->calvar[e][d][pre];
208     };
209     };
210     } else {
211     printf(" CALORIMETER - ERROR: problems finding a good calibration in this file! \n\n ");
212     fcrfile->Close();
213     return(-111);
214     };
215     fcrfile->Close();
216     delete ce;
217     return(0);
218     }
219    
220     void PamVMCCaloDig::DigitizeCALOCALIB() {
221    
222    
223     // Header of the four sections
224    
225     fSecCalo[0] = 0xAA00; // XE
226     fSecCalo[1] = 0xB100; // XO
227     fSecCalo[2] = 0xB600; // YE
228     fSecCalo[3] = 0xAD00; // YO
229    
230     // length of the data is 0x1215
231    
232     fSecCaloLength[0] = 0x1215; // XE
233     fSecCaloLength[1] = 0x1215; // XO
234     fSecCaloLength[2] = 0x1215; // YE
235     fSecCaloLength[3] = 0x1215; // YO
236    
237     Short_t CRC;
238    
239     for(Int_t sec=0; sec < 4; sec++){
240     //
241     // sec = 0 -> XE 1 -> XO 2-> YE 3 -> YO
242     //
243     fData.clear();
244    
245     fData.push_back(fSecCalo[sec]);
246     fData.push_back(fSecCaloLength[sec]);
247    
248     if ( sec==1 ){
249     FillCalPedBack(&fData,sec);
250     FillCalThrBack(&fData,sec);
251     FillCalRmsBack(&fData,sec);
252     FillCalBaseVarBack(&fData,sec);
253     } else {
254     FillCalPedNorm(&fData,sec);
255     FillCalThrNorm(&fData,sec);
256     FillCalRmsNorm(&fData,sec);
257     FillCalBaseVarNorm(&fData,sec);
258     }
259     //CRC
260    
261     fData.push_back(0x0000);
262     CRC = 0;
263    
264    
265     USBuffer::const_iterator p = fData.begin();
266     while( p!= fData.end() ){
267     CRC=crc(CRC,*p);
268     p++;
269     }
270     fData.push_back(CRC);
271    
272     //Writing calibrations
273     DigitizePSCU(fData.size()*2,0x18);
274     fraw->WritePSCU(&fDataPSCU);
275     fraw->CopyUShortToBuff(&fData);
276    
277     //No padding for now...
278     }
279     }
280    
281     void PamVMCCaloDig::FillCalPedNorm(USBuffer *buff, Int_t sec){
282    
283     Int_t chksum = 0;
284     for (Int_t plane=0; plane < 11; plane++){
285     for (Int_t strip=0; strip < 96; strip++){
286     chksum += (Int_t)fcalped[sec][plane][strip];
287     buff->push_back((Int_t)fcalped[sec][plane][strip]);
288     buff->push_back((Int_t)fcalgood[sec][plane][strip]);
289     }
290    
291     }//3-2114
292    
293     buff->push_back(((UShort_t)chksum));
294     buff->push_back(0x0000);
295     buff->push_back(((UShort_t)((Int_t)(chksum >> 16))));
296     }
297    
298     void PamVMCCaloDig::FillCalPedBack(USBuffer *buff, Int_t sec){
299    
300     Int_t chksum = 0;
301     for (Int_t plane=0; plane < 11; plane++){
302     for (Int_t strip=95; strip >= 0; strip--){
303     chksum += (Int_t)fcalped[sec][plane][strip];
304     buff->push_back((Int_t)fcalped[sec][plane][strip]);
305     buff->push_back((Int_t)fcalgood[sec][plane][strip]);
306     }
307     }//3-2114
308    
309     buff->push_back(((UShort_t)chksum));
310     buff->push_back(0x0000);
311     buff->push_back(((UShort_t)((Int_t)(chksum >> 16))));
312     }
313    
314    
315     void PamVMCCaloDig::FillCalThrNorm(USBuffer *buff, Int_t sec){
316    
317     Int_t chksum = 0;
318     for (Int_t plane=0; plane < 11; plane++){
319     for (Int_t strip=0; strip < 6; strip++){
320     chksum += (Int_t)fcalthr[sec][plane][strip];
321     buff->push_back(0x0000);
322     buff->push_back((Int_t)fcalthr[sec][plane][strip]);
323     }
324    
325     }//2118-2249
326    
327     buff->push_back(0x0000);
328     buff->push_back(((UShort_t)chksum));
329     buff->push_back(0x0000);
330     buff->push_back(((UShort_t)((Int_t)(chksum >> 16))));
331     }
332    
333     void PamVMCCaloDig::FillCalThrBack(USBuffer *buff, Int_t sec){
334    
335     Int_t chksum = 0;
336     for (Int_t plane=0; plane < 11; plane++){
337     for (Int_t strip=5; strip >= 0; strip--){
338     chksum += (Int_t)fcalthr[sec][plane][strip];
339     buff->push_back(0x0000);
340     buff->push_back((Int_t)fcalthr[sec][plane][strip]);
341     }
342     }//2118-2249
343    
344     buff->push_back(0x0000);
345     buff->push_back(((UShort_t)chksum));
346     buff->push_back(0x0000);
347     buff->push_back(((UShort_t)((Int_t)(chksum >> 16))));
348     }
349    
350     void PamVMCCaloDig::FillCalRmsNorm(USBuffer *buff, Int_t sec){
351    
352     for (Int_t plane=0; plane < 11; plane++){
353     for (Int_t strip=0; strip < 96; strip++){
354     buff->push_back(0x0000);
355     buff->push_back((Int_t)fcalrms[sec][plane][strip]);
356     }
357     }//2254-4365
358     }
359    
360     void PamVMCCaloDig::FillCalRmsBack(USBuffer *buff, Int_t sec){
361    
362     for (Int_t plane=0; plane < 11; plane++){
363     for (Int_t strip=95; strip >= 0; strip--){
364     buff->push_back(0x0000);
365     buff->push_back((Int_t)fcalrms[sec][plane][strip]);
366     }
367     }//2254-4365
368     }
369    
370     void PamVMCCaloDig::FillCalBaseVarNorm(USBuffer *buff, Int_t sec){
371    
372     for (Int_t plane=0; plane < 11; plane++){
373     for (Int_t strip=0; strip < 6; strip++){
374     buff->push_back(0x0000);
375     buff->push_back(((Int_t)fcalbase[sec][plane][strip]));
376     buff->push_back(0x0000);
377     buff->push_back(((Int_t)fcalvar[sec][plane][strip]));
378     }
379     }//4366-4629
380     }
381    
382     void PamVMCCaloDig::FillCalBaseVarBack(USBuffer *buff, Int_t sec){
383    
384     for (Int_t plane=0; plane < 11; plane++){
385     for (Int_t strip=5; strip >= 0; strip--){
386     buff->push_back(0x0000);
387     buff->push_back((Int_t)fcalbase[sec][plane][strip]);
388     buff->push_back(0x0000);
389     buff->push_back((Int_t)fcalvar[sec][plane][strip]);
390     }
391     }//4366-4629
392     }
393    
394     void PamVMCCaloDig::Digitize(){
395    
396     #ifdef DIG_DEBUG
397     cout<<"Digitizing CALO..."<<endl;
398     #endif
399    
400     Int_t ModCalo = 1; // 0 is RAW, 1 is COMPRESS, 2 is FULL
401     //####@@@@ should be given as input par @@@@####
402    
403    
404     // call different routines depending on the acq mode you want to simulate
405    
406     switch ( ModCalo ){
407     case 0:
408     this->DigitizeCaloRaw();
409     break;
410     case 1:
411     this->DigitizeCaloCompress();
412     break;
413     case 2:
414     this->DigitizeCaloFull();
415     break;
416     };
417     }
418    
419     void PamVMCCaloDig::DigitizeCaloRaw(){
420    
421    
422     // some variables
423     //
424     Float_t ens = 0.;
425     UInt_t adcsig = 0;
426     UInt_t adcbase = 0;
427     UInt_t adc = 0;
428     Int_t pre = 0;
429     UInt_t l = 0;
430     UInt_t lpl = 0;
431     UInt_t tstrip = 0;
432     UInt_t fSecPointer = 0;
433     Int_t fCALOlength;
434     Double_t pedenoise;
435     Float_t rms = 0.;
436     Float_t pedestal = 0.;
437     //
438     // clean the data structure
439     //
440    
441     UShort_t DataCALO[4264];
442     fData.clear();
443    
444     memset(DataCALO,0,sizeof(UShort_t)*4264);
445    
446     static const Float_t CALOGeV2MIPratio = 0.0001059994;
447    
448     // Header of the four sections
449    
450     fSecCalo[0] = 0xEA08; // XE
451     fSecCalo[1] = 0xF108; // XO
452     fSecCalo[2] = 0xF608; // YE
453     fSecCalo[3] = 0xED08; // YO
454    
455     // length of the data is 0x0428 in RAW mode
456    
457     fSecCaloLength[0] = 0x0428; // XE
458     fSecCaloLength[1] = 0x0428; // XO
459     fSecCaloLength[2] = 0x0428; // YE
460     fSecCaloLength[3] = 0x0428; // YO
461     // let's start
462     //
463     fCALOlength = 0;
464     //
465     for (Int_t sec=0; sec < 4; sec++){
466     //
467     // sec = 0 -> XE 1 -> XO 2-> YE 3 -> YO
468     //
469     l = 0; // XE and XO are Y planes
470     if ( sec < 2 ) l = 1; // while YE and YO are X planes
471     //
472     fSecPointer = fCALOlength;
473     //
474     // First of all we have section header and packet length
475     //
476     DataCALO[fCALOlength] = fSecCalo[sec]; //1
477     fCALOlength++;
478     DataCALO[fCALOlength] = fSecCaloLength[sec]; //2
479     fCALOlength++;
480     //
481     // selftrigger coincidences - in the future we should add here some code
482     // to simulate timing response of pre-amplifiers
483     //
484     for (Int_t autoplane=0; autoplane < 7; autoplane++){
485     DataCALO[fCALOlength] = 0x0000;
486     fCALOlength++;
487     }; //3-8
488     //
489     //
490     // here comes data
491     //
492     //
493     // Section XO is read in the opposite direction respect to the others
494     //
495     if ( sec == 1 ){
496     tstrip = 96*11 + fCALOlength; //tstrip = 1064
497     } else {
498     tstrip = 0;
499     };
500     //
501     pre = -1;
502     //
503     for (Int_t strip=0; strip < 96; strip++){
504     //
505     // which is the pre for this strip?
506     //
507     if (strip%16 == 0) {
508     pre++;
509     };
510     //
511     if ( sec == 1 ) tstrip -= 11;
512     //
513     for (Int_t plane=0; plane < 11; plane++){
514    
515     if ( sec == 0 || sec == 3 ) lpl = plane * 2;
516     if ( sec == 1 || sec == 2 ) lpl = (plane * 2) + 1;
517     //
518     // get the energy in GeV from the simulation for that strip
519     //
520     ens = this->GetCaloErel(sec,plane,strip);
521     //
522     // convert it into ADC channels
523     //
524     adcsig = int(ens*fCalomip[l][lpl][strip]/CALOGeV2MIPratio);
525     //
526     // sum baselines
527     //
528     adcbase = (UInt_t)fcalbase[sec][plane][pre];
529     //
530     // add noise and pedestals
531     //
532     pedestal = fcalped[sec][plane][strip];
533     rms = fcalrms[sec][plane][strip]/4.;
534     //
535     // Add random gaussian noise of RMS rms and Centered in the pedestal
536     //
537     pedenoise = frandom->Gaus((Double_t)pedestal,(Double_t)rms);
538     //
539     // Sum all contribution
540     //
541     adc = adcsig + adcbase + (Int_t)round(pedenoise);
542     //
543     // Signal saturation
544     //
545     if ( adc > 0x7FFF ) adc = 0x7FFF;
546     //
547     // save value
548     //
549     if ( sec == 1 ){
550     DataCALO[tstrip] = adc;
551     tstrip++;
552     } else {
553     DataCALO[fCALOlength] = adc;
554     };
555     fCALOlength++;
556     //
557     };
558     //
559     if ( sec == 1 ) tstrip -= 11;
560     //
561     };
562     //
563     // here we calculate and save the CRC
564     //
565     Short_t CRC = 0;
566     for (UInt_t i=0; i<(fCALOlength-fSecPointer); i++){
567     CRC=crc(CRC,DataCALO[i+fSecPointer]);
568     };
569     DataCALO[fCALOlength] = (UShort_t)CRC;
570     fCALOlength++;
571     };
572    
573     for (Int_t i = 0; i<fCALOlength; i++) fData.push_back(DataCALO[i]);
574     //reverse(fData.begin(),fData.end());
575     //for (Int_t i = 0; i<fCALOlength; i++) fData.push_back(DataCALO[fCALOlength-i-1]);
576     }
577    
578    
579     Float_t PamVMCCaloDig::GetCaloErel(Int_t sec, Int_t plane, Int_t strip){
580    
581     // determine plane and strip
582    
583     Int_t mplane = 0;
584    
585     switch (sec){
586     case 0: //3
587     mplane = plane * 4 ; // it must be 0, 4, 8, ... (+1) from plane = 0, 11 YO
588     break;
589     case 1: //2
590     mplane = plane * 4 + 2; // it must be 2, 6, 10, ... (+1) from plane = 0, 11 YE
591     break;
592     case 2: //1
593     mplane = plane * 4 + 3 ; // it must be 3, 7, 11, ... (+1) from plane = 0, 11 XO
594     break;
595     case 3: //0
596     mplane = plane * 4 + 1 ; // it must be 1, 5, 9, ... (+1) from plane = 0, 11 XE
597     break;
598     }
599    
600    
601     Float_t Erel = 0.;
602    
603     // search energy release in VMC output
604     if (fhc){
605     PamVMCDetectorHit * hit = (PamVMCDetectorHit*)fhc->At(mplane*96+strip);
606     if (hit) Erel=hit->GetEREL();
607     // if(Erel) cout<<"sec:"<<sec<<" plane:"<<plane<<" mplane:"<<mplane<<" AT:"<<mplane*96+strip<<" POS:"<<hit->GetPOS()<<endl;
608     } //else cout<<"CALO HIT Collection pointer not found"<<endl;
609    
610     // if ((sec == 2) && (plane == 0)) Erel = 0.01; else Erel = 0.;
611    
612     return Erel;
613     }
614    
615    
616    
617     void PamVMCCaloDig::DigitizeCaloCompress(){
618     //
619     // CompressMode implemented by C.Pizzolotto october 2009
620     //
621     // some variables
622     //
623     Float_t ens = 0.;
624     UInt_t adcsig = 0;
625     UInt_t adcbase = 0;
626     UInt_t adc[16];
627     Float_t rms = 0.;
628     Double_t pedenoise=0.;
629     Float_t pedestal = 0.;
630     UInt_t pedround[16];
631     Float_t thres[16];
632     Float_t goodflag[16];
633     UInt_t min_adc = 0x7FFF;
634     UInt_t min_adc_ch = 0;
635     UInt_t l = 0;
636     UInt_t lpl = 0;
637     Int_t plane = 0;
638     Int_t pre;
639     Int_t npre = 0; // number of pre between 0-5
640     UInt_t strip = 0;
641     UInt_t remainder;
642     Float_t basesum=0.;
643     Float_t basenof=0.;
644     UInt_t baseline=0;
645     UInt_t fSecPointer = 0;
646     UInt_t fNofTStripsPointer = 0;
647     UInt_t NofTransmittedStrips = 0 ;
648     Int_t fCALOlength;
649     UShort_t DataCALO[9040]; //TOO LONG? 4264ma non e' vero che e' cosi' lungo...... CECI CECI CECI
650     fData.clear();
651     static const Float_t CALOGeV2MIPratio = 0.0001059994;
652     //
653     // clean the data structure
654     //
655     memset(adc, 0,sizeof(adc));
656     memset(pedround, 0,sizeof(pedround));
657     memset(thres, 0,sizeof(thres));
658     memset(goodflag, 0,sizeof(goodflag));
659     //
660     memset(DataCALO,0,sizeof(UShort_t)*9040);
661     //
662     // Header of the four sections
663     //
664     fSecCalo[0] = 0xEA00; // XE
665     fSecCalo[1] = 0xF100; // XO
666     fSecCalo[2] = 0xF600; // YE
667     fSecCalo[3] = 0xED00; // YO
668     //
669     // here comes raw data
670     //
671     fCALOlength = 0;
672     //
673     for (Int_t sec=0; sec < 4; sec++){ //
674     //
675     // sec = 0 -> XE 1 -> XO 2-> YE 3 -> YO
676     //
677     l = 0; // XE and XO are Y planes
678     if ( sec < 2 ) l = 1; // while YE and YO are X planes
679     //
680     fSecPointer = fCALOlength;
681     //
682     // First of all we have section header and packet length
683     //
684     DataCALO[fCALOlength] = fSecCalo[sec];
685     fCALOlength++;
686     DataCALO[fCALOlength] = 0; // Unknown: length must be calculated on fly
687     fCALOlength++;
688     //
689     // selftrigger coincidences - in the future we should add here some code to simulate timing response of pre-amplifiers
690     //
691     for (Int_t autoplane=0; autoplane < 7; autoplane++){
692     DataCALO[fCALOlength] = 0x0000;
693     fCALOlength++;
694     };
695     //
696     // second level trigger
697     //
698     DataCALO[fCALOlength] = 0x0000;
699     fCALOlength++;
700     //
701     // Nof strips transmetted: must be calculated on fly
702     //
703     fNofTStripsPointer = fCALOlength;
704     DataCALO[fCALOlength] = 0x0000;
705     fCALOlength++;
706     NofTransmittedStrips=0;
707     //
708     // Identifier of calo data
709     //
710     DataCALO[fCALOlength] = 0xCA50;
711     fCALOlength++;
712     DataCALO[fCALOlength] = 0xCA50;
713     fCALOlength++;
714     DataCALO[fCALOlength] = 0xFFFF; // compresso
715     fCALOlength++;
716     //
717     // Pedestal threashold table checksum
718     //
719     DataCALO[fCALOlength] = 0x0000;
720     fCALOlength++;
721     //
722     // Calorimeter event counter
723     //
724     DataCALO[fCALOlength] = Getevtcalo() ;
725     fCALOlength++;
726     //cout<<" evtcalo?"<<Getevtcalo()<<endl;
727     //
728     // Start here with data
729     //
730     plane=-1;
731     npre =-1;
732     for (Int_t ipre=0; ipre< 66; ipre++){ // (11 planes*6 preampl)
733     //
734     // which plane
735     if ( (ipre % 6) == 0) {
736     plane++;
737     }
738     //
739     pre=ipre;
740     //
741     // Adjust counter for plane X0
742     if (sec==1) // conto invertito
743     {
744     remainder = pre % 6 ;
745     pre = ((plane+1)*6) - remainder ;
746     }
747     //
748     if ( sec == 0 || sec == 3 ) lpl = plane * 2;
749     if ( sec == 1 || sec == 2 ) lpl = (plane * 2) + 1;
750     //
751     // initialize min_adc
752     min_adc = 0x7FFF;
753     for (Int_t ch=0; ch <16; ch++){ // 16 channels each pre
754     //
755     // strip number
756     //
757     strip=((pre-(6*plane))*16)+ch;
758     if(sec==1) strip = ((pre-(6*plane))*16)+(15-ch)-16;
759     //
760     // calculate npre a number between 0-5
761     //
762     if( sec==1) {
763     if ( ((95-strip) % 16) == 0) {
764     npre++;
765     if(npre>5) npre=0;
766     }
767     } else {
768     if ( (strip % 16) == 0) {
769     npre++;
770     if(npre>5) npre=0;
771     }
772     }
773     //
774     ens = this->GetCaloErel(sec,plane,strip);
775     //
776     // convert it into ADC channels
777     //
778     adcsig = int(ens*fCalomip[l][lpl][strip]/CALOGeV2MIPratio);
779     //
780     // sum baselines
781     //
782     adcbase = (UInt_t)fcalbase[sec][plane][npre];
783     //
784     // add noise and pedestals
785     //
786     pedestal = fcalped[sec][plane][strip];
787     rms = fcalrms[sec][plane][strip]/4.;
788     //
789     // Add random gaussian noise of RMS rms and Centered in the pedestal
790     //
791     pedenoise = frandom->Gaus((Double_t)pedestal,(Double_t)rms);
792     //
793     // Sum all contribution
794     //
795     adc[ch] = adcsig + adcbase + (Int_t)round(pedenoise);
796     //
797     // Signal saturation
798     //
799     if ( adc[ch] > 0x7FFF ) adc[ch] = 0x7FFF;
800     //
801     // save infos
802     //
803     pedround[ch] = (Int_t)round(pedestal) ;
804     thres[ch] = ( fcalthr[sec][plane][npre] );
805     goodflag[ch] = ( fcalgood[sec][plane][strip] ); // if bad should be 255
806     //
807     // Find minimum adc in this preamp
808     //
809     if ( goodflag[ch]==0 && (adc[ch]-pedround[ch])<min_adc )
810     {
811     min_adc = ( adc[ch]-pedround[ch] ) ;
812     min_adc_ch = ch ;
813     }
814     }; // close channel loop ch
815     //
816     // Find how many channels are below threshold
817     //
818     Int_t nof_chs_below = 0;
819     for (Int_t ch=0; ch <16; ch++){ // 16 channels each pre
820     if ( goodflag[ch]==0 && ((adc[ch]-pedround[ch]) < (min_adc+thres[min_adc_ch])) )
821     nof_chs_below++;
822     };
823     //
824     // Transmit data: CASE nof_chs_below<9
825     //
826     if(nof_chs_below<9)
827     {
828     if(sec==1) {
829     DataCALO[fCALOlength] = 0x1000 + ipre ;
830     } else {
831     DataCALO[fCALOlength] = 0x1000 + pre ;
832     }
833     fCALOlength++;
834     for (Int_t ch=0; ch <16; ch++)
835     {
836     DataCALO[fCALOlength] = adc[ch];
837     fCALOlength++;
838     NofTransmittedStrips++;
839     };
840     }
841     else
842     //
843     // Transmit data: CASE nof_chs_below>=9
844     //
845     {
846     if(sec==1) {
847     DataCALO[fCALOlength] = 0x0800 + ipre ;
848     } else {
849     DataCALO[fCALOlength] = 0x0800 + pre;
850     }
851     fCALOlength++;
852     //
853     // calculate baseline and save it
854     //
855     basenof=0;
856     baseline=0;
857     basesum=0;
858     for (Int_t ch=0; ch <16; ch++){
859     if( goodflag[ch]==0 && ( (adc[ch]-pedround[ch])<(min_adc+thres[ch]) ) )
860     {
861     basesum = basesum + (adc[ch]-pedround[ch]) ;
862     basenof++;
863     }
864     };
865     baseline = (Int_t)round( basesum / basenof );
866     DataCALO[fCALOlength] = baseline;
867     fCALOlength++;
868     //
869     // Transmit only channels > (min_adc+thres[ch])
870     //
871     for (Int_t ch=0; ch <16; ch++){
872     if ( (adc[ch]-pedround[ch] )>(min_adc+thres[ch]) )
873     {
874     DataCALO[fCALOlength] = ch;
875     fCALOlength++;
876     DataCALO[fCALOlength] = adc[ch];
877     fCALOlength++;
878     NofTransmittedStrips++;
879     }
880     };
881     } // close if nof_chs_below
882     }; // close preampl loop
883     //
884     // Write the length
885     //
886     DataCALO[fSecPointer+1] = (fCALOlength-fSecPointer+1)-2 ;
887     // total length of the packet: -2: because the words with status and length are not included
888     DataCALO[fNofTStripsPointer] = NofTransmittedStrips ;
889     //
890     // here we calculate and save the CRC
891     //
892     Short_t CRC = 0;
893     for (UInt_t i=0; i<(fCALOlength-fSecPointer); i++){
894     CRC=crc(CRC,DataCALO[i+fSecPointer]);
895     };
896     DataCALO[fCALOlength] = (UShort_t)CRC;
897     fCALOlength++;
898     //
899     }; // close sec loop
900    
901     Incrementevtcalo();
902    
903     for (Int_t i = 0; i<fCALOlength; i++) fData.push_back(DataCALO[i]);
904    
905     return;
906    
907    
908    
909    
910     }
911    
912     void PamVMCCaloDig::DigitizeCaloFull(){
913    
914     cout<<"!!!WARNING PamVMCCaloDig FULL MODE STILL NOT IMPLEMENTED!!!"<<endl;
915    
916     this->DigitizeCaloRaw();
917     return;
918    
919     }

  ViewVC Help
Powered by ViewVC 1.1.23