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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.12 - (show annotations) (download)
Fri Jun 17 11:56:42 2011 UTC (13 years, 6 months ago) by mocchiut
Branch: MAIN
Changes since 1.11: +48 -4 lines
Bugs fixed and new methods implemented

1 #include <CaloBragg.h>
2
3
4 ClassImp(CaloBragg);
5 //--------------------------------------
6 /*
7 * Default constructor
8 */
9 CaloBragg::CaloBragg(){
10 Clear();
11 };
12
13 CaloBragg::CaloBragg(PamLevel2 *l2p){
14 //
15 Clear();
16 LoadParam();
17 //
18 L2 = l2p;
19 //
20 if ( !L2->IsORB() ) printf(" WARNING: OrbitalInfo Tree is needed, the plugin could not work properly without it \n");
21 //
22 OBT = 0;
23 PKT = 0;
24 atime = 0;
25 //
26 debug = false;
27 usetrack = false;
28 usepl18x = false;
29 newchi2 = false;
30 usenewBB = false;
31 fzeta = -1.;
32 //
33 };
34
35 void CaloBragg::Clear(){
36 //
37 ndf = 0;
38 tr = 0;
39 sntr = 0;
40 // qtchi2 = 0.;
41 // qtz = 0.;
42 // qtetot = 0.;
43 // qtpskip = 0.;
44 lpchi2 = 0.;
45 lpz = 0.;
46 lpetot = 0.;
47 lppskip = 0.;
48 memset(calorimetro,0,44*2*sizeof(Float_t));
49 memset(spessore,0,4*sizeof(Float_t));
50 memset(estremi,0,2*2*sizeof(Float_t));
51 Integrale=0.;
52
53 for(Int_t l=0;l<44;l++){
54 calorimetro[l][0]=-1.;
55 }
56
57 };
58
59 void CaloBragg::Print(){
60 //
61
62 if(!debug) Process();
63 //
64 printf("========================================================================\n");
65 printf(" OBT: %u PKT: %u ATIME: %u Track %i Use track %i \n",OBT,PKT,atime,tr,usetrack);
66 printf(" first plane: %f \n", estremi[0][0]);
67 printf(" last plane: %f \n", estremi[1][0]);
68 // printf(" chi 2 from truncated mean: %f \n", qtchi2);
69 // printf(" Z from truncated mean %f: \n", qtz);
70 // printf(" energy from truncated mean %f: \n", qtetot);
71 // printf(" plane not used for truncated mean %f: \n", qtpskip);
72 printf(" chi 2 from loop %f: \n", lpchi2);
73 printf(" Z from loop %f: \n", lpz);
74 printf(" energy from loop %f: \n", lpetot);
75 printf(" plane not used for loop %f: \n", lppskip);
76 printf(" ndf: %i \n",ndf);
77 printf("========================================================================\n");
78 //
79 };
80
81 void CaloBragg::Delete(){
82 Clear();
83 //delete this;
84 };
85
86
87 void CaloBragg::Process(){
88 Process(-1);
89 };
90
91
92 void CaloBragg::CleanPlanes(Float_t epiano[22][2]){
93 // return;
94 Int_t hitplanes = 0;
95 for (Int_t i = 0; i<22; i++){
96 for (Int_t j = 1; j>=0; j--){
97 if ( epiano[i][j] > 0.7 ) hitplanes++;
98 };
99 };
100 Float_t lowlim = 0.85;
101 Float_t dedxone = 0.;
102 Float_t step1 = 0.8*L2->GetCaloLevel2()->qtot/(Float_t)hitplanes;
103 while ( dedxone < step1 ){
104 for (Int_t i = 0; i<22; i++){
105 for (Int_t j = 1; j>=0; j--){
106 if ( epiano[i][j] >= step1 && dedxone < 0.7 ) dedxone = epiano[i][j];
107 };
108 };
109 }
110 if ( dedxone < 0.7 ){
111 for (Int_t i = 0; i<22; i++){
112 for (Int_t j = 1; j>=0; j--){
113 if ( epiano[i][j] > 0. && dedxone < 0.7 ) dedxone = epiano[i][j];
114 };
115 };
116 }
117 //
118 // printf(" dedxone = %f step1 %f \n",dedxone,step1);
119 Bool_t revulsera = false;
120 Bool_t nullius = false;
121 Int_t nulliferus = 0;
122 for (Int_t i = 0; i<22; i++){
123 for (Int_t j = 1; j>=0; j--){
124 if ( epiano[i][j] < dedxone*lowlim ){
125 // printf(" %i %i epiano %f limit %f nulliferus %i nullius %i \n",i,j,epiano[i][j],dedxone*lowlim,nulliferus,nullius);
126 epiano[i][j] = 0.;
127 } else {
128 //x printf(" %i %i epiano %f limit %f nulliferus %i nullius %i \n",i,j,epiano[i][j],dedxone*lowlim,nulliferus,nullius);
129 nulliferus = 0;
130 revulsera = true;
131 };
132 if ( epiano[i][j] < 0.7 && revulsera ) nulliferus++;
133 if ( nulliferus > 10 ) nullius = true;
134 if ( nullius ) epiano[i][j] = 0.;
135 };
136 };
137
138 }
139
140
141 void CaloBragg::Process(Int_t ntr){
142 //
143 if ( !L2 ){
144 printf(" ERROR: cannot find PamLevel2 object, use the correct constructor or check your program!\n");
145 printf(" ERROR: CaloBragg variables not filled \n");
146 return;
147 };
148 //
149 Bool_t newentry = false;
150 //
151 if ( L2->IsORB() ){
152 if ( L2->GetOrbitalInfo()->pkt_num != PKT || L2->GetOrbitalInfo()->OBT != OBT || L2->GetOrbitalInfo()->absTime != atime || ntr != sntr ){
153 newentry = true;
154 OBT = L2->GetOrbitalInfo()->OBT;
155 PKT = L2->GetOrbitalInfo()->pkt_num;
156 atime = L2->GetOrbitalInfo()->absTime;
157 sntr = ntr;
158 };
159 } else {
160 newentry = true;
161 };
162 //
163 if ( !newentry ) return;
164 //
165 tr = ntr;
166 //
167 if ( debug ) printf(" Processing event at OBT %u PKT %u time %u \n",OBT,PKT,atime);
168 //
169 Clear();
170
171 //
172 //
173 //
174 Int_t view = 0;
175 Int_t plane = 0;
176 Int_t strip = 0;
177 Float_t mip = 0.;
178 Float_t epiano[22][2];
179 memset(epiano,0,22*2*sizeof(Float_t));
180 for ( Int_t i=0; i<L2->GetCaloLevel1()->istrip; i++ ){
181 //
182 mip = L2->GetCaloLevel1()->DecodeEstrip(i,view,plane,strip);
183 //
184 if ( !usepl18x && view==0 && plane==18 ) mip = 0.;
185 //
186 epiano[plane][view]+=mip;
187 //
188 //
189 };
190 //
191 this->CleanPlanes(*&epiano);
192 //
193 PamTrack *ptrack = 0;
194 CaloTrkVar *track = 0;
195 //
196 if ( usetrack ){
197 if ( ntr >= 0 ){
198 ptrack = L2->GetTrack(ntr);
199 if ( ptrack ) track = ptrack->GetCaloTrack();
200 } else {
201 track = L2->GetCaloStoredTrack(ntr); //al momento e' vera solo questa riga
202 };
203 //
204 if ( !track && ntr >= 0 ){
205 printf(" ERROR: cannot find any track!\n");
206 printf(" ERROR: CaloBragg variables not completely filled \n");
207 return;
208 };
209 } else {
210 if ( ntr >= 0 ){
211 if ( debug ) printf(" ERROR: you asked not to use a track but you are looking for track number %i !\n",ntr);
212 if ( debug ) printf(" ERROR: CaloBragg variables not completely filled \n");
213 return;
214 };
215 };
216 //
217 if(L2->GetCaloLevel2()->npcfit[0]==0 && L2->GetCaloLevel2()->npcfit[1]==0 && L2->GetCaloLevel2()->npcfit[2]==0 && L2->GetCaloLevel2()->npcfit[3]==0) return;// controllo sulla traccia nel calorimetro
218
219 //
220 for(Int_t p=0; p<22; p++){
221 for(Int_t v=0; v<2; v++){
222 /*per usare traccia non del calo camboare cibar*/
223 calorimetro[(2*p)+1-v][0] = L2->GetCaloLevel2()->cibar[p][v];//strip attraversata
224 calorimetro[(2*p)+1-v][1] = epiano[p][v]; //energia del piano //(epiano[p][v])/0.89
225 };
226 };
227
228 /*per ogni evento calcolo la conversione mip e w attraversato in equivalente Si*/
229 conversione(); // out: 1) g/cm2 Si , 2) spessoreW equivalente in Si, 3)Mip corretta per inclinazione
230
231 /*settaggio della soglia per il loop sulla determinazione del piano di partenza */
232 Float_t ordplane[44];//mi serve per la media troncata
233 memset(ordplane,0,44*sizeof(Float_t));
234
235 for(Int_t ipla=0; ipla< 2*NPLA; ipla++) ordplane[ipla]=calorimetro[ipla][1]; //energia del piano
236
237
238 //ordino tutte le energie dei piani in ordine crescente
239
240 Long64_t work[200];
241 Int_t ind = 0;
242 //Int_t l = 0;
243 Int_t RN = 0;
244 Float_t sum4 = 0.;
245 Float_t qm = 0.;
246 while ( RN < 4 && ind < 44 ){
247 qm = TMath::KOrdStat((Long64_t)44,ordplane,(Long64_t)ind,work);
248 if (qm >= 0.7 ){
249 if ( RN < 4 ){
250 sum4 += qm;
251 RN++;
252 };
253 };
254 ind++;
255 };
256 //
257 //sum4 /= (Float_t)RN;
258 Float_t Zmean = (sqrt((sum4*MIP)/(((Float_t)RN)*spessore[2])));
259 if(Zmean ==0.) Zmean=1.;
260 if ( Zmean < 1. ) Zmean = 1.;
261
262
263 /*trova primo e ultimo piano attraversati*/
264 Int_t p = 0;//contatore piani
265 //per il primo parte da 0 e va in giu'
266 while( estremi[0][1] <= 0. && p<(2*NPLA) ){ // era ==0 ma ricorda i problemi con Float == !!!!!
267 // if( (calorimetro[p][0] != -1) && (calorimetro[p][1] >50.)){
268 // if( (calorimetro[p][0] >0) && (calorimetro[p][1]*MIP >0.3)){ //0.7 mip = 70MeV soglia minima
269 if( (calorimetro[p][0] >0) && (calorimetro[p][1]*MIP >Zmean*0.7)){ // 70% della MIP
270 estremi[0][0]=p;
271 estremi[0][1]=calorimetro[p][1] *MIP; //energia in MeV
272 };
273 p++;
274 };
275
276 //ultimo parte da 44 e sale
277 p=43;
278 while( (estremi[1][1] <= 0.) && (p>(int)estremi[0][0]) ){
279 if( (calorimetro[p][0] != -1) && (calorimetro[p][1] >0.7)){
280 estremi[1][0]=p;//
281 estremi[1][1]=calorimetro[p][1] *MIP;//energia in MeV
282 };
283 p = p-1;
284 };
285 //
286
287 Float_t lastok = 0.;
288 // Bool_t goback = false;
289 for ( int o = 0; o < estremi[1][0]; o++ ){
290 //
291 if ( calorimetro[o][1] > 0.7 ) lastok = calorimetro[o][1];
292 if ( calorimetro[o][1] < 0.7 && lastok > 0. ) calorimetro[o][1] = lastok;
293 // if ( calorimetro[o][1] < 0.7 ) goback = true;
294 //
295 };
296 lastok = 0.;
297 // if ( goback ){
298 for ( int o = estremi[1][0]; o >= 0; o-- ){
299 //
300 // printf(" goback1: o %i calo %f lastok %f \n",o,calorimetro[o][1],lastok);
301 if ( o < estremi[1][0] && calorimetro[o][1] > calorimetro[o+1][1]*1.2 && lastok > 0. ) calorimetro[o][1] = lastok;
302 if ( calorimetro[o][1] > 0.7 ) lastok = calorimetro[o][1];
303 if ( calorimetro[o][1] < 0.7 && lastok > 0. ) calorimetro[o][1] = lastok;
304 // printf(" goback2: o %i calo %f lastok %f \n",o,calorimetro[o][1],lastok);
305 //
306 };
307 // };
308
309 if ( startZero ) {
310 estremi[0][0] = 0.;
311 // estremi[0][1] = 0.;
312 }
313
314 /*integrale: energia totale rilasciata nel calo (aggiungendo quella 'teorica' nel W )*/
315 for(Int_t pl=0; pl<(2*NPLA); pl++){
316 // printf(" integrale: calorimetro %f \n",calorimetro[pl][1]);
317 //calcolo intergale in unita di spessori di silicio
318 Integrale += calorimetro[pl][1] * MIP;//piano di silicio
319 // se non e'il 1o dopo l'Y (tutti i pari) c'e' il W
320 if(pl%2!=0){ //equival W in Si
321 Integrale+= 0.5*((calorimetro[pl-1][1] * MIP)+(calorimetro[pl][1] * MIP))*(spessore[1]);
322 };
323 };
324 //Integrale=24000;//Integrale*1000;
325 Integrale *= 1000.;
326
327 /*z ed energia con media troncata*/
328 // mediatroncata(); // out: 1)chi2, 2)z, 3)Etot, 4)Pskip
329
330 /*z ed energia con loop*/
331 if ( debug ) printf(" call Zdaloop with integrale %f \n",Integrale);
332 Zdaloop(); // out: 1)chi2, 2)z, 3)Etot, 4)Pskip
333
334
335 if ( debug ) this->Print();
336 if ( debug ) printf(" fine evento \n");
337 //
338 };
339
340
341 Float_t CaloBragg::Integral(){
342 Process();
343
344 Float_t dEpianiloop[44];
345 Int_t tz1=(Int_t)lpz;
346 Enetrack(&tz1, &lpetot, &estremi[0][0],&estremi[1][0], dEpianiloop);//calcola rilascio energetico sui piani da loop
347
348
349 Float_t integ = 0.;
350 for(Int_t i=0;i<=estremi[1][0];i++){
351 // integ += dEplan[i];
352 //printf(" step %i integ %f deplan %f \n",i,integ,dEplan[i]);
353 integ += dEpianiloop[i];
354 // printf(" step %i integ %f deplan %f \n",i,integ,dEpianiloop[i]);
355 }
356 return integ;
357 }
358
359 Float_t CaloBragg::LastIntegral(){
360 Process();
361
362 Float_t integ = 0.;
363 for(Int_t i=0;i<=estremi[1][0];i++){
364 integ += dEplan[i];
365 //printf(" step %i integ %f deplan %f \n",i,integ,dEplan[i]);
366 }
367 return integ;
368 }
369
370
371 void CaloBragg::Draw(){
372
373 Process();
374
375 this->Draw(0.,0.);
376
377 }
378
379 void CaloBragg::Draw(Int_t Z, Float_t enetot){
380
381 // Float_t dEpianimean[44];
382 Float_t dEpianiloop[44];
383 Float_t Depth[44];
384 // Int_t tz=(Int_t)qtz;
385 Int_t tz1= Z;
386 Float_t enet = enetot;
387 // Float_t enet = lpetot;
388
389 if ( Z > 0. && enetot > 0. ){
390 estremi[0][0] = 0;
391 estremi[1][0] = 43;
392
393
394 Float_t ytgx = 0.;
395 Float_t ytgy = 0.;
396
397 //lunghezza effettiva di silicio attraversata (mm)
398 Float_t SiCross = sqrt(SQ(ySi) + SQ(ytgx) + SQ(ytgy));
399
400 spessore[0] = (SiCross/10.) * rhoSi; //spessore silicio in g/cm2
401
402 /*tungsteno*/
403
404 //rapporto tra rilasci energetici nei due materiali
405 Float_t WCross = sqrt((yW*yW) + (ytgx*ytgx) + (ytgy*ytgy));//mm* rapporto lunghezze rad
406 //gcm2W = WCross/10. * rhoW;
407
408 // (g/cm2W)/(g/cm2Si)
409 spessore[3] = (WCross/10.) * rhoW;
410 Float_t a=(WCross/SiCross)*(rhoW/rhoSi)*(1.145/1.664); //(gcm2W)/(SiCross/10. * rhoSi)* (1.145/1.664);
411 spessore[1] = a;
412 //riscala mip allo spessore attraversato
413 spessore[2] = MIP*(SiCross/ySi);
414
415 } else {
416 tz1=(Int_t)lpz;
417 enet = lpetot;
418 // Enetrack(&tz, &qtetot, &estremi[0][0],&estremi[1][0], dEpianimean);//calcola rilascio energetico sui piani da media troncata
419
420 }
421 Enetrack(&tz1, &enet, &estremi[0][0],&estremi[1][0], dEpianiloop);//calcola rilascio energetico sui piani da loop
422
423 Float_t sp= spessore[0]*spessore[1];
424 for(Int_t i=0;i<44;i++)Depth[i]=i*sp;
425 //
426 gStyle->SetLabelSize(0.04);
427 gStyle->SetNdivisions(510,"XY");
428 //
429 TString hid = Form("cCaloBragg");
430 TCanvas *tc = dynamic_cast<TCanvas*>(gDirectory->FindObject(hid));
431 if ( tc ){
432 // tc->Clear();
433 } else {
434 tc = new TCanvas(hid,hid);
435 // tc->Divide(1,2);
436 };
437 //
438 // TString thid = Form("hCaloBragg");
439 // TH2F *th = dynamic_cast<TH2F*>(gDirectory->FindObject(thid));
440 // if ( th ) th->Delete();
441 // th->Clear();
442 // th->Reset();
443 // } else {
444 // th = new TH2F(thid,thid,300,-0.5,300.,1000,0.,150.);
445 // th->SetMarkerStyle(20);
446 // };
447 //
448 tc->cd();
449 TString thid2 = Form("hCaloBragg2");
450 TH2F *th2 = dynamic_cast<TH2F*>(gDirectory->FindObject(thid2));
451 if ( th2 ) th2->Delete();
452 th2 = new TH2F(thid2,thid2,300,-0.5,300.,1000,0.,150.);
453 th2->SetMarkerStyle(20);
454 th2->SetMarkerColor(kRed);
455 //
456 TString thid3 = Form("hCaloBragg3");
457 TH2F *th3 = dynamic_cast<TH2F*>(gDirectory->FindObject(thid3));
458 if ( th3 ) th3->Delete();
459 th3 = new TH2F(thid3,thid3,300,-0.5,300.,1000,0.,150.);
460 th3->SetMarkerStyle(20);
461 th3->SetMarkerColor(kBlue);
462
463
464 // tc->cd(1);
465 //
466 // for(Int_t i=0;i<=estremi[1][0];i++)th->Fill(Depth[i],dEpianimean[i]);
467 for(Int_t i=0;i<=estremi[1][0];i++)th2->Fill(Depth[i],calorimetro[i][1]*MIP);
468 // th->Draw();
469 th2->Draw("same");
470
471 // tc->cd(2);
472 tc->cd();
473 //
474 for(Int_t i=0;i<=estremi[1][0];i++){
475 th3->Fill(Depth[i],dEpianiloop[i]);
476 // printf(" i %i Depth %f depianiloop %f \n",i,Depth[i],dEpianiloop[i]);
477 }
478 th3->Draw();
479 th2->Draw("same");
480
481 tc->Modified();
482 tc->Update();
483
484 //
485 gStyle->SetLabelSize(0);
486 gStyle->SetNdivisions(1,"XY");
487 //
488 };
489
490
491
492 void CaloBragg::LoadParam(){
493
494 //
495 elem[0] = 1.00794; //H 1
496 elem[1] = 4.002602; //He 2
497 elem[2] = 6.941; //Li 3
498 elem[3] = 9.012182;//Be 4
499 elem[4] = 10.811; //B 5
500 elem[5] = 12.0107; //C 6
501 elem[6] = 14.00674;//N 7
502 elem[7] = 15.9994; //O 8
503 elem[8] = 18.9984032; //F 9
504 elem[9] = 20.1797; //Ne 10
505 elem[10] = 22.98977;//Na 11
506 elem[11] = 24.3050; //Mg 12
507 elem[12] = 26.9815; //Al 13
508 elem[13] = 28.0855; //Si 14
509 elem[14] = 30.974; //P 15
510 elem[15] = 32.066; //S 16
511 elem[16] = 35.4527; //Cl 17
512 elem[17] = 39.948; //Ar 18
513 elem[18] = 39.0983; //K 19
514 elem[19] = 40.078; //Ca 20
515 elem[20] = 44.95591;//Sc 21
516 elem[21] = 47.867; //Ti 22
517 elem[22] = 50.9415; //V 23
518 elem[23] = 51.9961; //Cr 24
519 elem[24] = 54.938049;//Mn 25
520 elem[25] = 55.845; //Fe 26
521 elem[26] = 58.9332; //Co 27
522 elem[27] = 58.6934; //Ni 28
523 elem[28] = 63.546; //Cu 29
524 elem[29] = 65.39; //Zn 30
525 elem[30] = 69.723; //Ga 31
526 elem[31] = 72.61; //Ge 32
527
528
529 //parametri calorimetro
530 NPLA = 22;
531 NCHA = 96;
532 nView = 2;
533
534 AA = 0.96;//mm larghezza strip
535 ADIST = 80.5;//mm distanza tra pad
536 PIANO = 8.59;//mm distanza
537
538 ySi = 0.38;//mm spessore silicio
539 yW = 2.66;//mm spessore tungsteno
540 rhoSi = 2.33;//g/cm3 densita' silicio
541 rhoW = 19.3;//g/cm3 densita' tugsteno
542 MIP = 0.106;//Mev g/cm2 energia al minimo nel silicio per 0.38 mm
543
544 emin = 0.;
545
546 //parametri bethe-bloch
547 pigr = 3.1415;
548 Na = 6.02e-23;
549 ZA = 0.49; /*Z/A per Si*/
550 // ISi =182e-06; /*MeV*/
551 ISi = 171e-06; /*MeV*/
552 IW = 735e-06; /*MeV*/
553 // ISi =0.0001059994; /*GeV!!*/ no era giusto!!
554 Me = 0.511; /* MeV*/
555 MassP = 931.27;/*MeV*/
556 r2 = 7.95e-26; /*ro*ro in cm */
557
558 };
559
560
561
562 //
563 void CaloBragg::conversione(){
564
565 // calcolo spessore Si attraverato in funzione dell'inclinazione
566 // e conversione dello spessore di W in Si e correzione del valore
567 // della Mip pe lo spessore effettivo
568 //
569 // in : evento
570 //
571 // out: out[0] = gcm2Si = spessore silicio attraversato nel piano
572 // out[1] = WinSi = spessore equivalente in Si del W attraversato
573 // out[2] = Mip = fattore conversione energia riscalato allo spessore attrversatonel piano
574
575 Float_t SiCross=0.;
576 Float_t WCross = 0.;
577 Float_t ytgx = 0;
578 Float_t ytgy = 0;
579 Float_t a = 0.;
580
581 /*silicio*/
582 ytgx = ySi * L2->GetCaloLevel2()->tanx[0];
583 ytgy = ySi * L2->GetCaloLevel2()->tany[0];
584
585 //lunghezza effettiva di silicio attraversata (mm)
586 SiCross = sqrt(SQ(ySi) + SQ(ytgx) + SQ(ytgy));
587
588 spessore[0] = (SiCross/10.) * rhoSi; //spessore silicio in g/cm2
589
590 /*tungsteno*/
591 ytgx = yW * L2->GetCaloLevel2()->tanx[0];
592 ytgy = yW * L2->GetCaloLevel2()->tany[0];
593
594 //rapporto tra rilasci energetici nei due materiali
595 WCross = sqrt((yW*yW) + (ytgx*ytgx) + (ytgy*ytgy));//mm* rapporto lunghezze rad
596 //gcm2W = WCross/10. * rhoW;
597
598 // (g/cm2W)/(g/cm2Si)
599 spessore[3] = (WCross/10.) * rhoW;
600 a=(WCross/SiCross)*(rhoW/rhoSi)*(1.145/1.664); //(gcm2W)/(SiCross/10. * rhoSi)* (1.145/1.664);
601 spessore[1] = a;
602 //riscala mip allo spessore attraversato
603 spessore[2] = MIP*(SiCross/ySi);
604 };//end conversione
605
606
607
608
609
610 void CaloBragg::BetheBloch(Float_t *x, Float_t *z, Float_t *Mass, Float_t *gam, Float_t *Bet, Float_t *out, Float_t II){
611
612 //rilascio energetico con bethe bloch con correzioni
613 //in: x: g/cm2
614 // z: carica
615 // Mass: Massa uma
616 // Ene: energia particella MeV//tolta
617 // gam: (etot/massa)
618 // Bet: rad((g2-1)/g2)
619 //
620 //out: energia rilasciata MeV
621
622
623 Float_t eta =0.;
624 Float_t Wmax =0.;
625 Float_t lg =0.;
626 Float_t Energia=0.;
627 Float_t C=0.;
628 Float_t INo = ISi;
629
630 if ( usenewBB ) INo = II;
631
632 eta = (*gam)*(*Bet);
633
634 //Bet=3/gam; SQ(*gam) * SQ(*Bet)
635 Wmax = 2.* Me * SQ(eta) / (1. + 2.*(*gam)*Me/(*Mass) + SQ(Me)/SQ(*Mass));
636
637 lg = 2.* Me * SQ(eta) * Wmax / SQ(INo);
638 // Energia = x* 2 * pigr * Na * r2 * Me * rhoSi *ZA* SQ(z)/SQ(Bet) * lg;
639 C=(0.42237*pow(eta,-2.) + 0.0304*pow(eta,-4.) - 0.00038*pow(eta,-6.))*pow(10.,-6.)* pow(INo,2.) +
640 (3.858*pow(eta,-2.) - 0.1668*pow(eta,-4.) + 0.00158*pow(eta,-6.))*pow(10.,-9.)*pow(INo,3.);
641
642 if(eta <= 0.13) C= C * log(eta/0.0653)/log(0.13/0.0653);
643
644 Energia = (*x) * 0.307/28.09 * 14. *SQ(*z)/SQ(*Bet)*(0.5*log(lg) - SQ(*Bet) - C/14.);
645
646 *out =Energia;//out
647
648 };//end Bethebloch
649
650
651
652
653 void CaloBragg::ELOSS(Float_t *dx, Int_t *Z, Float_t *Etot, Float_t *out, Float_t II){
654
655 /*perdita di energia per ioni pesanti (come da routine geant)*/
656 // in : dx => spessore g/cm2
657 // Z => carica
658 // Etot => energia perticella
659 //
660 // out: energia persa
661
662
663 Float_t Q=0.;
664 Float_t v=0.;
665 Float_t gam=0.;
666 Float_t Bet=0.;
667 Float_t dEP=0.;
668
669 // gamma // Mass = A * MassP; /*in Mev/c2*/
670 gam = (*Etot)/(elem[*Z-1]*MassP); // E = gamma M c2
671
672
673 Bet = sqrt((SQ(gam) -1.)/SQ(gam));
674
675 // v= 121.4139*(Bet/pow((*Z),(2./3.))) + 0.0378*sin(190.7165*(Bet/pow((*Z),(2./3.))));
676 v= 121.4139*(Bet*pow((*Z),(2./3.))) + 0.0378*sin(190.7165*(Bet*pow((*Z),(2./3.)))); // EMI AAAAGGH!!
677
678 //carica effettiva
679 Q= (*Z)*(1- (1.034 - 0.1777*exp(-0.08114*(*Z)))*exp(-v));
680
681 //perdita energia per un protone
682 Float_t protone =1.;
683 // Float_t Mass=(elem[*Z-1]*MassP); //EMI
684 // BetheBloch(dx, &protone, &Mass, &gam, &Bet, &dEP);//ene non serve..go gamma.. BetheBloch(dx, 1, MassP, Etot/A, gam, Bet, &dEP);
685
686 BetheBloch(dx, &protone, &MassP, &gam, &Bet, &dEP, II);//ene non serve..go gamma.. BetheBloch(dx, 1, MassP, Etot/A, gam, Bet, &dEP); //EMI
687
688 *out= (SQ(Q)*(dEP));//*dx;
689
690
691 };//end ELOSS
692
693
694
695
696 void CaloBragg::Enetrack(Int_t* Z, Float_t* E0, Float_t* primo,Float_t* ultimo, Float_t out[]){
697
698 //calcola energia rilasciata sulla traccia (usa ELOSS)
699 // in : Z =>carica
700 // E0 =>energia
701 // spess2[3] => conversione spessore Si, Si in W, mip
702 // primo => posizione primo piano attraversato
703 //
704 // out: array[44] =>rilasci energetici calcolati per ogni piano[44] dopo il primo(estremi[0][0])
705
706
707
708 Float_t dE=0.; //energia rilasciata
709 Float_t Ezero= *E0;//energia iniziale
710
711 //azzero energia rilasciata sui piani
712 memset(out, 0, 2*NPLA*sizeof(Float_t));
713
714 Float_t Massa = (elem[(*Z)-1] * MassP);
715
716 for( Int_t ipla=((int)(*primo)); ipla<= ((int)(*ultimo)); ipla++){
717 dE=0.;
718 //spessore silicio corretto x inclinazione, z, energia, out:rilascio
719 ELOSS(&spessore[0], Z, &Ezero, &dE, ISi);//spessore in g/cm2!!
720 if((Ezero-dE) <= Massa){//se l'energia depositata e' maggiore dell'energia della perticella stop
721 out[ipla] = Ezero - Massa; //MeV
722 return;
723
724 }else{
725 out[ipla] = dE; //MeV
726 Ezero = Ezero - dE;//energia residua
727 if ( debug ) printf(" zompa %i out %f dE %f ezero %f \n",ipla,out[ipla],dE,Ezero);
728 };
729 //se sono su un piano Y (tutti i pari) dopo c'e' il tungsteno
730 if(ipla%2 == 0){
731 /*tungsteno*/
732 dE=0.;
733 Float_t sp = 0.;
734 Float_t II = ISi;
735 if ( usenewBB ){
736 sp = spessore[3];
737 II = IW;
738 } else {
739 sp = spessore[0]*spessore[1]; //((gcm2Si)*(WinSi))//spessore attraversato in g/cm2
740 }
741 // printf(" sp %f II %f \n",sp,II);
742 ELOSS(&sp, Z, &Ezero, &dE,II);
743 if((Ezero-dE) <= Massa){//se l'energia depositata e' maggiore dell'energia della perticella stop
744 return;
745 }else{
746 Ezero = Ezero -dE;//energia residua
747 };
748 };
749
750 };//fine loop piani
751
752
753 };//end Enetrack
754
755
756
757 void CaloBragg::chiquadro(Float_t dE[], Float_t out[]){
758
759 // calcola chi2 tra energia calcolata e misurata
760 // in : dE[44] =>energia calcolata
761 // calo3[44][2]=> [0]strip attraversata [1]energia misurata per ogni piano
762 // estr2 => array con primo[0][0] e ultimo[1][0] piano attraversati ed energie[][1]
763 //
764 // out: array[3]=> (chi2; piani scartati consecutivi(79= >3 quindi frammentato); piani scartati totale)
765
766
767 Float_t sum = 0.;
768 Float_t PianoPrecedente=0.;
769 Float_t badplane=0.;
770 Float_t badplanetot=0.;
771 Float_t w,wi;
772 //
773 if ( newchi2 ){
774 ndf = 0;
775 sum = 0.;
776 for( Int_t ipla=((int)(estremi[0][0])); ipla<= ((int)(estremi[1][0])); ipla++){
777 sum += pow((dE[ipla] - (calorimetro[ipla][1] * spessore[2]))/(0.05*dE[ipla]),2.);
778 // printf(" quiqui: dE %f calor %f spessore[2] %f \n",dE[ipla],spessore[2]*calorimetro[ipla][1],spessore[2]);
779 ndf++;
780 }
781 ndf -= 2;
782 if ( ndf > 0 ) sum /= (float)ndf;
783 out[0] = sum;
784 out[1] = 0.;
785 out[2] = (int)(estremi[1][0])-ndf;
786 // printf(" sum %f ndf %i \n ",sum,ndf);
787 } else {
788 for(Int_t ipla=0; ipla<2*NPLA; ipla++){
789 //tutti i piani attraversati dalla traiettoria
790 if(calorimetro[ipla][0] != -1.){ //
791 w=0.; //normalizzazione;
792 wi=1.;//peso
793
794 //tolgo piani attraversati dalla traccia ma precedenti il piano individuato come ingresso
795 if (ipla<estremi[0][0]) wi=0.;
796
797 //tolgo piani attraversati da traccia ma successivi all'ultimo se sono diversi da 0
798 //if((ipla>estremi[1][0]) && (calorimetro[ipla][1] >0.) ) wi=0.;
799 if((ipla>estremi[1][0])) wi=0.;
800
801 //normalizzazione
802 if (calorimetro[ipla][1] != 0.) w=1./(calorimetro[ipla][1]* MIP); //
803
804 //tolgo piani con rilasci inferiori al 30% del precedente
805 if(calorimetro[ipla][1] < (0.7*PianoPrecedente)){ // cosi' i piani senza rilascio non vengono considerati nel calcolo del chi2
806 wi=0.;
807 //se sono piani intermedi (non si e' fermta) li considero non buoni
808 if( (ipla <= estremi[1][0]) && (calorimetro[ipla][1] !=0.)){//
809 badplane+=1.;
810 badplanetot+=1.;
811 };
812 };
813
814 //meno peso ai piani con rilasci maggiori di 1000 MIP
815 // if(calorimetro[ipla][1] > 1000) wi=0.5;
816 if(calorimetro[ipla][1] > 1200.) wi=0.5;
817
818 Float_t arg = w*wi*(dE[ipla] - (calorimetro[ipla][1] * MIP));
819
820 sum += SQ(arg); // w*wi*(dEpiani[p][v]-(eplane[p][v]*MIP))));//( dEpiani[p][v] - (eplane[p][v]*MIP));
821 if(debug){
822 printf("dedx calcolata %f e reale %f \n",dE[ipla],(calorimetro[ipla][1] * MIP));
823 }
824 //se trovo piano non buono (tolto quindi wi=0) non modifico il piano precedente
825 if(wi != 0.){//
826 PianoPrecedente= calorimetro[ipla][1];//tengo piano precedente
827 badplane = 0.;//azzero contatore piani scartati consecutivi
828 };
829 };
830
831 //da Emi
832 if(badplane > 2){
833 // printf(" AAAAAAAAAAAAAAAAAAAAAAAAGGGGGGGGGGGGGGGGGGGGGGGGGGGGGG\n");
834 out[1] =79.;
835 break;
836 };
837
838 };//fine loop piani
839 //chi2,frammentato,pskip
840 out[0]=sum;
841 out[2]=badplanetot;
842 }
843 };//end chiquadro
844
845
846
847 void CaloBragg::loopze( Float_t step, Float_t E0,Float_t Zstart, Float_t Zlimite, Int_t nostep = 1000){
848 //
849 //loop su z ed energie per trovare miglior z (ed energia)
850 //in: nloop => energia massima da provare (nloop x E0)
851 // E0 => energia iniziale (intergale)
852 // Zstart => minimo z da cui patire
853 // Zlimite => z a cui fermarsi (z al minimo di ionizz sul 1o piano)
854 //
855 //out: array[4]=> chi2,Zbest,Ebest,piani saltati nel chi2
856 //
857
858
859 memset(dEplan,0,2*NPLA*sizeof(Float_t));
860
861 Int_t Z = 0;// z iniziale
862
863 Float_t Massa = 0.;
864
865 Float_t Stepint =(step)/(Float_t)nostep;//passo per il calcolo di energia
866
867 Float_t energia =0.;//energia del loop
868
869 Float_t chi2[3] = {0,0,0};//out dal calcolo chi2: chi2, piani consecutivi saltati, piani totali saltati
870
871 Int_t zmin = (int)Zstart;
872 Int_t max=32;//max z di cui so la massa :P
873 if((Zlimite)<=31) max=(int)(Zlimite) + 1;
874
875 if ( fzeta > 0. ){
876 zmin = fzeta;
877 max = fzeta+1;
878 }
879
880 Int_t colmax=32;
881 Int_t rowmax=3000;
882
883 Float_t matrixchi2[colmax][rowmax][3];
884 memset(matrixchi2, 0, colmax*rowmax*3*sizeof(Float_t));
885
886 Int_t imin = 1-nostep/2;
887 Int_t imax = nostep/2;
888
889 //loop elementi
890 for(Int_t inucl=zmin; inucl<max; inucl++){
891
892 Z= inucl;
893
894 Massa = elem[inucl-1]*MassP;
895
896 //loop energia
897 Int_t iene2 = 0;
898 // for(Int_t iene= 0; iene<1000; iene++){// da non cambiare in base a Stepint altrimenti cambia la matrice bestchi2!!!cosi' non raggiungo mai integrale!!!!! mettere <=??
899 for(Int_t iene= imin; iene<imax; iene++){// da non cambiare in base a Stepint altrimenti cambia la matrice bestchi2!!!cosi' non raggiungo mai integrale!!!!! mettere <=??
900
901 iene2++;
902 energia= Massa + (E0)+ iene*Stepint;//gli do un'energia totale (momento) massa+energia cinetica, aumentando la cinetica..
903
904 Enetrack(&Z, &energia, &estremi[0][0],&estremi[1][0], dEplan);//calcola rilascio energetico sui piani
905 //calcolo chi2
906 chiquadro(dEplan,chi2);
907
908 // printf(" last deplan from: Z = %i iene %i energia %f chi2 %f \n",inucl,iene,energia,chi2[0]);
909
910 if( (chi2[1] != 79.) ){//salto quelli che frammentano
911 matrixchi2[inucl][iene2][0]=chi2[0];//valore chi2 per questo z a questa energia
912 matrixchi2[inucl][iene2][1]=energia;//energia per questo chi2
913 matrixchi2[inucl][iene2][2]=chi2[2];//piani saltati nel chi2
914 } else {
915 matrixchi2[inucl][iene2][0]=1000.;//valore chi2 per questo z a questa energia
916 matrixchi2[inucl][iene2][1]=1000.;//energia per questo chi2
917 matrixchi2[inucl][iene2][2]=1000.;//piani saltati nel chi2
918 break;
919 }
920 }//fine loop energia
921
922
923 };//fine loop z
924
925
926 //Emi
927 for (Int_t nu=zmin; nu<max; nu++){
928 for (Int_t en=0; en<nostep; en++){
929 if((matrixchi2[nu][en][0]<bestchi2[0]) && (matrixchi2[nu][en][0] >0.)){
930 bestchi2[0]= matrixchi2[nu][en][0];// chi2
931 bestchi2[1]= (Float_t)nu; // z
932 bestchi2[2]= matrixchi2[nu][en][1];//energia;
933 bestchi2[3]= matrixchi2[nu][en][2];// totale piani saltati
934 }
935 }
936 }
937
938
939 };//endloopze
940
941
942
943
944
945 // void CaloBragg::mediatroncata(){
946 // //calcolo Z con media troncata e utilizzo questo Z per trovare l'energia migliore
947 // //in: ordplane[44] => array con energia dei piani
948 // // spess[3] => conversioni spessore di silicio, w, mip
949 // // estr[2][2] => primo[0][0] e ultimo[1][0] piano attraversati ed energie[][1]
950 // // calo[44][2]=> energia[][1] e strip[][0] passaggio su ogni piano
951 // // integrale => energia totale nel calorimetro considerando il W
952 // //
953 // // out[4] chi2,z,Etot,Pskip
954
955 // Float_t ordplane[44];//mi serve per la media troncata
956 // memset(ordplane,0,44*sizeof(Float_t));
957
958 // for(Int_t ipla=0; ipla< 2*NPLA; ipla++) ordplane[ipla]=calorimetro[ipla][1]; //energia del piano
959
960
961 // //ordino tutte le energie dei piani in ordine crescente
962
963 // Long64_t work[200];
964 // Int_t ind = 0;
965 // //Int_t l = 0;
966 // Int_t RN = 0;
967 // Float_t sum4 = 0.;
968 // Float_t qm = 0.;
969 // //
970 // //Float_t qmt = ethr*0.8; // *0.9
971 // //
972 // //Int_t uplim = TMath::Max(3,N);
973 // //
974 // while ( RN < 4 && ind < 44 ){
975 // qm = TMath::KOrdStat(44,ordplane,ind,work);
976 // if (qm >= 0.7 ){
977 // if ( RN < 4 ){
978 // sum4 += qm;
979 // RN++;
980 // };
981 // // l++;
982 // // if ( debug ) printf(" value no %i qm %f sum4 %f \n",l,qm,sum4);
983 // };
984 // ind++;
985 // };
986 // //
987 // sum4 /= (Float_t)RN;
988 // Float_t Zmean = (sqrt((sum4*MIP)/(((Float_t)RN)*spessore[2])));//ma non e'/1??
989 // if(Zmean ==0.) Zmean=1.;
990 // if ( Zmean < 1. ) Zmean = 1.;
991
992
993 // // Zmean =round(Zmean);
994 // // if(Zmean <1.) Zmean=1.;
995
996 // // if(Zmean >0.)Zmean =round(Zmean);
997
998 // //======== per i nuclei=======
999 // if (Zmean >=2.){
1000 // ind = 0;
1001 // RN = 0;
1002 // sum4 = 0.;
1003 // qm = 0.;
1004 // while ( RN < 4 && ind < 44 ){
1005 // qm = TMath::KOrdStat(44,ordplane,ind,work);
1006 // if (qm >= (Zmean*Zmean)-Zmean*Zmean*0.2 ){
1007 // if ( RN < 4 ){
1008 // sum4 += qm;
1009 // RN++;
1010 // };
1011 // };
1012 // ind++;
1013 // };
1014 // //
1015 // sum4 /= (Float_t)RN;
1016 // Zmean = (sqrt((sum4*MIP)/(4.*spessore[2])));//ma non e' /1??
1017 // }
1018
1019
1020 // //calcolo energia migliore per Z trovato con media troncata
1021 // // Float_t zmin=Zmean;
1022 // Float_t zmin=round(Zmean);
1023
1024 // bestchi2[0]=10000.;
1025 // bestchi2[1]=0.;
1026 // bestchi2[2]=0.;
1027 // bestchi2[3]=0.;
1028 // Float_t zero=0.;
1029
1030 // // step energia zstart zstop
1031 // loopze(Integrale,zero,zmin,zmin);
1032
1033
1034 // qtchi2=bestchi2[0];
1035 // qtz=bestchi2[1];
1036 // qtetot=bestchi2[2];
1037 // qtpskip=bestchi2[3];
1038 // };//end mediatroncata
1039
1040
1041
1042 void CaloBragg::Zdaloop(){
1043 //calcolo Z con un loop su tutti i possibli Z ed energie
1044 //in: ordplane[44]=> array con energia dei piani
1045 // spess1[3]=> conversioni spessore di silicio, w e mip
1046 // estr3[2][2]=> primo[0][0] e ultimo[1][0] piano ed energie
1047 // calo1[44][2]=> energia[][1] e strip[][0] passaggio su ogni piano
1048 // integrale=> energia totale nel calorimetro considerando il W
1049 //
1050 // out[4] chi2,z,Etot,Pskip
1051
1052
1053 /*z se particella fosse al minimo*/ //energia1piano/mip corretta
1054 // Float_t zmax = round(sqrt(estremi[0][1]/spessore[2]));
1055 // if(zmax<31)zmax=zmax+1;
1056
1057 /*calcolo Z ed E con loop sui vari elementi ed energie*/
1058
1059 Float_t zmin=1.;
1060 Float_t zmax=32.;
1061 Float_t bestchitemp[4] = {0,0,0,0};
1062
1063 bestchi2[0]=10000.;
1064 bestchi2[1]=0.;
1065 bestchi2[2]=0.;
1066 bestchi2[3]=0.;
1067 Float_t zero=0.;
1068 //------------primo loop ----------------------
1069 // energia ezero, zstart zstop
1070 // loopze(Integrale,zero,zmin,zmax);
1071
1072 //-> loopze(Integrale*1.2/500.,Integrale/1000.,zmin,zmax,50);
1073 loopze(Integrale*1.2/500.,Integrale/1000.,zmin,zmax,200);
1074
1075 // loopze(Integrale*2.,Integrale/100.,zmin,zmax);
1076 if ( debug ) printf(" Integrale %f , outene %f \n",Integrale,bestchi2[2]);
1077
1078 //------------secondo loop ----------------------
1079 for(Int_t i=0;i<4;i++) bestchitemp[i]=bestchi2[i];
1080 bestchi2[0] = 10000.;
1081 bestchi2[1] = 0.;
1082 bestchi2[2] = 0.;
1083 bestchi2[3] = 0.;//riazzero
1084
1085 Float_t step = bestchitemp[2];//
1086 zero=0.; // qualsiasi altro valore peggiora le cose
1087 // zmin=zmax=bestchitemp[1];
1088 zmin=bestchitemp[1]-1;
1089 zmax=bestchitemp[1]+1;
1090 // loopze(step,zero,zmin,zmax); //
1091
1092 //-> loopze(step,step/2.,zmin,zmax,200); //
1093 loopze(step,step/2.,zmin,zmax,500); //
1094
1095 if ( debug ) printf(" Integrale2 %f , outene %f step %f \n",Integrale,bestchi2[2],step);
1096
1097
1098 //chi2,z,Etot,Pskip
1099 lpchi2=bestchi2[0];
1100 lpz=bestchi2[1];
1101 lpetot=bestchi2[2];
1102 lppskip=bestchi2[3];
1103
1104 };//endZdaloop
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116

  ViewVC Help
Powered by ViewVC 1.1.23