#include "machdefs.h" #include #include #include #include "histos.h" #include "perrors.h" #include "poly.h" #include "strutil.h" namespace SOPHYA { /*! \class Histo \ingroup HiStats Classe d'histogrammes 1D */ /********* Methode *********/ /*! Constructeur par defaut */ Histo::Histo() : mData(NULL), mErr2(NULL) , mUnder(0.), mOver(0.), nHist(0.), nEntries(0) , mBins(0), mMin(0.), mMax(0.), binWidth(0.) { } /********* Methode *********/ /*! Constructeur d'un histo de nBin bins allant de xMin a xMax */ Histo::Histo(r_8 xMin, r_8 xMax, int_4 nBin) : mData(NULL), mErr2(NULL) { CreateOrResize(xMin,xMax,nBin); } /********* Methode *********/ /*! Constructeur d'un histo de nBin bins allant de xMin a xMax */ Histo::Histo(r_4 xMin, r_4 xMax, int_4 nBin) : mData(NULL), mErr2(NULL) { CreateOrResize((r_8)xMin,(r_8)xMax,nBin); } /********* Methode *********/ /*! Constructeur par copie */ Histo::Histo(const Histo& H) : mData((H.mBins>0)? new r_8[H.mBins] : NULL), mErr2((H.mBins>0 && H.mErr2!=NULL)? new r_8[H.mBins]: NULL), mUnder(H.mUnder), mOver(H.mOver), nHist(H.nHist), nEntries(H.nEntries), mBins(H.mBins), mMin(H.mMin), mMax(H.mMax), binWidth(H.binWidth) { if(mBins<=0) return; memcpy(mData,H.mData,mBins*sizeof(r_8)); if(H.mErr2) memcpy(mErr2, H.mErr2, mBins*sizeof(r_8)); } /********* Methode *********/ /*! Destructeur */ Histo::~Histo() { Delete(); } /********* Methode *********/ /*! Redimensionnement d'un histo */ void Histo::ReSize(r_8 xMin, r_8 xMax, int_4 nBin) { CreateOrResize(xMin,xMax,nBin); } /********* Methode *********/ /*! Gestion de l'allocation */ void Histo::CreateOrResize(r_8 xMin, r_8 xMax, int_4 nBin) { //cout<<"Histo::CreateOrResize()"<0 && mData==NULL) mData = new r_8[nBin]; if(mData) memset(mData,0,nBin*sizeof(r_8)); mBins = nBin; mMin = xMin; mMax = xMax; binWidth = (nBin>0) ? (mMax-mMin)/nBin : 0.; mUnder = mOver = nHist = 0.; nEntries = 0; } /********* Methode *********/ /*! Gestion de la des-allocation */ void Histo::Delete() { //cout<<"Histo::Delete()"< &v) const { v.Realloc(mBins); for(int_4 i=0;i &v) const { v.Realloc(mBins); for(int_4 i=0;i &v) const { v.Realloc(mBins); if(!mErr2) {for(int_4 i=0;i &v) const { v.Realloc(mBins); if(!mErr2) {for(int_4 i=0;i &v, int_4 ierr) { //if(v.NElts()<(uint_4) mBins) throw SzMismatchError(PExcLongMessage("")); uint_4 n = (v.NElts()<(uint_4) mBins) ? v.NElts(): (uint_4) mBins; if(n>0) for(uint_4 i=0;i &v, int_4 ierr) { //if(v.NElts()<(uint_4) mBins) throw SzMismatchError(PExcLongMessage("")); uint_4 n = (v.NElts()<(uint_4) mBins) ? v.NElts(): (uint_4) mBins; if(n>0) for(uint_4 i=0;i &v) { //if(v.NElts()<(uint_4) mBins) throw SzMismatchError(PExcLongMessage("")); uint_4 n = (v.NElts()<(uint_4) mBins) ? v.NElts(): (uint_4) mBins; if(n>0) { if(!mErr2) Errors(); for(uint_4 i=0;i &v) { //if(v.NElts()<(uint_4) mBins) throw SzMismatchError(PExcLongMessage("")); uint_4 n = (v.NElts()<(uint_4) mBins) ? v.NElts(): (uint_4) mBins; if(n>0) { if(!mErr2) Errors(); for(uint_4 i=0;i0.) mErr2[i] += v(i); } return; } /*! Remplissage des erreurs de l'histo avec les valeurs d'un vecteur */ void Histo::PutError(TVector &v) { //if(v.NElts()<(uint_4) mBins) throw SzMismatchError(PExcLongMessage("")); uint_4 n = (v.NElts()<(uint_4) mBins) ? v.NElts(): (uint_4) mBins; if(n>0) { if(!mErr2) Errors(); for(uint_4 i=0;i0.) mErr2[i]=v(i)*v(i); else mErr2[i]=-v(i)*v(i); } return; } /********* Methode *********/ /*! Addition du contenu de l'histo pour abscisse x poids w */ void Histo::Add(r_8 x, r_8 w) { int_4 numBin = FindBin(x); if (numBin<0) mUnder += w; else if (numBin>=mBins) mOver += w; else { mData[numBin] += w; if(mErr2!=NULL) mErr2[numBin] += w*w; nHist += w; nEntries++; } } /********* Methode *********/ /*! Addition du contenu de l'histo bin numBin poids w */ void Histo::AddBin(int_4 numBin, r_8 w) { if (numBin<0) mUnder += w; else if (numBin>=mBins) mOver += w; else { mData[numBin] += w; if(mErr2!=NULL) mErr2[numBin] += w*w; nHist += w; nEntries++; } } /********* Methode *********/ /*! Remplissage du contenu de l'histo pour abscisse x poids w */ void Histo::SetBin(r_8 x, r_8 w) { int_4 numBin = FindBin(x); SetBin(numBin,w); } /********* Methode *********/ /*! Remplissage du contenu de l'histo pour numBin poids w */ void Histo::SetBin(int_4 numBin, r_8 w) { if (numBin<0) mUnder = w; else if (numBin>=mBins) mOver = w; else { nHist -= mData[numBin]; mData[numBin] = w; nHist += w; } } /********* Methode *********/ /*! Remplissage des erreurs au carre pour abscisse x */ void Histo::SetErr2(r_8 x, r_8 e2) { int_4 numBin = FindBin(x); SetErr2(numBin,e2); } /********* Methode *********/ /*! Remplissage des erreurs au carre pour numBin poids */ void Histo::SetErr2(int_4 numBin, r_8 e2) { if( mErr2==NULL) return; if ( numBin<0 || numBin>=mBins ) return; mErr2[numBin] = e2; } /********* Methode *********/ /*! Remplissage des erreurs pour abscisse x */ void Histo::SetErr(r_8 x, r_8 e) { SetErr2(x, e*e); } /********* Methode *********/ /*! Remplissage des erreurs pour numBin */ void Histo::SetErr(int_4 numBin, r_8 e) { SetErr2(numBin, e*e); } /********* Methode *********/ /*! Methode virtuelle de mise a jour - Ne fait rien pour Histo - Voir HProf */ void Histo::UpdateHisto(bool force) const { return; } /********* Methode *********/ /*! Numero du bin ayant le contenu maximum */ int_4 Histo::IMax() const { if(mBins<=0) return 0; int_4 imx=0; if(mBins==1) return imx; r_8 mx=mData[0]; for (int_4 i=1; imx) {imx = i; mx=mData[i];} return imx; } /********* Methode *********/ /*! Numero du bin ayant le contenu minimum */ int_4 Histo::IMin() const { if(mBins<=0) return 0; int_4 imx=0; if(mBins==1) return imx; r_8 mx=mData[0]; for (int_4 i=1; imx) mx=mData[i]; return mx; } /********* Methode *********/ /*! Valeur le contenu minimum */ r_8 Histo::VMin() const { if(mBins<=0) return 0.; r_8 mx=mData[0]; if(mBins==1) return mx; for (int_4 i=1;i=0.) ? mData[i] : -mData[i]; n += v; sx += BinCenter(i)*v; } if(n>0.) return sx/n; else return mMin; } /********* Methode *********/ /*! Valeur du sigma */ r_8 Histo::Sigma() const { if(mBins<=0) return 0.; r_8 n = 0; r_8 sx = 0; r_8 sx2= 0; for (int_4 i=0; i=0.) ? mData[i] : -mData[i]; n += v; sx += BinCenter(i)*v; sx2+= BinCenter(i)*BinCenter(i)*v; } // attention a l'erreur d'arrondi si un seul bin rempli!! if( n == 0 ) return 0.; sx2 = sx2/n - (sx/n)*(sx/n); if( sx2>0. ) return sqrt( sx2 ); else return 0.; } /********* Methode *********/ /*! Valeur de la moyenne calculee entre les bin il et ih */ r_8 Histo::MeanLH(int_4 il,int_4 ih) const { if(mBins<=0) return 0.; if( ih < il ) { il = 0; ih = mBins-1; } if( il < 0 ) il = 0; if( ih >= mBins ) ih = mBins-1; r_8 n = 0; r_8 sx = 0; for (int_4 i=il; i<=ih; i++) { r_8 v = (mData[i]>=0.) ? mData[i] : -mData[i]; n += v; sx += BinCenter(i)*v; } if(n>0.) return sx/n; else return BinLowEdge(il); } /********* Methode *********/ /*! Valeur de la moyenne calculee entre les bin il et ih */ r_8 Histo::SigmaLH(int_4 il,int_4 ih) const { if(mBins<=0) return 0.; if( ih < il ) { il = 0; ih = mBins - 1; } if( il < 0 ) il = 0; if( ih >= mBins ) ih = mBins - 1; r_8 n = 0; r_8 sx = 0; r_8 sx2= 0; for (int_4 i=il; i<=ih; i++) { r_8 v = (mData[i]>=0.) ? mData[i] : -mData[i]; n += v; sx += BinCenter(i)*v; sx2+= BinCenter(i)*BinCenter(i)*v; } if( n == 0 ) return 0.; sx2 = sx2/n - (sx/n)*(sx/n); if( sx2>0. ) return sqrt( sx2 ); else return 0.; } /********* Methode *********/ /*! Valeur de la moyenne calculee entre x0-dx et x0+dx */ r_8 Histo::Mean(r_8 x0, r_8 dx) const { if(mBins<=0) return 0.; r_8 sdata = 0; r_8 sx = 0; int_4 iMin = FindBin(x0-dx); if (iMin<0) iMin = 0; int_4 iMax = FindBin(x0+dx); if (iMax>mBins) iMax=mBins; for (int_4 i=iMin; i=0.) ? mData[i] : -mData[i]; sx += BinCenter(i)*v; sdata += v; } if(sdata>0.) return sx/sdata; else return mMin; } /********* Methode *********/ /*! Valeur du sigma calcule entre x0-dx et x0+dx */ r_8 Histo::Sigma(r_8 x0, r_8 dx) const { if(mBins<=0) return 0.; r_8 sx = 0; r_8 sx2= 0; r_8 sdata = 0; int_4 iMin = FindBin(x0-dx); if (iMin<0) iMin = 0; int_4 iMax = FindBin(x0+dx); if (iMax>mBins) iMax=mBins; for (int_4 i=iMin; i=0.) ? mData[i] : -mData[i]; sx += BinCenter(i)*v; sx2+= BinCenter(i)*BinCenter(i)*v; sdata += v; } if(sdata>0.) return sqrt( sx2/sdata - (sx/sdata)*(sx/sdata) ); else return 0.; } /********* Methode *********/ /*! Valeur de la moyenne et du sigma nettoyes */ r_8 Histo::CleanedMean(r_8& sigma) const { if(mBins<=0) return 0.; if (!nHist) return 0; // on fait ca de facon bete, on pourra raffiner apres... r_8 x0 = Mean(); r_8 s = Sigma()+binWidth; for (int_4 i=0; i<3; i++) { r_8 xx0 = Mean(x0,5*s); s = Sigma(x0,5*s)+binWidth; x0 = xx0; } sigma = s; return x0; } /********* Methode *********/ /*! Valeur de la moyenne nettoyee */ r_8 Histo::CleanedMean() const { if(mBins<=0) return 0.; if (!nHist) return 0; r_8 s=0; return CleanedMean(s); } /********* Methode *********/ /*! Retourne le nombre de bins non-nul */ int_4 Histo::BinNonNul() const { if(mBins<=0) return -1; int_4 non=0; for (int_4 i=0;i= n ) break; } if( i == mBins ) i = mBins-1; return i; } /********* Methode *********/ /*! Renvoie les numeros de bins imin,imax (0==mBins ) return -2; r_8 N1 = 0.; for(int_4 i=0; i<=I; i++) N1 += mData[i]; N1 *= per; r_8 N2 = 0.; {for(int_4 i=I; i=0; imin--) { n += mData[imin]; if(n>=N1) break; } if( imin<0 ) imin = 0; // cout<<"I="<1) cout<<"Max histo i="<= degree+1 ) break; iLow--; iHigh++; if(iLow<0 && iHigh>=NBins()) { if (debug>1) cout<<"Mode histogramme = "<= NBins()) iHigh = NBins()-1; } TVector xFit(nLowHigh); TVector yFit(nLowHigh); TVector e2Fit(nLowHigh); TVector errcoef(1); int_4 ii = 0; for (int_4 j=iLow; j<=iHigh; j++) { if ((*this)(j)>0) { if(!mErr2) { xFit(ii) = BinCenter(j)-xCenter; yFit(ii) = (*this)(j); e2Fit(ii) = yFit(ii); ii++; } else if(Error2(j)>0.) { xFit(ii) = BinCenter(j)-xCenter; yFit(ii) = (*this)(j); e2Fit(ii) = Error2(j); ii++; } } } if(debug>4) { int_4 k; for(k=0;k1) cout << "resultat fit = " << pol << endl; pol.Derivate(); int_4 DPolDeg = pol.Degre(); r_8 fd = 0.; if (DPolDeg == 0) { // on est dans le cas d'un fit de droite if( pol[0] > 0. ) { // on a fitte une droite de pente >0. fd = xFit(nLowHigh-1) + binWidth/2. + xCenter; } else if( pol[0] < 0. ) { // on a fitte une droite de pente <0. fd = xFit(0) - binWidth/2. + xCenter; } else { // on a fitte une droite de pente =0. (creneau) fd = (xFit(0)+xFit(nLowHigh-1))/2. + xCenter; } } else if (DPolDeg == 1) { // on est dans le cas d'un fit de parabole r_8 r=0; if (pol.Root1(r)==0) throw ParmError(PExcLongMessage("")); fd = r + xCenter; } else if (DPolDeg == 2) { // on est dans le cas d'un fit de cubique r_8 r1=0; r_8 r2=0; if (pol.Root2(r1,r2) == 0) throw ParmError(PExcLongMessage("")); pol.Derivate(); fd = (pol(r1)<0) ? r1 + xCenter : r2 + xCenter; } else { // on est dans un cas non prevu throw ParmError(PExcLongMessage("")); } if(fd>mMax) fd = mMax; if(fd1) cout << "Mode histogramme = " << fd << " (DerPol.degre =" << DPolDeg << " pol.coeff[0] =" << pol[0] << ")" << endl; return fd; } /********* Methode *********/ /*! Calcul de la largeur a frac pourcent du maximum autour du bin du maximum. L'histo est suppose etre remplit de valeurs positives */ r_8 Histo::FindWidth(r_8 frac, int_4 debug) const { if(mBins<=0) return 0; r_8 xmax = BinCenter( IMax() ); return FindWidth(xmax,frac,debug); } /********* Methode *********/ /*! Calcul de la largeur a frac pourcent de la valeur du bin xmax. L'histo est suppose etre remplit de valeurs positives */ r_8 Histo::FindWidth(r_8 xmax,r_8 frac, int_4 debug) const { if(mBins<=0) return 0; if (frac <= 0 || frac >= 1.) frac = 0.5; if (debug > 1) cout << "Histo::FindWidth a " << frac << " de xmax= " << xmax << " , ndata= " << NData() << " , nent= " << NEntries() << " , nbin= " << NBins() << endl; if (NEntries() < 1) throw ParmError(PExcLongMessage("")); if (NBins() < 3) throw ParmError(PExcLongMessage("")); int_4 iMax = FindBin(xmax); if (iMax<0 || iMax>=NBins()) throw ParmError(PExcLongMessage("")); r_8 hmax = mData[iMax]; r_8 limit = frac*hmax; if (debug > 1) cout << " Max histo[" << iMax << "] = " << hmax << ", limite " << limit << endl; int_4 iLow = iMax; while (iLow>=0 && mData[iLow]>limit) iLow--; if( iLow < 0 ) iLow = 0; int_4 iHigh = iMax; while (iHighlimit) iHigh++; if( iHigh >=NBins() ) iHigh = NBins()-1; r_8 xlow = BinCenter(iLow); r_8 ylow = mData[iLow]; r_8 xlow1=xlow, ylow1=ylow; if(iLow+1=0) { xhigh1 = BinCenter(iHigh-1); yhigh1 = mData[iHigh-1]; } r_8 xlpred,xhpred,wd; if(ylow1>ylow) { xlpred = xlow + (xlow1-xlow)/(ylow1-ylow)*(limit-ylow); if(xlpred < xlow) xlpred = xlow; } else xlpred = xlow; if(yhigh1>yhigh) { xhpred = xhigh + (xhigh1-xhigh)/(yhigh1-yhigh)*(limit-yhigh); if(xhpred > xhigh) xhpred = xhigh; } else xhpred = xhigh; wd = xhpred - xlpred; if (debug > 1) { cout << "Limites histo: " << " Width " << wd << endl; cout << " Low: [" << iLow << "]=" << ylow << " pred-> " << xlpred << " - High: [" << iHigh << "]=" << yhigh << " pred-> " << xhpred << endl; } return wd; } /********* Methode *********/ /*! Cf suivant mais im est le bin du maximum de l'histo */ int_4 Histo::EstimeMax(r_8& xm,int_4 SzPav) const { int_4 im = IMax(); return EstimeMax(im,xm,SzPav); } /********* Methode *********/ /*! Determine l'abscisses du maximum donne par im en moyennant dans un pave SzPav autour du maximum \verbatim Return: 0 = si fit maximum reussi avec SzPav pixels 1 = si fit maximum reussi avec moins que SzPav pixels 2 = si fit maximum echoue et renvoit BinCenter() -1 = si echec: SzPav <= 0 ou im hors limites \endverbatim */ int_4 Histo::EstimeMax(int_4& im,r_8& xm,int_4 SzPav) const { if(mBins<=0) return -1; xm = 0; if( SzPav <= 0 ) return -1; if( im < 0 || im >= mBins ) return -1; if( SzPav%2 == 0 ) SzPav++; SzPav = (SzPav-1)/2; int_4 rc = 0; r_8 dxm = 0, wx = 0; for(int_4 i=im-SzPav;i<=im+SzPav;i++) { if( i<0 || i>= mBins ) {rc=1; continue;} dxm += BinCenter(i) * (*this)(i); wx += (*this)(i); } if( wx > 0. ) { xm = dxm/wx; return rc; } else { xm = BinCenter(im); return 2; } } /********* Methode *********/ /*! Determine la largeur a frac% du maximum a gauche (widthG) et a droite (widthD) */ void Histo::EstimeWidthS(r_8 frac,r_8& widthG,r_8& widthD) const { if(mBins<=0) {frac=0.; return;} int_4 i; widthG = widthD = -1.; if( mBins<=1 || frac<=0. || frac>=1. ) return; int_4 imax = 0; r_8 maxi = mData[0]; for(i=1;imaxi) {imax=i; maxi=mData[i];} r_8 xmax = BinCenter(imax); r_8 maxf = maxi * frac; // recherche du sigma a gauche widthG = 0.; for(i=imax;i>=0;i--) if( mData[i] <= maxf ) break; if(i<0) i=0; if(i=mBins) i=mBins-1; if(i>imax) { if( mData[i] != mData[i-1] ) { widthD = BinCenter(i) - binWidth * (maxf-mData[i])/(mData[i-1]-mData[i]); widthD -= xmax; } else widthD = BinCenter(i) - xmax; } } /********* Methode *********/ /*! Impression de l'histogramme \verbatim hdyn = nombre de colonnes pour imprimer les etoiles si =0 alors defaut(100), si <0 pas de print histo seulement infos hmin = minimum de la dynamique hmax = maximum de la dynamique si hmax<=hmin : hmin=VMin() hmax=VMax() si hmax<=hmin et hmin=0 : hmin=0 hmax=VMax() sinon : hmin hmax pflag < 0 : on imprime les informations (nbin,min,...) sans l'histogramme = 0 : on imprime "BinCenter(i) mData[i]" (note "... ...") bit 0 on : (v=1): numero_bin "... ..." bit 1 on : (v=2): "... ..." erreur \endverbatim */ void Histo::Show(ostream & os) const { os << " Histo::Show " << " nHist=" << nHist << " nEntries=" << nEntries << " mUnder=" << mUnder << " mOver=" << mOver << endl; os << " mBins=" << mBins << " min=" << mMin << " mMax=" << mMax << " binWidth=" << binWidth << endl; os << " mean=" << Mean() << " r.m.s=" << Sigma() << " Errors="< ih ) { il = 0; ih = mBins-1; } if( il < 0 ) il = 0; if( ih >= mBins ) ih = mBins-1; r_8 dhmin = hmin; r_8 dhmax = hmax; r_8 hb,hbmn,hbmx; if(hdyn==0) hdyn = 100; Show(cout); if(hdyn<0 || pflag<0 ) return; if(dhmax<=dhmin) { if(hmin != 0.) dhmin = VMin(); else dhmin=0.; dhmax = VMax(); } if(dhmin>dhmax) return; if(dhmin==dhmax) {dhmin -= 1.; dhmax += 1.;} r_8 wb = (dhmax-dhmin) / hdyn; // determination de la position de la valeur zero int_4 i0 = (int_4)(-dhmin/wb); // si le zero est dans la dynamique, // il doit imperativement etre une limite de bin if( 0 <= i0 && i0 < hdyn ) { hbmn = dhmin + i0*wb; hbmx = hbmn + wb; if( hbmn != 0. ) { hbmn *= -1.; if( hbmn < hbmx ) { // le zero est + pres du bord negatif du bin dhmin += hbmn; dhmax += hbmn; } else { // le zero est + pres du bord positif du bin dhmin -= hbmx; dhmax -= hbmx; } wb = (dhmax-dhmin) / hdyn; i0 = (int_4)(-dhmin/wb); } } cout <<" bin minimum="<=hdyn) ii = hdyn-1; // limite du bin hbmn = dhmin + ii*wb; hbmx = hbmn + wb; // remplissage de s[] en tenant compte des courbes +/- if(i0<0) { // courbe entierement positive for(int_4 k=0;k<=ii;k++) s[k] = 'X'; } else if(i0>=hdyn) { // courbe entierement negative for(int_4 k=hdyn-1;k>=ii;k--) s[k] = 'X'; } else { // courbe positive et negative s[i0] = '|'; if(ii>i0) for(int_4 k=i0+1;k<=ii;k++) s[k] = 'X'; else if(ii0.) ib = (int_4)( 10.*(hb-hbmn)/(hbmx-hbmn) ); else if(hb<0.) ib = (int_4)( 10.*(hbmx-hb)/(hbmx-hbmn) ); else ib = -1; if(ib==-1) s[ii] = '|'; else if(ib==0) s[ii] = '.'; else if(ib>0 && ib<10) s[ii] = (char)((int_4) '0' + ib); // traitement des saturations +/- if( hb < dhmin ) s[0] = '*'; else if( hb > dhmax ) s[hdyn-1] = '*'; if(pflag&1) printf("%3d ",i); printf("%10.4g %10.4g ",BinCenter(i),hb); if(pflag&2 && mErr2!=NULL) printf("%10.4g ",Error(i)); printf("= %s\n",s); }} // derniere ligne for(int_4 i=0;i=3)) const int_4 ndig = 7; char sn[2*ndig+10]; hb = ( fabs(dhmin) > fabs(dhmax) ) ? fabs(dhmin) : fabs(dhmax); int_4 n; if( hb > 0. ) n = (int_4)(log10(hb*1.00000001)); else n = 1; r_8 scaledig = pow(10.,(r_8) ndig-2); r_8 expo = scaledig/pow(10.,(r_8) n); // cout <<"n="< , pour SGI-CC en particulier */ void ObjFileIO::ReadSelf(PInPersist& is) { char strg[256]; if(dobj==NULL) dobj=new Histo; else dobj->Delete(); // Lecture entete is.GetLine(strg, 255); bool nentries_en_8_bytes=false; if(strg[0]=='V' && strg[1]=='_' && strg[2]=='2') nentries_en_8_bytes=true; is.GetLine(strg, 255); is.GetLine(strg, 255); // Lecture des valeurs is.Get(dobj->mBins); if(nentries_en_8_bytes) is.Get(dobj->nEntries); else {int_4 dum; is.Get(dum); dobj->nEntries=(uint_8)dum;} int_4 errok; is.Get(errok); is.Get(dobj->binWidth); is.Get(dobj->mMin); is.Get(dobj->mMax); is.Get(dobj->mUnder); is.Get(dobj->mOver); is.Get(dobj->nHist); // Lecture des donnees Histo 1D is.GetLine(strg, 255); if(dobj->mBins>0) dobj->mData = new r_8[dobj->mBins]; is.Get(dobj->mData, dobj->mBins); // Lecture des erreurs if(errok) { is.GetLine(strg, 255); if(dobj->mBins>0) dobj->mErr2 = new r_8[dobj->mBins]; is.Get(dobj->mErr2, dobj->mBins); } return; } DECL_TEMP_SPEC /* equivalent a template <> , pour SGI-CC en particulier */ void ObjFileIO::WriteSelf(POutPersist& os) const { if (dobj == NULL) return; char strg[256]; // Que faut-il ecrire? int_4 errok = (dobj->mErr2) ? 1 : 0; // Ecriture entete pour identifier facilement sprintf(strg,"V_2 mBins=%6d NEnt=%15llu errok=%1d",dobj->mBins,dobj->nEntries,errok); os.PutLine(strg); sprintf(strg,"binw=%g mMin=%g mMax=%g",dobj->binWidth,dobj->mMin,dobj->mMax); os.PutLine(strg); sprintf(strg, "mUnder=%g mOver=%g nHist=%g",dobj->mUnder,dobj->mOver,dobj->nHist); os.PutLine(strg); // Ecriture des valeurs os.Put(dobj->mBins); os.Put(dobj->nEntries); os.Put(errok); os.Put(dobj->binWidth); os.Put(dobj->mMin); os.Put(dobj->mMax); os.Put(dobj->mUnder); os.Put(dobj->mOver); os.Put(dobj->nHist); // Ecriture des donnees Histo 1D sprintf(strg,"Histo: Tableau des donnees %d",dobj->mBins); os.PutLine(strg); os.Put(dobj->mData, dobj->mBins); // Ecriture des erreurs if(errok) { sprintf(strg,"Histo: Tableau des erreurs %d",dobj->mBins); os.PutLine(strg); os.Put(dobj->mErr2, dobj->mBins); } return; } #ifdef __CXX_PRAGMA_TEMPLATES__ #pragma define_template ObjFileIO #endif #if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES) template class ObjFileIO; #endif } // FIN namespace SOPHYA