// // $Id: histos.cc,v 1.6 2000-06-30 13:21:52 ansari Exp $ // #include "machdefs.h" #include #include #include #include "histos.h" #include "perrors.h" #include "poly.h" #include "strutil.h" #include "generalfit.h" /*! \class SOPHYA::Histo \ingroup HiStats Classe d'histogrammes 1D */ /********* Methode *********/ /*! Constructeur par defaut */ Histo::Histo() : data(NULL), err2(NULL), under(0), over(0), nHist(0), nEntries(0), bins(0), min(0), max(0), binWidth(0) { END_CONSTRUCTOR } /********* Methode *********/ /*! Constructeur d'un histo de nBin bins allant de xMin a xMax */ Histo::Histo(float xMin, float xMax, int nBin) : data(new float[nBin]), err2(NULL), under(0), over(0), nHist(0), nEntries(0), bins(nBin), min(xMin), max(xMax), binWidth((max - min)/nBin) { Zero(); END_CONSTRUCTOR } /********* Methode *********/ /*! Constructeur par copie */ Histo::Histo(const Histo& H) : data(new float[H.bins]), err2(NULL), under(H.under), over(H.over), nHist(H.nHist), nEntries(H.nEntries), bins(H.bins), min(H.min), max(H.max), binWidth(H.binWidth) { memcpy(data, H.data, bins*sizeof(float)); if( H.err2 != NULL ) { err2 = new double[bins]; memcpy(err2, H.err2, bins*sizeof(double)); } END_CONSTRUCTOR } /********* Methode *********/ /*! Gestion de la des-allocation */ void Histo::Delete() { if( data != NULL ) { delete[] data; data = NULL;} if( err2 != NULL ) { delete[] err2; err2 = NULL;} under = over = min = max = binWidth= 0.; nHist = 0.; nEntries = bins = 0; } /********* Methode *********/ /*! Destructeur */ Histo::~Histo() { Delete(); } /********* Methode *********/ /*! Remise a zero */ void Histo::Zero() { memset(data, 0, bins*sizeof(float)); under = over = 0; nHist = 0; nEntries = 0; if( err2 != NULL ) memset(err2, 0, bins*sizeof(double)); } /********* Methode *********/ /*! Pour avoir le calcul des erreurs */ void Histo::Errors() { if( bins > 0 ) { if(err2==NULL) err2 = new double[bins]; memset(err2, 0, bins*sizeof(double)); } } /********* Methode *********/ /*! Operateur egal */ Histo& Histo::operator = (const Histo& h) { if(this == &h) return *this; if( h.bins > bins ) Delete(); if(!data) data = new float[h.bins]; if( !h.err2 && err2 ) { delete [] err2; err2=NULL;} if( h.err2 && !err2 ) err2 = new double[h.bins]; under = h.under; over = h.over; nHist = h.nHist; nEntries = h.nEntries; bins = h.bins; min = h.min; max = h.max; binWidth = h.binWidth; memcpy(data, h.data, bins*sizeof(float)); if(err2) memcpy(err2, h.err2, bins*sizeof(double)); return *this; } /********* Methode *********/ /*! Operateur de multiplication par une constante */ Histo& Histo::operator *= (double b) { double b2 = b*b; for(int i=0;i &v) { v.Realloc(bins); for(int i=0;i &v) { v.Realloc(bins); for(int i=0;i &v) { v.Realloc(bins); if(!err2) {for(int i=0;i &v) { v.Realloc(bins); if(!err2) {for(int i=0;i &v, int ierr) { if(v.NElts()<(uint_4) bins) THROW(sizeMismatchErr); for(int i=0;i &v, int ierr) { if(v.NElts()<(uint_4) bins) THROW(sizeMismatchErr); for(int i=0;i &v) { if(v.NElts()<(uint_4) bins) THROW(sizeMismatchErr); if(!err2) Errors(); for(int i=0;i &v) { if(v.NElts()<(uint_4) bins) THROW(sizeMismatchErr); if(!err2) Errors(); for(int i=0;i0.) err2[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) bins) THROW(sizeMismatchErr); if(!err2) Errors(); for(int i=0;i0.) err2[i]=v(i)*v(i); else err2[i]=-v(i)*v(i); return; } /********* Methode *********/ /*! Addition du contenu de l'histo pour abscisse x poids w */ void Histo::Add(float x, float w) { int numBin = FindBin(x); if (numBin<0) under += w; else if (numBin>=bins) over += w; else { data[numBin] += w; if(err2!=NULL) err2[numBin] += w*w; nHist += w; nEntries++; } } /********* Methode *********/ /*! Addition du contenu de l'histo bin numBin poids w */ void Histo::AddBin(int numBin, float w) { if (numBin<0) under += w; else if (numBin>=bins) over += w; else { data[numBin] += w; if(err2!=NULL) err2[numBin] += w*w; nHist += w; nEntries++; } } /********* Methode *********/ /*! Remplissage du contenu de l'histo pour abscisse x poids w */ void Histo::SetBin(float x, float w) { int numBin = FindBin(x); SetBin(numBin,w); } /********* Methode *********/ /*! Remplissage du contenu de l'histo pour numBin poids w */ void Histo::SetBin(int numBin, float w) { if (numBin<0) under = w; else if (numBin>=bins) over = w; else { nHist -= data[numBin]; data[numBin] = w; nHist += w; } } /********* Methode *********/ /*! Remplissage des erreurs au carre pour abscisse x */ void Histo::SetErr2(float x, double e2) { int numBin = FindBin(x); SetErr2(numBin,e2); } /********* Methode *********/ /*! Remplissage des erreurs au carre pour numBin poids */ void Histo::SetErr2(int numBin, double e2) { if( err2==NULL) return; if ( numBin<0 || numBin>=bins ) return; err2[numBin] = e2; } /********* Methode *********/ /*! Remplissage des erreurs pour abscisse x */ void Histo::SetErr(float x, float e) { SetErr2(x, (double) e*e); } /********* Methode *********/ /*! Remplissage des erreurs pour numBin */ void Histo::SetErr(int numBin, float e) { SetErr2(numBin, (double) e*e); } /********* Methode *********/ /*! Numero du bin ayant le contenu maximum */ int Histo::IMax() const { int imx=0; if(bins==1) return imx; float mx=data[0]; for (int i=1; imx) {imx = i; mx=data[i];} return imx; } /********* Methode *********/ /*! Numero du bin ayant le contenu minimum */ int Histo::IMin() const { int imx=0; if(bins==1) return imx; float mx=data[0]; for (int i=1; imx) mx=data[i]; return mx; } /********* Methode *********/ /*! Valeur le contenu minimum */ float Histo::VMin() const { float mx=data[0]; if(bins==1) return mx; for (int i=1;i=0.) ? data[i] : -data[i]; n += v; sx += BinCenter(i)*v; } if(n>0.) return sx/n; else return min; } /********* Methode *********/ /*! Valeur du sigma */ float Histo::Sigma() const { double n = 0; double sx = 0; double sx2= 0; for (int i=0; i=0.) ? data[i] : -data[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 */ float Histo::MeanLH(int il,int ih) const { if( ih < il ) { il = 0; ih = bins-1; } if( il < 0 ) il = 0; if( ih >= bins ) ih = bins-1; double n = 0; double sx = 0; for (int i=il; i<=ih; i++) { double v = (data[i]>=0.) ? data[i] : -data[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 */ float Histo::SigmaLH(int il,int ih) const { if( ih < il ) { il = 0; ih = bins - 1; } if( il < 0 ) il = 0; if( ih >= bins ) ih = bins - 1; double n = 0; double sx = 0; double sx2= 0; for (int i=il; i<=ih; i++) { double v = (data[i]>=0.) ? data[i] : -data[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 */ float Histo::Mean(float x0, float dx) const { double sdata = 0; double sx = 0; int iMin = FindBin(x0-dx); if (iMin<0) iMin = 0; int iMax = FindBin(x0+dx); if (iMax>bins) iMax=bins; for (int i=iMin; i=0.) ? data[i] : -data[i]; sx += BinCenter(i)*v; sdata += v; } if(sdata>0.) return sx/sdata; else return min; } /********* Methode *********/ /*! Valeur du sigma calcule entre x0-dx et x0+dx */ float Histo::Sigma(float x0, float dx) const { double sx = 0; double sx2= 0; double sdata = 0; int iMin = FindBin(x0-dx); if (iMin<0) iMin = 0; int iMax = FindBin(x0+dx); if (iMax>bins) iMax=bins; for (int i=iMin; i=0.) ? data[i] : -data[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 */ float Histo::CleanedMean(float& sigma) const { if (!nHist) return 0; // on fait ca de facon bete, on pourra raffiner apres... float x0 = Mean(); float s = Sigma()+binWidth; for (int i=0; i<3; i++) { float xx0 = Mean(x0,5*s); s = Sigma(x0,5*s)+binWidth; x0 = xx0; } sigma = s; return x0; } /********* Methode *********/ /*! Valeur de la moyenne nettoyee */ float Histo::CleanedMean() const { if (!nHist) return 0; float s=0; return CleanedMean(s); } /********* Methode *********/ /*! Retourne le nombre de bins non-nul */ int Histo::BinNonNul() const { int non=0; for (int i=0;i= n ) break; } if( i == bins ) i = bins-1; return i; } /********* Methode *********/ /*! Renvoie les numeros de bins imin,imax (0==bins ) return -2; double N1 = 0.; for(int i=0; i<=I; i++) N1 += data[i]; N1 *= per; double N2 = 0.; {for(int i=I; i=0; imin--) { n += data[imin]; if(n>=N1) break; } if( imin<0 ) imin = 0; // cout<<"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 ii = 0; for (int j=iLow; j<=iHigh; j++) { if ((*this)(j)>0) { if(!err2) { 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 k; for(k=0;k1) cout << "resultat fit = " << pol << endl; pol.Derivate(); } CATCHALL { THROW_SAME; } ENDTRY int DPolDeg = pol.Degre(); float 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 double r=0; if (pol.Root1(r)==0) THROW(inconsistentErr); fd = r + xCenter; } else if (DPolDeg == 2) { // on est dans le cas d'un fit de cubique double r1=0; double r2=0; if (pol.Root2(r1,r2) == 0) THROW(inconsistentErr); pol.Derivate(); fd = (pol(r1)<0) ? r1 + xCenter : r2 + xCenter; } else { // on est dans un cas non prevu THROW(inconsistentErr); } if(fd>max) fd = max; 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 */ float Histo::FindWidth(float frac, int debug) const { float 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 */ float Histo::FindWidth(float xmax,float frac, int debug) const { 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(inconsistentErr); if (NBins() < 3) THROW(inconsistentErr); int iMax = FindBin(xmax); if (iMax<0 || iMax>=NBins()) THROW(inconsistentErr); float hmax = data[iMax]; float limit = frac*hmax; if (debug > 1) cout << " Max histo[" << iMax << "] = " << hmax << ", limite " << limit << endl; int iLow = iMax; while (iLow>=0 && data[iLow]>limit) iLow--; if( iLow < 0 ) iLow = 0; int iHigh = iMax; while (iHighlimit) iHigh++; if( iHigh >=NBins() ) iHigh = NBins()-1; float xlow = BinCenter(iLow); float ylow = data[iLow]; float xlow1=xlow, ylow1=ylow; if(iLow+1=0) { xhigh1 = BinCenter(iHigh-1); yhigh1 = data[iHigh-1]; } float 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 Histo::EstimeMax(float& xm,int SzPav) { int 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 Histo::EstimeMax(int& im,float& xm,int SzPav) { xm = 0; if( SzPav <= 0 ) return -1; if( im < 0 || im >= bins ) return -1; if( SzPav%2 == 0 ) SzPav++; SzPav = (SzPav-1)/2; int rc = 0; double dxm = 0, wx = 0; for(int i=im-SzPav;i<=im+SzPav;i++) { if( i<0 || i>= bins ) {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(float frac,float& widthG,float& widthD) { int i; widthG = widthD = -1.; if( bins<=1 || frac<=0. || frac>=1. ) return; int imax = 0; float maxi = data[0]; for(i=1;imaxi) {imax=i; maxi=data[i];} float xmax = BinCenter(imax); float maxf = maxi * frac; // recherche du sigma a gauche widthG = 0.; for(i=imax;i>=0;i--) if( data[i] <= maxf ) break; if(i<0) i=0; if(i=bins) i=bins-1; if(i>imax) { if( data[i] != data[i-1] ) { widthD = BinCenter(i) - binWidth * (maxf-data[i])/(data[i-1]-data[i]); widthD -= xmax; } else widthD = BinCenter(i) - xmax; } } ////////////////////////////////////////////////////////// /*! Fit de l'histogramme par ``gfit''. \verbatim typ_err = 0 : - erreur attachee au bin si elle existe - sinon 1 typ_err = 1 : - erreur attachee au bin si elle existe - sinon max( sqrt(abs(bin) ,1 ) typ_err = 2 : - erreur forcee a 1 typ_err = 3 : - erreur forcee a max( sqrt(abs(bin) ,1 ) typ_err = 4 : - erreur forcee a 1, nulle si bin a zero. typ_err = 5 : - erreur forcee a max( sqrt(abs(bin) ,1 ), nulle si bin a zero. \endverbatim */ int Histo::Fit(GeneralFit& gfit,unsigned short typ_err) { if(NBins()<=0) return -1000; if(typ_err>5) typ_err=0; GeneralFitData mydata(1,NBins()); for(int i=0;i par = gfit.GetParm(); Histo h(*this); for(int i=0;iValue(&x,par.Data()); } return h; } /*! Retourne une classe contenant la fonction du fit ``gfit''. */ Histo Histo::FitFunction(GeneralFit& gfit) { if(NBins()<=0) throw(SzMismatchError("Histo::FitFunction: size mismatch\n")); GeneralFunction* f = gfit.GetFunction(); if(f==NULL) throw(NullPtrError("Histo::FitFunction: NULL pointer\n")); TVector par = gfit.GetParm(); Histo h(*this); for(int i=0;iValue(&x,par.Data()); } return h; } /********* Methode *********/ /*! Impression de l'histogramme dans le fichier fp \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) data[i]" (note "... ...") bit 0 on : (v=1): numero_bin "... ..." bit 1 on : (v=2): "... ..." erreur \endverbatim */ void Histo::PrintF(FILE * fp, int hdyn,float hmin, float hmax,int pflag, int il, int ih) { if( il > ih ) { il = 0; ih = bins-1; } if( il < 0 ) il = 0; if( ih >= bins ) ih = bins-1; double dhmin = (double) hmin; double dhmax = (double) hmax; double hb,hbmn,hbmx; if(hdyn==0) hdyn = 100; cout << "~Histo::Print " << " nHist=" << nHist << " nEntries=" << nEntries << " under=" << under << " over=" << over << endl; cout << " bins=" << bins << " min=" << min << " max=" << max << " binWidth=" << binWidth << endl; cout << " mean=" << Mean() << " r.m.s=" << Sigma() << endl; if(hdyn<0 || pflag<0 ) return; if(dhmax<=dhmin) { if(hmin != 0.) dhmin = (double) VMin(); else dhmin=0.; dhmax = (double) VMax(); } if(dhmin>dhmax) return; if(dhmin==dhmax) {dhmin -= 1.; dhmax += 1.;} double wb = (dhmax-dhmin) / (double) hdyn; // determination de la position de la valeur zero int i0 = (int)(-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) / (double) hdyn; i0 = (int)(-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 k=0;k<=ii;k++) s[k] = 'X'; } else if(i0>=hdyn) { // courbe entierement negative for(int k=hdyn-1;k>=ii;k--) s[k] = 'X'; } else { // courbe positive et negative s[i0] = '|'; if(ii>i0) for(int k=i0+1;k<=ii;k++) s[k] = 'X'; else if(ii0.) ib = (int)( 10.*(hb-hbmn)/(hbmx-hbmn) ); else if(hb<0.) ib = (int)( 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) '0' + ib); // traitement des saturations +/- if( hb < dhmin ) s[0] = '*'; else if( hb > dhmax ) s[hdyn-1] = '*'; if(pflag&1) fprintf( fp, "%3d ",i); fprintf( fp, "%10.4g %10.4g ",BinCenter(i),hb); if(pflag&2 && err2!=NULL) fprintf( fp, "%10.4g ",Error(i)); fprintf( fp, "= %s\n",s); }} // derniere ligne for(int i=0;i=3)) const int ndig = 7; char sn[2*ndig+10]; hb = ( fabs(dhmin) > fabs(dhmax) ) ? fabs(dhmin) : fabs(dhmax); int n; if( hb > 0. ) n = (int)(log10(hb*1.00000001)); else n = 1; double scaledig = pow(10.,(double) ndig-2); double expo = scaledig/pow(10.,(double) n); // cout <<"n="<::ReadSelf(PInPersist& is) { char strg[256]; if(dobj==NULL) dobj=new Histo; else dobj->Delete(); // Lecture entete is.GetLine(strg, 255); is.GetLine(strg, 255); is.GetLine(strg, 255); // Lecture des valeurs is.Get(dobj->bins); is.Get(dobj->nEntries); int_4 errok; is.Get(errok); is.Get(dobj->binWidth); is.Get(dobj->min); is.Get(dobj->max); is.Get(dobj->under); is.Get(dobj->over); is.Get(dobj->nHist); // Lecture des donnees Histo 1D dobj->data = new float[dobj->bins]; is.GetLine(strg, 255); is.Get(dobj->data, dobj->bins); // Lecture des erreurs if(errok) { is.GetLine(strg, 255); dobj->err2 = new double[dobj->bins]; is.Get(dobj->err2, dobj->bins); } return; } void ObjFileIO::WriteSelf(POutPersist& os) const { if (dobj == NULL) return; char strg[256]; // Que faut-il ecrire? int_4 errok = (dobj->err2) ? 1 : 0; // Ecriture entete pour identifier facilement sprintf(strg,"bins=%6d NEnt=%15d errok=%1d",dobj->bins,dobj->nEntries,errok); os.PutLine(strg); sprintf(strg,"binw=%g min=%g max=%g",dobj->binWidth,dobj->min,dobj->max); os.PutLine(strg); sprintf(strg, "under=%g over=%g nHist=%g",dobj->under,dobj->over,dobj->nHist); os.PutLine(strg); // Ecriture des valeurs os.Put(dobj->bins); os.Put(dobj->nEntries); os.Put(errok); os.Put(dobj->binWidth); os.Put(dobj->min); os.Put(dobj->max); os.Put(dobj->under); os.Put(dobj->over); os.Put(dobj->nHist); // Ecriture des donnees Histo 1D sprintf(strg,"Histo: Tableau des donnees %d",dobj->bins); os.PutLine(strg); os.Put(dobj->data, dobj->bins); // Ecriture des erreurs if(errok) { sprintf(strg,"Histo: Tableau des erreurs %d",dobj->bins); os.PutLine(strg); os.Put(dobj->err2, dobj->bins); } return; } #ifdef __CXX_PRAGMA_TEMPLATES__ #pragma define_template ObjFileIO #endif #if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES) template class ObjFileIO; #endif