/[PAMELA software]/DarthVader/CalorimeterLevel2/src/CaloLevel0.cpp
ViewVC logotype

Annotation of /DarthVader/CalorimeterLevel2/src/CaloLevel0.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.2 - (hide annotations) (download)
Mon Mar 26 14:02:06 2007 UTC (17 years, 8 months ago) by mocchiut
Branch: MAIN
CVS Tags: v3r04, v3r03
Changes since 1.1: +40 -13 lines
Bug in CaloStrip/Level0 fixed + new alignement parameters added

1 mocchiut 1.1 /**
2     * \file src/CaloLevel0.cpp
3     * \author Emiliano Mocchiutti
4     **/
5     //
6     // C/C++ headers
7     //
8     #include <sstream>
9     #include <fstream>
10     //
11     // ROOT headers
12     //
13     #include <TTree.h>
14     #include <TBranch.h>
15     #include <TFile.h>
16     #include <TObject.h>
17     //
18     // YODA headers
19     //
20     #include <PamelaRun.h>
21     #include <physics/calorimeter/CalorimeterEvent.h>
22     #include <CalibCalPedEvent.h>
23     //
24     //
25     //
26     #include <GLTables.h>
27     //
28     // this package headers
29     //
30     #include <delay.h>
31     #include <CaloLevel0.h>
32     //
33     //
34     // Declaration of the core fortran routines
35     //
36     #define calol2cm calol2cm_
37     extern "C" int calol2cm();
38     #define calol2tr calol2tr_
39     extern "C" int calol2tr();
40     //
41     using namespace std;
42     //
43     //
44     // Public methods
45     //
46    
47     CaloLevel0::~CaloLevel0(){
48     if ( de ) delete de;
49     delete this;
50     }
51    
52     CaloLevel0::CaloLevel0(){
53     //
54     extern struct FlCaLevel1 clevel1_;
55     extern struct FlCaLevel2 clevel2_;
56     clevel1 = &clevel1_;
57     clevel2 = &clevel2_;
58     //
59     trkseqno = 0;
60     ClearStructs();
61     //
62     memset(dexy, 0, 2*22*96*sizeof(Float_t));
63     memset(dexyc, 0, 2*22*96*sizeof(Float_t));
64     memset(mip, 0, 2*22*96*sizeof(Float_t));
65     memset(base, 0, 2*22*6*sizeof(Float_t));
66     memset(sbase, 0, 2*22*6*sizeof(Float_t));
67     calopar1 = true;
68     calopar2 = true;
69     calopar3 = true;
70     crosst = true;
71     ftcalopar1 = 0;
72     ttcalopar1 = 0;
73     ftcalopar2 = 0;
74     ttcalopar2 = 0;
75     ftcalopar3 = 0;
76     ttcalopar3 = 0;
77     }
78    
79     void CaloLevel0::SetCrossTalk(Bool_t ct){
80     crosst = ct;
81     };
82    
83     void CaloLevel0::SetVerbose(Bool_t ct){
84     verbose = ct;
85     };
86    
87     /**
88     * Initialize CaloLevel0 object
89     **/
90     void CaloLevel0::ProcessingInit(TSQLServer *dbc, UInt_t hs, Int_t &sgnl, TTree *l0tree, Bool_t isdeb, Bool_t isverb){
91     //
92     debug = isdeb;
93     verbose = isverb;
94     //
95     l0tr=(TTree*)l0tree;
96     de = new pamela::calorimeter::CalorimeterEvent();
97     l0calo = (TBranch*)l0tr->GetBranch("Calorimeter");
98     l0tr->SetBranchAddress("Calorimeter", &de);
99     //
100     trkseqno = 0;
101     ClearStructs();
102     //
103     GL_CALO_CALIB *glcalo = new GL_CALO_CALIB();
104     //
105     sgnl = 0;
106     UInt_t uptime = 0;
107     //
108     for (Int_t s = 0; s < 4; s++){
109     idcalib[s] = 0;
110     fromtime[s] = 0;
111     totime[s] = 0;
112     calibno[s] = 0;
113     ClearCalibVals(s);
114     //
115     sgnl = glcalo->Query_GL_CALO_CALIB(hs,uptime,s,dbc);
116     if ( sgnl < 0 ){
117     if ( verbose ) printf(" CALORIMETER - ERROR: error from GLTables\n");
118     return;
119     };
120     //
121     idcalib[s] = glcalo->ID_ROOT_L0;
122     fromtime[s] = glcalo->FROM_TIME;
123     if ( glcalo->TO_TIME < hs ){ // calibration is corrupted and we are using the one that preceed the good one
124     totime[s] = uptime;
125     } else {
126     totime[s] = glcalo->TO_TIME;
127     };
128     calibno[s] = glcalo->EV_ROOT;
129     //
130     if ( totime[s] == 0 ){
131     if ( verbose ) printf(" CALORIMETER - WARNING: data with no associated calibration\n");
132     ClearCalibVals(s);
133     sgnl = 100;
134     };
135     };
136     //
137     // determine path and name and entry of the calibration file
138     //
139     GL_ROOT *glroot = new GL_ROOT();
140     if ( verbose ) printf("\n");
141     for (Int_t s = 0; s < 4; s++){
142     if ( verbose ) printf(" ** SECTION %i **\n",s);
143     if ( totime[s] > 0 ){
144     //
145     sgnl = glroot->Query_GL_ROOT(idcalib[s],dbc);
146     if ( sgnl < 0 ){
147     if ( verbose ) printf(" CALORIMETER - ERROR: error from GLTables\n");
148     return;
149     };
150     //
151     stringstream name;
152     name.str("");
153     name << glroot->PATH.Data() << "/";
154     name << glroot->NAME.Data();
155     //
156     fcalname[s] = (TString)name.str().c_str();
157     if ( verbose ) printf(" - runheader at time %u. From time %u to time %u \n use file %s \n calibration at entry %i \n\n",hs,fromtime[s],totime[s],fcalname[s].Data(),calibno[s]);
158     } else {
159     if ( verbose ) printf(" - runheader at time %u. NO CALIBRATION INCLUDE THE RUNHEADER! ",hs);
160     };
161     sgnl = LoadCalib(s);
162     if ( sgnl ) break;
163     };
164     //
165     delete glcalo;
166     delete glroot;
167     //
168     return;
169     //
170     }
171    
172     Int_t CaloLevel0::ChkCalib(TSQLServer *dbc, UInt_t atime){
173     Int_t sgnl = 0;
174     for ( Int_t s = 0; s < 4; s++){
175     if ( atime > totime[s] ){
176     if ( !dbc->IsConnected() ) throw -116;
177     sgnl = Update(dbc,atime,s);
178     if ( sgnl < 0 ) return(sgnl);
179     };
180     };
181     return(sgnl);
182     }
183    
184 mocchiut 1.2 Int_t CaloLevel0::ChkParam(TSQLServer *dbc, UInt_t runheader, Bool_t mechal){
185 mocchiut 1.1 stringstream calfile;
186     stringstream bmfile;
187     stringstream aligfile;
188     Int_t error = 0;
189     FILE *f = 0;
190     ifstream badfile;
191     GL_PARAM *glparam = new GL_PARAM();
192     //
193     if ( calopar1 || ( ttcalopar1 != 0 && ttcalopar1 < runheader ) ){
194 mocchiut 1.2 //
195     //
196     //
197     if ( debug ) printf(" calopar1 %i ftcalopar1 %u ttcalopar1 %u runheader %u \n",calopar1,ftcalopar1,ttcalopar1,runheader);
198     //
199 mocchiut 1.1 calopar1 = false;
200     //
201     // determine where I can find calorimeter ADC to MIP conversion file
202     //
203     if ( verbose ) printf(" Querying DB for calorimeter parameters files...\n");
204     //
205     error = 0;
206     error = glparam->Query_GL_PARAM(runheader,101,dbc);
207     if ( error < 0 ) return(error);
208     //
209     calfile.str("");
210     calfile << glparam->PATH.Data() << "/";
211     calfile << glparam->NAME.Data();
212     ftcalopar1 = glparam->FROM_TIME;
213     ttcalopar1 = glparam->TO_TIME;
214     //
215     if ( verbose ) printf("\n Using ADC to MIP conversion file: \n %s \n",calfile.str().c_str());
216     f = fopen(calfile.str().c_str(),"rb");
217     if ( !f ){
218     if ( verbose ) printf(" CALORIMETER - ERROR: no ADC to MIP file!\n");
219     return(-105);
220     };
221     //
222     for (Int_t m = 0; m < 2 ; m++ ){
223     for (Int_t k = 0; k < 22; k++ ){
224     for (Int_t l = 0; l < 96; l++ ){
225     fread(&mip[m][k][l],sizeof(mip[m][k][l]),1,f);
226     };
227     };
228     };
229     fclose(f);
230     };
231     //
232     if ( calopar2 || ( ttcalopar2 != 0 && ttcalopar2 < runheader ) ){
233 mocchiut 1.2 if ( debug ) printf(" calopar2 %i ftcalopar2 %u ttcalopar2 %u runheader %u \n",calopar2,ftcalopar2,ttcalopar2,runheader);
234 mocchiut 1.1 calopar2 = false;
235     //
236     // determine where I can find calorimeter alignment file
237     //
238     //
239     error = 0;
240     error = glparam->Query_GL_PARAM(runheader,102,dbc);
241     if ( error < 0 ) return(error);
242     //
243     aligfile.str("");
244     aligfile << glparam->PATH.Data() << "/";
245     aligfile << glparam->NAME.Data();
246     ftcalopar2 = glparam->FROM_TIME;
247     ttcalopar2 = glparam->TO_TIME;
248     //
249 mocchiut 1.2 if ( verbose ) printf("\n Using parameter file: \n %s \n",aligfile.str().c_str());
250 mocchiut 1.1 f = fopen(aligfile.str().c_str(),"rb");
251     if ( !f ){
252 mocchiut 1.2 if ( verbose ) printf(" CALORIMETER - ERROR: no parameter file!\n");
253 mocchiut 1.1 return(-106);
254     };
255     //
256 mocchiut 1.2 if ( !mechal ){
257     //
258     fread(&clevel1->xalig,sizeof(clevel1->xalig),1,f);
259     if ( debug ) printf(" xalig = %f \n",clevel1->xalig);
260     fread(&clevel1->yalig,sizeof(clevel1->yalig),1,f);
261     if ( debug ) printf(" yalig = %f \n",clevel1->yalig);
262     fread(&clevel1->zalig,sizeof(clevel1->zalig),1,f);
263     if ( debug ) printf(" zalig = %f \n",clevel1->zalig);
264     } else {
265     if ( verbose ) printf("\n Using MECHANICAL alignement parameters \n");
266     //
267     CaloStrip cs = CaloStrip();
268     cs.UseMechanicalAlig();
269     clevel1->xalig = cs.GetXalig();
270     if ( debug ) printf(" xalig = %f \n",clevel1->xalig);
271     clevel1->yalig = cs.GetYalig();
272     if ( debug ) printf(" yalig = %f \n",clevel1->yalig);
273     clevel1->zalig = cs.GetZalig();
274     if ( debug ) printf(" zalig = %f \n",clevel1->zalig);
275     //
276     Float_t tmp = 0;
277     fread(&tmp,sizeof(clevel1->xalig),1,f);
278     fread(&tmp,sizeof(clevel1->yalig),1,f);
279     fread(&tmp,sizeof(clevel1->zalig),1,f);
280     //
281     };
282 mocchiut 1.1 fread(&clevel1->emin,sizeof(clevel1->emin),1,f);
283     if ( debug ) printf(" signal threshold = %f \n",clevel1->emin);
284     //
285     fclose(f);
286     };
287     //
288     // Load offline bad strip mask
289     //
290     if ( calopar3 || ( ttcalopar3 != 0 && ttcalopar3 < runheader ) ){
291 mocchiut 1.2 if ( debug ) printf(" calopar3 %i ftcalopar3 %u ttcalopar3 %u runheader %u \n",calopar3,ftcalopar3,ttcalopar3,runheader);
292 mocchiut 1.1 calopar3 = false;
293     //
294     // determine where I can find calorimeter alignment file
295     //
296     //
297     error = 0;
298     error = glparam->Query_GL_PARAM(runheader,103,dbc);
299     if ( error < 0 ) return(error);
300     //
301     bmfile.str("");
302     bmfile << glparam->PATH.Data() << "/";
303     bmfile << glparam->NAME.Data();
304     ftcalopar3 = glparam->FROM_TIME;
305     ttcalopar3 = glparam->TO_TIME;
306     //
307     if ( verbose ) printf("\n Using bad strip offline mask file: \n %s \n\n",bmfile.str().c_str());
308     badfile.open(bmfile.str().c_str());
309     if ( !badfile ){
310     if ( verbose ) printf(" CALORIMETER - ERROR: no bad strip offline mask file!\n");
311     return(-115);
312     };
313     //
314     Bool_t isdone = false;
315     Int_t bad = 0;
316     Int_t view = 1;
317     Int_t strip = 0;
318     Int_t plane = 21;
319     while ( !isdone ) {
320     badfile >> bad;
321     obadmask[view][plane][strip] = bad;
322     if ( debug && bad ) printf(" SETTING view %i plane %i strip %i BAD = %i \n",view,plane,strip,bad);
323     strip++;
324     if ( strip > 95 ){
325     strip = 0;
326     plane--;
327     if ( plane < 0 ){
328     plane = 21;
329     view--;
330     };
331     if ( view < 0 ) isdone = true;
332     };
333     };
334     //
335     badfile.close();
336     };
337     //
338     delete glparam;
339     //
340     return(0);
341     }
342    
343    
344    
345     void CaloLevel0::FindBaseRaw(Int_t l, Int_t m, Int_t pre){
346     Float_t minstrip = 100000.;
347     Float_t rms = 0.;
348     base[l][m][pre] = 0.;
349     for (Int_t e = pre*16; e < (pre+1)*16 ; e++){
350     if ( calgood[l][m][e] == 0. && obadmask[l][m][e] == 0 && dexy[l][m][e]-calped[l][m][e] < minstrip && dexy[l][m][e] > 0.) {
351     minstrip = dexy[l][m][e]-calped[l][m][e];
352     rms = calthr[l][m][pre];
353     };
354     };
355     if ( minstrip != 100000. ) {
356     Float_t strip6s = 0.;
357     for (Int_t e = pre*16; e < (pre+1)*16 ; e++){
358     if ( (dexy[l][m][e]-calped[l][m][e]) >= minstrip && (dexy[l][m][e]-calped[l][m][e]) <= (minstrip+rms) ) {
359     strip6s += 1.;
360     base[l][m][pre] += (dexy[l][m][e] - calped[l][m][e]);
361     };
362     //
363     // compression
364     //
365     if ( abs((int)(dexy[l][m][e]-calped[l][m][e])) <= (minstrip+rms) ) {
366     dexyc[l][m][e] = 0.;
367     } else {
368     dexyc[l][m][e] = dexy[l][m][e];
369     };
370     };
371     if ( strip6s >= 9. ){
372     Double_t arro = base[l][m][pre]/strip6s;
373     Float_t deci = 1000.*((float)arro - float(int(arro)));
374     if ( deci < 500. ) {
375     arro = double(int(arro));
376     } else {
377     arro = 1. + double(int(arro));
378     };
379     base[l][m][pre] = arro;
380     } else {
381     base[l][m][pre] = 31000.;
382     for (Int_t e = pre*16; e < (pre+1)*16 ; e++){
383     dexyc[l][m][e] = dexy[l][m][e];
384     };
385     };
386     } else {
387     base[l][m][pre] = 31000.;
388     };
389     }
390    
391     Int_t CaloLevel0::Calibrate(Int_t ei){
392     //
393     // get entry ei
394     //
395     l0calo->GetEntry(ei);
396     //
397     // if it was not a selftrigger event, could it ever been a selftrigger event? if so trigty = 3.
398     //
399     Int_t val = 0;
400     Int_t del = 1100;
401     if ( clevel2->trigty != 2. ){
402     Bool_t ck = false;
403     for (Int_t sec = 0; sec < 4; sec++){
404     val = (Int_t)de->calselftrig[sec][6];
405     del = delay(val);
406     if ( del < 1100 ){
407     clevel2->wartrig = 0.;
408     clevel2->trigty = 3.;
409     ck = true;
410     break;
411     };
412     };
413     if ( !ck ) clevel2->wartrig = 100.;
414     } else {
415     Bool_t ck = false;
416     for (Int_t sec = 0; sec < 4; sec++){
417     val = (Int_t)de->calselftrig[sec][6];
418     del = delay(val);
419     if ( del < 1100 ){
420     clevel2->wartrig = 0.;
421     ck = true;
422     };
423     };
424     if ( !ck ) clevel2->wartrig = 100.;
425     };
426     //
427     Int_t se = 5;
428     Int_t done = 0;
429     Int_t pre = -1;
430     Bool_t isCOMP = false;
431     Bool_t isFULL = false;
432     Bool_t isRAW = false;
433     Float_t ener;
434     Int_t doneb = 0;
435     Int_t donec = 0;
436     Int_t ck = 0;
437     Int_t ipre = 0;
438     Int_t ip[3] = {0};
439     Float_t base0, base1, base2;
440     base0 = 0.;
441     base1 = 0.;
442     base2 = 0.;
443     Float_t qpre[6] = {0.,0.,0.,0.,0.,0.};
444     Float_t ene[96];
445     Int_t chdone[4] = {0,0,0,0};
446     Int_t pe = 0;
447     //
448     Float_t ener0 = 0.;
449     Float_t cbase0 = 0.;
450     Bool_t pproblem = false;
451     //
452     Float_t tim = 0.;
453     Int_t plo = 0;
454     Int_t fbi = 0;
455     Int_t cle = 0;
456     //
457     // run over views and planes
458     //
459     for (Int_t l = 0; l < 2; l++){
460     for (Int_t m = 0; m < 22; m++){
461     //
462 mocchiut 1.2 // determine the section number
463 mocchiut 1.1 //
464     se = 5;
465 mocchiut 1.2 if (l == 0 && m%2 == 0) se = 3;
466 mocchiut 1.1 if (l == 0 && m%2 != 0) se = 2;
467 mocchiut 1.2 if (l == 1 && m%2 != 0) se = 1;
468     if (l == 1 && m%2 == 0) se = 0;
469 mocchiut 1.1 //
470     // determine what kind of event we are going to analyze
471     //
472     isCOMP = false;
473     isFULL = false;
474     isRAW = false;
475     if ( de->stwerr[se] & (1 << 16) ) isCOMP = true;
476     if ( de->stwerr[se] & (1 << 17) ) isFULL = true;
477     if ( de->stwerr[se] & (1 << 3) ) isRAW = true;
478     if ( !chdone[se] ){
479     //
480     // check for any error in the event
481     //
482     clevel2->crc[se] = 0;
483     if ( de->perror[se] == 132 ){
484     clevel2->crc[se] = 1;
485     pe++;
486     };
487     clevel2->perr[se] = 0;
488     if ( de->perror[se] != 0 ){
489     clevel2->perr[se] = 1;
490     pe++;
491     };
492     clevel2->swerr[se] = 0;
493     for (Int_t j = 0; j < 7 ; j++){
494     if ( (j != 3) && (de->stwerr[se] & (1 << j)) ){
495     clevel2->swerr[se] = 1;
496     pe++;
497     };
498     };
499     chdone[se] = 1;
500     };
501     if ( clevel2->crc[se] == 0 && (clevel1->good2 == 1 || clevel2->trigty >= 2) ){
502     pre = -1;
503     //
504     for (Int_t nn = 0; nn < 96; nn++){
505     ene[nn] = 0.;
506     dexy[l][m][nn] = de->dexy[l][m][nn] ;
507     dexyc[l][m][nn] = de->dexyc[l][m][nn] ;
508     };
509     //
510     // run over preamplifiers
511     //
512     pre = -1;
513     cbase0 = 0.;
514     for (Int_t i = 0; i < 3; i++){
515     for (Int_t j = 0; j < 2; j++){
516     pre = j + i*2;
517     //
518     // baseline check and calculation
519     //
520     if ( !isRAW ) {
521     base[l][m][pre] = de->base[l][m][pre] ;
522     cbase0 += base[l][m][pre];
523     } else {
524     //
525     // if it is a raw event and we haven't checked
526     // yet, calculate the baseline.
527     //
528     FindBaseRaw(l,m,pre);
529     cbase0 += base[l][m][pre];
530     };
531     };
532     };
533     //
534     // run over strips
535     //
536     pre = -1;
537     ener0 = 0.;
538     for (Int_t i = 0 ; i < 3 ; i++){
539     ip[i] = 0;
540     for (Int_t n = i*32 ; n < (i+1)*32 ; n++){
541     if (n%16 == 0) {
542     ck = 0;
543     done = 0;
544     doneb = 0;
545     donec = 0;
546     pre++;
547     qpre[pre] = 0.;
548     };
549     //
550     // baseline check and calculation
551     //
552     // no suitable new baseline, use old ones!
553     //
554     if ( !done ){
555     if ( (base[l][m][pre] == 31000. || base[l][m][pre] == 0.) ){
556     ck = 1;
557     if (pre%2 == 0) {
558     ip[i] = pre + 1;
559     } else {
560     ip[i] = pre - 1;
561     };
562     if ( (base[l][m][ip[i]] == 31000. || base[l][m][ip[i]] == 0.) ){
563     //
564     ck = 2;
565     if ( sbase[l][m][pre] == 31000. || sbase[l][m][pre] == 0. ) {
566     ck = 3;
567     };
568     };
569     done = 1;
570     };
571     };
572     //
573     // CALIBRATION ALGORITHM
574     //
575     if ( !doneb ){
576     switch (ck) {
577     case 0:
578     base0 = base[l][m][pre];
579     base2 = calbase[l][m][pre];
580     break;
581     case 1:
582     base0 = base[l][m][ip[i]];
583     base2 = calbase[l][m][ip[i]];
584     break;
585     case 2:
586     base0 = sbase[l][m][pre];
587     base2 = calbase[l][m][pre];
588     break;
589     case 3:
590     base0 = calbase[l][m][pre];
591     base2 = calbase[l][m][pre];
592     break;
593     };
594     base1 = calbase[l][m][pre];
595     doneb = 1;
596     };
597     ener = dexyc[l][m][n];
598     ener0 += ener;
599     clevel1->estrip[n][m][l] = 0.;
600     if ( base0>0 && base0 < 30000. ){
601     if ( !donec && (base0 - base1 + base2) != 0. ){
602     sbase[l][m][pre] = base0 - base1 + base2;
603     donec = 1;
604     };
605     if ( ener > 0. ){
606     clevel1->estrip[n][m][l] = (ener - calped[l][m][n] - base0 - base1 + base2)/mip[l][m][n] ;
607     //
608     // OK, now in estrip we have the energy deposit in MIP of all the strips for this event (at the end of loops of course)
609     //
610     qpre[pre] += clevel1->estrip[n][m][l];
611     };
612     };
613     };
614     if ( crosst ){
615     if (ck == 1){
616     if (ip[i]%2 == 0) {
617     ipre = ip[i] + 1;
618     } else {
619     ipre = ip[i] - 1;
620     };
621     for (Int_t j = ipre*16 ; j < (ipre+1)*16 ; j++){
622     clevel1->estrip[j][m][l] += (qpre[ipre] - qpre[ip[i]]) * 0.00478;
623     };
624     };
625     if (ck == 2){
626     for (Int_t j = i*32 ; j < (i+1)*32 ; j++){
627     ipre = j/16 + 1;
628     clevel1->estrip[j][m][l] += qpre[ipre] * 0.00478;
629     };
630     };
631     };
632     };
633     //
634     if ( ener0 == 0. && cbase0 == 0. && !pproblem && clevel2->perr[se] == 0){
635     if ( verbose ) printf(" L0 entry %i : calorimeter power problems! event marked as bad perr %f swerr %X view %i plane %i \n",ei,de->perror[se],de->stwerr[se],l,m);
636     pproblem = true;
637     pe++;
638     };
639     //
640     Int_t j4 = -4;
641     Int_t jjj = -3;
642     Int_t jj = -2;
643     for (Int_t j = 0 ; j < 100 ; j++){
644     jj++;
645     jjj++;
646     j4++;
647     if ( j < 96 ) ene[j] = clevel1->estrip[j][m][l];
648     if ( crosst ){
649     if ( jj >= 0 && jj < 96 ){
650     if ( jj != 0 && jj != 32 && jj != 64 ) ene[jj-1] += -clevel1->estrip[jj][m][l] * 0.01581;
651     if ( jj != 31 && jj != 63 && jj != 95 ) ene[jj+1] += -clevel1->estrip[jj][m][l] * 0.01581;
652     };
653     if ( jjj >= 0 && jjj < 96 ){
654     if ( jjj != 0 && jjj != 32 && jjj != 64 ) clevel1->estrip[jjj-1][m][l] += -ene[jjj] * 0.01581;
655     if ( jjj != 31 && jjj != 63 && jjj != 95 ) clevel1->estrip[jjj+1][m][l] += -ene[jjj] * 0.01581;
656     };
657     };
658     if ( j4 >= 0 && j4 < 96 ){
659     //
660     // NOTICE: THE FOLLOWING LINE EXCLUDE ALL STRIPS FOR WHICH THE RMS*4 IS GREATER THAN 26 !!! <=============== IMPORTANT! =================>
661     //
662     if ( obadmask[l][m][j4] == 1 || clevel1->estrip[j4][m][l] <= clevel1->emin || calrms[l][m][j4] > 26 ){
663     clevel1->estrip[j4][m][l] = 0.;
664     };
665     //
666     // code and save the energy for each strip in svstrip
667     //
668     if ( clevel1->estrip[j4][m][l] > clevel1->emin ){
669     //
670     tim = 100000.;
671     plo = m;
672     fbi = 0;
673     if ( clevel1->estrip[j4][m][l] > 0.99995 ){
674     tim = 10000.;
675     plo = m;
676     fbi = 1;
677     };
678     if ( clevel1->estrip[j4][m][l] > 9.9995 ){
679     tim = 1000.;
680     plo = 22 + m;
681     fbi = 1;
682     };
683     if ( clevel1->estrip[j4][m][l] > 99.995 ){
684     tim = 100.;
685     plo = 22 + m;
686     fbi = 0;
687     };
688     if ( clevel1->estrip[j4][m][l] > 999.95 ){
689     tim = 10.;
690     plo = 44 + m;
691     fbi = 0;
692     };
693     if ( clevel1->estrip[j4][m][l] > 9999.5 ){
694     tim = 1.;
695     plo = 66 + m;
696     fbi = 0;
697     };
698     //
699     cle = (Int_t)lroundf(tim*clevel1->estrip[j4][m][l]);
700     //
701     if ( l == 0 ){
702     //
703     // +-PPSSmmmm.mmmm
704     //
705     svstrip[istrip] = fbi*1000000000 + plo*10000000 + j4*100000 + cle;
706     } else {
707     svstrip[istrip] = -(fbi*1000000000 + plo*10000000 + j4*100000 + cle);
708     };
709     //
710     // if ( ei >= -770 ) printf(" j %i l %i m %i estrip %f \n",j4,l,m,clevel1->estrip[j4][m][l]);
711     // if ( ei >= -770 ) printf(" num lim %i fbi %i tim %f plo %i cle %i \n",numeric_limits<Int_t>::max(),fbi,tim,plo,cle);
712     // if ( ei >= -770 ) printf(" svstrip %i \n",svstrip[istrip]);
713     //
714     istrip++;
715     };
716     };
717     };
718     //
719     } else {
720     for (Int_t nn = 0; nn < 96; nn++){
721     clevel1->estrip[nn][m][l] = 0.;
722     };
723     };
724     };
725     };
726     if ( !pe ){
727     clevel2->good = 1;
728     } else {
729     clevel2->good = 0;
730     };
731     return(0);
732     }
733    
734     void CaloLevel0::GetTrkVar(){
735     calol2tr();
736     }
737    
738     void CaloLevel0::FillTrkVar(CaloLevel2 *ca, Int_t nutrk){
739     //
740     CaloTrkVar *t_ca = new CaloTrkVar();
741     //
742     t_ca->trkseqno = trkseqno;
743     t_ca->ncore = (Int_t)clevel2->ncore;
744     t_ca->qcore = clevel2->qcore;
745     t_ca->noint = (Int_t)clevel2->noint;
746     t_ca->ncyl = (Int_t)clevel2->ncyl;
747     t_ca->qcyl = clevel2->qcyl;
748     t_ca->qtrack = clevel2->qtrack;
749     t_ca->qtrackx = clevel2->qtrackx;
750     t_ca->qtracky = clevel2->qtracky;
751     t_ca->dxtrack = clevel2->dxtrack;
752     t_ca->dytrack = clevel2->dytrack;
753     t_ca->qlast = clevel2->qlast;
754     t_ca->nlast = (Int_t)clevel2->nlast;
755     t_ca->qpre = clevel2->qpre;
756     t_ca->npre = (Int_t)clevel2->npre;
757     t_ca->qpresh = clevel2->qpresh;
758     t_ca->npresh = (Int_t)clevel2->npresh;
759     t_ca->qtr = clevel2->qtr;
760     t_ca->ntr = (Int_t)clevel2->ntr;
761     t_ca->planetot = (Int_t)clevel2->planetot;
762     t_ca->qmean = clevel2->qmean;
763     t_ca->dX0l = clevel2->dX0l;
764     t_ca->qlow = clevel2->qlow;
765     t_ca->nlow = (Int_t)clevel2->nlow;
766     //
767     if ( trkseqno == -1 ){
768     // ca->impx = clevel2->impx;
769     // ca->impy = clevel2->impy;
770     ca->tanx[1] = clevel2->tanx;
771     ca->tany[1] = clevel2->tany;
772     ca->elen = clevel2->elen;
773     ca->selen = clevel2->selen;
774     // memcpy(ca->cibar,clevel2->cibar,sizeof(clevel2->cibar));
775     // memcpy(ca->cbar,clevel2->cbar,sizeof(clevel2->cbar));
776     memcpy(t_ca->tibar,clevel2->cibar,sizeof(clevel2->cibar));
777     memcpy(t_ca->tbar,clevel2->cbar,sizeof(clevel2->cbar));
778     memcpy(ca->planemax,clevel2->planemax,sizeof(clevel2->planemax));
779     ca->varcfit[2] = clevel2->varcfit[0];
780     ca->varcfit[3] = clevel2->varcfit[1];
781     ca->npcfit[2] = clevel2->npcfit[0];
782     ca->npcfit[3] = clevel2->npcfit[1];
783     // memcpy(ca->varcfit,clevel2->varcfit,sizeof(clevel2->varcfit));
784     // memcpy(ca->npcfit,clevel2->npcfit,sizeof(clevel2->npcfit));
785     } else {
786     memcpy(t_ca->tibar,clevel2->tibar,sizeof(clevel2->tibar));
787     memcpy(t_ca->tbar,clevel2->tbar,sizeof(clevel2->tbar));
788     };
789     //
790     //
791     if ( !(ca->CaloTrk) ) ca->CaloTrk = new TClonesArray("CaloTrkVar",1); //ELENA
792     TClonesArray &t = *ca->CaloTrk;
793     new(t[nutrk]) CaloTrkVar(*t_ca);
794     //
795     delete t_ca;
796     //
797     ClearTrkVar();
798     }
799    
800     void CaloLevel0::GetCommonVar(){
801     calol2cm();
802     }
803    
804     void CaloLevel0::FillCommonVar(CaloLevel1 *c1, CaloLevel2 *ca){
805     //
806     ca->good = clevel2->good;
807     if ( clevel2->trigty == 2. ){
808     ca->selftrigger = 1;
809     } else {
810     ca->selftrigger = 0;
811     };
812     //
813     ca->selftrigger += (Int_t)clevel2->wartrig;
814     //
815     memcpy(ca->perr,clevel2->perr,sizeof(clevel2->perr));
816     memcpy(ca->swerr,clevel2->swerr,sizeof(clevel2->swerr));
817     memcpy(ca->crc,clevel2->crc,sizeof(clevel2->crc));
818     ca->nstrip = (Int_t)clevel2->nstrip;
819     ca->qtot = clevel2->qtot;
820     // ca->impx = clevel2->impx;
821     // ca->impy = clevel2->impy;
822     ca->tanx[0] = clevel2->tanx;
823     ca->tany[0] = clevel2->tany;
824     ca->nx22 = (Int_t)clevel2->nx22;
825     ca->qx22 = clevel2->qx22;
826     ca->qmax = clevel2->qmax;
827     ca->elen = clevel2->elen;
828     ca->selen = clevel2->selen;
829     memcpy(ca->qq,clevel2->qq,sizeof(clevel2->qq));
830     memcpy(ca->planemax,clevel2->planemax,sizeof(clevel2->planemax));
831     ca->varcfit[0] = clevel2->varcfit[0];
832     ca->varcfit[1] = clevel2->varcfit[1];
833     ca->npcfit[0] = clevel2->npcfit[0];
834     ca->npcfit[1] = clevel2->npcfit[1];
835     ca->fitmode[0] = clevel2->fmode[0];
836     ca->fitmode[1] = clevel2->fmode[1];
837     // memcpy(ca->varcfit,clevel2->varcfit,sizeof(clevel2->varcfit));
838     // memcpy(ca->npcfit,clevel2->npcfit,sizeof(clevel2->npcfit));
839     memcpy(ca->cibar,clevel2->cibar,sizeof(clevel2->cibar));
840     memcpy(ca->cbar,clevel2->cbar,sizeof(clevel2->cbar));
841     //
842     if ( c1 ){
843     c1->istrip = istrip;
844     c1->estrip = TArrayI(istrip,svstrip);
845     };
846     //
847     }
848    
849     void CaloLevel0::ClearStructs(){
850     ClearTrkVar();
851     ClearCommonVar();
852     }
853    
854     void CaloLevel0::RunClose(){
855     l0tr->Delete();
856     ClearStructs();
857     //
858     memset(dexy, 0, 2*22*96*sizeof(Float_t));
859     memset(dexyc, 0, 2*22*96*sizeof(Float_t));
860     memset(base, 0, 2*22*6*sizeof(Float_t));
861     memset(sbase, 0, 2*22*6*sizeof(Float_t));
862     //
863     }
864    
865     //
866     // Private methods
867     //
868    
869     void CaloLevel0::ClearTrkVar(){
870     clevel2->ncore = 0;
871     clevel2->qcore = 0.;
872     clevel2->noint = 0.;
873     clevel2->ncyl = 0.;
874     clevel2->qcyl = 0.;
875     clevel2->qtrack = 0.;
876     clevel2->qtrackx = 0.;
877     clevel2->qtracky = 0.;
878     clevel2->dxtrack = 0.;
879     clevel2->dytrack = 0.;
880     clevel2->qlast = 0.;
881     clevel2->nlast = 0.;
882     clevel2->qpre = 0.;
883     clevel2->npre = 0.;
884     clevel2->qpresh = 0.;
885     clevel2->npresh = 0.;
886     clevel2->qlow = 0.;
887     clevel2->nlow = 0.;
888     clevel2->qtr = 0.;
889     clevel2->ntr = 0.;
890     clevel2->planetot = 0.;
891     clevel2->qmean = 0.;
892     clevel2->dX0l = 0.;
893     clevel2->elen = 0.;
894     clevel2->selen = 0.;
895     memset(clevel1->al_p, 0, 5*2*sizeof(Double_t));
896     memset(clevel2->tibar, 0, 2*22*sizeof(Int_t));
897     memset(clevel2->tbar, 0, 2*22*sizeof(Float_t));
898     }
899    
900     void CaloLevel0::ClearCommonVar(){
901     istrip = 0;
902     clevel2->trigty = -1.;
903     clevel2->wartrig = 0.;
904     clevel2->good = 0;
905     clevel2->nstrip = 0.;
906     clevel2->qtot = 0.;
907     // clevel2->impx = 0.;
908     // clevel2->impy = 0.;
909     clevel2->tanx = 0.; // this is correct since it refers to the fortran structure
910     clevel2->tany = 0.; // this is correct since it refers to the fortran structure
911     clevel2->qmax = 0.;
912     clevel2->nx22 = 0.;
913     clevel2->qx22 = 0.;
914     memset(clevel2->perr, 0, 4*sizeof(Int_t));
915     memset(clevel2->swerr, 0, 4*sizeof(Int_t));
916     memset(clevel2->crc, 0, 4*sizeof(Int_t));
917     memset(clevel2->qq, 0, 4*sizeof(Int_t));
918     memset(clevel2->varcfit, 0, 4*sizeof(Float_t));
919     memset(clevel2->npcfit, 0, 4*sizeof(Int_t));
920     memset(clevel2->planemax, 0, 2*sizeof(Int_t));
921     memset(clevel2->fmode, 0, 2*sizeof(Int_t));
922     memset(clevel2->cibar, 0, 2*22*sizeof(Int_t));
923     memset(clevel2->cbar, 0, 2*22*sizeof(Float_t));
924     }
925    
926     void CaloLevel0::ClearCalibVals(Int_t s){
927     //
928     for ( Int_t d=0 ; d<11 ;d++ ){
929     Int_t pre = -1;
930     for ( Int_t j=0; j<96 ;j++){
931     if ( j%16 == 0 ) pre++;
932     if ( s == 2 ){
933     calped[0][2*d+1][j] = 0.;
934     cstwerr[3] = 0.;
935     cperror[3] = 0.;
936     calgood[0][2*d+1][j] = 0.;
937     calthr[0][2*d+1][pre] = 0.;
938     calrms[0][2*d+1][j] = 0.;
939     calbase[0][2*d+1][pre] = 0.;
940     calvar[0][2*d+1][pre] = 0.;
941     };
942     if ( s == 3 ){
943     calped[0][2*d][j] = 0.;
944     cstwerr[1] = 0.;
945     cperror[1] = 0.;
946     calgood[0][2*d][j] = 0.;
947     calthr[0][2*d][pre] = 0.;
948     calrms[0][2*d][j] = 0.;
949     calbase[0][2*d][pre] = 0.;
950     calvar[0][2*d][pre] = 0.;
951     };
952     if ( s == 0 ){
953     calped[1][2*d][j] = 0.;
954     cstwerr[0] = 0.;
955     cperror[0] = 0.;
956     calgood[1][2*d][j] = 0.;
957     calthr[1][2*d][pre] = 0.;
958     calrms[1][2*d][j] = 0.;
959     calbase[1][2*d][pre] = 0.;
960     calvar[1][2*d][pre] = 0.;
961     };
962     if ( s == 1 ){
963     calped[1][2*d+1][j] = 0.;
964     cstwerr[2] = 0.;
965     cperror[2] = 0.;
966     calgood[1][2*d+1][j] = 0.;
967     calthr[1][2*d+1][pre] = 0.;
968     calrms[1][2*d+1][j] = 0.;
969     calbase[1][2*d+1][pre] = 0.;
970     calvar[1][2*d+1][pre] = 0.;
971     };
972     };
973     };
974     return;
975     }
976    
977     Int_t CaloLevel0::Update(TSQLServer *dbc, UInt_t atime, Int_t s){
978     //
979     Int_t sgnl = 0;
980     //
981     GL_CALO_CALIB *glcalo = new GL_CALO_CALIB();
982     //
983     sgnl = 0;
984     //
985     idcalib[s] = 0;
986     fromtime[s] = 0;
987     totime[s] = 0;
988     calibno[s] = 0;
989     ClearCalibVals(s);
990     //
991     UInt_t uptime = 0;
992     //
993     sgnl = glcalo->Query_GL_CALO_CALIB(atime,uptime,s,dbc);
994     if ( sgnl < 0 ){
995     if ( verbose ) printf(" CALORIMETER - ERROR: error from GLTables\n");
996     return(sgnl);
997     };
998     //
999     idcalib[s] = glcalo->ID_ROOT_L0;
1000     fromtime[s] = glcalo->FROM_TIME;
1001     if ( glcalo->TO_TIME < atime ){ // calibration is corrupted and we are using the one that preceed the good one
1002     totime[s] = uptime;
1003     } else {
1004     totime[s] = glcalo->TO_TIME;
1005     };
1006     // totime[s] = glcalo->TO_TIME;
1007     calibno[s] = glcalo->EV_ROOT;
1008     //
1009     if ( totime[s] == 0 ){
1010     if ( verbose ) printf(" CALORIMETER - WARNING: data with no associated calibration\n");
1011     ClearCalibVals(s);
1012     sgnl = 100;
1013     };
1014     //
1015     // determine path and name and entry of the calibration file
1016     //
1017     GL_ROOT *glroot = new GL_ROOT();
1018     if ( verbose ) printf("\n");
1019     if ( verbose ) printf(" ** SECTION %i **\n",s);
1020     //
1021     sgnl = glroot->Query_GL_ROOT(idcalib[s],dbc);
1022     if ( sgnl < 0 ){
1023     if ( verbose ) printf(" CALORIMETER - ERROR: error from GLTables\n");
1024     return(sgnl);
1025     };
1026     //
1027     stringstream name;
1028     name.str("");
1029     name << glroot->PATH.Data() << "/";
1030     name << glroot->NAME.Data();
1031     //
1032     fcalname[s] = (TString)name.str().c_str();
1033     if ( verbose ) printf(" - event at time %u. From time %u to time %u \n use file %s \n calibration at entry %i \n\n",atime,fromtime[s],totime[s],fcalname[s].Data(),calibno[s]);
1034     //
1035     sgnl = LoadCalib(s);
1036     //
1037     if ( sgnl != 0 ) return(sgnl);
1038     delete glcalo;
1039     delete glroot;
1040     //
1041     return(0);
1042     //
1043     }
1044    
1045     Int_t CaloLevel0::LoadCalib(Int_t s){
1046     //
1047     ifstream myfile;
1048     myfile.open(fcalname[s].Data());
1049     if ( !myfile ){
1050     return(-107);
1051     };
1052     myfile.close();
1053     //
1054     TFile *File = new TFile(fcalname[s].Data());
1055     if ( !File ) return(-108);
1056     TTree *tr = (TTree*)File->Get("CalibCalPed");
1057     if ( !tr ) return(-109);
1058     //
1059     TBranch *calo = tr->GetBranch("CalibCalPed");
1060     //
1061     pamela::CalibCalPedEvent *ce = 0;
1062     tr->SetBranchAddress("CalibCalPed", &ce);
1063     //
1064     Long64_t ncalibs = calo->GetEntries();
1065     //
1066     if ( !ncalibs ) return(-110);
1067     //
1068     calo->GetEntry(calibno[s]);
1069     //
1070     if (ce->cstwerr[s] != 0 && ce->cperror[s] == 0 ) {
1071     for ( Int_t d=0 ; d<11 ;d++ ){
1072     Int_t pre = -1;
1073     for ( Int_t j=0; j<96 ;j++){
1074     if ( j%16 == 0 ) pre++;
1075     if ( s == 2 ){
1076     calped[0][2*d+1][j] = ce->calped[3][d][j];
1077     cstwerr[3] = ce->cstwerr[3];
1078     cperror[3] = ce->cperror[3];
1079     calgood[0][2*d+1][j] = ce->calgood[3][d][j];
1080     calthr[0][2*d+1][pre] = ce->calthr[3][d][pre];
1081     calrms[0][2*d+1][j] = ce->calrms[3][d][j];
1082     calbase[0][2*d+1][pre] = ce->calbase[3][d][pre];
1083     calvar[0][2*d+1][pre] = ce->calvar[3][d][pre];
1084     };
1085     if ( s == 3 ){
1086     calped[0][2*d][j] = ce->calped[1][d][j];
1087     cstwerr[1] = ce->cstwerr[1];
1088     cperror[1] = ce->cperror[1];
1089     calgood[0][2*d][j] = ce->calgood[1][d][j];
1090     calthr[0][2*d][pre] = ce->calthr[1][d][pre];
1091     calrms[0][2*d][j] = ce->calrms[1][d][j];
1092     calbase[0][2*d][pre] = ce->calbase[1][d][pre];
1093     calvar[0][2*d][pre] = ce->calvar[1][d][pre];
1094     };
1095     if ( s == 0 ){
1096     calped[1][2*d][j] = ce->calped[0][d][j];
1097     cstwerr[0] = ce->cstwerr[0];
1098     cperror[0] = ce->cperror[0];
1099     calgood[1][2*d][j] = ce->calgood[0][d][j];
1100     calthr[1][2*d][pre] = ce->calthr[0][d][pre];
1101     calrms[1][2*d][j] = ce->calrms[0][d][j];
1102     calbase[1][2*d][pre] = ce->calbase[0][d][pre];
1103     calvar[1][2*d][pre] = ce->calvar[0][d][pre];
1104     };
1105     if ( s == 1 ){
1106     calped[1][2*d+1][j] = ce->calped[2][d][j];
1107     cstwerr[2] = ce->cstwerr[2];
1108     cperror[2] = ce->cperror[2];
1109     calgood[1][2*d+1][j] = ce->calgood[2][d][j];
1110     calthr[1][2*d+1][pre] = ce->calthr[2][d][pre];
1111     calrms[1][2*d+1][j] = ce->calrms[2][d][j];
1112     calbase[1][2*d+1][pre] = ce->calbase[2][d][pre];
1113     calvar[1][2*d+1][pre] = ce->calvar[2][d][pre];
1114     };
1115     };
1116     };
1117     } else {
1118     if ( verbose ) printf(" CALORIMETER - ERROR: problems finding a good calibration in this file! \n\n ");
1119     return(-111);
1120     };
1121     File->Close();
1122     return(0);
1123     }

  ViewVC Help
Powered by ViewVC 1.1.23