source: Sophya/trunk/SophyaLib/HiStats/histos.cc@ 1155

Last change on this file since 1155 was 1144, checked in by ansari, 25 years ago

passage dans histo.cc d'argument int en int_4 pour compatibilite mac.

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