/[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.4 - (hide annotations) (download)
Wed Sep 6 11:03:31 2006 UTC (18 years, 3 months ago) by mocchiut
Branch: MAIN
Changes since 1.3: +33 -40 lines
Adapted to the new profiler

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

  ViewVC Help
Powered by ViewVC 1.1.23