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