1 | //
|
---|
2 | // $Id: histos.cc,v 1.7 1999-12-10 16:55:51 ansari Exp $
|
---|
3 | //
|
---|
4 |
|
---|
5 | #include "machdefs.h"
|
---|
6 | #include <string.h>
|
---|
7 | #include <stdio.h>
|
---|
8 | #include <math.h>
|
---|
9 | #include "histos.h"
|
---|
10 | #include "perrors.h"
|
---|
11 | #include "poly.h"
|
---|
12 | #include "strutil.h"
|
---|
13 | #include "generalfit.h"
|
---|
14 |
|
---|
15 | //++
|
---|
16 | // Class Histo
|
---|
17 | // Lib Outils++
|
---|
18 | // include histos.h
|
---|
19 | //
|
---|
20 | // Classe d'histogrammes 1D
|
---|
21 | //--
|
---|
22 |
|
---|
23 | //++
|
---|
24 | // Titre Constructeurs
|
---|
25 | //--
|
---|
26 |
|
---|
27 | /********* Methode *********/
|
---|
28 | //++
|
---|
29 | Histo::Histo()
|
---|
30 | //
|
---|
31 | //--
|
---|
32 | : data(NULL), err2(NULL),
|
---|
33 | under(0), over(0), nHist(0), nEntries(0),
|
---|
34 | bins(0), min(0), max(0),
|
---|
35 | binWidth(0)
|
---|
36 | {
|
---|
37 | END_CONSTRUCTOR
|
---|
38 | }
|
---|
39 |
|
---|
40 | /********* Methode *********/
|
---|
41 | //++
|
---|
42 | Histo::Histo(float xMin, float xMax, int nBin)
|
---|
43 | //
|
---|
44 | //--
|
---|
45 | : data(new float[nBin]), err2(NULL),
|
---|
46 | under(0), over(0), nHist(0), nEntries(0),
|
---|
47 | bins(nBin), min(xMin), max(xMax),
|
---|
48 | binWidth((max - min)/nBin)
|
---|
49 | {
|
---|
50 | Zero();
|
---|
51 | END_CONSTRUCTOR
|
---|
52 | }
|
---|
53 |
|
---|
54 | /********* Methode *********/
|
---|
55 | //++
|
---|
56 | Histo::Histo(const Histo& H)
|
---|
57 | //
|
---|
58 | //--
|
---|
59 | : data(new float[H.bins]), err2(NULL),
|
---|
60 | under(H.under), over(H.over), nHist(H.nHist), nEntries(H.nEntries),
|
---|
61 | bins(H.bins), min(H.min), max(H.max),
|
---|
62 | binWidth(H.binWidth)
|
---|
63 | {
|
---|
64 | memcpy(data, H.data, bins*sizeof(float));
|
---|
65 | if( H.err2 != NULL ) {
|
---|
66 | err2 = new double[bins];
|
---|
67 | memcpy(err2, H.err2, bins*sizeof(double));
|
---|
68 | }
|
---|
69 | END_CONSTRUCTOR
|
---|
70 | }
|
---|
71 |
|
---|
72 | /********* Methode *********/
|
---|
73 | void Histo::Delete()
|
---|
74 | {
|
---|
75 | if( data != NULL ) { delete[] data; data = NULL;}
|
---|
76 | if( err2 != NULL ) { delete[] err2; err2 = NULL;}
|
---|
77 | under = over = min = max = binWidth= 0.;
|
---|
78 | nHist = 0.;
|
---|
79 | nEntries = bins = 0;
|
---|
80 | }
|
---|
81 |
|
---|
82 | /********* Methode *********/
|
---|
83 | Histo::~Histo()
|
---|
84 | {
|
---|
85 | Delete();
|
---|
86 | }
|
---|
87 |
|
---|
88 | //++
|
---|
89 | // Titre Methodes
|
---|
90 | //--
|
---|
91 |
|
---|
92 | /********* Methode *********/
|
---|
93 | //++
|
---|
94 | void Histo::Zero()
|
---|
95 | //
|
---|
96 | // Remise a zero
|
---|
97 | //--
|
---|
98 | {
|
---|
99 | memset(data, 0, bins*sizeof(float));
|
---|
100 | under = over = 0;
|
---|
101 | nHist = 0;
|
---|
102 | nEntries = 0;
|
---|
103 | if( err2 != NULL ) memset(err2, 0, bins*sizeof(double));
|
---|
104 | }
|
---|
105 |
|
---|
106 | /********* Methode *********/
|
---|
107 | //++
|
---|
108 | void Histo::Errors()
|
---|
109 | //
|
---|
110 | // Pour avoir le calcul des erreurs
|
---|
111 | //--
|
---|
112 | {
|
---|
113 | if( bins > 0 ) {
|
---|
114 | if(err2==NULL) err2 = new double[bins];
|
---|
115 | memset(err2, 0, bins*sizeof(double));
|
---|
116 | }
|
---|
117 | }
|
---|
118 |
|
---|
119 | /********* Methode *********/
|
---|
120 | //++
|
---|
121 | Histo& Histo::operator = (const Histo& h)
|
---|
122 | //
|
---|
123 | //--
|
---|
124 | {
|
---|
125 | if(this == &h) return *this;
|
---|
126 | if( h.bins > bins ) Delete();
|
---|
127 | if(!data) data = new float[h.bins];
|
---|
128 | if( !h.err2 && err2 ) { delete [] err2; err2=NULL;}
|
---|
129 | if( h.err2 && !err2 ) err2 = new double[h.bins];
|
---|
130 |
|
---|
131 | under = h.under;
|
---|
132 | over = h.over;
|
---|
133 | nHist = h.nHist;
|
---|
134 | nEntries = h.nEntries;
|
---|
135 | bins = h.bins;
|
---|
136 | min = h.min;
|
---|
137 | max = h.max;
|
---|
138 | binWidth = h.binWidth;
|
---|
139 |
|
---|
140 | memcpy(data, h.data, bins*sizeof(float));
|
---|
141 | if(err2) memcpy(err2, h.err2, bins*sizeof(double));
|
---|
142 |
|
---|
143 | return *this;
|
---|
144 | }
|
---|
145 |
|
---|
146 | /********* Methode *********/
|
---|
147 | //++
|
---|
148 | Histo& Histo::operator *= (double b)
|
---|
149 | //
|
---|
150 | //--
|
---|
151 | {
|
---|
152 | double b2 = b*b;
|
---|
153 | for(int i=0;i<bins;i++) {
|
---|
154 | data[i] *= b;
|
---|
155 | if(err2) err2[i] *= b2;
|
---|
156 | }
|
---|
157 | under *= b;
|
---|
158 | over *= b;
|
---|
159 | nHist *= b;
|
---|
160 |
|
---|
161 | return *this;
|
---|
162 | }
|
---|
163 |
|
---|
164 | //++
|
---|
165 | Histo& Histo::operator /= (double b)
|
---|
166 | //
|
---|
167 | //--
|
---|
168 | {
|
---|
169 | if (b==0.) THROW(inconsistentErr);
|
---|
170 | double b2 = b*b;
|
---|
171 | for(int i=0;i<bins;i++) {
|
---|
172 | data[i] /= b;
|
---|
173 | if(err2) err2[i] /= b2;
|
---|
174 | }
|
---|
175 | under /= b;
|
---|
176 | over /= b;
|
---|
177 | nHist /= b;
|
---|
178 |
|
---|
179 | return *this;
|
---|
180 | }
|
---|
181 |
|
---|
182 | //++
|
---|
183 | Histo& Histo::operator += (double b)
|
---|
184 | //
|
---|
185 | //--
|
---|
186 | {
|
---|
187 | for(int i=0;i<bins;i++) data[i] += b;
|
---|
188 | under += b;
|
---|
189 | over += b;
|
---|
190 | nHist += bins*b;
|
---|
191 |
|
---|
192 | return *this;
|
---|
193 | }
|
---|
194 |
|
---|
195 | //++
|
---|
196 | Histo& Histo::operator -= (double b)
|
---|
197 | //
|
---|
198 | //--
|
---|
199 | {
|
---|
200 | for(int i=0;i<bins;i++) data[i] -= b;
|
---|
201 | under -= b;
|
---|
202 | over -= b;
|
---|
203 | nHist -= bins*b;
|
---|
204 |
|
---|
205 | return *this;
|
---|
206 | }
|
---|
207 |
|
---|
208 | /********* Methode *********/
|
---|
209 | //++
|
---|
210 | Histo operator * (const Histo& a, double b)
|
---|
211 | //
|
---|
212 | //--
|
---|
213 | {
|
---|
214 | Histo result(a);
|
---|
215 | return (result *= b);
|
---|
216 | }
|
---|
217 |
|
---|
218 | //++
|
---|
219 | Histo operator * (double b, const Histo& a)
|
---|
220 | //
|
---|
221 | //--
|
---|
222 | {
|
---|
223 | Histo result(a);
|
---|
224 | return (result *= b);
|
---|
225 | }
|
---|
226 |
|
---|
227 | //++
|
---|
228 | Histo operator / (const Histo& a, double b)
|
---|
229 | //
|
---|
230 | //--
|
---|
231 | {
|
---|
232 | Histo result(a);
|
---|
233 | return (result /= b);
|
---|
234 | }
|
---|
235 |
|
---|
236 | //++
|
---|
237 | Histo operator + (const Histo& a, double b)
|
---|
238 | //
|
---|
239 | //--
|
---|
240 | {
|
---|
241 | Histo result(a);
|
---|
242 | return (result += b);
|
---|
243 | }
|
---|
244 |
|
---|
245 | //++
|
---|
246 | Histo operator + (double b, const Histo& a)
|
---|
247 | //
|
---|
248 | //--
|
---|
249 | {
|
---|
250 | Histo result(a);
|
---|
251 | return (result += b);
|
---|
252 | }
|
---|
253 |
|
---|
254 | //++
|
---|
255 | Histo operator - (const Histo& a, double b)
|
---|
256 | //
|
---|
257 | //--
|
---|
258 | {
|
---|
259 | Histo result(a);
|
---|
260 | return (result -= b);
|
---|
261 | }
|
---|
262 |
|
---|
263 | //++
|
---|
264 | Histo operator - (double b, const Histo& a)
|
---|
265 | //
|
---|
266 | //--
|
---|
267 | {
|
---|
268 | Histo result(a);
|
---|
269 | result *= -1.;
|
---|
270 | return (result += b);
|
---|
271 | }
|
---|
272 |
|
---|
273 | /********* Methode *********/
|
---|
274 | //++
|
---|
275 | Histo& Histo::operator += (const Histo& a)
|
---|
276 | //
|
---|
277 | //--
|
---|
278 | {
|
---|
279 | if(bins!=a.bins) THROW(sizeMismatchErr);
|
---|
280 | for(int i=0;i<bins;i++) {
|
---|
281 | data[i] += a(i);
|
---|
282 | if(err2 && a.err2) err2[i] += a.Error2(i);
|
---|
283 | else err2[i] = 0.;
|
---|
284 | }
|
---|
285 | under += a.under;
|
---|
286 | over += a.over;
|
---|
287 | nHist += a.nHist;
|
---|
288 | nEntries += a.nEntries;
|
---|
289 |
|
---|
290 | return *this;
|
---|
291 | }
|
---|
292 |
|
---|
293 | //++
|
---|
294 | Histo& Histo::operator -= (const Histo& a)
|
---|
295 | //
|
---|
296 | //--
|
---|
297 | {
|
---|
298 | if(bins!=a.bins) THROW(sizeMismatchErr);
|
---|
299 | for(int i=0;i<bins;i++) {
|
---|
300 | data[i] -= a(i);
|
---|
301 | if(err2 && a.err2) err2[i] += a.Error2(i);
|
---|
302 | else err2[i] = 0.;
|
---|
303 | }
|
---|
304 | under -= a.under;
|
---|
305 | over -= a.over;
|
---|
306 | nHist -= a.nHist;
|
---|
307 | nEntries += a.nEntries;
|
---|
308 |
|
---|
309 | return *this;
|
---|
310 | }
|
---|
311 |
|
---|
312 | //++
|
---|
313 | Histo& Histo::operator *= (const Histo& a)
|
---|
314 | //
|
---|
315 | //--
|
---|
316 | {
|
---|
317 | if(bins!=a.bins) THROW(sizeMismatchErr);
|
---|
318 | nHist = 0.;
|
---|
319 | for(int i=0;i<bins;i++) {
|
---|
320 | if(err2 && a.err2)
|
---|
321 | err2[i] = a(i)*a(i)*err2[i] + data[i]*data[i]*a.Error2(i);
|
---|
322 | else err2[i] = 0.;
|
---|
323 | data[i] *= a(i);
|
---|
324 | nHist += data[i];
|
---|
325 | }
|
---|
326 | under *= a.under;
|
---|
327 | over *= a.over;
|
---|
328 | nEntries += a.nEntries;
|
---|
329 |
|
---|
330 | return *this;
|
---|
331 | }
|
---|
332 |
|
---|
333 | //++
|
---|
334 | Histo& Histo::operator /= (const Histo& a)
|
---|
335 | //
|
---|
336 | //--
|
---|
337 | {
|
---|
338 | if(bins!=a.bins) THROW(sizeMismatchErr);
|
---|
339 | nHist = 0.;
|
---|
340 | for(int i=0;i<bins;i++) {
|
---|
341 | if(a(i)==0.) {
|
---|
342 | data[i]=0.;
|
---|
343 | if(err2) err2[i]=0.;
|
---|
344 | continue;
|
---|
345 | }
|
---|
346 | if(err2 && a.err2)
|
---|
347 | err2[i] = (err2[i] + data[i]/a(i)*data[i]/a(i)*a.Error2(i))
|
---|
348 | /(a(i)*a(i));
|
---|
349 | else err2[i] = 0.;
|
---|
350 | data[i] /= a(i);
|
---|
351 | nHist += data[i];
|
---|
352 | }
|
---|
353 | if(a.under!=0.) under /= a.under; else under = 0.;
|
---|
354 | if(a.over!=0.) over /= a.over; else over = 0.;
|
---|
355 | nEntries += a.nEntries;
|
---|
356 |
|
---|
357 | return *this;
|
---|
358 | }
|
---|
359 |
|
---|
360 | /********* Methode *********/
|
---|
361 | //++
|
---|
362 | Histo operator + (const Histo& a, const Histo& b)
|
---|
363 | //
|
---|
364 | //--
|
---|
365 | {
|
---|
366 | if (b.NBins()!=a.NBins()) THROW(sizeMismatchErr);
|
---|
367 | Histo c(a);
|
---|
368 | return (c += b);
|
---|
369 | }
|
---|
370 |
|
---|
371 | //++
|
---|
372 | Histo operator - (const Histo& a, const Histo& b)
|
---|
373 | //
|
---|
374 | //--
|
---|
375 | {
|
---|
376 | if (b.NBins()!=a.NBins()) THROW(sizeMismatchErr);
|
---|
377 | Histo c(a);
|
---|
378 | return (c -= b);
|
---|
379 | }
|
---|
380 |
|
---|
381 | //++
|
---|
382 | Histo operator * (const Histo& a, const Histo& b)
|
---|
383 | //
|
---|
384 | //--
|
---|
385 | {
|
---|
386 | if (b.NBins()!=a.NBins()) THROW(sizeMismatchErr);
|
---|
387 | Histo c(a);
|
---|
388 | return (c *= b);
|
---|
389 | }
|
---|
390 |
|
---|
391 | //++
|
---|
392 | Histo operator / (const Histo& a, const Histo& b)
|
---|
393 | //
|
---|
394 | //--
|
---|
395 | {
|
---|
396 | if (b.NBins()!=a.NBins()) THROW(sizeMismatchErr);
|
---|
397 | Histo c(a);
|
---|
398 | return (c /= b);
|
---|
399 | }
|
---|
400 |
|
---|
401 | /********* Methode *********/
|
---|
402 | //++
|
---|
403 | void Histo::GetAbsc(Vector &v)
|
---|
404 | //
|
---|
405 | // Remplissage d'un tableau avec la valeur des abscisses
|
---|
406 | //--
|
---|
407 | {
|
---|
408 | v.Realloc(bins);
|
---|
409 | for(int i=0;i<bins;i++) v(i) = BinLowEdge(i);
|
---|
410 | return;
|
---|
411 | }
|
---|
412 |
|
---|
413 | //++
|
---|
414 | void Histo::GetValue(Vector &v)
|
---|
415 | //
|
---|
416 | // Remplissage d'un tableau avec la valeur du contenu
|
---|
417 | //--
|
---|
418 | {
|
---|
419 | v.Realloc(bins);
|
---|
420 | for(int i=0;i<bins;i++) v(i) = data[i];
|
---|
421 | return;
|
---|
422 | }
|
---|
423 |
|
---|
424 | //++
|
---|
425 | void Histo::GetError2(Vector &v)
|
---|
426 | //
|
---|
427 | // Remplissage d'un tableau avec la valeur des erreurs au carre
|
---|
428 | //--
|
---|
429 | {
|
---|
430 | v.Realloc(bins);
|
---|
431 | if(!err2) {for(int i=0;i<bins;i++) v(i) = 0.; return;}
|
---|
432 | for(int i=0;i<bins;i++) v(i) = err2[i];
|
---|
433 | return;
|
---|
434 | }
|
---|
435 |
|
---|
436 | //++
|
---|
437 | void Histo::GetError(Vector &v)
|
---|
438 | //
|
---|
439 | // Remplissage d'un tableau avec la valeur des erreurs
|
---|
440 | //--
|
---|
441 | {
|
---|
442 | v.Realloc(bins);
|
---|
443 | if(!err2) {for(int i=0;i<bins;i++) v(i) = 0.; return;}
|
---|
444 | for(int i=0;i<bins;i++) v(i) = Error(i);
|
---|
445 | return;
|
---|
446 | }
|
---|
447 |
|
---|
448 | /********* Methode *********/
|
---|
449 | //++
|
---|
450 | void Histo::PutValue(Vector &v, int ierr)
|
---|
451 | //
|
---|
452 | // Remplissage du contenu de l'histo avec les valeurs d'un vecteur
|
---|
453 | //--
|
---|
454 | {
|
---|
455 | if(v.NElts()<bins) THROW(sizeMismatchErr);
|
---|
456 | for(int i=0;i<bins;i++) {
|
---|
457 | data[i] = v(i);
|
---|
458 | if(err2&&ierr) err2[i] = fabs(v(i));
|
---|
459 | }
|
---|
460 | return;
|
---|
461 | }
|
---|
462 |
|
---|
463 | //++
|
---|
464 | void Histo::PutValueAdd(Vector &v, int ierr)
|
---|
465 | //
|
---|
466 | // Addition du contenu de l'histo avec les valeurs d'un vecteur
|
---|
467 | //--
|
---|
468 | {
|
---|
469 | if(v.NElts()<bins) THROW(sizeMismatchErr);
|
---|
470 | for(int i=0;i<bins;i++) {
|
---|
471 | data[i] += v(i);
|
---|
472 | if(err2&&ierr) err2[i] += fabs(v(i));
|
---|
473 | }
|
---|
474 | return;
|
---|
475 | }
|
---|
476 |
|
---|
477 | //++
|
---|
478 | void Histo::PutError2(Vector &v)
|
---|
479 | //
|
---|
480 | // Remplissage des erreurs au carre de l'histo avec les valeurs d'un vecteur
|
---|
481 | //--
|
---|
482 | {
|
---|
483 | if(v.NElts()<bins) THROW(sizeMismatchErr);
|
---|
484 | if(!err2) Errors();
|
---|
485 | for(int i=0;i<bins;i++) err2[i] = v(i);
|
---|
486 | return;
|
---|
487 | }
|
---|
488 |
|
---|
489 | //++
|
---|
490 | void Histo::PutError2Add(Vector &v)
|
---|
491 | //
|
---|
492 | // Addition des erreurs au carre de l'histo avec les valeurs d'un vecteur
|
---|
493 | //--
|
---|
494 | {
|
---|
495 | if(v.NElts()<bins) THROW(sizeMismatchErr);
|
---|
496 | if(!err2) Errors();
|
---|
497 | for(int i=0;i<bins;i++) if(v(i)>0.) err2[i] += v(i);
|
---|
498 | return;
|
---|
499 | }
|
---|
500 |
|
---|
501 | //++
|
---|
502 | void Histo::PutError(Vector &v)
|
---|
503 | //
|
---|
504 | // Remplissage des erreurs de l'histo avec les valeurs d'un vecteur
|
---|
505 | //--
|
---|
506 | {
|
---|
507 | if(v.NElts()<bins) THROW(sizeMismatchErr);
|
---|
508 | if(!err2) Errors();
|
---|
509 | for(int i=0;i<bins;i++)
|
---|
510 | if(v(i)>0.) err2[i]=v(i)*v(i); else err2[i]=-v(i)*v(i);
|
---|
511 | return;
|
---|
512 | }
|
---|
513 |
|
---|
514 | /********* Methode *********/
|
---|
515 | //++
|
---|
516 | void Histo::Add(float x, float w)
|
---|
517 | //
|
---|
518 | // Addition du contenu de l'histo pour abscisse x poids w
|
---|
519 | //--
|
---|
520 | {
|
---|
521 | int numBin = FindBin(x);
|
---|
522 | if (numBin<0) under += w;
|
---|
523 | else if (numBin>=bins) over += w;
|
---|
524 | else {
|
---|
525 | data[numBin] += w;
|
---|
526 | if(err2!=NULL) err2[numBin] += w*w;
|
---|
527 | nHist += w;
|
---|
528 | nEntries++;
|
---|
529 | }
|
---|
530 | }
|
---|
531 |
|
---|
532 | /********* Methode *********/
|
---|
533 | //++
|
---|
534 | void Histo::AddBin(int numBin, float w)
|
---|
535 | //
|
---|
536 | // Addition du contenu de l'histo bin numBin poids w
|
---|
537 | //--
|
---|
538 | {
|
---|
539 | if (numBin<0) under += w;
|
---|
540 | else if (numBin>=bins) over += w;
|
---|
541 | else {
|
---|
542 | data[numBin] += w;
|
---|
543 | if(err2!=NULL) err2[numBin] += w*w;
|
---|
544 | nHist += w;
|
---|
545 | nEntries++;
|
---|
546 | }
|
---|
547 | }
|
---|
548 |
|
---|
549 | /********* Methode *********/
|
---|
550 | //++
|
---|
551 | void Histo::SetBin(float x, float w)
|
---|
552 | //
|
---|
553 | // Remplissage du contenu de l'histo pour abscisse x poids w
|
---|
554 | //--
|
---|
555 | {
|
---|
556 | int numBin = FindBin(x);
|
---|
557 | SetBin(numBin,w);
|
---|
558 | }
|
---|
559 |
|
---|
560 | /********* Methode *********/
|
---|
561 | //++
|
---|
562 | void Histo::SetBin(int numBin, float w)
|
---|
563 | //
|
---|
564 | // Remplissage du contenu de l'histo pour numBin poids w
|
---|
565 | //--
|
---|
566 | {
|
---|
567 | if (numBin<0) under = w;
|
---|
568 | else if (numBin>=bins) over = w;
|
---|
569 | else {
|
---|
570 | nHist -= data[numBin];
|
---|
571 | data[numBin] = w;
|
---|
572 | nHist += w;
|
---|
573 | }
|
---|
574 | }
|
---|
575 |
|
---|
576 | /********* Methode *********/
|
---|
577 | //++
|
---|
578 | void Histo::SetErr2(float x, double e2)
|
---|
579 | //
|
---|
580 | // Remplissage des erreurs au carre pour abscisse x
|
---|
581 | //--
|
---|
582 | {
|
---|
583 | int numBin = FindBin(x);
|
---|
584 | SetErr2(numBin,e2);
|
---|
585 | }
|
---|
586 |
|
---|
587 | /********* Methode *********/
|
---|
588 | //++
|
---|
589 | void Histo::SetErr2(int numBin, double e2)
|
---|
590 | //
|
---|
591 | // Remplissage des erreurs au carre pour numBin poids
|
---|
592 | //--
|
---|
593 | {
|
---|
594 | if( err2==NULL) return;
|
---|
595 | if ( numBin<0 || numBin>=bins ) return;
|
---|
596 | err2[numBin] = e2;
|
---|
597 | }
|
---|
598 |
|
---|
599 | /********* Methode *********/
|
---|
600 | //++
|
---|
601 | void Histo::SetErr(float x, float e)
|
---|
602 | //
|
---|
603 | // Remplissage des erreurs pour abscisse x
|
---|
604 | //--
|
---|
605 | {
|
---|
606 | SetErr2(x, (double) e*e);
|
---|
607 | }
|
---|
608 |
|
---|
609 | /********* Methode *********/
|
---|
610 | //++
|
---|
611 | void Histo::SetErr(int numBin, float e)
|
---|
612 | //
|
---|
613 | // Remplissage des erreurs pour numBin
|
---|
614 | //--
|
---|
615 | {
|
---|
616 | SetErr2(numBin, (double) e*e);
|
---|
617 | }
|
---|
618 |
|
---|
619 | /********* Methode *********/
|
---|
620 | //++
|
---|
621 | int Histo::IMax() const
|
---|
622 | //
|
---|
623 | // Numero du bin ayant le contenu maximum
|
---|
624 | //--
|
---|
625 | {
|
---|
626 | int imx=0;
|
---|
627 | if(bins==1) return imx;
|
---|
628 | float mx=data[0];
|
---|
629 | for (int i=1; i<bins; i++)
|
---|
630 | if (data[i]>mx) {imx = i; mx=data[i];}
|
---|
631 | return imx;
|
---|
632 | }
|
---|
633 |
|
---|
634 | /********* Methode *********/
|
---|
635 | //++
|
---|
636 | int Histo::IMin() const
|
---|
637 | //
|
---|
638 | // Numero du bin ayant le contenu minimum
|
---|
639 | //--
|
---|
640 | {
|
---|
641 | int imx=0;
|
---|
642 | if(bins==1) return imx;
|
---|
643 | float mx=data[0];
|
---|
644 | for (int i=1; i<bins; i++)
|
---|
645 | if (data[i]<mx) {imx = i; mx=data[i];}
|
---|
646 | return imx;
|
---|
647 | }
|
---|
648 |
|
---|
649 | /********* Methode *********/
|
---|
650 | //++
|
---|
651 | float Histo::VMax() const
|
---|
652 | //
|
---|
653 | // Valeur le contenu maximum
|
---|
654 | //--
|
---|
655 | {
|
---|
656 | float mx=data[0];
|
---|
657 | if(bins==1) return mx;
|
---|
658 | for (int i=1;i<bins;i++) if(data[i]>mx) mx=data[i];
|
---|
659 | return mx;
|
---|
660 | }
|
---|
661 |
|
---|
662 | /********* Methode *********/
|
---|
663 | //++
|
---|
664 | float Histo::VMin() const
|
---|
665 | //
|
---|
666 | // Valeur le contenu minimum
|
---|
667 | //--
|
---|
668 | {
|
---|
669 | float mx=data[0];
|
---|
670 | if(bins==1) return mx;
|
---|
671 | for (int i=1;i<bins;i++) if(data[i]<mx) mx=data[i];
|
---|
672 | return mx;
|
---|
673 | }
|
---|
674 |
|
---|
675 | /********* Methode *********/
|
---|
676 | //++
|
---|
677 | float Histo::Mean() const
|
---|
678 | //
|
---|
679 | // Valeur moyenne
|
---|
680 | //--
|
---|
681 | {
|
---|
682 | double n = 0;
|
---|
683 | double sx = 0;
|
---|
684 | for (int i=0; i<bins; i++) {
|
---|
685 | double v = (data[i]>=0.) ? data[i] : -data[i];
|
---|
686 | n += v;
|
---|
687 | sx += BinCenter(i)*v;
|
---|
688 | }
|
---|
689 | if(n>0.) return sx/n; else return min;
|
---|
690 | }
|
---|
691 |
|
---|
692 | /********* Methode *********/
|
---|
693 | //++
|
---|
694 | float Histo::Sigma() const
|
---|
695 | //
|
---|
696 | // Valeur du sigma
|
---|
697 | //--
|
---|
698 | {
|
---|
699 | double n = 0;
|
---|
700 | double sx = 0;
|
---|
701 | double sx2= 0;
|
---|
702 | for (int i=0; i<bins; i++) {
|
---|
703 | double v = (data[i]>=0.) ? data[i] : -data[i];
|
---|
704 | n += v;
|
---|
705 | sx += BinCenter(i)*v;
|
---|
706 | sx2+= BinCenter(i)*BinCenter(i)*v;
|
---|
707 | }
|
---|
708 | // attention a l'erreur d'arrondi si un seul bin rempli!!
|
---|
709 | if( n == 0 ) return 0.;
|
---|
710 | sx2 = sx2/n - (sx/n)*(sx/n);
|
---|
711 | if( sx2>0. ) return sqrt( sx2 );
|
---|
712 | else return 0.;
|
---|
713 | }
|
---|
714 |
|
---|
715 | /********* Methode *********/
|
---|
716 | //++
|
---|
717 | float Histo::MeanLH(int il,int ih) const
|
---|
718 | //
|
---|
719 | // Valeur de la moyenne calculee entre les bin il et ih
|
---|
720 | //--
|
---|
721 | {
|
---|
722 | if( ih < il ) { il = 0; ih = bins-1; }
|
---|
723 | if( il < 0 ) il = 0;
|
---|
724 | if( ih >= bins ) ih = bins-1;
|
---|
725 | double n = 0;
|
---|
726 | double sx = 0;
|
---|
727 | for (int i=il; i<=ih; i++) {
|
---|
728 | double v = (data[i]>=0.) ? data[i] : -data[i];
|
---|
729 | n += v;
|
---|
730 | sx += BinCenter(i)*v;
|
---|
731 | }
|
---|
732 | if(n>0.) return sx/n; else return BinLowEdge(il);
|
---|
733 | }
|
---|
734 |
|
---|
735 | /********* Methode *********/
|
---|
736 | //++
|
---|
737 | float Histo::SigmaLH(int il,int ih) const
|
---|
738 | //
|
---|
739 | // Valeur de la moyenne calculee entre les bin il et ih
|
---|
740 | //--
|
---|
741 | {
|
---|
742 | if( ih < il ) { il = 0; ih = bins - 1; }
|
---|
743 | if( il < 0 ) il = 0;
|
---|
744 | if( ih >= bins ) ih = bins - 1;
|
---|
745 | double n = 0;
|
---|
746 | double sx = 0;
|
---|
747 | double sx2= 0;
|
---|
748 | for (int i=il; i<=ih; i++) {
|
---|
749 | double v = (data[i]>=0.) ? data[i] : -data[i];
|
---|
750 | n += v;
|
---|
751 | sx += BinCenter(i)*v;
|
---|
752 | sx2+= BinCenter(i)*BinCenter(i)*v;
|
---|
753 | }
|
---|
754 | if( n == 0 ) return 0.;
|
---|
755 | sx2 = sx2/n - (sx/n)*(sx/n);
|
---|
756 | if( sx2>0. ) return sqrt( sx2 );
|
---|
757 | else return 0.;
|
---|
758 | }
|
---|
759 |
|
---|
760 | /********* Methode *********/
|
---|
761 | //++
|
---|
762 | float Histo::Mean(float x0, float dx) const
|
---|
763 | //
|
---|
764 | // Valeur de la moyenne calculee entre x0-dx et x0+dx
|
---|
765 | //--
|
---|
766 | {
|
---|
767 | double sdata = 0;
|
---|
768 | double sx = 0;
|
---|
769 | int iMin = FindBin(x0-dx);
|
---|
770 | if (iMin<0) iMin = 0;
|
---|
771 | int iMax = FindBin(x0+dx);
|
---|
772 | if (iMax>bins) iMax=bins;
|
---|
773 | for (int i=iMin; i<iMax; i++) {
|
---|
774 | double v = (data[i]>=0.) ? data[i] : -data[i];
|
---|
775 | sx += BinCenter(i)*v;
|
---|
776 | sdata += v;
|
---|
777 | }
|
---|
778 | if(sdata>0.) return sx/sdata; else return min;
|
---|
779 | }
|
---|
780 |
|
---|
781 | /********* Methode *********/
|
---|
782 | //++
|
---|
783 | float Histo::Sigma(float x0, float dx) const
|
---|
784 | //
|
---|
785 | // Valeur du sigma calcule entre x0-dx et x0+dx
|
---|
786 | //--
|
---|
787 | {
|
---|
788 | double sx = 0;
|
---|
789 | double sx2= 0;
|
---|
790 | double sdata = 0;
|
---|
791 | int iMin = FindBin(x0-dx);
|
---|
792 | if (iMin<0) iMin = 0;
|
---|
793 | int iMax = FindBin(x0+dx);
|
---|
794 | if (iMax>bins) iMax=bins;
|
---|
795 | for (int i=iMin; i<iMax; i++) {
|
---|
796 | double v = (data[i]>=0.) ? data[i] : -data[i];
|
---|
797 | sx += BinCenter(i)*v;
|
---|
798 | sx2+= BinCenter(i)*BinCenter(i)*v;
|
---|
799 | sdata += v;
|
---|
800 | }
|
---|
801 | if(sdata>0.) return sqrt( sx2/sdata - (sx/sdata)*(sx/sdata) );
|
---|
802 | else return 0.;
|
---|
803 | }
|
---|
804 |
|
---|
805 | /********* Methode *********/
|
---|
806 | //++
|
---|
807 | float Histo::CleanedMean(float& sigma) const
|
---|
808 | //
|
---|
809 | // Valeur de la moyenne et du sigma nettoyes
|
---|
810 | //--
|
---|
811 | {
|
---|
812 | if (!nHist) return 0;
|
---|
813 | // on fait ca de facon bete, on pourra raffiner apres...
|
---|
814 | float x0 = Mean();
|
---|
815 | float s = Sigma()+binWidth;
|
---|
816 |
|
---|
817 | for (int i=0; i<3; i++) {
|
---|
818 | float xx0 = Mean(x0,5*s);
|
---|
819 | s = Sigma(x0,5*s)+binWidth;
|
---|
820 | x0 = xx0;
|
---|
821 | }
|
---|
822 | sigma = s;
|
---|
823 | return x0;
|
---|
824 | }
|
---|
825 |
|
---|
826 | /********* Methode *********/
|
---|
827 | //++
|
---|
828 | float Histo::CleanedMean() const
|
---|
829 | //
|
---|
830 | // Valeur de la moyenne nettoyee
|
---|
831 | //--
|
---|
832 | {
|
---|
833 | if (!nHist) return 0;
|
---|
834 | float s=0;
|
---|
835 | return CleanedMean(s);
|
---|
836 | }
|
---|
837 |
|
---|
838 | /********* Methode *********/
|
---|
839 | //++
|
---|
840 | int Histo::BinNonNul() const
|
---|
841 | //
|
---|
842 | // Retourne le nombre de bins non-nul
|
---|
843 | //--
|
---|
844 | {
|
---|
845 | int non=0;
|
---|
846 | for (int i=0;i<bins;i++) if( data[i] != 0. ) non++;
|
---|
847 | return non;
|
---|
848 | }
|
---|
849 |
|
---|
850 | /********* Methode *********/
|
---|
851 | //++
|
---|
852 | int Histo::ErrNonNul() const
|
---|
853 | //
|
---|
854 | // Retourne le nombre de bins ayant une erreur non-nulle
|
---|
855 | //--
|
---|
856 | {
|
---|
857 | if(err2==NULL) return -1;
|
---|
858 | int non=0;
|
---|
859 | for (int i=0;i<bins;i++) if( err2[i] != 0. ) non++;
|
---|
860 | return non;
|
---|
861 | }
|
---|
862 |
|
---|
863 | /********* Methode *********/
|
---|
864 | //++
|
---|
865 | int Histo::BinPercent(float per) const
|
---|
866 | //
|
---|
867 | // Renvoie le numero de bin tel que il y ait "per" pourcent entrees
|
---|
868 | // entre le bin 0 et ce bin (y compris le contenu de ce bin).
|
---|
869 | //--
|
---|
870 | {
|
---|
871 | double n = nHist*per;
|
---|
872 | double s = 0.;
|
---|
873 | int i;
|
---|
874 | for(i=0; i<bins; i++ ) {
|
---|
875 | s += data[i];
|
---|
876 | if( s >= n ) break;
|
---|
877 | }
|
---|
878 | if( i == bins ) i = bins-1;
|
---|
879 | return i;
|
---|
880 | }
|
---|
881 |
|
---|
882 | /********* Methode *********/
|
---|
883 | //++
|
---|
884 | int Histo::BinPercent(float x,float per,int& imin,int& imax)
|
---|
885 | //
|
---|
886 | // Renvoie les numeros de bins imin,imax (0=<i<bins) tels que:
|
---|
887 | //| notons I = bin contenant l'abscisse x
|
---|
888 | //| N1 = nombre d'entrees entre bin 0 et I compris
|
---|
889 | //| N2 = nombre d'entrees entre bin I et bins-1 compris
|
---|
890 | //| imin = bin tel que nombre d'entrees entre imin et I = N1 * per
|
---|
891 | //| imax = bin tel que nombre d'entrees entre I et imax = N2 * per
|
---|
892 | //| Return: <0 si echec
|
---|
893 | //| min(I-imin,imax-I) si ok
|
---|
894 | //--
|
---|
895 | {
|
---|
896 | imin = imax = -1;
|
---|
897 | if( per <= 0. ) return -1;
|
---|
898 |
|
---|
899 | int I = FindBin(x);
|
---|
900 | if( I<0 || I>=bins ) return -2;
|
---|
901 |
|
---|
902 | double N1 = 0.; for(int i=0; i<=I; i++) N1 += data[i]; N1 *= per;
|
---|
903 | double N2 = 0.; {for(int i=I; i<bins; i++) N2 += data[i]; N2 *= per;}
|
---|
904 |
|
---|
905 | double n = 0.;
|
---|
906 | for(imin=I; imin>=0; imin--) { n += data[imin]; if(n>=N1) break; }
|
---|
907 | if( imin<0 ) imin = 0;
|
---|
908 | // cout<<"I="<<I<<" imin="<<imin<<" n="<<n<<" N1="<<N1<<endl;
|
---|
909 |
|
---|
910 | n = 0.;
|
---|
911 | for(imax=I; imax<bins; imax++) { n += data[imax]; if(n>=N2) break; }
|
---|
912 | if( imax>=bins ) imax = bins-1;
|
---|
913 | // cout<<"I="<<I<<" imax="<<imax<<" n="<<n<<" N2="<<N2<<endl;
|
---|
914 |
|
---|
915 | return ( imax-I > I-imin ) ? I-imin : imax-I ;
|
---|
916 | }
|
---|
917 |
|
---|
918 | /********* Methode *********/
|
---|
919 | //++
|
---|
920 | int Histo::BinPercent(float x,float per,float& xmin,float& xmax)
|
---|
921 | //
|
---|
922 | // Idem precedent mais renvoie xmin et xmax
|
---|
923 | //--
|
---|
924 | {
|
---|
925 | xmin = xmax = 0.;
|
---|
926 | int imin,imax;
|
---|
927 | int rc = BinPercent(x,per,imin,imax);
|
---|
928 | if( rc >= 0 ) { xmin = BinLowEdge(imin); xmax = BinHighEdge(imax); }
|
---|
929 | return rc;
|
---|
930 | }
|
---|
931 |
|
---|
932 | /********* Methode *********/
|
---|
933 | //++
|
---|
934 | void Histo::HInteg(float norm)
|
---|
935 | //
|
---|
936 | // Remplace l'histogramme par son integrale normalise a norm:
|
---|
937 | // si norm <= 0 : pas de normalisation, integration seule
|
---|
938 | //--
|
---|
939 | {
|
---|
940 | if(bins>1)
|
---|
941 | for(int i=1; i<bins; i++) {
|
---|
942 | data[i] += data[i-1];
|
---|
943 | if(err2!=NULL) err2[i] += err2[i-1];
|
---|
944 | }
|
---|
945 | if( norm <= 0. ) return;
|
---|
946 | norm /= data[bins-1];
|
---|
947 | for(int i=0; i<bins; i++) {
|
---|
948 | data[i] *= norm;
|
---|
949 | if(err2!=NULL) err2[i] *= norm*norm;
|
---|
950 | }
|
---|
951 | }
|
---|
952 |
|
---|
953 | /********* Methode *********/
|
---|
954 | //++
|
---|
955 | void Histo::HDeriv()
|
---|
956 | //
|
---|
957 | // Remplace l'histogramme par sa derivee
|
---|
958 | //--
|
---|
959 | {
|
---|
960 | if( bins <= 1 ) return;
|
---|
961 | float *temp = new float[bins];
|
---|
962 | memcpy(temp, data, bins*sizeof(float));
|
---|
963 | if(bins>=3) for(int i=1; i<bins-1; i++)
|
---|
964 | temp[i] = (data[i+1]-data[i-1])/(2.*binWidth);
|
---|
965 | temp[0] = (data[1]-data[0])/binWidth;
|
---|
966 | temp[bins-1] = (data[bins-1]-data[bins-2])/binWidth;
|
---|
967 | memcpy(data, temp, bins*sizeof(float));
|
---|
968 | delete [] temp;
|
---|
969 | }
|
---|
970 |
|
---|
971 | /********* Methode *********/
|
---|
972 | //++
|
---|
973 | void Histo::HRebin(int nbinew)
|
---|
974 | //
|
---|
975 | // Pour rebinner l'histogramme sur nbinew bins
|
---|
976 | //--
|
---|
977 | {
|
---|
978 | if( nbinew <= 0 ) return;
|
---|
979 |
|
---|
980 | // memoraisation de l'histogramme original
|
---|
981 | Histo H(*this);
|
---|
982 |
|
---|
983 | // Le nombre de bins est il plus grand ??
|
---|
984 | if( nbinew > bins ) {
|
---|
985 | delete [] data; data = NULL;
|
---|
986 | data = new float[nbinew];
|
---|
987 | if( err2 != NULL ) {
|
---|
988 | delete [] err2; err2 = NULL;
|
---|
989 | err2 = new double[nbinew];
|
---|
990 | }
|
---|
991 | }
|
---|
992 |
|
---|
993 | // remise en forme de this
|
---|
994 | bins = nbinew;
|
---|
995 | this->Zero();
|
---|
996 | binWidth = (max-min)/bins;
|
---|
997 | // On recopie les parametres qui n'ont pas changes
|
---|
998 | under = H.under;
|
---|
999 | over = H.over;
|
---|
1000 |
|
---|
1001 |
|
---|
1002 | // Remplissage
|
---|
1003 | for(int i=0;i<bins;i++) {
|
---|
1004 | float xmi = BinLowEdge(i);
|
---|
1005 | float xma = BinHighEdge(i);
|
---|
1006 | int imi = H.FindBin(xmi);
|
---|
1007 | if( imi < 0 ) imi = 0;
|
---|
1008 | int ima = H.FindBin(xma);
|
---|
1009 | if( ima >= H.bins ) ima = H.bins-1;
|
---|
1010 | double w = 0.;
|
---|
1011 | if( ima == imi ) {
|
---|
1012 | w = H(imi) * binWidth/H.binWidth;
|
---|
1013 | } else {
|
---|
1014 | w += H(imi) * (H.BinHighEdge(imi)-xmi)/H.binWidth;
|
---|
1015 | w += H(ima) * (xma -H.BinLowEdge(ima))/H.binWidth;
|
---|
1016 | if( ima > imi+1 )
|
---|
1017 | for(int ii=imi+1;ii<ima;ii++) w += H(ii);
|
---|
1018 | }
|
---|
1019 | AddBin(i,(float) w);
|
---|
1020 | }
|
---|
1021 |
|
---|
1022 | }
|
---|
1023 |
|
---|
1024 | /********* Methode *********/
|
---|
1025 | //++
|
---|
1026 | int Histo::MaxiLocal(float& maxi,int& imax,float& maxn,int& imaxn)
|
---|
1027 | //
|
---|
1028 | // Retourne le nombre de maxima locaux, la valeur du maximum (maxi)
|
---|
1029 | // et sa position (imax), ainsi que la valeur du second maximum
|
---|
1030 | // local (maxn) et sa position (imaxn).
|
---|
1031 | // Attention: si un maximum a la meme valeur sur plusieurs bins
|
---|
1032 | // consecutifs, le bin le plus a droite est pris.
|
---|
1033 | //--
|
---|
1034 | {
|
---|
1035 | int nml = 0;
|
---|
1036 | imax = imaxn = -1;
|
---|
1037 | maxi = maxn = -1.f;
|
---|
1038 |
|
---|
1039 | bool up = true;
|
---|
1040 | bool down = false;
|
---|
1041 | for(int i=0;i<bins;i++) {
|
---|
1042 | if( !up ) if( data[i] > data[i-1] ) up = true;
|
---|
1043 | if( up && !down ) {
|
---|
1044 | if( i == bins-1 ) down=true;
|
---|
1045 | else if( data[i] > data[i+1] ) down=true;
|
---|
1046 | }
|
---|
1047 |
|
---|
1048 | if( up && down ) {
|
---|
1049 | nml++;
|
---|
1050 | up = down = false;
|
---|
1051 | if( imax < 0 ) {
|
---|
1052 | imax = i; maxi = data[i];
|
---|
1053 | } else if( data[i] >= maxi ) {
|
---|
1054 | imaxn = imax; maxn = maxi;
|
---|
1055 | imax = i; maxi = data[i];
|
---|
1056 | } else {
|
---|
1057 | if( imaxn < 0 || data[i] >= maxn ) { imaxn = i; maxn = data[i]; }
|
---|
1058 | }
|
---|
1059 | }
|
---|
1060 |
|
---|
1061 | }
|
---|
1062 | return nml;
|
---|
1063 | }
|
---|
1064 |
|
---|
1065 | /********* Methode *********/
|
---|
1066 | //++
|
---|
1067 | float Histo::FitMax(int degree, float frac, int debug) const
|
---|
1068 | //
|
---|
1069 | // Fit de la position du maximum de l'histo par un polynome
|
---|
1070 | // de degre `degree' a `frac' pourcent du maximum.
|
---|
1071 | // L'histo est suppose etre remplit de valeurs positives
|
---|
1072 | //--
|
---|
1073 | {
|
---|
1074 | if (degree < 2 || degree > 3) degree = 3;
|
---|
1075 | if (frac <= 0. || frac >= 1.) frac = 0.5;
|
---|
1076 |
|
---|
1077 | if (debug > 1)
|
---|
1078 | cout<<"Histo::FitMax : Nb Entrees histo ="<<NEntries()<<endl;
|
---|
1079 |
|
---|
1080 | if (NEntries() < 1) THROW(inconsistentErr);
|
---|
1081 |
|
---|
1082 | int iMax = IMax();
|
---|
1083 | float hmax = (*this)(iMax);
|
---|
1084 | float xCenter = BinCenter(iMax);
|
---|
1085 |
|
---|
1086 | if(debug>1)
|
---|
1087 | cout<<"Max histo i="<<iMax<<" x="<<xCenter<<" v="<<hmax<<endl;
|
---|
1088 |
|
---|
1089 | /* find limits at frac*hmax */
|
---|
1090 |
|
---|
1091 | float limit = frac*hmax;
|
---|
1092 |
|
---|
1093 | volatile int iLow = iMax;
|
---|
1094 | while (iLow>0 && (*this)(iLow)>limit) iLow--;
|
---|
1095 |
|
---|
1096 | volatile int iHigh = iMax;
|
---|
1097 | while (iHigh<NBins()-1 && (*this)(iHigh)>limit) iHigh++;
|
---|
1098 |
|
---|
1099 | int nLowHigh;
|
---|
1100 | for(;;) {
|
---|
1101 | nLowHigh = 0;
|
---|
1102 | for (int i=iLow; i<=iHigh; i++) if((*this)(i)>0) {
|
---|
1103 | if(!err2) nLowHigh++;
|
---|
1104 | else if(Error2(i)>0.) nLowHigh++;
|
---|
1105 | }
|
---|
1106 | if (debug > 1) cout <<"Limites histo "<<iLow<<" - "<<iHigh
|
---|
1107 | <<" ("<<nLowHigh<<" non nuls)"<<endl;
|
---|
1108 | if( nLowHigh >= degree+1 ) break;
|
---|
1109 | iLow--; iHigh++;
|
---|
1110 | if(iLow<0 && iHigh>=NBins()) {
|
---|
1111 | if (debug>1)
|
---|
1112 | cout<<"Mode histogramme = "<<xCenter
|
---|
1113 | <<" BinCenter("<<iMax<<")"<<endl;
|
---|
1114 | return xCenter;
|
---|
1115 | }
|
---|
1116 | if(iLow < 0 ) iLow = 0;
|
---|
1117 | if(iHigh >= NBins()) iHigh = NBins()-1;
|
---|
1118 | }
|
---|
1119 |
|
---|
1120 | Vector xFit(nLowHigh);
|
---|
1121 | Vector yFit(nLowHigh);
|
---|
1122 | Vector e2Fit(nLowHigh);
|
---|
1123 | Vector errcoef(1);
|
---|
1124 | int ii = 0;
|
---|
1125 | for (int j=iLow; j<=iHigh; j++) {
|
---|
1126 | if ((*this)(j)>0) {
|
---|
1127 | if(!err2) {
|
---|
1128 | xFit(ii) = BinCenter(j)-xCenter;
|
---|
1129 | yFit(ii) = (*this)(j);
|
---|
1130 | e2Fit(ii) = yFit(ii);
|
---|
1131 | ii++;
|
---|
1132 | } else if(Error2(j)>0.) {
|
---|
1133 | xFit(ii) = BinCenter(j)-xCenter;
|
---|
1134 | yFit(ii) = (*this)(j);
|
---|
1135 | e2Fit(ii) = Error2(j);
|
---|
1136 | ii++;
|
---|
1137 | }
|
---|
1138 | }
|
---|
1139 | }
|
---|
1140 | if(debug>4) {
|
---|
1141 | int k;
|
---|
1142 | for(k=0;k<nLowHigh;k++) cout<<" "<<xFit(k); cout<<endl;
|
---|
1143 | for(k=0;k<nLowHigh;k++) cout<<" "<<yFit(k); cout<<endl;
|
---|
1144 | for(k=0;k<nLowHigh;k++) cout<<" "<<e2Fit(k); cout<<endl;
|
---|
1145 | }
|
---|
1146 | if( ii != nLowHigh ) THROW(inconsistentErr);
|
---|
1147 | Poly pol(degree);
|
---|
1148 | TRY {
|
---|
1149 | pol.Fit(xFit, yFit, e2Fit, degree, errcoef);
|
---|
1150 | if (debug>1) cout << "resultat fit = " << pol << endl;
|
---|
1151 | pol.Derivate();
|
---|
1152 | } CATCHALL {
|
---|
1153 | THROW_SAME;
|
---|
1154 | } ENDTRY
|
---|
1155 |
|
---|
1156 | int DPolDeg = pol.Degre();
|
---|
1157 | float fd = 0.;
|
---|
1158 |
|
---|
1159 | if (DPolDeg == 0) {
|
---|
1160 | // on est dans le cas d'un fit de droite
|
---|
1161 | if( pol[0] > 0. ) {
|
---|
1162 | // on a fitte une droite de pente >0.
|
---|
1163 | fd = xFit(nLowHigh-1) + binWidth/2. + xCenter;
|
---|
1164 | } else if( pol[0] < 0. ) {
|
---|
1165 | // on a fitte une droite de pente <0.
|
---|
1166 | fd = xFit(0) - binWidth/2. + xCenter;
|
---|
1167 | } else {
|
---|
1168 | // on a fitte une droite de pente =0. (creneau)
|
---|
1169 | fd = (xFit(0)+xFit(nLowHigh-1))/2. + xCenter;
|
---|
1170 | }
|
---|
1171 | } else if (DPolDeg == 1) {
|
---|
1172 | // on est dans le cas d'un fit de parabole
|
---|
1173 | double r=0;
|
---|
1174 | if (pol.Root1(r)==0) THROW(inconsistentErr);
|
---|
1175 | fd = r + xCenter;
|
---|
1176 | } else if (DPolDeg == 2) {
|
---|
1177 | // on est dans le cas d'un fit de cubique
|
---|
1178 | double r1=0;
|
---|
1179 | double r2=0;
|
---|
1180 | if (pol.Root2(r1,r2) == 0) THROW(inconsistentErr);
|
---|
1181 | pol.Derivate();
|
---|
1182 | fd = (pol(r1)<0) ? r1 + xCenter : r2 + xCenter;
|
---|
1183 | } else {
|
---|
1184 | // on est dans un cas non prevu
|
---|
1185 | THROW(inconsistentErr);
|
---|
1186 | }
|
---|
1187 |
|
---|
1188 | if(fd>max) fd = max;
|
---|
1189 | if(fd<min) fd = min;
|
---|
1190 |
|
---|
1191 | if (debug>1)
|
---|
1192 | cout << "Mode histogramme = " << fd
|
---|
1193 | << " (DerPol.degre =" << DPolDeg
|
---|
1194 | << " pol.coeff[0] =" << pol[0]
|
---|
1195 | << ")" << endl;
|
---|
1196 |
|
---|
1197 | return fd;
|
---|
1198 | }
|
---|
1199 |
|
---|
1200 |
|
---|
1201 | /********* Methode *********/
|
---|
1202 | //++
|
---|
1203 | float Histo::FindWidth(float frac, int debug) const
|
---|
1204 | //
|
---|
1205 | // Calcul de la largeur a frac pourcent du maximum
|
---|
1206 | // autour du bin du maximum.
|
---|
1207 | // L'histo est suppose etre remplit de valeurs positives
|
---|
1208 | //--
|
---|
1209 | {
|
---|
1210 | float xmax = BinCenter( IMax() );
|
---|
1211 | return FindWidth(xmax,frac,debug);
|
---|
1212 | }
|
---|
1213 |
|
---|
1214 | /********* Methode *********/
|
---|
1215 | //++
|
---|
1216 | float Histo::FindWidth(float xmax,float frac, int debug) const
|
---|
1217 | //
|
---|
1218 | // Calcul de la largeur a frac pourcent de la valeur du bin xmax.
|
---|
1219 | // L'histo est suppose etre remplit de valeurs positives
|
---|
1220 | //--
|
---|
1221 | {
|
---|
1222 | if (frac <= 0 || frac >= 1.) frac = 0.5;
|
---|
1223 |
|
---|
1224 | if (debug > 1)
|
---|
1225 | cout << "Histo::FindWidth a " << frac
|
---|
1226 | << " de xmax= " << xmax
|
---|
1227 | << " , ndata= " << NData()
|
---|
1228 | << " , nent= " << NEntries()
|
---|
1229 | << " , nbin= " << NBins() << endl;
|
---|
1230 |
|
---|
1231 | if (NEntries() < 1) THROW(inconsistentErr);
|
---|
1232 | if (NBins() < 3) THROW(inconsistentErr);
|
---|
1233 |
|
---|
1234 | int iMax = FindBin(xmax);
|
---|
1235 | if (iMax<0 || iMax>=NBins()) THROW(inconsistentErr);
|
---|
1236 | float hmax = data[iMax];
|
---|
1237 | float limit = frac*hmax;
|
---|
1238 | if (debug > 1)
|
---|
1239 | cout << " Max histo[" << iMax << "] = " << hmax
|
---|
1240 | << ", limite " << limit << endl;
|
---|
1241 |
|
---|
1242 | int iLow = iMax;
|
---|
1243 | while (iLow>=0 && data[iLow]>limit) iLow--;
|
---|
1244 | if( iLow < 0 ) iLow = 0;
|
---|
1245 |
|
---|
1246 | int iHigh = iMax;
|
---|
1247 | while (iHigh<NBins() && data[iHigh]>limit) iHigh++;
|
---|
1248 | if( iHigh >=NBins() ) iHigh = NBins()-1;
|
---|
1249 |
|
---|
1250 | float xlow = BinCenter(iLow);
|
---|
1251 | float ylow = data[iLow];
|
---|
1252 |
|
---|
1253 | float xlow1=xlow, ylow1=ylow;
|
---|
1254 | if(iLow+1<NBins()) {
|
---|
1255 | xlow1 = BinCenter(iLow+1);
|
---|
1256 | ylow1 = data[iLow+1];
|
---|
1257 | }
|
---|
1258 |
|
---|
1259 | float xhigh = BinCenter(iHigh);
|
---|
1260 | float yhigh = data[iHigh];
|
---|
1261 |
|
---|
1262 | float xhigh1=xhigh, yhigh1=yhigh;
|
---|
1263 | if(iHigh-1>=0) {
|
---|
1264 | xhigh1 = BinCenter(iHigh-1);
|
---|
1265 | yhigh1 = data[iHigh-1];
|
---|
1266 | }
|
---|
1267 |
|
---|
1268 | float xlpred,xhpred,wd;
|
---|
1269 |
|
---|
1270 | if(ylow1>ylow) {
|
---|
1271 | xlpred = xlow + (xlow1-xlow)/(ylow1-ylow)*(limit-ylow);
|
---|
1272 | if(xlpred < xlow) xlpred = xlow;
|
---|
1273 | } else xlpred = xlow;
|
---|
1274 |
|
---|
1275 | if(yhigh1>yhigh) {
|
---|
1276 | xhpred = xhigh + (xhigh1-xhigh)/(yhigh1-yhigh)*(limit-yhigh);
|
---|
1277 | if(xhpred > xhigh) xhpred = xhigh;
|
---|
1278 | } else xhpred = xhigh;
|
---|
1279 |
|
---|
1280 | wd = xhpred - xlpred;
|
---|
1281 |
|
---|
1282 | if (debug > 1) {
|
---|
1283 | cout << "Limites histo: " << " Width " << wd << endl;
|
---|
1284 | cout << " Low: [" << iLow << "]=" << ylow << " pred-> " << xlpred
|
---|
1285 | << " - High: [" << iHigh << "]=" << yhigh << " pred-> " << xhpred
|
---|
1286 | << endl;
|
---|
1287 | }
|
---|
1288 |
|
---|
1289 | return wd;
|
---|
1290 | }
|
---|
1291 |
|
---|
1292 |
|
---|
1293 | /********* Methode *********/
|
---|
1294 | //++
|
---|
1295 | int Histo::EstimeMax(float& xm,int SzPav)
|
---|
1296 | //
|
---|
1297 | // Cf suivant mais im est le bin du maximum de l'histo
|
---|
1298 | //--
|
---|
1299 | {
|
---|
1300 | int im = IMax();
|
---|
1301 | return EstimeMax(im,xm,SzPav);
|
---|
1302 | }
|
---|
1303 |
|
---|
1304 | /********* Methode *********/
|
---|
1305 | //++
|
---|
1306 | int Histo::EstimeMax(int& im,float& xm,int SzPav)
|
---|
1307 | //
|
---|
1308 | // Determine l'abscisses du maximum donne par im
|
---|
1309 | // en moyennant dans un pave SzPav autour du maximum
|
---|
1310 | //| Return:
|
---|
1311 | //| 0 = si fit maximum reussi avec SzPav pixels
|
---|
1312 | //| 1 = si fit maximum reussi avec moins que SzPav pixels
|
---|
1313 | //| 2 = si fit maximum echoue et renvoit BinCenter()
|
---|
1314 | //| -1 = si echec: SzPav <= 0 ou im hors limites
|
---|
1315 | //--
|
---|
1316 | {
|
---|
1317 | xm = 0;
|
---|
1318 | if( SzPav <= 0 ) return -1;
|
---|
1319 | if( im < 0 || im >= bins ) return -1;
|
---|
1320 |
|
---|
1321 | if( SzPav%2 == 0 ) SzPav++;
|
---|
1322 | SzPav = (SzPav-1)/2;
|
---|
1323 |
|
---|
1324 | int rc = 0;
|
---|
1325 | double dxm = 0, wx = 0;
|
---|
1326 | for(int i=im-SzPav;i<=im+SzPav;i++) {
|
---|
1327 | if( i<0 || i>= bins ) {rc=1; continue;}
|
---|
1328 | dxm += BinCenter(i) * (*this)(i);
|
---|
1329 | wx += (*this)(i);
|
---|
1330 | }
|
---|
1331 |
|
---|
1332 | if( wx > 0. ) {
|
---|
1333 | xm = dxm/wx;
|
---|
1334 | return rc;
|
---|
1335 | } else {
|
---|
1336 | xm = BinCenter(im);
|
---|
1337 | return 2;
|
---|
1338 | }
|
---|
1339 |
|
---|
1340 | }
|
---|
1341 |
|
---|
1342 | /********* Methode *********/
|
---|
1343 | //++
|
---|
1344 | void Histo::EstimeWidthS(float frac,float& widthG,float& widthD)
|
---|
1345 | //
|
---|
1346 | // Determine la largeur a frac% du maximum a gauche (widthG)
|
---|
1347 | // et a droite (widthD)
|
---|
1348 | //--
|
---|
1349 | {
|
---|
1350 | int i;
|
---|
1351 | widthG = widthD = -1.;
|
---|
1352 | if( bins<=1 || frac<=0. || frac>=1. ) return;
|
---|
1353 |
|
---|
1354 | int imax = 0;
|
---|
1355 | float maxi = data[0];
|
---|
1356 | for(i=1;i<bins;i++) if(data[i]>maxi) {imax=i; maxi=data[i];}
|
---|
1357 | float xmax = BinCenter(imax);
|
---|
1358 | float maxf = maxi * frac;
|
---|
1359 |
|
---|
1360 | // recherche du sigma a gauche
|
---|
1361 | widthG = 0.;
|
---|
1362 | for(i=imax;i>=0;i--) if( data[i] <= maxf ) break;
|
---|
1363 | if(i<0) i=0;
|
---|
1364 | if(i<imax) {
|
---|
1365 | if( data[i+1] != data[i] ) {
|
---|
1366 | widthG = BinCenter(i) + binWidth
|
---|
1367 | * (maxf-data[i])/(data[i+1]-data[i]);
|
---|
1368 | widthG = xmax - widthG;
|
---|
1369 | } else widthG = xmax - BinCenter(i);
|
---|
1370 | }
|
---|
1371 |
|
---|
1372 | // recherche du sigma a droite
|
---|
1373 | widthD = 0.;
|
---|
1374 | for(i=imax;i<bins;i++) if( data[i] <= maxf ) break;
|
---|
1375 | if(i>=bins) i=bins-1;
|
---|
1376 | if(i>imax) {
|
---|
1377 | if( data[i] != data[i-1] ) {
|
---|
1378 | widthD = BinCenter(i) - binWidth
|
---|
1379 | * (maxf-data[i])/(data[i-1]-data[i]);
|
---|
1380 | widthD -= xmax;
|
---|
1381 | } else widthD = BinCenter(i) - xmax;
|
---|
1382 | }
|
---|
1383 |
|
---|
1384 | }
|
---|
1385 |
|
---|
1386 | //////////////////////////////////////////////////////////
|
---|
1387 | //++
|
---|
1388 | int Histo::Fit(GeneralFit& gfit,unsigned short typ_err)
|
---|
1389 | //
|
---|
1390 | // Fit de l'histogramme par ``gfit''.
|
---|
1391 | //| typ_err = 0 :
|
---|
1392 | //| - erreur attachee au bin si elle existe
|
---|
1393 | //| - sinon 1
|
---|
1394 | //| typ_err = 1 :
|
---|
1395 | //| - erreur attachee au bin si elle existe
|
---|
1396 | //| - sinon max( sqrt(abs(bin) ,1 )
|
---|
1397 | //| typ_err = 2 :
|
---|
1398 | //| - erreur forcee a 1
|
---|
1399 | //| typ_err = 3 :
|
---|
1400 | //| - erreur forcee a max( sqrt(abs(bin) ,1 )
|
---|
1401 | //| typ_err = 4 :
|
---|
1402 | //| - erreur forcee a 1, nulle si bin a zero.
|
---|
1403 | //| typ_err = 5 :
|
---|
1404 | //| - erreur forcee a max( sqrt(abs(bin) ,1 ),
|
---|
1405 | //| nulle si bin a zero.
|
---|
1406 | //--
|
---|
1407 | {
|
---|
1408 | if(NBins()<=0) return -1000;
|
---|
1409 | if(typ_err>5) typ_err=0;
|
---|
1410 |
|
---|
1411 | GeneralFitData mydata(1,NBins());
|
---|
1412 |
|
---|
1413 | for(int i=0;i<NBins();i++) {
|
---|
1414 | double x = (double) BinCenter(i);
|
---|
1415 | double f = (double) (*this)(i);
|
---|
1416 | double saf = sqrt(fabs((double) f)); if(saf<1.) saf=1.;
|
---|
1417 | double e;
|
---|
1418 | if(typ_err==0) {if(HasErrors()) e=Error(i); else e=1.;}
|
---|
1419 | else if(typ_err==1) {if(HasErrors()) e=Error(i); else e=saf;}
|
---|
1420 | else if(typ_err==2) e=1.;
|
---|
1421 | else if(typ_err==3) e=saf;
|
---|
1422 | else if(typ_err==4) e=(f==0.)?0.:1.;
|
---|
1423 | else if(typ_err==5) e=(f==0.)?0.:saf;
|
---|
1424 | mydata.AddData1(x,f,e);
|
---|
1425 | }
|
---|
1426 |
|
---|
1427 | gfit.SetData(&mydata);
|
---|
1428 |
|
---|
1429 | return gfit.Fit();
|
---|
1430 | }
|
---|
1431 |
|
---|
1432 | //++
|
---|
1433 | Histo Histo::FitResidus(GeneralFit& gfit)
|
---|
1434 | //
|
---|
1435 | // Retourne une classe contenant les residus du fit ``gfit''.
|
---|
1436 | //--
|
---|
1437 | {
|
---|
1438 | if(NBins()<=0)
|
---|
1439 | throw(SzMismatchError("Histo::FitResidus: size mismatch\n"));
|
---|
1440 | GeneralFunction* f = gfit.GetFunction();
|
---|
1441 | if(f==NULL)
|
---|
1442 | throw(NullPtrError("Histo::FitResidus: NULL pointer\n"));
|
---|
1443 | Vector par = gfit.GetParm();
|
---|
1444 | Histo h(*this);
|
---|
1445 | for(int i=0;i<NBins();i++) {
|
---|
1446 | double x = (double) BinCenter(i);
|
---|
1447 | h(i) -= (float) f->Value(&x,par.Data());
|
---|
1448 | }
|
---|
1449 | return h;
|
---|
1450 | }
|
---|
1451 |
|
---|
1452 | //++
|
---|
1453 | Histo Histo::FitFunction(GeneralFit& gfit)
|
---|
1454 | //
|
---|
1455 | // Retourne une classe contenant la fonction du fit ``gfit''.
|
---|
1456 | //--
|
---|
1457 | {
|
---|
1458 | if(NBins()<=0)
|
---|
1459 | throw(SzMismatchError("Histo::FitFunction: size mismatch\n"));
|
---|
1460 | GeneralFunction* f = gfit.GetFunction();
|
---|
1461 | if(f==NULL)
|
---|
1462 | throw(NullPtrError("Histo::FitFunction: NULL pointer\n"));
|
---|
1463 | Vector par = gfit.GetParm();
|
---|
1464 | Histo h(*this);
|
---|
1465 | for(int i=0;i<NBins();i++) {
|
---|
1466 | double x = (double) BinCenter(i);
|
---|
1467 | h(i) = (float) f->Value(&x,par.Data());
|
---|
1468 | }
|
---|
1469 | return h;
|
---|
1470 | }
|
---|
1471 |
|
---|
1472 | /********* Methode *********/
|
---|
1473 | //++
|
---|
1474 | void Histo::PrintF(FILE * fp, int hdyn,float hmin, float hmax,int pflag,
|
---|
1475 | int il, int ih)
|
---|
1476 | //
|
---|
1477 | // Impression de l'histogramme dans le fichier fp
|
---|
1478 | //| hdyn = nombre de colonnes pour imprimer les etoiles
|
---|
1479 | //| si =0 alors defaut(100),
|
---|
1480 | //| si <0 pas de print histo seulement infos
|
---|
1481 | //| hmin = minimum de la dynamique
|
---|
1482 | //| hmax = maximum de la dynamique
|
---|
1483 | //| si hmax<=hmin : hmin=VMin() hmax=VMax()
|
---|
1484 | //| si hmax<=hmin et hmin=0 : hmin=0 hmax=VMax()
|
---|
1485 | //| sinon : hmin hmax
|
---|
1486 | //| pflag < 0 : on imprime les informations (nbin,min,...) sans l'histogramme
|
---|
1487 | //| = 0 : on imprime "BinCenter(i) data[i]" (note "... ...")
|
---|
1488 | //| bit 0 on : (v=1): numero_bin "... ..."
|
---|
1489 | //| bit 1 on : (v=2): "... ..." erreur
|
---|
1490 | //--
|
---|
1491 | {
|
---|
1492 |
|
---|
1493 |
|
---|
1494 | if( il > ih ) { il = 0; ih = bins-1; }
|
---|
1495 | if( il < 0 ) il = 0;
|
---|
1496 | if( ih >= bins ) ih = bins-1;
|
---|
1497 |
|
---|
1498 | double dhmin = (double) hmin;
|
---|
1499 | double dhmax = (double) hmax;
|
---|
1500 | double hb,hbmn,hbmx;
|
---|
1501 |
|
---|
1502 | if(hdyn==0) hdyn = 100;
|
---|
1503 |
|
---|
1504 | cout << "~Histo::Print "
|
---|
1505 | << " nHist=" << nHist << " nEntries=" << nEntries
|
---|
1506 | << " under=" << under << " over=" << over << endl;
|
---|
1507 | cout << " bins=" << bins
|
---|
1508 | << " min=" << min << " max=" << max
|
---|
1509 | << " binWidth=" << binWidth << endl;
|
---|
1510 | cout << " mean=" << Mean() << " r.m.s=" << Sigma() << endl;
|
---|
1511 |
|
---|
1512 | if(hdyn<0 || pflag<0 ) return;
|
---|
1513 |
|
---|
1514 | if(dhmax<=dhmin) { if(hmin != 0.) dhmin = (double) VMin(); else dhmin=0.;
|
---|
1515 | dhmax = (double) VMax(); }
|
---|
1516 | if(dhmin>dhmax) return;
|
---|
1517 | if(dhmin==dhmax) {dhmin -= 1.; dhmax += 1.;}
|
---|
1518 | double wb = (dhmax-dhmin) / (double) hdyn;
|
---|
1519 |
|
---|
1520 | // determination de la position de la valeur zero
|
---|
1521 | int i0 = (int)(-dhmin/wb);
|
---|
1522 |
|
---|
1523 | // si le zero est dans la dynamique,
|
---|
1524 | // il doit imperativement etre une limite de bin
|
---|
1525 | if( 0 <= i0 && i0 < hdyn ) {
|
---|
1526 | hbmn = dhmin + i0*wb;
|
---|
1527 | hbmx = hbmn + wb;
|
---|
1528 | if( hbmn != 0. ) {
|
---|
1529 | hbmn *= -1.;
|
---|
1530 | if( hbmn < hbmx ) {
|
---|
1531 | // le zero est + pres du bord negatif du bin
|
---|
1532 | dhmin += hbmn;
|
---|
1533 | dhmax += hbmn;
|
---|
1534 | } else {
|
---|
1535 | // le zero est + pres du bord positif du bin
|
---|
1536 | dhmin -= hbmx;
|
---|
1537 | dhmax -= hbmx;
|
---|
1538 | }
|
---|
1539 | wb = (dhmax-dhmin) / (double) hdyn;
|
---|
1540 | i0 = (int)(-dhmin/wb);
|
---|
1541 | }
|
---|
1542 | }
|
---|
1543 |
|
---|
1544 | cout <<" bin minimum="<<dhmin<<" bin maximum="<<dhmax
|
---|
1545 | <<" binw="<<wb<< endl;
|
---|
1546 |
|
---|
1547 | char* s = new char[hdyn+1];
|
---|
1548 | s[hdyn] = '\0';
|
---|
1549 |
|
---|
1550 | // premiere ligne
|
---|
1551 | {for(int i=0;i<hdyn;i++) s[i] = '=';}
|
---|
1552 | if( 0 <= i0 && i0 < hdyn ) s[i0] = '0';
|
---|
1553 | if(pflag&1) fprintf( fp, "====");
|
---|
1554 | fprintf( fp, "======================");
|
---|
1555 | if(pflag&2 && err2!=NULL) fprintf( fp, "===========");
|
---|
1556 | fprintf( fp, "==%s\n",s);
|
---|
1557 |
|
---|
1558 | // histogramme
|
---|
1559 | {for(int i=il;i<=ih;i++) {
|
---|
1560 | for(int k=0;k<hdyn;k++) s[k] = ' ';
|
---|
1561 | hb = (*this)(i);
|
---|
1562 |
|
---|
1563 | //choix du bin (hdyn bin entre [dhmin,dhmax[
|
---|
1564 | int ii = (int)( (hb-dhmin)/wb );
|
---|
1565 | if(ii<0) ii = 0; else if(ii>=hdyn) ii = hdyn-1;
|
---|
1566 |
|
---|
1567 | // limite du bin
|
---|
1568 | hbmn = dhmin + ii*wb;
|
---|
1569 | hbmx = hbmn + wb;
|
---|
1570 |
|
---|
1571 | // remplissage de s[] en tenant compte des courbes +/-
|
---|
1572 | if(i0<0) {
|
---|
1573 | // courbe entierement positive
|
---|
1574 | for(int k=0;k<=ii;k++) s[k] = 'X';
|
---|
1575 | } else if(i0>=hdyn) {
|
---|
1576 | // courbe entierement negative
|
---|
1577 | for(int k=hdyn-1;k>=ii;k--) s[k] = 'X';
|
---|
1578 | } else {
|
---|
1579 | // courbe positive et negative
|
---|
1580 | s[i0] = '|';
|
---|
1581 | if(ii>i0) for(int k=i0+1;k<=ii;k++) s[k] = 'X';
|
---|
1582 | else if(ii<i0) for(int k=ii;k<i0;k++) s[k] = 'X';
|
---|
1583 | else s[ii] = 'X';
|
---|
1584 | }
|
---|
1585 |
|
---|
1586 | // valeur a mettre dans le bin le plus haut/bas
|
---|
1587 | int ib;
|
---|
1588 | if(hb>0.) ib = (int)( 10.*(hb-hbmn)/(hbmx-hbmn) );
|
---|
1589 | else if(hb<0.) ib = (int)( 10.*(hbmx-hb)/(hbmx-hbmn) );
|
---|
1590 | else ib = -1;
|
---|
1591 | if(ib==-1) s[ii] = '|';
|
---|
1592 | else if(ib==0) s[ii] = '.';
|
---|
1593 | else if(ib>0 && ib<10) s[ii] = (char)((int) '0' + ib);
|
---|
1594 |
|
---|
1595 | // traitement des saturations +/-
|
---|
1596 | if( hb < dhmin ) s[0] = '*';
|
---|
1597 | else if( hb > dhmax ) s[hdyn-1] = '*';
|
---|
1598 |
|
---|
1599 | if(pflag&1) fprintf( fp, "%3d ",i);
|
---|
1600 | fprintf( fp, "%10.4g %10.4g ",BinCenter(i),hb);
|
---|
1601 | if(pflag&2 && err2!=NULL) fprintf( fp, "%10.4g ",Error(i));
|
---|
1602 | fprintf( fp, "= %s\n",s);
|
---|
1603 | }}
|
---|
1604 |
|
---|
1605 | // derniere ligne
|
---|
1606 | for(int i=0;i<hdyn;i++) s[i] = '=';
|
---|
1607 | if( 0 <= i0 && i0 < hdyn ) s[i0] = '0';
|
---|
1608 | if(pflag&1) fprintf( fp, "====");
|
---|
1609 | fprintf( fp, "======================");
|
---|
1610 | if(pflag&2 && err2!=NULL) fprintf( fp, "===========");
|
---|
1611 | fprintf( fp, "==%s\n",s);
|
---|
1612 |
|
---|
1613 | // valeur basse des bins (sur ["ndig-1" digits + signe] = ndig char (>=3))
|
---|
1614 | const int ndig = 7;
|
---|
1615 | char sn[2*ndig+10];
|
---|
1616 | hb = ( fabs(dhmin) > fabs(dhmax) ) ? fabs(dhmin) : fabs(dhmax);
|
---|
1617 | int n;
|
---|
1618 | if( hb > 0. ) n = (int)(log10(hb*1.00000001)); else n = 1;
|
---|
1619 | double scaledig = pow(10.,(double) ndig-2);
|
---|
1620 | double expo = scaledig/pow(10.,(double) n);
|
---|
1621 | // cout <<"n="<<n<<" hb="<<hb<<" scaledig="<<scaledig<<" expo="<<expo<<endl;
|
---|
1622 | for(int k=0;k<ndig;k++) {
|
---|
1623 | for(int i=0;i<hdyn;i++) s[i] = ' ';
|
---|
1624 | {for(int i=0;i<hdyn;i++) {
|
---|
1625 | n = (int)( (dhmin + i*wb)*expo );
|
---|
1626 | for(int j=0;j<2*ndig+10;j++) sn[j] = ' ';
|
---|
1627 | sprintf(sn,"%d%c",n,'\0');
|
---|
1628 | strip(sn,'B',' ');
|
---|
1629 | // cout <<"n="<<n<<" sn=("<<sn<<") l="<<strlen(sn)<<" k="<<k;
|
---|
1630 | if( (int) strlen(sn) > k ) s[i] = sn[k];
|
---|
1631 | // cout <<" s=("<<s<<")"<<endl;
|
---|
1632 | }}
|
---|
1633 | if(pflag&1) fprintf( fp, " ");
|
---|
1634 | fprintf( fp, " ");
|
---|
1635 | if(pflag&2 && err2!=NULL) fprintf( fp, " ");
|
---|
1636 | fprintf( fp, " %s\n",s);
|
---|
1637 | }
|
---|
1638 | fprintf( fp, " (valeurs a multiplier par %.0e)\n",1./expo);
|
---|
1639 |
|
---|
1640 | delete[] s;
|
---|
1641 | }
|
---|
1642 |
|
---|
1643 | /********* Methode *********/
|
---|
1644 | //++
|
---|
1645 | void Histo::Print(int hdyn,float hmin, float hmax,int pflag,
|
---|
1646 | int il, int ih)
|
---|
1647 | //
|
---|
1648 | // Impression de l'histogramme sur stdout
|
---|
1649 | //--
|
---|
1650 | {
|
---|
1651 | Histo::PrintF(stdout, hdyn, hmin, hmax, pflag, il, ih);
|
---|
1652 | }
|
---|
1653 |
|
---|
1654 |
|
---|
1655 | // Rappel des inlines functions pour commentaires
|
---|
1656 | //++
|
---|
1657 | // inline float XMin() const
|
---|
1658 | // Retourne l'abscisse minimum
|
---|
1659 | //--
|
---|
1660 | //++
|
---|
1661 | // inline float XMax() const
|
---|
1662 | // Retourne l'abscisse maximum
|
---|
1663 | //--
|
---|
1664 | //++
|
---|
1665 | // inline int NBins() const
|
---|
1666 | // Retourne le nombre de bins
|
---|
1667 | //--
|
---|
1668 | //++
|
---|
1669 | // inline float BinWidth() const
|
---|
1670 | // Retourne la largeur du bin
|
---|
1671 | //--
|
---|
1672 | //++
|
---|
1673 | // inline float* Bins() const
|
---|
1674 | // Retourne le pointeur sur le tableaux des contenus
|
---|
1675 | //--
|
---|
1676 | //++
|
---|
1677 | // inline float operator()(int i) const
|
---|
1678 | // Retourne le contenu du bin i
|
---|
1679 | //--
|
---|
1680 | //++
|
---|
1681 | // inline float& operator()(int i)
|
---|
1682 | // Remplit le contenu du bin i
|
---|
1683 | //--
|
---|
1684 | //++
|
---|
1685 | // inline float Error(int i) const
|
---|
1686 | // Retourne l'erreur du bin i
|
---|
1687 | //--
|
---|
1688 | //++
|
---|
1689 | // inline double Error2(int i) const
|
---|
1690 | // Retourne l'erreur au carre du bin i
|
---|
1691 | //--
|
---|
1692 | //++
|
---|
1693 | // inline double& Error2(int i)
|
---|
1694 | // Remplit l'erreur au carre du bin i
|
---|
1695 | //--
|
---|
1696 | //++
|
---|
1697 | // inline float NData() const
|
---|
1698 | // Retourne la somme ponderee
|
---|
1699 | //--
|
---|
1700 | //++
|
---|
1701 | // inline float NEntries()
|
---|
1702 | // Retourne le nombre d'entrees
|
---|
1703 | //--
|
---|
1704 | //++
|
---|
1705 | // inline float NOver() const
|
---|
1706 | // Retourne le nombre d'overflow
|
---|
1707 | //--
|
---|
1708 | //++
|
---|
1709 | // inline float NUnder() const
|
---|
1710 | // Retourne le nombre d'underflow
|
---|
1711 | //--
|
---|
1712 | //++
|
---|
1713 | // inline float BinLowEdge(int i) const
|
---|
1714 | // Retourne l'abscisse du bord inferieur du bin i
|
---|
1715 | //--
|
---|
1716 | //++
|
---|
1717 | // inline float BinCenter(int i) const
|
---|
1718 | // Retourne l'abscisse du centre du bin i
|
---|
1719 | //--
|
---|
1720 | //++
|
---|
1721 | // inline float BinHighEdge(int i) const
|
---|
1722 | // Retourne l'abscisse du bord superieur du bin i
|
---|
1723 | //--
|
---|
1724 | //++
|
---|
1725 | // inline int FindBin(float x) const
|
---|
1726 | // Retourne le numero du bin contenant l'abscisse x
|
---|
1727 | //--
|
---|
1728 |
|
---|
1729 |
|
---|
1730 | ///////////////////////////////////////////////////////////
|
---|
1731 | // --------------------------------------------------------
|
---|
1732 | // Les objets delegues pour la gestion de persistance
|
---|
1733 | // --------------------------------------------------------
|
---|
1734 | ///////////////////////////////////////////////////////////
|
---|
1735 |
|
---|
1736 | void ObjFileIO<Histo>::ReadSelf(PInPersist& is)
|
---|
1737 | {
|
---|
1738 | char strg[256];
|
---|
1739 |
|
---|
1740 | if(dobj==NULL) dobj=new Histo;
|
---|
1741 | else dobj->Delete();
|
---|
1742 |
|
---|
1743 | // Lecture entete
|
---|
1744 | is.GetLine(strg, 255);
|
---|
1745 | is.GetLine(strg, 255);
|
---|
1746 | is.GetLine(strg, 255);
|
---|
1747 |
|
---|
1748 | // Lecture des valeurs
|
---|
1749 | is.Get(dobj->bins);
|
---|
1750 | is.Get(dobj->nEntries);
|
---|
1751 | int_4 errok;
|
---|
1752 | is.Get(errok);
|
---|
1753 |
|
---|
1754 | is.Get(dobj->binWidth);
|
---|
1755 | is.Get(dobj->min);
|
---|
1756 | is.Get(dobj->max);
|
---|
1757 | is.Get(dobj->under);
|
---|
1758 | is.Get(dobj->over);
|
---|
1759 |
|
---|
1760 | is.Get(dobj->nHist);
|
---|
1761 |
|
---|
1762 | // Lecture des donnees Histo 1D
|
---|
1763 | dobj->data = new float[dobj->bins];
|
---|
1764 | is.GetLine(strg, 255);
|
---|
1765 | is.Get(dobj->data, dobj->bins);
|
---|
1766 |
|
---|
1767 | // Lecture des erreurs
|
---|
1768 | if(errok) {
|
---|
1769 | is.GetLine(strg, 255);
|
---|
1770 | dobj->err2 = new double[dobj->bins];
|
---|
1771 | is.Get(dobj->err2, dobj->bins);
|
---|
1772 | }
|
---|
1773 |
|
---|
1774 | return;
|
---|
1775 | }
|
---|
1776 |
|
---|
1777 | void ObjFileIO<Histo>::WriteSelf(POutPersist& os) const
|
---|
1778 | {
|
---|
1779 | if (dobj == NULL) return;
|
---|
1780 | char strg[256];
|
---|
1781 |
|
---|
1782 | // Que faut-il ecrire?
|
---|
1783 | int_4 errok = (dobj->err2) ? 1 : 0;
|
---|
1784 |
|
---|
1785 | // Ecriture entete pour identifier facilement
|
---|
1786 | sprintf(strg,"bins=%6d NEnt=%15d errok=%1d",dobj->bins,dobj->nEntries,errok);
|
---|
1787 | os.PutLine(strg);
|
---|
1788 | sprintf(strg,"binw=%g min=%g max=%g",dobj->binWidth,dobj->min,dobj->max);
|
---|
1789 | os.PutLine(strg);
|
---|
1790 | sprintf(strg, "under=%g over=%g nHist=%g",dobj->under,dobj->over,dobj->nHist);
|
---|
1791 | os.PutLine(strg);
|
---|
1792 |
|
---|
1793 | // Ecriture des valeurs
|
---|
1794 | os.Put(dobj->bins);
|
---|
1795 | os.Put(dobj->nEntries);
|
---|
1796 | os.Put(errok);
|
---|
1797 |
|
---|
1798 | os.Put(dobj->binWidth);
|
---|
1799 | os.Put(dobj->min);
|
---|
1800 | os.Put(dobj->max);
|
---|
1801 | os.Put(dobj->under);
|
---|
1802 | os.Put(dobj->over);
|
---|
1803 |
|
---|
1804 | os.Put(dobj->nHist);
|
---|
1805 |
|
---|
1806 | // Ecriture des donnees Histo 1D
|
---|
1807 | sprintf(strg,"Histo: Tableau des donnees %d",dobj->bins);
|
---|
1808 | os.PutLine(strg);
|
---|
1809 | os.Put(dobj->data, dobj->bins);
|
---|
1810 |
|
---|
1811 | // Ecriture des erreurs
|
---|
1812 | if(errok) {
|
---|
1813 | sprintf(strg,"Histo: Tableau des erreurs %d",dobj->bins);
|
---|
1814 | os.PutLine(strg);
|
---|
1815 | os.Put(dobj->err2, dobj->bins);
|
---|
1816 | }
|
---|
1817 |
|
---|
1818 | return;
|
---|
1819 | }
|
---|
1820 |
|
---|
1821 | #ifdef __CXX_PRAGMA_TEMPLATES__
|
---|
1822 | #pragma define_template ObjFileIO<Histo>
|
---|
1823 | #endif
|
---|
1824 |
|
---|
1825 | #if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES)
|
---|
1826 | template class ObjFileIO<Histo>;
|
---|
1827 | #endif
|
---|