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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.9 - (show annotations) (download)
Thu Jan 11 09:32:54 2007 UTC (17 years, 11 months ago) by pam-fi
Branch: MAIN
Changes since 1.8: +1 -0 lines
memory-leak bug fixed

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

  ViewVC Help
Powered by ViewVC 1.1.23