/[PAMELA software]/calo/flight/CaloProfile/src/CaloProfile.cpp
ViewVC logotype

Contents of /calo/flight/CaloProfile/src/CaloProfile.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.15 - (show annotations) (download)
Mon Sep 14 13:56:50 2009 UTC (15 years, 5 months ago) by mocchiut
Branch: MAIN
Changes since 1.14: +32 -6 lines
Bug in Longitudinal fit fixed

1 #include <CaloProfile.h>
2 //
3 ClassImp(CaloLat);
4 ClassImp(CaloLong);
5 ClassImp(Calo2D);
6
7
8 int isempty(struct stack *s){
9 return (s->top == EMPTY) ? 1 : 0;
10 }
11
12 void emptystack(struct stack* s){
13 s->top=EMPTY;
14 strcpy(s->data," ");
15 }
16
17 void push(struct stack* s,int item){
18 if(s->top == (MAX-1)){
19 printf("\n ERROR! STACK FULL (too many digits in SetLower/UpperLimit)\n");
20 } else {
21 s->data[++s->top]=(char)item;
22 };
23 }
24
25 char pop(struct stack* s){
26 char ret=(char)EMPTY;
27 if(!isempty(s)){
28 ret= s->data[s->top--];
29 };
30 return ret;
31 }
32
33
34 int ipop(struct stack* s){
35 int ret=EMPTY;
36 if(s->top == EMPTY)
37 printf("\n ERROR! STACK EMPTY (too few digits in SetLower/UpperLimit)\n");
38 else {
39 ret= s->data[s->top--];
40 };
41 return ret;
42 }
43
44
45 void display(struct stack s){
46 int ss = s.top;
47 while(s.top != EMPTY){
48 printf("\n%d",s.data[s.top--]);
49 };
50 printf(" s.top %i \n",ss);
51 }
52
53 int isoperator(char e){
54 if(e == '+' || e == '-' || e == '*' || e == '/' || e == '%')
55 return 1;
56 else
57 return 0;
58 }
59
60
61 int priority(char e){
62 int pri = 0;
63 if(e == '*' || e == '/' || e =='%')
64 pri = 2;
65 else {
66 if(e == '+' || e == '-') pri = 1;
67 };
68 return pri;
69 }
70
71 void infix2postfix(char* infix, char * postfix, int insertspace){
72 char *i,*p;
73 struct stack X;
74 char n1;
75 emptystack(&X);
76 i = &infix[0];
77 p = &postfix[0];
78 while(*i){
79 while(*i == ' ' || *i == '\t' ){
80 i++;
81 };
82 TString c=i;
83 if( isdigit(*i) || isalpha(*i) || c.BeginsWith(".") ){
84 while( isdigit(*i) || isalpha(*i) || c.BeginsWith(".") ){
85 *p = *i;
86 p++;
87 i++;
88 c=i;
89 };
90 /*SPACE CODE*/
91 if(insertspace){
92 *p = ' ';
93 p++;
94 };
95 /*END SPACE CODE*/
96 };
97 if( *i == '(' ){
98 push(&X,*i);
99 i++;
100 };
101 if( *i == ')'){
102 n1 = pop(&X);
103 while( n1 != '(' ){
104 *p = n1;
105 p++;
106 /*SPACE CODE*/
107 if(insertspace){
108 *p = ' ';
109 p++;
110 };
111 /*END SPACE CODE*/
112 n1 = pop(&X);
113 };
114 i++;
115 };
116 if( isoperator(*i) ){
117 if(isempty(&X)){
118 push(&X,*i);
119 } else {
120 n1 = pop(&X);
121 while( priority(n1) >= priority(*i) ){
122 *p = n1;
123 p++;
124 /*SPACE CODE*/
125 if(insertspace){
126 *p = ' ';
127 p++;
128 };
129 /*END SPACE CODE*/
130 n1 = pop(&X);
131 };
132 push(&X,n1);
133 push(&X,*i);
134 };
135 i++;
136 };
137 };
138 while(!isempty(&X)){
139 n1 = pop(&X);
140 *p = n1;
141 p++;
142 /*SPACE CODE*/
143 if(insertspace){
144 *p = ' ';
145 p++;
146 };
147 /*END SPACE CODE*/
148 };
149 *p = '\0';
150 }
151
152
153 Float_t evaluate(char *postfix){
154 //
155 Float_t op1 = 0.;
156 Float_t op2 = 0.;
157 Float_t result = 0.;
158 //
159 TString e = postfix;
160 Float_t st[50];
161 memset(st,0,50*sizeof(Float_t));
162 TObjArray *ae = e.Tokenize(" ");
163 Int_t o = 0;
164 Int_t a = 0;
165 Int_t ap = 0;
166 //
167 while ( o < ae->GetEntries() ){
168 if ( ((TString)(((TObjString*)ae->At(o))->GetString())).IsFloat() ){
169 st[a] = (Float_t)((TString)(((TObjString*)ae->At(o))->GetString())).Atof();;
170 a++;
171 } else {
172 ap = a-1;
173 op1 = st[ap--];
174 op2 = st[ap];
175 const char *p=((TString)(((TObjString*)ae->At(o))->GetString())).Data();
176 switch(*p){
177 case '+':
178 result = op2 + op1;
179 break;
180 case '-':
181 result = op2 - op1;
182 break;
183 case '/':
184 result = op2 / op1;
185 break;
186 case '*':
187 result = op2 * op1;
188 break;
189 case '%':
190 result = (Int_t)round(op2) % (Int_t)round(op1);
191 break;
192 default:
193 printf("\nInvalid Operator: %s \n",((TString)(((TObjString*)ae->At(o))->GetString())).Data());
194 return 0;
195 };
196 st[ap] = result;
197 a = ap+1;
198 //
199 };
200 o++;
201 };
202 return result;
203 //
204 }
205
206
207
208 ////////////////////////////////////////////////////////////////////////
209 /**
210 * 1-dimension function describing lateral distribution of the
211 * shower as viewed by calorimeter
212 * (projection of 3d function in one direction)
213 *
214 * xi[0] = x or y coordinate relative to shower axis
215 * parmin[0] = rt
216 * parmin[1] = p
217 * parmin[2] = rc
218 *
219 */
220 ////////////////////////////////////////////////////////////////////////
221 Double_t cfradx(Double_t *xi, Double_t *parmin) {
222
223 double fradxmin2,p,rt,rc,es,x,pig,norm,c;
224 x=*xi;
225 pig = acos(-1.);
226 rt=parmin[0];
227 p=parmin[1];
228 rc=parmin[2];
229 norm=parmin[3];
230 c=parmin[4];
231 x=x-c;
232 es=1.5;
233 fradxmin2=p*pig*pow(rc,2)/pow((pow(x,2)+pow(rc,2)),es);
234 fradxmin2=fradxmin2+(1-p)*pig*pow(rt,2)/pow((pow(x,2)+pow(rt,2)),es);
235 fradxmin2=norm*fradxmin2/(2*pig);
236 //cout<<"x,fradxmin2 "<< x<<" "<<fradxmin2 <<endl;
237 return fradxmin2;
238 }
239
240
241 Double_t cfunc(Double_t *tin,Double_t *par)
242 {
243
244 Double_t gaa,a,tm,et,value,t,t0;
245 t=*tin;
246 et=par[0];
247 a=par[1];
248 tm=par[2];
249 t0=par[3];
250 gaa=TMath::Gamma(a);
251
252 value=et*((a-1)/(tm*gaa))*pow(((a-1)*(t-t0)/tm),(a-1))*exp(-(a-1)*(t-t0)/tm);
253 return value;
254 }
255
256
257 //--------------------------------------
258 /**
259 * Default constructor
260 */
261 CaloLat::CaloLat(){
262 Clear();
263 };
264
265 /**
266 * Default constructor
267 */
268 Calo2D::Calo2D(){
269 Clear();
270 };
271
272 CaloLat::CaloLat(PamLevel2 *l2p){
273 //
274 Clear();
275 //
276 L2 = l2p;
277 //
278 if ( !L2->IsORB() ) printf(" WARNING: OrbitalInfo Tree is needed, the plugin could not work properly without it \n");
279 //
280 OBT = 0;
281 PKT = 0;
282 atime = 0;
283 //
284 suf = "";
285 debug = false;
286 //
287 };
288
289 Calo2D::Calo2D(PamLevel2 *l2p){
290 //
291 Clear();
292 //
293 L2 = l2p;
294 //
295 if ( !L2->IsORB() ) printf(" WARNING: OrbitalInfo Tree is needed, the plugin could not work properly without it \n");
296 //
297 OBT = 0;
298 PKT = 0;
299 atime = 0;
300 //
301 suf = "";
302 debug = false;
303 //
304 };
305
306 void CaloLat::Clear(){
307 //
308 };
309
310 void Calo2D::Clear(){
311 //
312 };
313
314 void CaloLat::Print(){
315 //
316 Process();
317 //
318 printf("==================== Calorimeter Lateral Profile =======================\n");
319 printf(" OBT: %u PKT: %u ATIME: %u \n",OBT,PKT,atime);
320 // printf(" nx [number of X combination]:.. %i\n",nx);
321 // printf(" ny [number of Y combination]:.. %i\n",ny);
322 printf("========================================================================\n");
323 //
324 };
325
326 void Calo2D::Print(){
327 //
328 Process();
329 //
330 printf("==================== Calorimeter 2D Profile =======================\n");
331 printf(" OBT: %u PKT: %u ATIME: %u \n",OBT,PKT,atime);
332 // printf(" nx [number of X combination]:.. %i\n",nx);
333 // printf(" ny [number of Y combination]:.. %i\n",ny);
334 printf("========================================================================\n");
335 //
336 };
337
338 void CaloLat::Draw(){
339 //
340 Process();
341 Draw(-1,-1);
342 };
343
344 void Calo2D::Draw(){
345 //
346 Process();
347 Draw(-1);
348 };
349
350 void CaloLat::Draw(Int_t view,Int_t plane){
351 //
352 Int_t minv = 0;
353 Int_t maxv = 0;
354 Int_t minp = 0;
355 Int_t maxp = 0;
356 //
357 if ( view == -1 ){
358 minv = 0;
359 maxv = 2;
360 } else {
361 minv = view;
362 maxv = view+1;
363 };
364
365 if ( plane == -1 ){
366 minp = 0;
367 maxp = 22;
368 } else {
369 minp = plane;
370 maxp = plane+1;
371 };
372 //
373 Process();
374 //
375 gStyle->SetLabelSize(0.04);
376 gStyle->SetNdivisions(510,"XY");
377 //
378 for (Int_t v=minv; v<maxv;v++){
379 for (Int_t p=minp; p<maxp;p++){
380 TString hid = Form("clatv%ip%i%s",v,p,suf.Data());
381 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
382 if ( tc ){
383 // tc->Clear();
384 } else {
385 tc = new TCanvas(hid,hid);
386 };
387 //
388 TString thid = Form("hlatv%ip%i%s",v,p,suf.Data());
389 TH1F *th = dynamic_cast<TH1F*>(gDirectory->FindObject(thid));
390 if ( th ) th->Delete();
391 // th->Clear();
392 // th->Reset();
393 // } else {
394 th = new TH1F(thid,thid,96,-0.5,95.5);
395 // };
396 tc->cd();
397 //
398 for (Int_t st=0;st<96;st++){
399 th->Fill(st,estrip[v][p][st]);
400 };
401 th->Draw();
402 tc->Modified();
403 tc->Update();
404 };
405 };
406 //
407 gStyle->SetLabelSize(0);
408 gStyle->SetNdivisions(1,"XY");
409 //
410 };
411
412 void Calo2D::Draw(Int_t plane){
413 //
414 Int_t minp = 0;
415 Int_t maxp = 0;
416 //
417 if ( plane == -1 ){
418 minp = 0;
419 maxp = 23;
420 } else {
421 minp = plane;
422 maxp = plane+1;
423 };
424 //
425 Process();
426 //
427 gStyle->SetLabelSize(0.04);
428 gStyle->SetNdivisions(510,"XY");
429 //
430 for (Int_t p=minp; p<maxp;p++){
431 TString hid = Form("c2dp%i%s",p,suf.Data());
432 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
433 if ( tc ){
434 // tc->Clear();
435 } else {
436 tc = new TCanvas(hid,hid);
437 };
438 //
439 TString thid = Form("h2dp%i%s",p,suf.Data());
440 TH2F *th = dynamic_cast<TH2F*>(gDirectory->FindObject(thid));
441 if ( th ) th->Delete();
442 // th->Clear();
443 // th->Reset();
444 // } else {
445 Int_t minx = smax[p] - 10;
446 if ( minx < 0 ) minx = 0;
447 Int_t maxx = minx + 20;
448 if ( maxx > 95 ){
449 maxx = 95;
450 minx = 75;
451 };
452 Int_t miny = smay[p] - 10;
453 if ( miny < 0 ) miny = 0;
454 Int_t maxy = miny + 20;
455 if ( maxy > 95 ){
456 maxy = 95;
457 miny = 75;
458 };
459 th = new TH2F(thid,thid,20,(Float_t)minx-0.5,(Float_t)maxx-0.5,20,(Float_t)miny-0.5,(Float_t)maxy-0.5);
460 // th = new TH2F(thid,thid,96,-0.5,95.5,96,-0.5,95.5);
461 // };
462 tc->cd();
463 //
464 for (Int_t stx=minx;stx<maxx+1;stx++){
465 for (Int_t sty=miny;sty<maxy+1;sty++){
466 th->Fill(stx,sty,estrip[p][stx][sty]);
467 };
468 };
469 gStyle->SetPalette(1);
470 // tc->SetLogz();
471 // th->Draw("colbox");
472 th->Draw("cont4");
473 tc->Modified();
474 tc->Update();
475 };
476 //
477 gStyle->SetLabelSize(0);
478 gStyle->SetNdivisions(1,"XY");
479 //
480 };
481
482 void CaloLat::Delete(){
483 Clear();
484 //delete this;
485 };
486
487 void Calo2D::Delete(){
488 Clear();
489 //delete this;
490 };
491
492
493 void CaloLat::Process(){
494 //
495 if ( !L2 ){
496 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
497 printf(" ERROR: CaloHough variables not filled \n");
498 return;
499 };
500 //
501 Bool_t newentry = false;
502 //
503 if ( L2->IsORB() ){
504 if ( L2->GetOrbitalInfo()->pkt_num != PKT || L2->GetOrbitalInfo()->OBT != OBT || L2->GetOrbitalInfo()->absTime != atime ){
505 newentry = true;
506 OBT = L2->GetOrbitalInfo()->OBT;
507 PKT = L2->GetOrbitalInfo()->pkt_num;
508 atime = L2->GetOrbitalInfo()->absTime;
509 };
510 } else {
511 newentry = true;
512 };
513 //
514 if ( !newentry ) return;
515 //
516 if ( debug ) printf(" Start processing event at OBT %u PKT %u time %u \n",OBT,PKT,atime);
517 //
518 Clear();
519 //
520 // let's start
521 //
522 memset(estrip,0, 4224*sizeof(Float_t));
523 Float_t mip1 = 0.;
524 Int_t view1 = 0;
525 Int_t plane1 = 0;
526 Int_t strip1 = 0;
527 //
528 for (Int_t i=0; i<L2->GetCaloLevel1()->istrip ; i++){
529 mip1 = L2->GetCaloLevel1()->DecodeEstrip(i,view1,plane1,strip1);
530 estrip[view1][plane1][strip1] = mip1;
531 };
532 //
533 if ( debug ) this->Print();
534 if ( debug ) printf(" exit \n");
535 //
536 };
537
538 void Calo2D::Process(){
539 //
540 if ( !L2 ){
541 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
542 printf(" ERROR: CaloHough variables not filled \n");
543 return;
544 };
545 //
546 Bool_t newentry = false;
547 //
548 if ( L2->IsORB() ){
549 if ( L2->GetOrbitalInfo()->pkt_num != PKT || L2->GetOrbitalInfo()->OBT != OBT || L2->GetOrbitalInfo()->absTime != atime ){
550 newentry = true;
551 OBT = L2->GetOrbitalInfo()->OBT;
552 PKT = L2->GetOrbitalInfo()->pkt_num;
553 atime = L2->GetOrbitalInfo()->absTime;
554 };
555 } else {
556 newentry = true;
557 };
558 //
559 if ( !newentry ) return;
560 //
561 if ( debug ) printf(" Start processing event at OBT %u PKT %u time %u \n",OBT,PKT,atime);
562 //
563 Clear();
564 //
565 // let's start
566 //
567 Float_t es[2][22][96];
568 memset(es,0, 4224*sizeof(Float_t));
569 memset(estrip,0, 4224*sizeof(Float_t));
570 Float_t mip1 = 0.;
571 Int_t view1 = 0;
572 Int_t plane1 = 0;
573 Int_t strip1 = 0;
574 //
575 for (Int_t i=0; i<L2->GetCaloLevel1()->istrip ; i++){
576 mip1 = L2->GetCaloLevel1()->DecodeEstrip(i,view1,plane1,strip1);
577 es[view1][plane1][strip1] = mip1;
578 };
579 //
580 Int_t plane2 = 0;
581 Float_t emax[23];
582 memset(emax,0,sizeof(Float_t)*23);
583 memset(smax,0,sizeof(Int_t)*23);
584 memset(smay,0,sizeof(Int_t)*23);
585 //
586 //
587 for (Int_t p=0; p < 23 ; p++){
588 //
589 plane1 = p-1;
590 plane2 = p;
591 //
592 if ( p == 0 ){
593 plane1 = -1;
594 plane2 = 0;
595 };
596 if ( p == 22 ){
597 plane1 = 21;
598 plane2 = -1;
599 };
600 //
601 for (Int_t s=0; s < 96 ; s++){ // x
602 for (Int_t ss=0; ss < 96 ; ss++){ // y
603 if ( p == 0 ){
604 estrip[p][s][ss] += es[1][plane2][ss];
605 if ( (es[1][plane2][ss]) > emax[p] ){
606 smax[p] = 45;
607 smay[p] = ss;
608 emax[p] = es[1][plane2][ss] ;
609 };
610 };
611 if ( p > 0 && p < 22 ){
612 estrip[p][s][ss] += es[0][plane1][s] + es[1][plane2][ss];
613 if ( (es[0][plane1][s] + es[1][plane2][ss]) > emax[p] ){
614 smax[p] = s;
615 smay[p] = ss;
616 emax[p] = es[0][plane1][s] + es[1][plane2][ss] ;
617 };
618 };
619 if ( p == 22 ){
620 estrip[p][s][ss] += es[0][plane1][s];
621 if ( (es[1][plane2][s]) > emax[p] ){
622 smax[p] = s;
623 smay[p] = 45;
624 emax[p] = es[1][plane2][s] ;
625 };
626 };
627 };
628 };
629 //
630 };
631 //
632 if ( debug ) this->Print();
633 if ( debug ) printf(" exit \n");
634 //
635 };
636
637
638 /**
639 * Default constructor
640 */
641 CaloLong::CaloLong(){
642 Clear();
643 };
644
645 CaloLong::CaloLong(PamLevel2 *l2p){
646 //
647 Clear();
648 //
649 L2 = l2p;
650 //
651 if ( !L2->IsORB() ) printf(" WARNING: OrbitalInfo Tree is needed, the plugin could not work properly without it \n");
652 //
653 OBT = 0;
654 PKT = 0;
655 atime = 0;
656 //
657 clp = L2->GetCaloLevel2();
658 //
659 sel = true;
660 cont = false;
661 N = 0;
662 NC = 22;
663 mask18b = -1;
664 //
665 no18x = true;
666 debug = false;
667 maskXE = false;
668 maskXO = false;
669 maskYE = false;
670 maskYO = false;
671 //
672 lmax = 0.;
673 umax = 100.;
674 slmax = "";
675 sumax = "";
676 suf = "";
677 xyaverage = true;
678 //
679 letmax = -1;
680 heavytail = false;
681 lmipth = 100.;
682 //
683 };
684
685 void CaloLong::MaskSection(TString sec){
686 sec.ToUpper();
687 if ( sec.Contains("XO") ) maskXO = true;
688 if ( sec.Contains("YO") ) maskYO = true;
689 if ( sec.Contains("XE") ) maskXE = true;
690 if ( sec.Contains("YE") ) maskYE = true;
691 }
692
693 void CaloLong::UnMaskSections(){
694 this->UnMaskSection("XEXOYEYO");
695 }
696
697 void CaloLong::UnMaskSection(TString sec){
698 sec.ToUpper();
699 if ( sec.Contains("XO") ) maskXO = false;
700 if ( sec.Contains("YO") ) maskYO = false;
701 if ( sec.Contains("XE") ) maskXE = false;
702 if ( sec.Contains("YE") ) maskYE = false;
703 }
704
705 void CaloLong::Clear(){
706 //
707 if ( debug ) printf(" Clear called \n");
708 //
709 memset(eplane,0, 2*22*sizeof(Float_t));
710 //
711 chi2 = 0.;
712 ndf = 0.;
713 E0 = 0.;
714 defE0 = 0.;
715 a = 0.;
716 b = 0.;
717 errE0 = 0.;
718 erra = 0.;
719 errb = 0.;
720 etmax = 0.;
721 asymm = 0.;
722 fitresult = 0;
723 //
724 X0pl = 0.76;
725 //
726 };
727
728 void CaloLong::Print(){
729 //
730 Fit();
731 //
732 printf("==================== Calorimeter Longitudinal Profile =======================\n");
733 printf(" OBT: %u PKT: %u ATIME: %u \n",OBT,PKT,atime);
734 printf(" fitresult:.. %i\n",fitresult);
735 printf(" chi2 :.. %f\n",chi2);
736 printf(" ndf :.. %f\n",ndf);
737 printf(" nchi2 :.. %f\n",chi2/ndf);
738 printf(" E0 :.. %f\n",E0);
739 printf(" E0/260. :.. %f\n",E0/260.);
740 printf(" defE0 :.. %f\n",defE0);
741 printf(" lower :.. %f\n",lmax);
742 printf(" upper :.. %f\n",umax);
743 printf(" s lower :.. %s\n",slmax.Data());
744 printf(" s upper :.. %s\n",sumax.Data());
745 printf(" a :.. %f\n",a);
746 printf(" b :.. %f\n",b);
747 printf(" errE0 :.. %f\n",errE0);
748 printf(" erra :.. %f\n",erra);
749 printf(" errb :.. %f\n",errb);
750 printf(" asymm :.. %f\n",asymm);
751 printf(" tmax :.. %f\n",((a-1.)/b));
752 printf(" etmax :.. %f\n",etmax);
753 printf(" X0pl :.. %f\n",X0pl);
754 printf("========================================================================\n");
755 //
756 };
757
758 void CaloLong::SetNoWpreSampler(Int_t n){
759 //
760 if ( NC+n <= 22 && NC+n >= 0 ){
761 N = n;
762 } else {
763 printf(" ERROR! Calorimeter is made of 22 W planes\n");
764 printf(" you are giving N presampler = %i and N calo = %i \n",n,NC);
765 printf(" WARNING: using default values NWpre = 0, NWcalo = 22\n");
766 NC = 22;
767 N = 0;
768 };
769 }
770
771 void CaloLong::SetNoWcalo(Int_t n){
772 if ( N+n <= 22 && N+n >= 0 ){
773 NC = n;
774 } else {
775 printf(" ERROR! Calorimeter is made of 22 W planes\n");
776 printf(" you are giving N W presampler = %i and N W calo = %i \n",N,n);
777 printf(" WARNING: using default values NWpre = 0, NWcalo = 22\n");
778 NC = 22;
779 N = 0;
780 };
781 }
782
783 void CaloLong::SplitInto(Int_t NoWpreSampler, Int_t NoWcalo){
784 this->SetNoWpreSampler(0);
785 this->SetNoWcalo(0);
786 if ( NoWpreSampler < NoWcalo ){
787 this->SetNoWpreSampler(NoWpreSampler);
788 this->SetNoWcalo(NoWcalo);
789 } else {
790 this->SetNoWcalo(NoWcalo);
791 this->SetNoWpreSampler(NoWpreSampler);
792 };
793 }
794
795 void CaloLong::Process(){
796 //
797 if ( !L2 ){
798 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
799 printf(" ERROR: CaloHough variables not filled \n");
800 return;
801 };
802 //
803 Bool_t newentry = false;
804 //
805 if ( L2->IsORB() ){
806 if ( L2->GetOrbitalInfo()->pkt_num != PKT || L2->GetOrbitalInfo()->OBT != OBT || L2->GetOrbitalInfo()->absTime != atime ){
807 newentry = true;
808 OBT = L2->GetOrbitalInfo()->OBT;
809 PKT = L2->GetOrbitalInfo()->pkt_num;
810 atime = L2->GetOrbitalInfo()->absTime;
811 };
812 } else {
813 newentry = true;
814 };
815 //
816 if ( !newentry ) return;
817 //
818 if ( debug ) printf(" Start processing event at OBT %u PKT %u time %u \n",OBT,PKT,atime);
819 //
820 Clear();
821 //
822 // let's start
823 //
824 if ( cont ){
825 for (Int_t i=0; i<22; i++){
826 if ( i == (18+N) ){
827 mask18b = 18 + N;
828 break;
829 };
830 };
831 };
832 //
833 if ( sel ){
834 for (Int_t i=0; i<22; i++){
835 if ( i == (18-N) ){
836 mask18b = 18 - N;
837 break;
838 };
839 };
840 };
841 //
842 // if ( mask18b == 18 ) mask18b = -1;
843 //
844 Int_t view = 0;
845 Int_t plane = 0;
846 Int_t strip = 0;
847 Float_t mip = 0.;
848 Bool_t gof = true;
849 for (Int_t i=0; i < L2->GetCaloLevel1()->istrip; i++){
850 mip = L2->GetCaloLevel1()->DecodeEstrip(i,view,plane,strip);
851 gof = true;
852 if ( maskXE && (plane%2)==0 && view==1 ) gof = false;
853 if ( maskXO && (plane%2)!=0 && view==1 ) gof = false;
854 if ( maskYE && (plane%2)!=0 && view==0 ) gof = false;
855 if ( maskYO && (plane%2)==0 && view==0 ) gof = false;
856 if ( gof ) eplane[view][plane] += mip;
857 };
858 //
859 // inclination factor (stolen from Daniele's code)
860 //
861 Float_t ytgx = 0.;
862 Float_t ytgy = 0.;
863 ytgx = 0.76 * clp->tanx[0];
864 ytgy = 0.76 * clp->tany[0];
865 X0pl = sqrt( pow(0.76,2.) + pow(ytgx,2.) + pow(ytgy,2.) );
866 //
867 // Find experimental plane of maximum
868 //
869 Int_t pmax = 0;
870 Int_t vmax = 0;
871 Float_t emax = 0.;
872 for (Int_t v=0; v<2; v++){
873 for (Int_t i=0; i<22; i++){
874 if ( eplane[v][i] > emax ){
875 emax = eplane[v][i];
876 vmax = v;
877 pmax = i;
878 };
879 };
880 };
881 //
882 //
883 //
884 if ( vmax == 0 ) pmax++;
885 etmax = pmax * X0pl;
886 //
887 if ( debug ) this->Print();
888 if ( debug ) printf(" exit \n");
889 //
890 };
891
892 void CaloLong::SetEnergies(Float_t myene[][22]){
893 //
894 if ( !L2 ){
895 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
896 printf(" ERROR: CaloHough variables not filled \n");
897 return;
898 };
899 //
900 Bool_t newentry = false;
901 //
902 if ( L2->IsORB() ){
903 if ( L2->GetOrbitalInfo()->pkt_num != PKT || L2->GetOrbitalInfo()->OBT != OBT || L2->GetOrbitalInfo()->absTime != atime ){
904 newentry = true;
905 OBT = L2->GetOrbitalInfo()->OBT;
906 PKT = L2->GetOrbitalInfo()->pkt_num;
907 atime = L2->GetOrbitalInfo()->absTime;
908 };
909 } else {
910 newentry = true;
911 };
912 //
913 if ( !newentry ) return;
914 //
915 if ( debug ) printf(" SET ENERGIES Start processing event at OBT %u PKT %u time %u \n",OBT,PKT,atime);
916 //
917 Clear();
918 //
919 // let's start
920 //
921 if ( cont ){
922 for (Int_t i=0; i<22; i++){
923 if ( i == (18+N) ){
924 mask18b = 18 + N;
925 break;
926 };
927 };
928 };
929 //
930 if ( sel ){
931 for (Int_t i=0; i<22; i++){
932 if ( i == (18-N) ){
933 mask18b = 18 - N;
934 break;
935 };
936 };
937 };
938 //
939 // if ( mask18b == 18 ) mask18b = -1;
940 //
941 Int_t view = 0;
942 Int_t plane = 0;
943 Bool_t gof = true;
944 for (view=0; view < 2; view++){
945 for (plane=0; plane < 22; plane++){
946 gof = true;
947 if ( maskXE && (plane%2)==0 && view==1 ) gof = false;
948 if ( maskXO && (plane%2)!=0 && view==1 ) gof = false;
949 if ( maskYE && (plane%2)!=0 && view==0 ) gof = false;
950 if ( maskYO && (plane%2)==0 && view==0 ) gof = false;
951 if ( gof ) eplane[view][plane] = myene[view][plane];
952 if ( debug ) printf(" XE %i XO %i YE %i YO %i eplane %i %i = %f ........... myene %i %i = %f\n", maskXE, maskXO, maskYE, maskYO,view,plane,eplane[view][plane],view,plane,myene[view][plane]);
953 };
954 };
955 //
956 if ( !clp && debug ) printf(" ERROR!! no CaloLevel2 obj!!\n");
957 //
958 // inclination factor (stolen from Daniele's code)
959 //
960 Float_t ytgx = 0.;
961 Float_t ytgy = 0.;
962 ytgx = 0.76 * clp->tanx[0];
963 ytgy = 0.76 * clp->tany[0];
964 X0pl = sqrt( pow(0.76,2.) + pow(ytgx,2.) + pow(ytgy,2.) );
965 if ( debug ) printf(" X0pl %f tanx %f tany %f \n",X0pl,ytgx,ytgy);
966 //
967 // Find experimental plane of maximum
968 //
969 Int_t pmax = 0;
970 Int_t vmax = 0;
971 Float_t emax = 0.;
972 for (Int_t v=0; v<2; v++){
973 for (Int_t i=0; i<22; i++){
974 if ( eplane[v][i] > emax ){
975 emax = eplane[v][i];
976 vmax = v;
977 pmax = i;
978 };
979 };
980 };
981 //
982 //
983 //
984 if ( vmax == 0 ) pmax++;
985 etmax = pmax * X0pl;
986 letmax = etmax;
987 //
988 if ( debug ) this->Print();
989 if ( debug ) printf(" exit \n");
990 //
991 };
992
993 Double_t ccurve(Double_t *ti,Double_t *par){
994 //
995 Double_t t = *ti;
996 Double_t cE0 = par[0];
997 Double_t ca = par[1];
998 Double_t cb = par[2];
999 Double_t gammaa = TMath::Gamma(ca);
1000 //
1001 Double_t value = cE0 * cb * ( pow((cb * t),(ca - 1)) * exp( -cb * t ) ) / gammaa;
1002 //
1003 return value;
1004 //
1005 }
1006
1007 void CaloLong::Fit(){
1008 this->Fit(false);
1009 };
1010
1011 Float_t CaloLong::Evaluate(TString s, Float_t tmax, Float_t X0pl){
1012 /* SAMPLE OUTPUT:
1013 Enter Infix Expression : A + B + C / (E - F)
1014 Postfix Expression is : A B + C E F - / +
1015 */
1016 if ( !s.Contains("tmax") && !s.Contains("X0pl") ){
1017 printf(" ERROR, the input formula must contain \"t\"\n");
1018 return 0.;
1019 };
1020 if ( tmax != tmax ){
1021 printf(" ERROR, tmax is nan! \n");
1022 return 0.;
1023 };
1024 TString g=Form("%f",tmax);
1025 TString h=Form("%f",X0pl);
1026 TString *ts= new TString("");
1027 ts->Prepend(s.Data());
1028 ts->ReplaceAll("tmax",4,g.Data(),g.Capacity());
1029 ts->ReplaceAll("X0pl",4,h.Data(),h.Capacity());
1030 ts->Prepend("(");
1031 ts->Append(")");
1032 if ( debug ) printf(" ts %s tssize %i capac %i s %s g %s \n",ts->Data(),ts->Sizeof(),ts->Capacity(),s.Data(),g.Data());
1033 char in[50],post[50];
1034 strcpy(&in[0]," ");
1035 strcpy(&in[0],ts->Data());
1036 strcpy(&post[0]," ");
1037 if ( debug ) printf("Infix Expression is : ##%s##\n",&in[0]);
1038 infix2postfix(&in[0],&post[0],1);
1039 if ( debug ) printf("Postfix Expression is : ##%s##\n",&post[0]);
1040 /* SAMPLE OUTPUT:
1041 Enter Postfix Expression : 3 5 + 2 /
1042 3 5 + 2 / EQUALS 4
1043 */
1044 Float_t res = evaluate(&post[0]);
1045 if ( debug ) printf("%s EQUALS %f\n",post,res);
1046 return res;
1047 }
1048
1049 void CaloLong::Fit(Bool_t draw){
1050 //
1051 Process();
1052 //
1053 //
1054 if ( !L2 ){
1055 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
1056 printf(" ERROR: CaloHough variables not filled \n");
1057 return;
1058 };
1059 //
1060 Bool_t newentry = false;
1061 //
1062 if ( L2->IsORB() ){
1063 if ( L2->GetOrbitalInfo()->pkt_num != fPKT || L2->GetOrbitalInfo()->OBT != fOBT || L2->GetOrbitalInfo()->absTime != fatime ){
1064 newentry = true;
1065 fOBT = L2->GetOrbitalInfo()->OBT;
1066 fPKT = L2->GetOrbitalInfo()->pkt_num;
1067 fatime = L2->GetOrbitalInfo()->absTime;
1068 };
1069 } else {
1070 newentry = true;
1071 };
1072 //
1073 if ( !newentry ) return;
1074 //
1075 if ( debug ) printf(" Start fitting event at OBT %u PKT %u time %u \n",fOBT,fPKT,fatime);
1076 //
1077 if ( draw ){
1078 gStyle->SetLabelSize(0.04);
1079 gStyle->SetNdivisions(510,"XY");
1080 };
1081 //
1082 TString hid = Form("clongfit%s",suf.Data());
1083 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
1084 // if ( tc ) tc->Delete();
1085 // if ( tc ) tc->Close();
1086 if ( !tc && draw ){
1087 tc = new TCanvas(hid,hid);
1088 } else {
1089 if ( tc ) tc->cd();
1090 };
1091 //
1092 TString thid = Form("hlongfit%s",suf.Data());
1093 TH2F *th = dynamic_cast<TH2F*>(gDirectory->FindObject(thid));
1094 if ( th ) th->Delete();
1095 //
1096 TString ghid = Form("glongfit%s",suf.Data());
1097 TGraphErrors *gh = dynamic_cast<TGraphErrors*>(gDirectory->FindObject(ghid));
1098 if ( gh ) gh->Delete();
1099 //
1100 Float_t xpos = 0.;
1101 Float_t enemip = 0.;
1102 Float_t xmax = NC * X0pl + 0.2;
1103 //
1104 Double_t xxx[50];
1105 Double_t exx[50];
1106 Double_t yyy[50];
1107 Double_t eyy[50];
1108 Int_t numpo = 0;
1109 memset(xxx,0,50*sizeof(Double_t));
1110 memset(exx,0,50*sizeof(Double_t));
1111 memset(yyy,0,50*sizeof(Double_t));
1112 memset(eyy,0,50*sizeof(Double_t));
1113 //
1114 // AGH, BUG!
1115 //
1116 Int_t mmin = 0;
1117 Int_t mmax = 0;
1118 if ( cont ){
1119 mmin = N;
1120 mmax = NC+N;
1121 } else {
1122 mmin = 0;
1123 mmax = NC;
1124 };
1125 //
1126 Float_t emax = 0.;
1127 Float_t qtotparz = 0.;
1128 for (Int_t st=mmin;st<mmax+1;st++){
1129 enemip = 0.;
1130 xpos = (st - mmin) * X0pl;
1131 //
1132 if ( xyaverage ){
1133 //
1134 if ( st > mmin && st < mmax ){
1135 if ( no18x && ( st == 18+1 || st == mask18b+1 )){
1136 if ( !maskYO ){
1137 enemip = 2. * eplane[1][st];
1138 } else {
1139 enemip = eplane[1][st];
1140 };
1141 } else {
1142 enemip = eplane[0][st-1] + eplane[1][st];
1143 };
1144 } else {
1145 if ( st == mmin ){
1146 if ( !maskYE ){
1147 enemip = 2. * eplane[1][st];
1148 } else {
1149 enemip = eplane[1][st];
1150 };
1151 };
1152 if ( st == mmax ){
1153 if ( !maskXE ){
1154 enemip = 2. * eplane[0][st-1];
1155 } else {
1156 enemip = eplane[0][st-1];
1157 };
1158 };
1159 };
1160 //
1161 qtotparz += enemip;
1162 if ( enemip > emax ) emax = enemip;
1163 if ( enemip > 0. ){
1164 xxx[numpo] = xpos;
1165 exx[numpo] = 0.1;
1166 yyy[numpo] = enemip;
1167 eyy[numpo] = sqrt(enemip*3.)+sqrt(5.);
1168 if ( xpos > letmax && enemip > lmipth && heavytail) eyy[numpo] = (sqrt(enemip*3.)+sqrt(5.))/numpo;
1169 numpo++;
1170 // th->Fill(xpos,enemip);
1171 if ( debug ) printf(" AVE Filling: st %i xpos %f energy %f error %f letmax %f lmpith %f \n",st,xpos,enemip,eyy[numpo-1],letmax,lmipth);
1172 };
1173 } else {
1174 for (Int_t jj=0; jj<2; jj++){
1175 if ( st > mmin && st < mmax ){
1176 if ( jj == 0 && no18x && ( st == 18+1 || st == mask18b+1 )){
1177 enemip = eplane[1][st];
1178 } else {
1179 if ( jj == 0 ){
1180 enemip = eplane[jj][st-1];
1181 } else {
1182 enemip = eplane[jj][st];
1183 };
1184 };
1185 } else {
1186 if ( st == mmin && jj == 1 ){
1187 enemip = eplane[jj][st];
1188 };
1189 if ( st == mmax && jj == 0){
1190 enemip = eplane[jj][st-1];
1191 };
1192 };
1193 //
1194 qtotparz += enemip;
1195 if ( enemip > emax ) emax = enemip;
1196 if ( enemip > 0. ){
1197 xxx[numpo] = xpos;
1198 exx[numpo] = 0.1;
1199 yyy[numpo] = enemip;
1200 eyy[numpo] = sqrt(enemip*3.)+sqrt(5.);
1201 if ( xpos > letmax && enemip > lmipth && heavytail ) eyy[numpo] = (sqrt(enemip*3.)+sqrt(5.))/numpo;
1202 // eyy[numpo] = sqrt(enemip)/(st*0.95);
1203 numpo++;
1204 // th->Fill(xpos,enemip);
1205 if ( debug ) printf(" XY Filling: st %i xpos %f energy %f \n",st,xpos,enemip);
1206 };
1207 };
1208 };
1209
1210 //
1211 // for (Int_t v=1; v>=0;v--)// {
1212 // //
1213 // if ( v == 1 ){
1214 // xpos = (st - N) * X0pl;
1215 // } else {
1216 // xpos = (st + 1 - N) * X0pl;
1217 // };
1218 // //
1219 // if ( no18x && st == 18 && v == 0 ){
1220 // // skip plane 18x
1221 // } else {
1222 // if ( v == 1 && st == mask18b ){
1223 // // emulate plane 18x
1224 // } else {
1225 // if ( eplane[v][st] > 0. ){
1226 // th->Fill(xpos,eplane[v][st]);
1227 // if ( debug ) printf(" Filling: st %i v %i xpos %f energy %f \n",st,v,xpos,eplane[v][st]);
1228 // };
1229 // };
1230 // };
1231 // //
1232 // };
1233 };
1234 //
1235 // th = new TH2F(thid,thid,int(NC*1.5),-0.2,xmax);
1236 th = new TH2F(thid,thid,1000,-0.2,xmax,1000,0.,emax*1.2);
1237 gh = new TGraphErrors(numpo,xxx,yyy,exx,eyy);
1238 TString fnam=Form("lfit%s",suf.Data());
1239 TF1 *lfit = dynamic_cast<TF1*>(gDirectory->FindObject(fnam));
1240 if ( lfit ) lfit->Delete();
1241 lfit = new TF1(fnam,ccurve,0.,xmax,3);
1242 // if ( !lfit ){
1243 // lfit = new TF1(fnam,ccurve,0.,xmax,3);
1244 // };
1245 // lfit->Clear();
1246 // lfit->SetName(fnam);
1247 // lfit->SetRange(0.,xmax);
1248 // lfit->Print();
1249 //
1250 if ( debug ) printf("qtot %f qtotparz %f \n",clp->qtot,qtotparz);
1251 E0 = qtotparz;
1252 // E0 = clp->qtot;
1253 a = 5.;
1254 b = 0.5;
1255 if ( debug ) printf(" STARTING PARAMETERS: E0 %f a %f b %f \n",E0,a,b);
1256 lfit->SetParameters(E0,a,b);
1257 // lfit->SetParLimits(0,0.,1000.);
1258 // lfit->SetParLimits(1,-1.,80.);
1259 // lfit->SetParLimits(2,-1.,10.);
1260 // TString optio = "ROW"; // "RO"
1261 TString optio = "RO"; // "RO"
1262 if ( debug ){
1263 optio += "NV";
1264 if ( draw ) optio += "V";
1265 } else {
1266 optio += "NQ";
1267 if ( draw ) optio += "Q";
1268 };
1269 //
1270 if ( debug ) printf(" OK, start the fitting procedure... options: %s \n",optio.Data());
1271 //
1272 // fitresult = th->Fit("lfit",optio);
1273 // lfit->Update();
1274 fitresult = gh->Fit("lfit",optio);
1275 // fitresult = gh->Fit("lfit","ROQW");
1276 //
1277 if ( debug ) printf(" the fit is done! result: %i \n",fitresult);
1278 //
1279 E0 = lfit->GetParameter(0);
1280 a = lfit->GetParameter(1);
1281 b = lfit->GetParameter(2);
1282 errE0 = lfit->GetParError(0);
1283 erra = lfit->GetParError(1);
1284 errb = lfit->GetParError(2);
1285 chi2 = lfit->GetChisquare();
1286 ndf = lfit->GetNDF();
1287 Float_t tmax = 0.;
1288 if ( debug ) printf(" Parameters are retrieved \n");
1289 if ( b != 0 ) tmax = (a - 1.)/b;
1290 //
1291 if ( fitresult != 0 ){
1292 if ( debug ) printf(" The fit failed, no integrals calculation and asymm is set to -1. \n");
1293 asymm = -1.;
1294 defE0 = -1.;
1295 } else {
1296 if ( slmax.MaybeRegexp() ) lmax = this->Evaluate(slmax,tmax,X0pl);
1297 if ( sumax.MaybeRegexp() ) umax = this->Evaluate(sumax,tmax,X0pl);
1298 Int_t npp = 1000;
1299 double *xpp=new double[npp];
1300 double *wpp=new double[npp];
1301 lfit->CalcGaussLegendreSamplingPoints(npp,xpp,wpp,1e-12);
1302 defE0 = lfit->IntegralFast(npp,xpp,wpp,lmax,umax);
1303 //
1304 double *xp=new double[npp];
1305 double *wp=new double[npp];
1306 lfit->CalcGaussLegendreSamplingPoints(npp,xp,wp,1e-12);
1307 Float_t imax = lfit->IntegralFast(npp,xp,wp,0.,tmax);
1308 // Float_t imax = lfit->Integral(0.,tmax);
1309 if ( debug ) printf(" Integral till maximum (%f): %f \n",tmax,imax);
1310 Int_t np = 1000;
1311 double *x=new double[np];
1312 double *w=new double[np];
1313 lfit->CalcGaussLegendreSamplingPoints(np,x,w,1e-12);
1314 Float_t i10max = lfit->IntegralFast(np,x,w,0.,10.*tmax);
1315 delete x;
1316 delete w;
1317 delete xp;
1318 delete wp;
1319 delete xpp;
1320 delete wpp;
1321 // Float_t i10max = lfit->Integral(0.,10.*tmax);
1322 if ( debug ) printf(" Integral: %f \n",i10max);
1323 //
1324 if ( i10max != imax ){
1325 asymm = imax / (i10max-imax);
1326 } else {
1327 if ( debug ) printf(" i10max == imax, asymm undefined\n");
1328 asymm = -2.;
1329 };
1330 if ( asymm != asymm ){
1331 if ( debug ) printf(" asymm is nan \n");
1332 asymm = -3.;
1333 };
1334 //lfit->Integral(0.,tmax)/(lfit->Integral(0.,10.*tmax)-lfit->Integral(0.,tmax));
1335 if ( debug ) printf(" Asymmetry has been calculated \n");
1336 };
1337 //
1338 if ( asymm < 0. || ndf <= 0. || chi2 < 0. || tmax < 0. ){
1339 if ( debug ) printf(" Funny asymm||ndf||chi2||tmax values, fit failed \n");
1340 fitresult = 100;
1341 };
1342 //
1343 if ( draw ){
1344 //
1345 tc->cd();
1346 // gStyle->SetOptStat(11111);
1347 tc->SetTitle();
1348 th->SetTitle("");
1349 th->SetName("");
1350 th->SetMarkerStyle(20);
1351 // axis titles
1352 th->SetXTitle("Depth [X0]");
1353 th->SetYTitle("Energy [MIP]");
1354 // th->DrawCopy("Perror");
1355 th->DrawCopy();
1356 gh->SetMarkerSize(1.);
1357 gh->SetMarkerStyle(21);
1358 // gh->SetMarkerColor(kRed);
1359 gh->Draw("Psame");
1360 lfit->Draw("same");
1361 tc->Modified();
1362 tc->Update();
1363 //
1364 gStyle->SetLabelSize(0);
1365 gStyle->SetNdivisions(1,"XY");
1366 //
1367 } else {
1368 if ( th ) th->Delete();
1369 };
1370 //
1371 // gMinuit->SetPrintLevel(1);
1372 // gMinuit->mnamin();
1373 // // gMinuit->mnfree(0);
1374 // gMinuit->mncler();
1375 // gMinuit->mnrset(1);
1376 // gMinuit->Delete();
1377 //ROOT::Math::Minimizer::Clear();
1378 // gMinuit->mnhelp("*");
1379 //
1380 // delete lfit;
1381 //
1382 };
1383
1384 void CaloLong::Draw(){
1385 //
1386 Process();
1387 Draw(-1);
1388 };
1389
1390 void CaloLong::Draw(Int_t view){
1391 //
1392 Int_t minv = 0;
1393 Int_t maxv = 0;
1394 //
1395 if ( view == -1 ){
1396 maxv = -1;
1397 } else {
1398 minv = view;
1399 maxv = view+1;
1400 };
1401 //
1402 Process();
1403 //
1404 gStyle->SetLabelSize(0.04);
1405 gStyle->SetNdivisions(510,"XY");
1406 //
1407 if ( maxv != -1 ){
1408 for (Int_t v=minv; v<maxv;v++){
1409 TString hid = Form("clongv%i%s",v,suf.Data());
1410 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
1411 if ( tc ){
1412 // tc->Clear();
1413 } else {
1414 tc = new TCanvas(hid,hid);
1415 };
1416 //
1417 TString thid = Form("hlongv%i%s",v,suf.Data());
1418 TH1F *th = dynamic_cast<TH1F*>(gDirectory->FindObject(thid));
1419 if ( th ) th->Delete();
1420 // th->Clear();
1421 // th->Reset();
1422 // } else {
1423 th = new TH1F(thid,thid,22,-0.5,21.5);
1424 // };
1425 tc->cd();
1426 //
1427 for (Int_t st=0;st<22;st++){
1428 th->Fill(st,eplane[v][st]);
1429 };
1430 th->Draw();
1431 tc->Modified();
1432 tc->Update();
1433 };
1434 } else {
1435 //
1436 TString hid = Form("clongvyvx%s",suf.Data());
1437 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
1438 if ( tc ){
1439 } else {
1440 tc = new TCanvas(hid,hid);
1441 };
1442 //
1443 TString thid = Form("hlongvyvx%s",suf.Data());
1444 TH1F *th = dynamic_cast<TH1F*>(gDirectory->FindObject(thid));
1445 if ( th ) th->Delete();
1446 th = new TH1F(thid,thid,44,-0.5,43.5);
1447 tc->cd();
1448 Int_t pp=0;
1449 for (Int_t st=0;st<22;st++){
1450 for (Int_t v=1; v>=0;v--){
1451 //
1452 th->Fill(pp,eplane[v][st]);
1453 //
1454 pp++;
1455 };
1456 };
1457 th->Draw();
1458 tc->Modified();
1459 tc->Update();
1460 };
1461 //
1462 gStyle->SetLabelSize(0);
1463 gStyle->SetNdivisions(1,"XY");
1464 //
1465 };
1466
1467 void CaloLong::Delete(){
1468 Clear();
1469 //delete this;
1470 };
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490

  ViewVC Help
Powered by ViewVC 1.1.23