/[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.13 - (show annotations) (download)
Wed Aug 19 07:25:34 2009 UTC (15 years, 6 months ago) by mocchiut
Branch: MAIN
Changes since 1.12: +17 -13 lines
Features added

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 };
680
681 void CaloLong::MaskSection(TString sec){
682 sec.ToUpper();
683 if ( sec.Contains("XO") ) maskXO = true;
684 if ( sec.Contains("YO") ) maskYO = true;
685 if ( sec.Contains("XE") ) maskXE = true;
686 if ( sec.Contains("YE") ) maskYE = true;
687 }
688
689 void CaloLong::UnMaskSections(){
690 this->UnMaskSection("XEXOYEYO");
691 }
692
693 void CaloLong::UnMaskSection(TString sec){
694 sec.ToUpper();
695 if ( sec.Contains("XO") ) maskXO = false;
696 if ( sec.Contains("YO") ) maskYO = false;
697 if ( sec.Contains("XE") ) maskXE = false;
698 if ( sec.Contains("YE") ) maskYE = false;
699 }
700
701 void CaloLong::Clear(){
702 //
703 memset(eplane,0, 2*22*sizeof(Float_t));
704 //
705 chi2 = 0.;
706 ndf = 0.;
707 E0 = 0.;
708 defE0 = 0.;
709 a = 0.;
710 b = 0.;
711 errE0 = 0.;
712 erra = 0.;
713 errb = 0.;
714 etmax = 0.;
715 asymm = 0.;
716 fitresult = 0;
717 //
718 X0pl = 0.76;
719 //
720 };
721
722 void CaloLong::Print(){
723 //
724 Fit();
725 //
726 printf("==================== Calorimeter Longitudinal Profile =======================\n");
727 printf(" OBT: %u PKT: %u ATIME: %u \n",OBT,PKT,atime);
728 printf(" fitresult:.. %i\n",fitresult);
729 printf(" chi2 :.. %f\n",chi2);
730 printf(" ndf :.. %f\n",ndf);
731 printf(" nchi2 :.. %f\n",chi2/ndf);
732 printf(" E0 :.. %f\n",E0);
733 printf(" E0/260. :.. %f\n",E0/260.);
734 printf(" defE0 :.. %f\n",defE0);
735 printf(" lower :.. %f\n",lmax);
736 printf(" upper :.. %f\n",umax);
737 printf(" s lower :.. %s\n",slmax.Data());
738 printf(" s upper :.. %s\n",sumax.Data());
739 printf(" a :.. %f\n",a);
740 printf(" b :.. %f\n",b);
741 printf(" errE0 :.. %f\n",errE0);
742 printf(" erra :.. %f\n",erra);
743 printf(" errb :.. %f\n",errb);
744 printf(" asymm :.. %f\n",asymm);
745 printf(" tmax :.. %f\n",((a-1.)/b));
746 printf(" etmax :.. %f\n",etmax);
747 printf(" X0pl :.. %f\n",X0pl);
748 printf("========================================================================\n");
749 //
750 };
751
752 void CaloLong::SetNoWpreSampler(Int_t n){
753 //
754 if ( NC+n <= 22 && NC+n >= 0 ){
755 N = n;
756 } else {
757 printf(" ERROR! Calorimeter is made of 22 W planes\n");
758 printf(" you are giving N presampler = %i and N calo = %i \n",n,NC);
759 printf(" WARNING: using default values NWpre = 0, NWcalo = 22\n");
760 NC = 22;
761 N = 0;
762 };
763 }
764
765 void CaloLong::SetNoWcalo(Int_t n){
766 if ( N+n <= 22 && N+n >= 0 ){
767 NC = n;
768 } else {
769 printf(" ERROR! Calorimeter is made of 22 W planes\n");
770 printf(" you are giving N W presampler = %i and N W calo = %i \n",N,n);
771 printf(" WARNING: using default values NWpre = 0, NWcalo = 22\n");
772 NC = 22;
773 N = 0;
774 };
775 }
776
777 void CaloLong::SplitInto(Int_t NoWpreSampler, Int_t NoWcalo){
778 this->SetNoWpreSampler(0);
779 this->SetNoWcalo(0);
780 if ( NoWpreSampler < NoWcalo ){
781 this->SetNoWpreSampler(NoWpreSampler);
782 this->SetNoWcalo(NoWcalo);
783 } else {
784 this->SetNoWcalo(NoWcalo);
785 this->SetNoWpreSampler(NoWpreSampler);
786 };
787 }
788
789 void CaloLong::Process(){
790 //
791 if ( !L2 ){
792 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
793 printf(" ERROR: CaloHough variables not filled \n");
794 return;
795 };
796 //
797 Bool_t newentry = false;
798 //
799 if ( L2->IsORB() ){
800 if ( L2->GetOrbitalInfo()->pkt_num != PKT || L2->GetOrbitalInfo()->OBT != OBT || L2->GetOrbitalInfo()->absTime != atime ){
801 newentry = true;
802 OBT = L2->GetOrbitalInfo()->OBT;
803 PKT = L2->GetOrbitalInfo()->pkt_num;
804 atime = L2->GetOrbitalInfo()->absTime;
805 };
806 } else {
807 newentry = true;
808 };
809 //
810 if ( !newentry ) return;
811 //
812 if ( debug ) printf(" Start processing event at OBT %u PKT %u time %u \n",OBT,PKT,atime);
813 //
814 Clear();
815 //
816 // let's start
817 //
818 if ( cont ){
819 for (Int_t i=0; i<22; i++){
820 if ( i == (18+N) ){
821 mask18b = 18 + N;
822 break;
823 };
824 };
825 };
826 //
827 if ( sel ){
828 for (Int_t i=0; i<22; i++){
829 if ( i == (18-N) ){
830 mask18b = 18 - N;
831 break;
832 };
833 };
834 };
835 //
836 // if ( mask18b == 18 ) mask18b = -1;
837 //
838 Int_t view = 0;
839 Int_t plane = 0;
840 Int_t strip = 0;
841 Float_t mip = 0.;
842 Bool_t gof = true;
843 for (Int_t i=0; i < L2->GetCaloLevel1()->istrip; i++){
844 mip = L2->GetCaloLevel1()->DecodeEstrip(i,view,plane,strip);
845 gof = true;
846 if ( maskXE && (plane%2)==0 && view==1 ) gof = false;
847 if ( maskXO && (plane%2)!=0 && view==1 ) gof = false;
848 if ( maskYE && (plane%2)!=0 && view==0 ) gof = false;
849 if ( maskYO && (plane%2)==0 && view==0 ) gof = false;
850 if ( gof ) eplane[view][plane] += mip;
851 };
852 //
853 // inclination factor (stolen from Daniele's code)
854 //
855 Float_t ytgx = 0;
856 Float_t ytgy = 0;
857 ytgx = 0.76 * clp->tanx[0];
858 ytgy = 0.76 * clp->tany[0];
859 X0pl = sqrt( pow(0.76,2.) + pow(ytgx,2.) + pow(ytgy,2.) );
860 //
861 // Find experimental plane of maximum
862 //
863 Int_t pmax = 0;
864 Int_t vmax = 0;
865 Float_t emax = 0.;
866 for (Int_t v=0; v<2; v++){
867 for (Int_t i=0; i<22; i++){
868 if ( eplane[v][i] > emax ){
869 emax = eplane[v][i];
870 vmax = v;
871 pmax = i;
872 };
873 };
874 };
875 //
876 //
877 //
878 if ( vmax == 0 ) pmax++;
879 etmax = pmax * X0pl;
880 //
881 if ( debug ) this->Print();
882 if ( debug ) printf(" exit \n");
883 //
884 };
885
886 void CaloLong::SetEnergies(Float_t myene[][22]){
887 //
888 if ( !L2 ){
889 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
890 printf(" ERROR: CaloHough variables not filled \n");
891 return;
892 };
893 //
894 Bool_t newentry = false;
895 //
896 if ( L2->IsORB() ){
897 if ( L2->GetOrbitalInfo()->pkt_num != PKT || L2->GetOrbitalInfo()->OBT != OBT || L2->GetOrbitalInfo()->absTime != atime ){
898 newentry = true;
899 OBT = L2->GetOrbitalInfo()->OBT;
900 PKT = L2->GetOrbitalInfo()->pkt_num;
901 atime = L2->GetOrbitalInfo()->absTime;
902 };
903 } else {
904 newentry = true;
905 };
906 //
907 if ( !newentry ) return;
908 //
909 if ( debug ) printf(" SET ENERGIES Start processing event at OBT %u PKT %u time %u \n",OBT,PKT,atime);
910 //
911 Clear();
912 //
913 // let's start
914 //
915 if ( cont ){
916 for (Int_t i=0; i<22; i++){
917 if ( i == (18+N) ){
918 mask18b = 18 + N;
919 break;
920 };
921 };
922 };
923 //
924 if ( sel ){
925 for (Int_t i=0; i<22; i++){
926 if ( i == (18-N) ){
927 mask18b = 18 - N;
928 break;
929 };
930 };
931 };
932 //
933 // if ( mask18b == 18 ) mask18b = -1;
934 //
935 Int_t view = 0;
936 Int_t plane = 0;
937 Bool_t gof = true;
938 for (view=0; view < 2; view++){
939 for (plane=0; plane < 22; plane++){
940 gof = true;
941 if ( maskXE && (plane%2)==0 && view==1 ) gof = false;
942 if ( maskXO && (plane%2)!=0 && view==1 ) gof = false;
943 if ( maskYE && (plane%2)!=0 && view==0 ) gof = false;
944 if ( maskYO && (plane%2)==0 && view==0 ) gof = false;
945 if ( gof ) eplane[view][plane] = myene[view][plane];
946 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]);
947 };
948 };
949 //
950 // inclination factor (stolen from Daniele's code)
951 //
952 Float_t ytgx = 0;
953 Float_t ytgy = 0;
954 ytgx = 0.76 * clp->tanx[0];
955 ytgy = 0.76 * clp->tany[0];
956 X0pl = sqrt( pow(0.76,2.) + pow(ytgx,2.) + pow(ytgy,2.) );
957 //
958 // Find experimental plane of maximum
959 //
960 Int_t pmax = 0;
961 Int_t vmax = 0;
962 Float_t emax = 0.;
963 for (Int_t v=0; v<2; v++){
964 for (Int_t i=0; i<22; i++){
965 if ( eplane[v][i] > emax ){
966 emax = eplane[v][i];
967 vmax = v;
968 pmax = i;
969 };
970 };
971 };
972 //
973 //
974 //
975 if ( vmax == 0 ) pmax++;
976 etmax = pmax * X0pl;
977 //
978 if ( debug ) this->Print();
979 if ( debug ) printf(" exit \n");
980 //
981 };
982
983 Double_t ccurve(Double_t *ti,Double_t *par){
984 //
985 Double_t t = *ti;
986 Double_t cE0 = par[0];
987 Double_t ca = par[1];
988 Double_t cb = par[2];
989 Double_t gammaa = TMath::Gamma(ca);
990 //
991 Double_t value = cE0 * cb * ( pow((cb * t),(ca - 1)) * exp( -cb * t ) ) / gammaa;
992 //
993 return value;
994 //
995 }
996
997 void CaloLong::Fit(){
998 this->Fit(false);
999 };
1000
1001 Float_t CaloLong::Evaluate(TString s, Float_t tmax, Float_t X0pl){
1002 /* SAMPLE OUTPUT:
1003 Enter Infix Expression : A + B + C / (E - F)
1004 Postfix Expression is : A B + C E F - / +
1005 */
1006 if ( !s.Contains("tmax") && !s.Contains("X0pl") ){
1007 printf(" ERROR, the input formula must contain \"t\"\n");
1008 return 0.;
1009 };
1010 if ( tmax != tmax ){
1011 printf(" ERROR, tmax is nan! \n");
1012 return 0.;
1013 };
1014 TString g=Form("%f",tmax);
1015 TString h=Form("%f",X0pl);
1016 TString *ts= new TString("");
1017 ts->Prepend(s.Data());
1018 ts->ReplaceAll("tmax",4,g.Data(),g.Capacity());
1019 ts->ReplaceAll("X0pl",4,h.Data(),h.Capacity());
1020 ts->Prepend("(");
1021 ts->Append(")");
1022 if ( debug ) printf(" ts %s tssize %i capac %i s %s g %s \n",ts->Data(),ts->Sizeof(),ts->Capacity(),s.Data(),g.Data());
1023 char in[50],post[50];
1024 strcpy(&in[0]," ");
1025 strcpy(&in[0],ts->Data());
1026 strcpy(&post[0]," ");
1027 if ( debug ) printf("Infix Expression is : ##%s##\n",&in[0]);
1028 infix2postfix(&in[0],&post[0],1);
1029 if ( debug ) printf("Postfix Expression is : ##%s##\n",&post[0]);
1030 /* SAMPLE OUTPUT:
1031 Enter Postfix Expression : 3 5 + 2 /
1032 3 5 + 2 / EQUALS 4
1033 */
1034 Float_t res = evaluate(&post[0]);
1035 if ( debug ) printf("%s EQUALS %f\n",post,res);
1036 return res;
1037 }
1038
1039 void CaloLong::Fit(Bool_t draw){
1040 //
1041 Process();
1042 //
1043 //
1044 if ( !L2 ){
1045 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
1046 printf(" ERROR: CaloHough variables not filled \n");
1047 return;
1048 };
1049 //
1050 Bool_t newentry = false;
1051 //
1052 if ( L2->IsORB() ){
1053 if ( L2->GetOrbitalInfo()->pkt_num != fPKT || L2->GetOrbitalInfo()->OBT != fOBT || L2->GetOrbitalInfo()->absTime != fatime ){
1054 newentry = true;
1055 fOBT = L2->GetOrbitalInfo()->OBT;
1056 fPKT = L2->GetOrbitalInfo()->pkt_num;
1057 fatime = L2->GetOrbitalInfo()->absTime;
1058 };
1059 } else {
1060 newentry = true;
1061 };
1062 //
1063 if ( !newentry ) return;
1064 //
1065 if ( debug ) printf(" Start fitting event at OBT %u PKT %u time %u \n",fOBT,fPKT,fatime);
1066 //
1067 if ( draw ){
1068 gStyle->SetLabelSize(0.04);
1069 gStyle->SetNdivisions(510,"XY");
1070 };
1071 //
1072 TString hid = Form("clongfit%s",suf.Data());
1073 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
1074 // if ( tc ) tc->Delete();
1075 // if ( tc ) tc->Close();
1076 if ( !tc && draw ){
1077 tc = new TCanvas(hid,hid);
1078 } else {
1079 if ( tc ) tc->cd();
1080 };
1081 //
1082 TString thid = Form("hlongfit%s",suf.Data());
1083 TH2F *th = dynamic_cast<TH2F*>(gDirectory->FindObject(thid));
1084 if ( th ) th->Delete();
1085 //
1086 TString ghid = Form("glongfit%s",suf.Data());
1087 TGraphErrors *gh = dynamic_cast<TGraphErrors*>(gDirectory->FindObject(ghid));
1088 if ( gh ) gh->Delete();
1089 //
1090 Float_t xpos = 0.;
1091 Float_t enemip = 0.;
1092 Float_t xmax = NC * X0pl + 0.2;
1093 //
1094 Double_t xxx[50];
1095 Double_t exx[50];
1096 Double_t yyy[50];
1097 Double_t eyy[50];
1098 Int_t numpo = 0;
1099 memset(xxx,0,50*sizeof(Double_t));
1100 memset(exx,0,50*sizeof(Double_t));
1101 memset(yyy,0,50*sizeof(Double_t));
1102 memset(eyy,0,50*sizeof(Double_t));
1103 //
1104 // AGH, BUG!
1105 //
1106 Int_t mmin = 0;
1107 Int_t mmax = 0;
1108 if ( cont ){
1109 mmin = N;
1110 mmax = NC+N;
1111 } else {
1112 mmin = 0;
1113 mmax = NC;
1114 };
1115 //
1116 Float_t emax = 0.;
1117 Float_t qtotparz = 0.;
1118 for (Int_t st=mmin;st<mmax+1;st++){
1119 enemip = 0.;
1120 xpos = (st - mmin) * X0pl;
1121 //
1122 if ( xyaverage ){
1123 //
1124 if ( st > mmin && st < mmax ){
1125 if ( no18x && ( st == 18+1 || st == mask18b+1 )){
1126 if ( !maskYO ){
1127 enemip = 2. * eplane[1][st];
1128 } else {
1129 enemip = eplane[1][st];
1130 };
1131 } else {
1132 enemip = eplane[0][st-1] + eplane[1][st];
1133 };
1134 } else {
1135 if ( st == mmin ){
1136 if ( !maskYE ){
1137 enemip = 2. * eplane[1][st];
1138 } else {
1139 enemip = eplane[1][st];
1140 };
1141 };
1142 if ( st == mmax ){
1143 if ( !maskXE ){
1144 enemip = 2. * eplane[0][st-1];
1145 } else {
1146 enemip = eplane[0][st-1];
1147 };
1148 };
1149 };
1150 //
1151 qtotparz += enemip;
1152 if ( enemip > emax ) emax = enemip;
1153 if ( enemip > 0. ){
1154 xxx[numpo] = xpos;
1155 exx[numpo] = 0.1;
1156 yyy[numpo] = enemip;
1157 eyy[numpo] = sqrt(enemip*3.)+sqrt(5.);
1158 numpo++;
1159 // th->Fill(xpos,enemip);
1160 if ( debug ) printf(" AVE Filling: st %i xpos %f energy %f \n",st,xpos,enemip);
1161 };
1162 } else {
1163 for (Int_t jj=0; jj<2; jj++){
1164 if ( st > mmin && st < mmax ){
1165 if ( jj == 0 && no18x && ( st == 18+1 || st == mask18b+1 )){
1166 enemip = eplane[1][st];
1167 } else {
1168 if ( jj == 0 ){
1169 enemip = eplane[jj][st-1];
1170 } else {
1171 enemip = eplane[jj][st];
1172 };
1173 };
1174 } else {
1175 if ( st == mmin && jj == 1 ){
1176 enemip = eplane[jj][st];
1177 };
1178 if ( st == mmax && jj == 0){
1179 enemip = eplane[jj][st-1];
1180 };
1181 };
1182 //
1183 qtotparz += enemip;
1184 if ( enemip > emax ) emax = enemip;
1185 if ( enemip > 0. ){
1186 xxx[numpo] = xpos;
1187 exx[numpo] = 0.1;
1188 yyy[numpo] = enemip;
1189 eyy[numpo] = sqrt(enemip*3.)+sqrt(5.);
1190 // eyy[numpo] = sqrt(enemip)/(st*0.95);
1191 numpo++;
1192 // th->Fill(xpos,enemip);
1193 if ( debug ) printf(" XY Filling: st %i xpos %f energy %f \n",st,xpos,enemip);
1194 };
1195 };
1196 };
1197
1198 //
1199 // for (Int_t v=1; v>=0;v--)// {
1200 // //
1201 // if ( v == 1 ){
1202 // xpos = (st - N) * X0pl;
1203 // } else {
1204 // xpos = (st + 1 - N) * X0pl;
1205 // };
1206 // //
1207 // if ( no18x && st == 18 && v == 0 ){
1208 // // skip plane 18x
1209 // } else {
1210 // if ( v == 1 && st == mask18b ){
1211 // // emulate plane 18x
1212 // } else {
1213 // if ( eplane[v][st] > 0. ){
1214 // th->Fill(xpos,eplane[v][st]);
1215 // if ( debug ) printf(" Filling: st %i v %i xpos %f energy %f \n",st,v,xpos,eplane[v][st]);
1216 // };
1217 // };
1218 // };
1219 // //
1220 // };
1221 };
1222 //
1223 // th = new TH2F(thid,thid,int(NC*1.5),-0.2,xmax);
1224 th = new TH2F(thid,thid,1000,-0.2,xmax,1000,0.,emax*1.2);
1225 gh = new TGraphErrors(numpo,xxx,yyy,exx,eyy);
1226 TString fnam=Form("lfit%s",suf.Data());
1227 TF1 *lfit = dynamic_cast<TF1*>(gDirectory->FindObject(fnam));
1228 if ( lfit ) lfit->Delete();
1229 lfit = new TF1(fnam,ccurve,0.,xmax,3);
1230 if ( debug ) printf("qtot %f qtotparz %f \n",clp->qtot,qtotparz);
1231 E0 = qtotparz;
1232 // E0 = clp->qtot;
1233 a = 5.;
1234 b = 0.5;
1235 if ( debug ) printf(" STARTING PARAMETERS: E0 %f a %f b %f \n",E0,a,b);
1236 lfit->SetParameters(E0,a,b);
1237 // lfit->SetParLimits(0,0.,1000.);
1238 // lfit->SetParLimits(1,-1.,80.);
1239 // lfit->SetParLimits(2,-1.,10.);
1240 // TString optio = "ROW"; // "RO"
1241 TString optio = "RO"; // "RO"
1242 if ( debug ){
1243 optio += "NV";
1244 if ( draw ) optio += "V";
1245 } else {
1246 optio += "NQ";
1247 if ( draw ) optio += "Q";
1248 };
1249 //
1250 if ( debug ) printf(" OK, start the fitting procedure...\n");
1251 //
1252 // fitresult = th->Fit("lfit",optio);
1253 fitresult = gh->Fit("lfit",optio);
1254 //
1255 if ( debug ) printf(" the fit is done! result: %i \n",fitresult);
1256 //
1257 E0 = lfit->GetParameter(0);
1258 a = lfit->GetParameter(1);
1259 b = lfit->GetParameter(2);
1260 errE0 = lfit->GetParError(0);
1261 erra = lfit->GetParError(1);
1262 errb = lfit->GetParError(2);
1263 chi2 = lfit->GetChisquare();
1264 ndf = lfit->GetNDF();
1265 Float_t tmax = 0.;
1266 if ( debug ) printf(" Parameters are retrieved \n");
1267 if ( b != 0 ) tmax = (a - 1.)/b;
1268 //
1269 if ( fitresult != 0 ){
1270 if ( debug ) printf(" The fit failed, no integrals calculation and asymm is set to -1. \n");
1271 asymm = -1.;
1272 defE0 = -1.;
1273 } else {
1274 if ( slmax.MaybeRegexp() ) lmax = this->Evaluate(slmax,tmax,X0pl);
1275 if ( sumax.MaybeRegexp() ) umax = this->Evaluate(sumax,tmax,X0pl);
1276 Int_t npp = 1000;
1277 double *xpp=new double[npp];
1278 double *wpp=new double[npp];
1279 lfit->CalcGaussLegendreSamplingPoints(npp,xpp,wpp,1e-12);
1280 defE0 = lfit->IntegralFast(npp,xpp,wpp,lmax,umax);
1281 //
1282 double *xp=new double[npp];
1283 double *wp=new double[npp];
1284 lfit->CalcGaussLegendreSamplingPoints(npp,xp,wp,1e-12);
1285 Float_t imax = lfit->IntegralFast(npp,xp,wp,0.,tmax);
1286 // Float_t imax = lfit->Integral(0.,tmax);
1287 if ( debug ) printf(" Integral till maximum (%f): %f \n",tmax,imax);
1288 Int_t np = 1000;
1289 double *x=new double[np];
1290 double *w=new double[np];
1291 lfit->CalcGaussLegendreSamplingPoints(np,x,w,1e-12);
1292 Float_t i10max = lfit->IntegralFast(np,x,w,0.,10.*tmax);
1293 delete x;
1294 delete w;
1295 delete xp;
1296 delete wp;
1297 // Float_t i10max = lfit->Integral(0.,10.*tmax);
1298 if ( debug ) printf(" Integral: %f \n",i10max);
1299 //
1300 if ( i10max != imax ){
1301 asymm = imax / (i10max-imax);
1302 } else {
1303 if ( debug ) printf(" i10max == imax, asymm undefined\n");
1304 asymm = -2.;
1305 };
1306 if ( asymm != asymm ){
1307 if ( debug ) printf(" asymm is nan \n");
1308 asymm = -3.;
1309 };
1310 //lfit->Integral(0.,tmax)/(lfit->Integral(0.,10.*tmax)-lfit->Integral(0.,tmax));
1311 if ( debug ) printf(" Asymmetry has been calculated \n");
1312 };
1313 //
1314 if ( asymm < 0. || ndf <= 0. || chi2 < 0. || tmax < 0. ){
1315 if ( debug ) printf(" Funny asymm||ndf||chi2||tmax values, fit failed \n");
1316 fitresult = 100;
1317 };
1318 //
1319 if ( draw ){
1320 //
1321 tc->cd();
1322 // gStyle->SetOptStat(11111);
1323 tc->SetTitle();
1324 th->SetTitle("");
1325 th->SetName("");
1326 th->SetMarkerStyle(20);
1327 // axis titles
1328 th->SetXTitle("Depth [X0]");
1329 th->SetYTitle("Energy [MIP]");
1330 // th->DrawCopy("Perror");
1331 th->DrawCopy();
1332 gh->SetMarkerSize(1.);
1333 gh->SetMarkerStyle(21);
1334 // gh->SetMarkerColor(kRed);
1335 gh->Draw("Psame");
1336 lfit->Draw("same");
1337 tc->Modified();
1338 tc->Update();
1339 //
1340 gStyle->SetLabelSize(0);
1341 gStyle->SetNdivisions(1,"XY");
1342 //
1343 } else {
1344 if ( th ) th->Delete();
1345 };
1346 //
1347 // delete lfit;
1348 //
1349 };
1350
1351 void CaloLong::Draw(){
1352 //
1353 Process();
1354 Draw(-1);
1355 };
1356
1357 void CaloLong::Draw(Int_t view){
1358 //
1359 Int_t minv = 0;
1360 Int_t maxv = 0;
1361 //
1362 if ( view == -1 ){
1363 maxv = -1;
1364 } else {
1365 minv = view;
1366 maxv = view+1;
1367 };
1368 //
1369 Process();
1370 //
1371 gStyle->SetLabelSize(0.04);
1372 gStyle->SetNdivisions(510,"XY");
1373 //
1374 if ( maxv != -1 ){
1375 for (Int_t v=minv; v<maxv;v++){
1376 TString hid = Form("clongv%i%s",v,suf.Data());
1377 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
1378 if ( tc ){
1379 // tc->Clear();
1380 } else {
1381 tc = new TCanvas(hid,hid);
1382 };
1383 //
1384 TString thid = Form("hlongv%i%s",v,suf.Data());
1385 TH1F *th = dynamic_cast<TH1F*>(gDirectory->FindObject(thid));
1386 if ( th ) th->Delete();
1387 // th->Clear();
1388 // th->Reset();
1389 // } else {
1390 th = new TH1F(thid,thid,22,-0.5,21.5);
1391 // };
1392 tc->cd();
1393 //
1394 for (Int_t st=0;st<22;st++){
1395 th->Fill(st,eplane[v][st]);
1396 };
1397 th->Draw();
1398 tc->Modified();
1399 tc->Update();
1400 };
1401 } else {
1402 //
1403 TString hid = Form("clongvyvx%s",suf.Data());
1404 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
1405 if ( tc ){
1406 } else {
1407 tc = new TCanvas(hid,hid);
1408 };
1409 //
1410 TString thid = Form("hlongvyvx%s",suf.Data());
1411 TH1F *th = dynamic_cast<TH1F*>(gDirectory->FindObject(thid));
1412 if ( th ) th->Delete();
1413 th = new TH1F(thid,thid,44,-0.5,43.5);
1414 tc->cd();
1415 Int_t pp=0;
1416 for (Int_t st=0;st<22;st++){
1417 for (Int_t v=1; v>=0;v--){
1418 //
1419 th->Fill(pp,eplane[v][st]);
1420 //
1421 pp++;
1422 };
1423 };
1424 th->Draw();
1425 tc->Modified();
1426 tc->Update();
1427 };
1428 //
1429 gStyle->SetLabelSize(0);
1430 gStyle->SetNdivisions(1,"XY");
1431 //
1432 };
1433
1434 void CaloLong::Delete(){
1435 Clear();
1436 //delete this;
1437 };
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457

  ViewVC Help
Powered by ViewVC 1.1.23