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