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

Last change on this file since 2822 was 2808, checked in by ansari, 20 years ago

MAJ documentation - Reza 14/6/2005

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