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