source: Sophya/trunk/SophyaLib/HiStats/ntuple.cc@ 3718

Last change on this file since 3718 was 3619, checked in by cmv, 16 years ago

add various #include<> for g++ 4.3 (jaunty 9.04), cmv 05/05/2009

File size: 19.9 KB
RevLine 
[763]1#include <stdio.h>
[3619]2#include <stdlib.h>
[763]3#include <string.h>
4
5#include "strutil.h"
6#include "perrors.h"
7#include "ntuple.h"
[3392]8#include "thsafeop.h"
[763]9
[3236]10namespace SOPHYA {
[763]11
12#define LENNAME 8
13#define LENNAME1 (LENNAME+1)
[2663]14#define BADVAL -1.e19
[763]15
[1371]16/*!
[3236]17 \class NTuple
[1371]18 \ingroup HiStats
[2667]19 Simple NTuple class (a Table or 2-D data set) with double or
20 single precision floating point value columns.
[3392]21 \warning
22 Thread safe fill operation can be activated using the method SetThreadSafe()
23 Default mode is NON thread-safe fill.
[1405]24 \sa SOPHYA::ObjFileIO<NTuple>
25
26 \code
27 #include "ntuple.h"
28 // ...
[3572]29 const char * names[3] = {"XPos", "YPos", "Val"};
[2667]30 // NTuple (Table) creation with 3 columns (double precision)
[1405]31 NTuple nt(3, names);
32 // Filling the NTuple
[2667]33 r_8 x[3];
[1405]34 for(int i=0; i<63; i++) {
35 x[0] = (i%9)-4.; x[1] = (i/9)-3.; x[2] = x[0]*x[0]+x[1]*x[1];
36 nt.Fill(x);
37 }
38 // Printing table info
39 cout << nt ;
40 // Saving object into a PPF file
41 POutPersist po("nt.ppf");
42 po << nt ;
43 \endcode
[1371]44*/
45
[763]46//++
47// Class NTuple
48// Lib Outils++
49// include ntuple.h
50//
51// Classe de ntuples
52//--
53//++
54// Links Parents
55// PPersist
56// NTupleInterface
57//--
58
59/* --Methode-- */
[2808]60//! Default constructor - To be used when reading in an NTuple.
[763]61//++
62NTuple::NTuple()
63//
64// Createur par defaut
65//--
66{
67mNVar = mNEnt = mBlk = mNBlk = 0;
68mVar = NULL;
69mVarD = NULL;
[2663]70mFgDouble = true;
[763]71mInfo = NULL;
[3392]72mThS = NULL;
[763]73}
74
75
[3572]76//! Constructor with specification of number of columns and column names
[1405]77/*!
78 \param nvar : Number of columns in the table
[3572]79 \param noms : Array of column names
[1405]80 \param blk : Optional argument specifying number of table lines
81 in a given data block
[2663]82 \param fgdouble : if \b true: internal data kept as double precision values (r_8),
83 simple precision (r_4) otherwise
[1405]84 */
[763]85//++
[2663]86NTuple::NTuple(int nvar, char** noms, int blk, bool fgdouble)
[763]87//
88// Createur d'un ntuple de `nvar' variables dont les
[3112]89// noms sont dans le tableau de chaines de caracteres `noms'
[763]90// avec `blk' d'evenements par blocks.
91//--
92{
[3572]93if (nvar <= 0) throw ParmError("NTuple::NTuple(nvar<=0) with char** noms");
94Initialize(nvar,blk,fgdouble);
95for(int i=0; i<nvar; i++) mNames.push_back(noms[i]);
96return;
97}
[2663]98
[3572]99
100//! Constructor with specification of number of columns and column names
101/*!
102 \param nvar : Number of columns in the table
103 \param noms : Array of column names
104 \param blk : Optional argument specifying number of table lines
105 in a given data block
106 \param fgdouble : if \b true: internal data kept as double precision values (r_8),
107 simple precision (r_4) otherwise
108 */
109//++
110NTuple::NTuple(int nvar, const char** noms, int blk, bool fgdouble)
111//
112// Createur d'un ntuple de `nvar' variables dont les
113// noms sont dans le tableau de chaines de caracteres `noms'
114// avec `blk' d'evenements par blocks.
115//--
116{
117if (nvar <= 0) throw ParmError("NTuple::NTuple(nvar<=0) with const char** noms");
118Initialize(nvar,blk,fgdouble);
[2663]119for(int i=0; i<nvar; i++) mNames.push_back(noms[i]);
[763]120return;
121}
122
[3572]123//! Constructor with specification of number of columns and column names as a string vector
[3112]124/*!
125 \param noms : Array of column names (length(name) < 8 characters)
126 \param blk : Optional argument specifying number of table lines
127 in a given data block
128 \param fgdouble : if \b true: internal data kept as double precision values (r_8),
129 simple precision (r_4) otherwise
130 */
131NTuple::NTuple(vector<string>& noms, int blk, bool fgdouble)
132{
[3572]133int nvar = noms.size();
134if (nvar <= 0) throw ParmError("NTuple::NTuple(nvar<=0) with vector<string>& noms");
135Initialize(nvar,blk,fgdouble);
136for(int i=0; i<nvar; i++) mNames.push_back(noms[i]);
137return;
138}
139
140
141/* --Methode-- */
142/* Initialisation pour Createurs (fonction privee) */
143void NTuple::Initialize(int nvar, int blk, bool fgdouble)
144{
[3112]145mNVar = mNEnt = mBlk = mNBlk = 0;
146mVar = NULL;
147mVarD = NULL;
148mInfo = NULL;
[3392]149mThS = NULL;
[3112]150mNVar = nvar;
151mVar = new r_4[nvar];
152mVarD = new r_8[nvar];
153if (blk < 10) blk = 10;
154mBlk = blk;
155
156if (fgdouble) {
157 r_8* pt = new r_8[nvar*blk];
158 mNBlk = 1;
159 mPtrD.push_back(pt);
160 mFgDouble = true;
161}
162else {
163 r_4* pt = new r_4[nvar*blk];
164 mNBlk = 1;
165 mPtr.push_back(pt);
166 mFgDouble = false;
167}
[3572]168
[3112]169return;
170}
171
[3392]172/*!
173 Copy constructor - Copies the table definition and associated data,
174 as well as the tread safety state
175*/
[763]176// cmv 8/10/99
177//++
178NTuple::NTuple(const NTuple& NT)
179//
180// Createur par copie (clone).
181//--
182: mNVar(0), mNEnt(0), mBlk(0), mNBlk(0)
[3392]183 , mVar(NULL), mVarD(NULL), mInfo(NULL),mThS(NULL)
[763]184{
185if(NT.mNVar<=0) return; // cas ou NT est cree par defaut
186mNVar = NT.mNVar;
187mBlk = NT.mBlk;
188mVar = new r_4[NT.mNVar];
189mVarD = new r_8[NT.mNVar];
190
[2663]191mNames = NT.mNames;
192
[763]193int i;
[2663]194mFgDouble = NT.mFgDouble;
195if (mFgDouble) {
196 r_8* pt = new r_8[mNVar*mBlk];
197 mNBlk = 1; mPtrD.push_back(pt);
198 if(NT.mNEnt > 0)
199 for(i=0;i<NT.mNEnt;i++) {pt=NT.GetVecD(i,NULL); Fill(pt);}
200}
201else {
202 r_4* pt = new r_4[mNVar*mBlk];
203 mNBlk = 1; mPtr.push_back(pt);
204 for(i=0;i<NT.mNEnt;i++) {pt=NT.GetVec(i,NULL); Fill(pt);}
205}
[763]206
[2663]207mNames = NT.mNames;
[763]208
209if(NT.mInfo!=NULL) {mInfo = new DVList; *mInfo = *(NT.mInfo);}
[3392]210mThS = NULL;
211if(NT.mThS!=NULL) SetThreadSafe(true);
[763]212return;
213}
214
215
216/* --Methode-- */
217NTuple::~NTuple()
218{
219Clean();
220}
[3572]221
[3392]222/* --Methode-- */
223/*!
224 \brief Activate or deactivate thread-safe \b Fill() operations
225 If fg==true, create an ThSafeOp object in order to insure atomic Fill()
226 operations. if fg==false, the ThSafeOp object (mThS) of the target NTuple
227 is destroyed.
228 \warning The default Fill() operation mode for NTuples is NOT thread-safe.
229 Please note also that the thread-safety state is NOt saved to PPF (or FITS) streams.
230*/
231void NTuple::SetThreadSafe(bool fg)
232{
233if (fg) {
234 if (mThS) return;
235 else mThS = new ThSafeOp();
236}
237else {
238 if (mThS) delete mThS;
239 mThS = NULL;
240}
241}
[763]242
243/* --Methode-- */
[1405]244//! Clear the data table definition and deletes the associated data
[763]245void NTuple::Clean()
246{
247if (mVar) delete[] mVar;
248if (mVarD) delete[] mVarD;
249if (mInfo) delete mInfo;
[3392]250if (mThS) delete mThS;
[763]251int i;
[2663]252if(mNBlk>0) {
253 if (mFgDouble) for(i=0; i<mNBlk; i++) delete[] mPtrD[i];
254 else for(i=0; i<mNBlk; i++) delete[] mPtr[i];
255}
256if (mFgDouble) mPtrD.erase(mPtrD.begin(), mPtrD.end());
257else mPtr.erase(mPtr.begin(), mPtr.end());
[763]258mNVar = mNEnt = mBlk = mNBlk = 0;
259mVar = NULL;
260mVarD = NULL;
261mInfo = NULL;
[3392]262mThS = NULL;
[763]263return;
264}
265
266/* --Methode-- cmv 08/10/99 */
[1405]267//! = operator, copies the data table definition and its contents
[763]268//++
269NTuple& NTuple::operator = (const NTuple& NT)
270//
271// Operateur egal (clone).
272//--
273{
274if(this == &NT) return *this;
275Clean();
276if(NT.mNVar<=0) return *this; // cas ou NT est cree par defaut
277mNVar = NT.mNVar;
278mBlk = NT.mBlk;
279mVar = new r_4[NT.mNVar];
280mVarD = new r_8[NT.mNVar];
281
[2663]282mNames = NT.mNames;
283
[763]284int i;
[2663]285mFgDouble = NT.mFgDouble;
286if (mFgDouble) {
287 r_8* pt = new r_8[mNVar*mBlk];
288 mNBlk = 1; mPtrD.push_back(pt);
289 if(NT.mNEnt > 0)
290 for(i=0;i<NT.mNEnt;i++) {pt=NT.GetVecD(i,NULL); Fill(pt);}
291}
292else {
293 r_4* pt = new r_4[mNVar*mBlk];
294 mNBlk = 1; mPtr.push_back(pt);
295 for(i=0;i<NT.mNEnt;i++) {pt=NT.GetVec(i,NULL); Fill(pt);}
296}
[763]297
298if(NT.mInfo!=NULL) {mInfo = new DVList; *mInfo = *(NT.mInfo);}
299
300return *this;
301}
302
303/* --Methode-- */
[1405]304//! Adds an entry (a line) to the table
305/*!
306 \param x : content of the line to be appended to the table
307 */
[763]308//++
309void NTuple::Fill(r_4* x)
310//
311// Remplit le ntuple avec le tableau cd reels `x'.
312//--
313{
[3392]314if (mThS) mThS->lock(); // Thread-safe operation
[763]315int numb = mNEnt/mBlk;
316if (numb >= mNBlk) {
[2663]317 if (mFgDouble) {
318 r_8* pt = new r_8[mNVar*mBlk];
319 mNBlk++;
320 mPtrD.push_back(pt);
321 }
322 else {
323 r_4* pt = new r_4[mNVar*mBlk];
324 mNBlk++;
325 mPtr.push_back(pt);
326 }
[763]327}
[2663]328
[763]329int offb = mNEnt-numb*mBlk;
[2663]330if (mFgDouble)
331 for(int i=0; i<mNVar; i++) (mPtrD[numb]+offb*mNVar)[i] = x[i];
332else memcpy((mPtr[numb]+offb*mNVar), x, mNVar*sizeof(r_4));
[763]333mNEnt++;
[3392]334if (mThS) mThS->unlock(); // Thread-safe operation
[763]335return;
336}
337
[2663]338/* --Methode-- */
339//! Adds an entry (a line) to the table (double precision)
340/*!
341 \param x : content of the line to be appended to the table
342 */
343//++
344void NTuple::Fill(r_8* x)
345//
346// Remplit le ntuple avec le tableau double precision `x'.
347//--
348{
[3392]349if (mThS) mThS->lock(); // Thread-safe operation
[2663]350int numb = mNEnt/mBlk;
351if (numb >= mNBlk) {
352 if (mFgDouble) {
353 r_8* pt = new r_8[mNVar*mBlk];
354 mNBlk++;
355 mPtrD.push_back(pt);
356 }
357 else {
358 r_4* pt = new r_4[mNVar*mBlk];
359 mNBlk++;
360 mPtr.push_back(pt);
361 }
362}
[763]363
[2663]364int offb = mNEnt-numb*mBlk;
365if (mFgDouble)
366 memcpy((mPtrD[numb]+offb*mNVar), x, mNVar*sizeof(r_8));
367else
368 for(int i=0; i<mNVar; i++) (mPtr[numb]+offb*mNVar)[i] = x[i];
369
370mNEnt++;
[3392]371if (mThS) mThS->unlock(); // Thread-safe operation
[2663]372return;
373}
374
375
[763]376/* --Methode-- */
377//++
378float NTuple::GetVal(int n, int k) const
379//
380// Retourne la valeur de la variable `k' de l'evenement `n'.
381//--
382{
383if (n >= mNEnt) return(BADVAL);
384if ( (k < 0) || (k >= mNVar) ) return(BADVAL);
385int numb = n/mBlk;
386int offb = n-numb*mBlk;
[2663]387if ( mFgDouble) return(*(mPtrD[numb]+offb*mNVar+k));
388else return(*(mPtr[numb]+offb*mNVar+k));
[763]389}
390
391
392/* --Methode-- */
393//++
394int NTuple::IndexNom(const char* nom) const
395//
396// Retourne le numero de la variable de nom `nom'.
397//--
398{
399int i;
[2663]400string snom = nom;
[763]401for(i=0; i<mNVar; i++)
[2663]402 if ( mNames[i] == snom) return(i);
[763]403return(-1);
404}
405
406
[2663]407static char nomretour[256];
[763]408/* --Methode-- */
409//++
410char* NTuple::NomIndex(int k) const
411//
412// Retourne le nom de la variable numero 'k'
413//--
414{
415nomretour[0] = '\0';
[2663]416if ((k >= 0) && (k < mNVar)) strncpy(nomretour, mNames[k].c_str(), 255);
[763]417return(nomretour);
418}
419
420
421/* --Methode-- */
422//++
423r_4* NTuple::GetVec(int n, r_4* ret) const
424//
425// Retourne l'evenement `n' dans le vecteur `ret'.
426//--
427{
428int i;
429if (ret == NULL) ret = mVar;
430if (n >= mNEnt) {
431 for(i=0; i<mNVar; i++) ret[i] = BADVAL;
432 return(ret);
433}
434
435int numb = n/mBlk;
436int offb = n-numb*mBlk;
[2663]437if (mFgDouble) for(i=0; i<mNVar; i++) ret[i] = (mPtrD[numb]+offb*mNVar)[i];
438else memcpy(ret, (mPtr[numb]+offb*mNVar), mNVar*sizeof(r_4));
[763]439return(ret);
440}
441
442/* --Methode-- */
443//++
444r_8* NTuple::GetVecD(int n, r_8* ret) const
445//
446// Retourne l'evenement `n' dans le vecteur `ret'.
447//--
448{
449int i;
450if (ret == NULL) ret = mVarD;
[2663]451if (n >= mNEnt) {
452 for(i=0; i<mNVar; i++) ret[i] = BADVAL;
453 return(ret);
454}
455
456int numb = n/mBlk;
457int offb = n-numb*mBlk;
[2667]458if (mFgDouble) memcpy(ret, (mPtrD[numb]+offb*mNVar), mNVar*sizeof(r_8));
459else for(i=0; i<mNVar; i++) ret[i] = (mPtr[numb]+offb*mNVar)[i];
[763]460return(ret);
461}
462
463
464
465/* --Methode-- */
466//++
467DVList& NTuple::Info()
468//
469// Renvoie une référence sur l'objet DVList Associé
470//--
471{
472if (mInfo == NULL) mInfo = new DVList;
473return(*mInfo);
474}
475
476/* --Methode-- */
477//++
478void NTuple::Print(int num, int nmax) const
479//
480// Imprime `nmax' evenements a partir du numero `num'.
481//--
482{
483int i,j;
484
485printf("Num ");
[2663]486for(i=0; i<mNVar; i++) printf("%8s ", mNames[i].c_str());
[763]487putchar('\n');
488
489if (nmax <= 0) nmax = 1;
490if (num < 0) num = 0;
491nmax += num;
492if (nmax > mNEnt) nmax = mNEnt;
493for(i=num; i<nmax; i++) {
494 GetVec(i, NULL);
495 printf("%6d ", i);
496 for(j=0; j<mNVar; j++) printf("%8g ", (float)mVar[j]);
497 putchar('\n');
498}
499return;
500}
501
502/* --Methode-- */
[1405]503//! Prints table definition and number of entries
[763]504//++
505void NTuple::Show(ostream& os) const
506//
507// Imprime l'information generale sur le ntuple.
508//--
509{
[3572]510const char * tt = "float";
[2663]511if (mFgDouble) tt = "double";
512os << "NTuple T=" << tt << " : NVar= " << mNVar << " NEnt=" << mNEnt
[763]513 << " (Blk Sz,Nb= " << mBlk << " ," << mNBlk << ")\n";
514os << " Variables Min Max \n";
515int i;
516double min, max;
517char buff[128];
518for(i=0; i<mNVar; i++) {
519 GetMinMax(i, min, max);
[2663]520 sprintf(buff, "%3d %16s %10g %10g \n", i, mNames[i].c_str(), min, max);
[763]521 os << (string)buff ;
522 }
523os << endl;
524}
525
526
527/* --Methode-- */
[1405]528//! Fills the table, reading lines from an ascii file
529/*!
530 \param fn : file name
531 \param defval : default value for empty cells in the ascii file
532 */
[763]533//++
534int NTuple::FillFromASCIIFile(string const& fn, float defval)
535//
536// Remplit le ntuple a partir d'un fichier ASCII.
537// Renvoie le nombre de lignes ajoutees.
538//--
539{
540if (NbColumns() < 1) {
541 cout << "NTuple::FillFromASCIIFile() Ntuple has " << NbColumns() << " columns" << endl;
542 return(-1);
543 }
544FILE * fip = fopen(fn.c_str(), "r");
545if (fip == NULL) {
546 cout << "NTuple::FillFromASCIIFile() Error opening file " << fn << endl;
547 return(-2);
548 }
549
550char lineb[4096];
551char *line;
552char* ccp;
[809]553int j,kk;
[763]554int postab, posb;
[2663]555double* xv = new double[NbColumns()];
[763]556
557int nlread = 0;
558int nvar = NbColumns();
559// On boucle sur toutes les lignes
560while (fgets(lineb,4096,fip) != NULL) {
561 lineb[4095] = '\0';
562 j = 0; line = lineb;
563// On enleve les espaces et tab de debut
564 while ( (line[j] != '\0') && ((line[j] == ' ') || (line[j] == '\t')) ) j++;
565 line = lineb+j;
566// Il faut que le premier caractere non-espace soit un digit, ou + ou - ou .
567 if (!( isdigit(line[0]) || (line[0] == '+') || (line[0] == '-') || (line[0] == '.') )) continue;
568 ccp = line;
569 for(kk=0; kk<nvar; kk++) xv[kk] = defval;
570 for(kk=0; kk<nvar; kk++) {
571// Les mots sont separes par des espaces ou des tab
572 postab = posc(ccp, '\t' );
573 posb = posc(ccp, ' ' );
574 if (postab >= 0) {
575 if (posb < 0) posb = postab;
576 else if (postab < posb) posb = postab;
577 }
578 if (posb >= 0) ccp[posb] = '\0';
579 if ( isdigit(line[0]) || (line[0] == '+') || (line[0] == '-') || (line[0] == '.') )
580 xv[kk] = atof(ccp);
581 if (posb < 0) break;
582 ccp += posb+1; j = 0;
583 while ( (ccp[j] != '\0') && ((ccp[j] == ' ') || (ccp[j] == '\t')) ) j++;
584 ccp += j;
585 }
586 Fill(xv);
587 nlread++;
588 }
589
590delete[] xv;
591cout << "NTuple::FillFromASCIIFile( " << fn << ") " << nlread << " fill from file " << endl;
592return(nlread);
593}
594
595
596// ------- Implementation de l interface NTuple ---------
597
598/* --Methode-- */
[2682]599sa_size_t NTuple::NbLines() const
[763]600{
601return(NEntry());
602}
603/* --Methode-- */
[2682]604sa_size_t NTuple::NbColumns() const
[763]605{
606return(NVar());
607}
608
609/* --Methode-- */
[2682]610r_8 * NTuple::GetLineD(sa_size_t n) const
[763]611{
612return(GetVecD(n));
613}
614
615/* --Methode-- */
[2682]616r_8 NTuple::GetCell(sa_size_t n, sa_size_t k) const
[763]617{
618return(GetVal(n, k));
619}
620
621/* --Methode-- */
[2682]622r_8 NTuple::GetCell(sa_size_t n, string const & nom) const
[763]623{
624return(GetVal(n, nom.c_str()));
625}
626
627/* --Methode-- */
628//++
[2682]629void NTuple::GetMinMax(sa_size_t k, double& min, double& max) const
[763]630//
631// Retourne le minimum et le maximum de la variable `k'.
632//--
633{
634min = 9.e19; max = -9.e19;
635if ( (k < 0) || (k >= mNVar) ) return;
636int jb,ib,i;
637double x;
638i=0;
639for(jb=0; jb< mNBlk; jb++)
640 for(ib=0; ib< mBlk; ib++) {
641 if (i >= mNEnt) break;
642 i++;
[2663]643 x = (mFgDouble) ? *(mPtrD[jb]+ib*mNVar+k) : *(mPtr[jb]+ib*mNVar+k);
[763]644 if(i==1) {min = x; max = x;}
645 if (x < min) min = x;
646 if (x > max) max = x;
647 }
648return;
649}
650
651/* --Methode-- */
652void NTuple::GetMinMax(string const & nom, double& min, double& max) const
653{
654GetMinMax(IndexNom(nom.c_str()), min, max);
655}
656
657/* --Methode-- */
[2682]658sa_size_t NTuple::ColumnIndex(string const & nom) const
[763]659{
660return(IndexNom(nom.c_str()));
661}
662
663/* --Methode-- */
[2682]664string NTuple::ColumnName(sa_size_t k) const
[763]665{
666return(NomIndex(k));
667}
668
669/* --Methode-- */
670//++
671string NTuple::VarList_C(const char* nomx) const
672//
673// Retourne une chaine de caracteres avec la declaration des noms de
674// variables. si "nomx!=NULL" , des instructions d'affectation
675// a partir d'un tableau "nomx[i]" sont ajoutees.
676//--
677{
678string rets="";
679int i;
680for(i=0; i<mNVar; i++) {
681 if ( (i%5 == 0) && (i > 0) ) rets += ";";
682 if (i%5 == 0) rets += "\ndouble ";
683 else rets += ",";
[2663]684 rets += mNames[i];
[763]685 }
686rets += "; \n";
687if (nomx) {
688 char buff[256];
689 for(i=0; i<mNVar; i++) {
[2663]690 sprintf(buff,"%s=%s[%d]; ", mNames[i].c_str(), nomx, i);
[763]691 rets += buff;
692 if ( (i%3 == 0) && (i > 0) ) rets += "\n";
693 }
694 }
695
696return(rets);
697}
698
699
700/* --Methode-- */
701//++
702string NTuple::LineHeaderToString() const
703// Retourne une chaine de caracteres avec la liste des noms de
704// variables, utilisables pour une impression
705//--
706{
707char buff[32];
708string rets=" Num ";
709for(int i=0; i<mNVar; i++) {
[2663]710 sprintf(buff, "%8s ", mNames[i].c_str());
[763]711 rets += buff;
712 }
713rets += '\n';
714return(rets);
715}
716
717/* --Methode-- */
718//++
[2682]719string NTuple::LineToString(sa_size_t n) const
[763]720// Retourne une chaine de caracteres avec le contenu de la ligne "n"
721// utilisable pour une impression
722//--
723{
724char buff[32];
725double* val;
726val = GetLineD(n);
727sprintf(buff,"%6d: ",n);
728string rets=buff;
729int i;
730for(i=0; i<mNVar; i++) {
731 sprintf(buff, "%8.3g ", val[i]);
732 rets += buff;
733 }
734rets += '\n';
735return(rets);
736}
737
[1405]738/*!
[3236]739 \class ObjFileIO<NTuple>
[1405]740 \ingroup HiStats
741 Persistence (serialisation) handler for class NTuple
742*/
[763]743
[846]744/* --Methode-- */
[763]745//++
[2341]746DECL_TEMP_SPEC /* equivalent a template <> , pour SGI-CC en particulier */
[763]747void ObjFileIO<NTuple>::WriteSelf(POutPersist& s) const
748//
749// Ecriture ppersist du ntuple.
750//--
751{
[1155]752if (dobj == NULL) return;
753
[2663]754// On ecrit cette chaine pour compatibilite avec les anciennes version (1,2)
[1155]755string strg = "NTuple";
756s.Put(strg);
757
758// On ecrit 3 uint_4 ....
[2663]759// [0]: Numero de version ;
760// [1] : bit1 non nul -> has info, bit 2 non nul mFgDouble=true
761// [2] : reserve
[1155]762uint_4 itab[3];
[2663]763itab[0] = 3; // Numero de version a 3
[1155]764itab[1] = itab[2] = 0;
765if (dobj->mInfo) itab[1] = 1;
[2663]766if (dobj->mFgDouble) itab[1] += 2;
[1155]767s.Put(itab, 3);
768
769s.Put(dobj->mNVar);
[2663]770for(int k=0; k<dobj->mNVar; k++) s.Put(dobj->mNames[k]);
[1155]771s.Put(dobj->mNEnt);
772s.Put(dobj->mBlk);
773s.Put(dobj->mNBlk);
[2663]774
[763]775if (dobj->mInfo) s << (*(dobj->mInfo));
776int jb;
[2663]777if (dobj->mFgDouble)
778 for(jb=0; jb<dobj->mNBlk; jb++)
779 s.Put(dobj->mPtrD[jb], dobj->mNVar*dobj->mBlk);
780else
781 for(jb=0; jb<dobj->mNBlk; jb++)
782 s.Put(dobj->mPtr[jb], dobj->mNVar*dobj->mBlk);
[763]783return;
784}
785
786/* --Methode-- */
787//++
[2341]788DECL_TEMP_SPEC /* equivalent a template <> , pour SGI-CC en particulier */
[763]789void ObjFileIO<NTuple>::ReadSelf(PInPersist& s)
790//
791// Lecture ppersist du ntuple.
792//--
793{
794
[1155]795if (dobj == NULL) dobj = new NTuple;
796else dobj->Clean();
[763]797
798bool hadinfo = false;
[1155]799string strg;
800s.Get(strg);
[2663]801uint_4 itab[3] = {0,0,0};
[1155]802if (strg == "NTuple") {
803 s.Get(itab, 3);
[2663]804 if ( ((itab[0] < 3) && (itab[1] != 0)) ||
805 ((itab[0] >= 3) && ((itab[1]&1) == 1)) ) hadinfo = true;
[1155]806}
807else {
808// Ancienne version de PPF NTuple - Pour savoir s'il y avait un DVList Info associe
809 char buff[256];
810 strcpy(buff, strg.c_str());
811 if (strncmp(buff+strlen(buff)-7, "HasInfo", 7) == 0) hadinfo = true;
812}
[2663]813if (itab[0] < 3) { // Lecture version anterieures V= 1 , 2
814 s.Get(dobj->mNVar);
815 dobj->mVar = new r_4[dobj->mNVar];
816 dobj->mVarD = new r_8[dobj->mNVar];
817
818 char * names = new char[dobj->mNVar*LENNAME1];
819 s.GetBytes(names, dobj->mNVar*LENNAME1);
820 for(int k=0; k<dobj->mNVar; k++) dobj->mNames.push_back(names+k*LENNAME1);
821 s.Get(dobj->mNEnt);
822 s.Get(dobj->mBlk);
823 s.Get(dobj->mNBlk);
824 if (hadinfo) { // Lecture eventuelle du DVList Info
825 if (dobj->mInfo == NULL) dobj->mInfo = new DVList;
826 s >> (*(dobj->mInfo));
[763]827 }
[2663]828 // Il n'y avait que des NTuple avec data float pour V < 3
829 dobj->mFgDouble = false;
830 int jb;
831 for(jb=0; jb<dobj->mNBlk; jb++) {
832 r_4* pt = new r_4[dobj->mNVar*dobj->mBlk];
833 dobj->mPtr.push_back(pt);
834 s.Get(dobj->mPtr[jb], dobj->mNVar*dobj->mBlk);
835 }
[763]836}
[2663]837else { // Lecture version V 3
838 s.Get(dobj->mNVar);
839 dobj->mVar = new r_4[dobj->mNVar];
840 dobj->mVarD = new r_8[dobj->mNVar];
841 string nom;
842 for(int k=0; k<dobj->mNVar; k++) {
843 s.Get(nom);
844 dobj->mNames.push_back(nom);
845 }
846 s.Get(dobj->mNEnt);
847 s.Get(dobj->mBlk);
848 s.Get(dobj->mNBlk);
849 if (hadinfo) { // Lecture eventuelle du DVList Info
850 if (dobj->mInfo == NULL) dobj->mInfo = new DVList;
851 s >> (*(dobj->mInfo));
852 }
853 // Il n'y avait que des NTuple avec data float pour V < 3
854 dobj->mFgDouble = ((itab[1]&2) == 2) ? true : false;
855 int jb;
856 if (dobj->mFgDouble) {
857 for(jb=0; jb<dobj->mNBlk; jb++) {
858 r_8* pt = new r_8[dobj->mNVar*dobj->mBlk];
859 dobj->mPtrD.push_back(pt);
860 s.Get(dobj->mPtrD[jb], dobj->mNVar*dobj->mBlk);
861 }
862 }
863 else {
864 for(jb=0; jb<dobj->mNBlk; jb++) {
865 r_4* pt = new r_4[dobj->mNVar*dobj->mBlk];
866 dobj->mPtr.push_back(pt);
867 s.Get(dobj->mPtr[jb], dobj->mNVar*dobj->mBlk);
868 }
869 }
870} // Fin lecture V3
[763]871
872}
873
[846]874
[763]875#ifdef __CXX_PRAGMA_TEMPLATES__
876#pragma define_template ObjFileIO<NTuple>
877#endif
878
879#if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES)
[3236]880template class ObjFileIO<NTuple>;
[763]881#endif
[3236]882
883} // FIN namespace SOPHYA
Note: See TracBrowser for help on using the repository browser.