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

Last change on this file since 1371 was 1371, checked in by ansari, 25 years ago

MAJ documentation, Makefile, ... - Reza 5/1/2001

File size: 12.7 KB
RevLine 
[763]1#include <stdio.h>
2#include <string.h>
3
4#include "strutil.h"
5#include "perrors.h"
6#include "ntuple.h"
7
8
9#define BADVAL -1.e19
10#define LENNAME 8
11#define LENNAME1 (LENNAME+1)
12
[1371]13/*!
14 \class SOPHYA::NTuple
15 \ingroup HiStats
16 Simple NTuple class (a Table or 2-D data set) with floating value
17 columns.
18*/
19
[763]20//++
21// Class NTuple
22// Lib Outils++
23// include ntuple.h
24//
25// Classe de ntuples
26//--
27//++
28// Links Parents
29// PPersist
30// NTupleInterface
31//--
32
33/* --Methode-- */
34//++
35NTuple::NTuple()
36//
37// Createur par defaut
38//--
39{
40mNVar = mNEnt = mBlk = mNBlk = 0;
41mVar = NULL;
42mVarD = NULL;
43mNames = NULL;
44mInfo = NULL;
45}
46
47
48//++
49NTuple::NTuple(int nvar, char** noms, int blk)
50//
51// Createur d'un ntuple de `nvar' variables dont les
52// noms sont dans le tableau de cahines de caracteres `noms'
53// avec `blk' d'evenements par blocks.
54//--
55{
56mNVar = mNEnt = mBlk = mNBlk = 0;
57mVar = NULL;
58mVarD = NULL;
59mNames = NULL;
60mInfo = NULL;
61if (nvar <= 0) THROW(sizeMismatchErr);
62mNVar = nvar;
63mVar = new r_4[nvar];
64mVarD = new r_8[nvar];
65if (blk < 10) blk = 10;
66mBlk = blk;
67// On prend des noms de LENNAME char pour le moment
68mNames = new char[nvar*LENNAME1];
69r_4* pt = new r_4[nvar*blk];
70mNBlk = 1;
71mPtr.push_back(pt);
72int i;
73for(i=0; i<nvar; i++)
74 { strncpy(mNames+i*LENNAME1, noms[i], LENNAME);
75 mNames[i*LENNAME1+LENNAME] = '\0'; }
76return;
77}
78
79// cmv 8/10/99
80//++
81NTuple::NTuple(const NTuple& NT)
82//
83// Createur par copie (clone).
84//--
85: mNVar(0), mNEnt(0), mBlk(0), mNBlk(0)
86, mVar(NULL), mVarD(NULL), mNames(NULL), mInfo(NULL)
87{
88if(NT.mNVar<=0) return; // cas ou NT est cree par defaut
89mNVar = NT.mNVar;
90mBlk = NT.mBlk;
91mVar = new r_4[NT.mNVar];
92mVarD = new r_8[NT.mNVar];
93mNames = new char[NT.mNVar*LENNAME1];
94
95int i;
96r_4* pt = new r_4[mNVar*mBlk];
97mNBlk = 1; mPtr.push_back(pt);
98
99for(i=0;i<mNVar;i++) strcpy(mNames+i*LENNAME1,NT.NomIndex(i));
100
101if(NT.mInfo!=NULL) {mInfo = new DVList; *mInfo = *(NT.mInfo);}
102
103if(NT.mNEnt<=0) return;
104for(i=0;i<NT.mNEnt;i++) {r_4* x=NT.GetVec(i,NULL); Fill(x);}
105
106return;
107}
108
109/* --Methode-- */
110//++
111NTuple::NTuple(char *flnm)
112//
113// Createur lecture fichier ppersist.
114//--
115{
116mNVar = mNEnt = mBlk = mNBlk = 0;
117mVar = NULL;
118mVarD = NULL;
119mNames = NULL;
120mInfo = NULL;
121PInPersist s(flnm);
122ObjFileIO<NTuple> fiont(this);
123fiont.Read(s);
124}
125
126/* --Methode-- */
127NTuple::~NTuple()
128{
129Clean();
130}
131
132/* --Methode-- */
133void NTuple::Clean()
134{
135if (mVar) delete[] mVar;
136if (mVarD) delete[] mVarD;
137if (mNames) delete[] mNames;
138if (mInfo) delete mInfo;
139int i;
140if(mNBlk>0) for(i=0; i<mNBlk; i++) delete[] mPtr[i];
141mPtr.erase(mPtr.begin(), mPtr.end());
142mNVar = mNEnt = mBlk = mNBlk = 0;
143mVar = NULL;
144mVarD = NULL;
145mNames = NULL;
146mInfo = NULL;
147return;
148}
149
150/* --Methode-- cmv 08/10/99 */
151//++
152NTuple& NTuple::operator = (const NTuple& NT)
153//
154// Operateur egal (clone).
155//--
156{
157if(this == &NT) return *this;
158Clean();
159if(NT.mNVar<=0) return *this; // cas ou NT est cree par defaut
160mNVar = NT.mNVar;
161mBlk = NT.mBlk;
162mVar = new r_4[NT.mNVar];
163mVarD = new r_8[NT.mNVar];
164mNames = new char[NT.mNVar*LENNAME1];
165
166int i;
167r_4* pt = new r_4[mNVar*mBlk];
168mNBlk = 1; mPtr.push_back(pt);
169
170for(i=0;i<mNVar;i++) strcpy(mNames+i*LENNAME1,NT.NomIndex(i));
171
172if(NT.mInfo!=NULL) {mInfo = new DVList; *mInfo = *(NT.mInfo);}
173
174if(NT.mNEnt<=0) return *this;
175for(i=0;i<NT.mNEnt;i++) {r_4* x=NT.GetVec(i,NULL); Fill(x);}
176
177// En fait il faudrait un createur par copie qui partage les donnees
178// quand l'objet est temporaire... trop complique A FAIRE ! cmv.
179return *this;
180}
181
182/* --Methode-- */
183//++
184void NTuple::Fill(r_4* x)
185//
186// Remplit le ntuple avec le tableau cd reels `x'.
187//--
188{
189int numb = mNEnt/mBlk;
190if (numb >= mNBlk) {
191 r_4* pt = new r_4[mNVar*mBlk];
192 mNBlk++;
193 mPtr.push_back(pt);
194}
195int offb = mNEnt-numb*mBlk;
196memcpy((mPtr[numb]+offb*mNVar), x, mNVar*sizeof(r_4));
197mNEnt++;
198return;
199}
200
201
202/* --Methode-- */
203//++
204float NTuple::GetVal(int n, int k) const
205//
206// Retourne la valeur de la variable `k' de l'evenement `n'.
207//--
208{
209if (n >= mNEnt) return(BADVAL);
210if ( (k < 0) || (k >= mNVar) ) return(BADVAL);
211int numb = n/mBlk;
212int offb = n-numb*mBlk;
213return(*(mPtr[numb]+offb*mNVar+k));
214}
215
216
217/* --Methode-- */
218//++
219int NTuple::IndexNom(const char* nom) const
220//
221// Retourne le numero de la variable de nom `nom'.
222//--
223{
224int i;
225for(i=0; i<mNVar; i++)
226 if ( strcmp(nom, mNames+i*LENNAME1) == 0) return(i);
227return(-1);
228}
229
230
231static char nomretour[2*LENNAME1];
232/* --Methode-- */
233//++
234char* NTuple::NomIndex(int k) const
235//
236// Retourne le nom de la variable numero 'k'
237//--
238{
239nomretour[0] = '\0';
240if ((k >= 0) && (k < mNVar)) strcpy(nomretour, mNames+k*LENNAME1);
241return(nomretour);
242}
243
244
245/* --Methode-- */
246//++
247r_4* NTuple::GetVec(int n, r_4* ret) const
248//
249// Retourne l'evenement `n' dans le vecteur `ret'.
250//--
251{
252int i;
253if (ret == NULL) ret = mVar;
254if (n >= mNEnt) {
255 for(i=0; i<mNVar; i++) ret[i] = BADVAL;
256 return(ret);
257}
258
259int numb = n/mBlk;
260int offb = n-numb*mBlk;
261memcpy(ret, (mPtr[numb]+offb*mNVar), mNVar*sizeof(r_4));
262return(ret);
263}
264
265/* --Methode-- */
266//++
267r_8* NTuple::GetVecD(int n, r_8* ret) const
268//
269// Retourne l'evenement `n' dans le vecteur `ret'.
270//--
271{
272int i;
273if (ret == NULL) ret = mVarD;
274float * fr = GetVec(n);
275for(i=0; i<mNVar; i++) ret[i] = fr[i];
276return(ret);
277}
278
279
280
281/* --Methode-- */
282//++
283DVList& NTuple::Info()
284//
285// Renvoie une référence sur l'objet DVList Associé
286//--
287{
288if (mInfo == NULL) mInfo = new DVList;
289return(*mInfo);
290}
291
292/* --Methode-- */
293//++
294void NTuple::Print(int num, int nmax) const
295//
296// Imprime `nmax' evenements a partir du numero `num'.
297//--
298{
299int i,j;
300
301printf("Num ");
302for(i=0; i<mNVar; i++) printf("%8s ", mNames+i*LENNAME1);
303putchar('\n');
304
305if (nmax <= 0) nmax = 1;
306if (num < 0) num = 0;
307nmax += num;
308if (nmax > mNEnt) nmax = mNEnt;
309for(i=num; i<nmax; i++) {
310 GetVec(i, NULL);
311 printf("%6d ", i);
312 for(j=0; j<mNVar; j++) printf("%8g ", (float)mVar[j]);
313 putchar('\n');
314}
315return;
316}
317
318/* --Methode-- */
319//++
320void NTuple::Show(ostream& os) const
321//
322// Imprime l'information generale sur le ntuple.
323//--
324{
325os << "NTuple: NVar= " << mNVar << " NEnt=" << mNEnt
326 << " (Blk Sz,Nb= " << mBlk << " ," << mNBlk << ")\n";
327os << " Variables Min Max \n";
328int i;
329double min, max;
330char buff[128];
331for(i=0; i<mNVar; i++) {
332 GetMinMax(i, min, max);
[1092]333 sprintf(buff, "%3d %16s %10g %10g \n", i, mNames+i*LENNAME1, min, max);
[763]334 os << (string)buff ;
335 }
336os << endl;
337}
338
339
340/* --Methode-- */
341//++
342int NTuple::FillFromASCIIFile(string const& fn, float defval)
343//
344// Remplit le ntuple a partir d'un fichier ASCII.
345// Renvoie le nombre de lignes ajoutees.
346//--
347{
348if (NbColumns() < 1) {
349 cout << "NTuple::FillFromASCIIFile() Ntuple has " << NbColumns() << " columns" << endl;
350 return(-1);
351 }
352FILE * fip = fopen(fn.c_str(), "r");
353if (fip == NULL) {
354 cout << "NTuple::FillFromASCIIFile() Error opening file " << fn << endl;
355 return(-2);
356 }
357
358char lineb[4096];
359char *line;
360char* ccp;
[809]361int j,kk;
[763]362int postab, posb;
363float* xv = new float[NbColumns()];
364
365int nlread = 0;
366int nvar = NbColumns();
367// On boucle sur toutes les lignes
368while (fgets(lineb,4096,fip) != NULL) {
369 lineb[4095] = '\0';
370 j = 0; line = lineb;
371// On enleve les espaces et tab de debut
372 while ( (line[j] != '\0') && ((line[j] == ' ') || (line[j] == '\t')) ) j++;
373 line = lineb+j;
374// Il faut que le premier caractere non-espace soit un digit, ou + ou - ou .
375 if (!( isdigit(line[0]) || (line[0] == '+') || (line[0] == '-') || (line[0] == '.') )) continue;
376 ccp = line;
377 for(kk=0; kk<nvar; kk++) xv[kk] = defval;
378 for(kk=0; kk<nvar; kk++) {
379// Les mots sont separes par des espaces ou des tab
380 postab = posc(ccp, '\t' );
381 posb = posc(ccp, ' ' );
382 if (postab >= 0) {
383 if (posb < 0) posb = postab;
384 else if (postab < posb) posb = postab;
385 }
386 if (posb >= 0) ccp[posb] = '\0';
387 if ( isdigit(line[0]) || (line[0] == '+') || (line[0] == '-') || (line[0] == '.') )
388 xv[kk] = atof(ccp);
389 if (posb < 0) break;
390 ccp += posb+1; j = 0;
391 while ( (ccp[j] != '\0') && ((ccp[j] == ' ') || (ccp[j] == '\t')) ) j++;
392 ccp += j;
393 }
394 Fill(xv);
395 nlread++;
396 }
397
398delete[] xv;
399cout << "NTuple::FillFromASCIIFile( " << fn << ") " << nlread << " fill from file " << endl;
400return(nlread);
401}
402
403
404// ------- Implementation de l interface NTuple ---------
405
406/* --Methode-- */
407uint_4 NTuple::NbLines() const
408{
409return(NEntry());
410}
411/* --Methode-- */
412uint_4 NTuple::NbColumns() const
413{
414return(NVar());
415}
416
417/* --Methode-- */
418r_8 * NTuple::GetLineD(int n) const
419{
420return(GetVecD(n));
421}
422
423/* --Methode-- */
424r_8 NTuple::GetCell(int n, int k) const
425{
426return(GetVal(n, k));
427}
428
429/* --Methode-- */
430r_8 NTuple::GetCell(int n, string const & nom) const
431{
432return(GetVal(n, nom.c_str()));
433}
434
435/* --Methode-- */
436//++
437void NTuple::GetMinMax(int k, double& min, double& max) const
438//
439// Retourne le minimum et le maximum de la variable `k'.
440//--
441{
442min = 9.e19; max = -9.e19;
443if ( (k < 0) || (k >= mNVar) ) return;
444int jb,ib,i;
445double x;
446i=0;
447for(jb=0; jb< mNBlk; jb++)
448 for(ib=0; ib< mBlk; ib++) {
449 if (i >= mNEnt) break;
450 i++;
451 x = *(mPtr[jb]+ib*mNVar+k);
452 if(i==1) {min = x; max = x;}
453 if (x < min) min = x;
454 if (x > max) max = x;
455 }
456return;
457}
458
459/* --Methode-- */
460void NTuple::GetMinMax(string const & nom, double& min, double& max) const
461{
462GetMinMax(IndexNom(nom.c_str()), min, max);
463}
464
465/* --Methode-- */
466int NTuple::ColumnIndex(string const & nom) const
467{
468return(IndexNom(nom.c_str()));
469}
470
471/* --Methode-- */
472string NTuple::ColumnName(int k) const
473{
474return(NomIndex(k));
475}
476
477/* --Methode-- */
478//++
479string NTuple::VarList_C(const char* nomx) const
480//
481// Retourne une chaine de caracteres avec la declaration des noms de
482// variables. si "nomx!=NULL" , des instructions d'affectation
483// a partir d'un tableau "nomx[i]" sont ajoutees.
484//--
485{
486string rets="";
487int i;
488for(i=0; i<mNVar; i++) {
489 if ( (i%5 == 0) && (i > 0) ) rets += ";";
490 if (i%5 == 0) rets += "\ndouble ";
491 else rets += ",";
492 rets += mNames+i*LENNAME1;
493 }
494rets += "; \n";
495if (nomx) {
496 char buff[256];
497 for(i=0; i<mNVar; i++) {
498 sprintf(buff,"%s=%s[%d]; ", mNames+i*LENNAME1, nomx, i);
499 rets += buff;
500 if ( (i%3 == 0) && (i > 0) ) rets += "\n";
501 }
502 }
503
504return(rets);
505}
506
507
508/* --Methode-- */
509//++
510string NTuple::LineHeaderToString() const
511// Retourne une chaine de caracteres avec la liste des noms de
512// variables, utilisables pour une impression
513//--
514{
515char buff[32];
516string rets=" Num ";
517for(int i=0; i<mNVar; i++) {
518 sprintf(buff, "%8s ", mNames+i*LENNAME1);
519 rets += buff;
520 }
521rets += '\n';
522return(rets);
523}
524
525/* --Methode-- */
526//++
527string NTuple::LineToString(int n) const
528// Retourne une chaine de caracteres avec le contenu de la ligne "n"
529// utilisable pour une impression
530//--
531{
532char buff[32];
533double* val;
534val = GetLineD(n);
535sprintf(buff,"%6d: ",n);
536string rets=buff;
537int i;
538for(i=0; i<mNVar; i++) {
539 sprintf(buff, "%8.3g ", val[i]);
540 rets += buff;
541 }
542rets += '\n';
543return(rets);
544}
545
546
[846]547/* --Methode-- */
[763]548//++
549void ObjFileIO<NTuple>::WriteSelf(POutPersist& s) const
550//
551// Ecriture ppersist du ntuple.
552//--
553{
[1155]554if (dobj == NULL) return;
555
556// On ecrit cette chaine pour compatibilite avec l'ancienne version
557string strg = "NTuple";
558s.Put(strg);
559
560// On ecrit 3 uint_4 ....
561// 0: Numero de version, 1 : non nul -> has info, 2 : reserve
562uint_4 itab[3];
563itab[0] = 2; // Numero de version a 1
564itab[1] = itab[2] = 0;
565if (dobj->mInfo) itab[1] = 1;
566s.Put(itab, 3);
567
568s.Put(dobj->mNVar);
[763]569s.PutBytes(dobj->mNames, dobj->mNVar*LENNAME1);
[1155]570s.Put(dobj->mNEnt);
571s.Put(dobj->mBlk);
572s.Put(dobj->mNBlk);
[763]573if (dobj->mInfo) s << (*(dobj->mInfo));
574int jb;
575for(jb=0; jb<dobj->mNBlk; jb++)
[1155]576 s.Put(dobj->mPtr[jb], dobj->mNVar*dobj->mBlk);
[763]577return;
578}
579
580/* --Methode-- */
581//++
582void ObjFileIO<NTuple>::ReadSelf(PInPersist& s)
583//
584// Lecture ppersist du ntuple.
585//--
586{
587
[1155]588if (dobj == NULL) dobj = new NTuple;
589else dobj->Clean();
[763]590
591bool hadinfo = false;
[1155]592string strg;
593s.Get(strg);
594if (strg == "NTuple") {
595 uint_4 itab[3];
596 s.Get(itab, 3);
597 if (itab[1] != 0) hadinfo = true;
598}
599else {
600// Ancienne version de PPF NTuple - Pour savoir s'il y avait un DVList Info associe
601 char buff[256];
602 strcpy(buff, strg.c_str());
603 if (strncmp(buff+strlen(buff)-7, "HasInfo", 7) == 0) hadinfo = true;
604}
605s.Get(dobj->mNVar);
[763]606dobj->mNames = new char[dobj->mNVar*LENNAME1];
607dobj->mVar = new r_4[dobj->mNVar];
608dobj->mVarD = new r_8[dobj->mNVar];
609s.GetBytes(dobj->mNames, dobj->mNVar*LENNAME1);
[1155]610s.Get(dobj->mNEnt);
611s.Get(dobj->mBlk);
612s.Get(dobj->mNBlk);
[763]613
614if (hadinfo) { // Lecture eventuelle du DVList Info
615 if (dobj->mInfo == NULL) dobj->mInfo = new DVList;
616 s >> (*(dobj->mInfo));
617 }
618
619int jb;
620for(jb=0; jb<dobj->mNBlk; jb++) {
621 r_4* pt = new r_4[dobj->mNVar*dobj->mBlk];
622 dobj->mPtr.push_back(pt);
[1155]623 s.Get(dobj->mPtr[jb], dobj->mNVar*dobj->mBlk);
[763]624}
625
626}
627
[846]628
[763]629#ifdef __CXX_PRAGMA_TEMPLATES__
630#pragma define_template ObjFileIO<NTuple>
631#endif
632
633#if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES)
634template class ObjFileIO<NTuple>;
635#endif
Note: See TracBrowser for help on using the repository browser.