/[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.2 - (show annotations) (download)
Wed Jun 18 15:24:12 2008 UTC (16 years, 6 months ago) by pamelats
Branch: MAIN
Changes since 1.1: +128 -123 lines
primo debug fatto

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

  ViewVC Help
Powered by ViewVC 1.1.23