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