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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.3 - (hide annotations) (download)
Fri Aug 4 10:31:27 2006 UTC (18 years, 4 months ago) by mocchiut
Branch: MAIN
Changes since 1.2: +32 -109 lines
Small memory leaks fixed

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

  ViewVC Help
Powered by ViewVC 1.1.23