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