/[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.17 - (show annotations) (download)
Mon Dec 14 14:49:16 2009 UTC (15 years, 2 months ago) by mocchiut
Branch: MAIN
Changes since 1.16: +19 -2 lines
Do not use plane 18x by default

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

  ViewVC Help
Powered by ViewVC 1.1.23