source: Sophya/trunk/SophyaLib/NTools/ntuple.cc@ 278

Last change on this file since 278 was 278, checked in by ansari, 26 years ago

Implementation PPersist par deleguation pour DVList Reza 28/04/99

File size: 6.6 KB
Line 
1#include <stdio.h>
2#include <string.h>
3
4#include "perrors.h"
5#include "ntuple.h"
6
7
8#define BADVAL -1.e19
9#define LENNAME 8
10#define LENNAME1 (LENNAME+1)
11
12//++
13// Class NTuple
14// Lib Outils++
15// include ntuple.h
16//
17// Classe de ntuples
18//--
19
20/* --Methode-- */
21//++
22NTuple::NTuple()
23//
24// Createur par defaut
25//--
26{
27mNVar = mNEnt = mBlk = mNBlk = 0;
28mVar = NULL;
29mNames = NULL;
30mInfo = NULL;
31}
32
33
34//++
35NTuple::NTuple(int nvar, char** noms, int blk)
36//
37// Createur d'un ntuple de `nvar' variables dont les
38// noms sont dans le tableau de cahines de caracteres `noms'
39// avec `blk' d'evenements par blocks.
40//--
41{
42mNVar = mNEnt = mBlk = mNBlk = 0;
43mVar = NULL;
44mNames = NULL;
45mInfo = NULL;
46if (nvar <= 0) THROW(sizeMismatchErr);
47mNVar = nvar;
48mVar = new r_4[nvar];
49if (blk < 10) blk = 10;
50mBlk = blk;
51// On prend des noms de LENNAME char pour le moment
52mNames = new char[nvar*LENNAME1];
53r_4* pt = new r_4[nvar*blk];
54mNBlk = 1;
55mPtr.push_back(pt);
56int i;
57for(i=0; i<nvar; i++)
58 { strncpy(mNames+i*LENNAME1, noms[i], LENNAME);
59 mNames[i*LENNAME1+LENNAME] = '\0'; }
60return;
61}
62
63/* --Methode-- */
64//++
65NTuple::NTuple(char *flnm)
66//
67// Createur lecture fichier ppersist.
68//--
69{
70mNVar = mNEnt = mBlk = mNBlk = 0;
71mVar = NULL;
72mNames = NULL;
73mInfo = NULL;
74PInPersist s(flnm);
75Read(s);
76}
77
78/* --Methode-- */
79NTuple::~NTuple()
80{
81Clean();
82}
83
84/* --Methode-- */
85void NTuple::Clean()
86{
87if (mVar) delete[] mVar;
88if (mNames) delete[] mNames;
89if (mInfo) delete mInfo;
90int i;
91for(i=0; i<mNBlk; i++) delete[] mPtr[i];
92mPtr.erase(mPtr.begin(), mPtr.end());
93mNVar = mNEnt = mBlk = mNBlk = 0;
94mVar = NULL;
95mNames = NULL;
96return;
97}
98
99/* --Methode-- */
100//++
101void NTuple::Fill(r_4* x)
102//
103// Remplit le ntuple avec le tableau cd reels `x'.
104//--
105{
106int numb = mNEnt/mBlk;
107if (numb >= mNBlk) {
108 r_4* pt = new r_4[mNVar*mBlk];
109 mNBlk++;
110 mPtr.push_back(pt);
111}
112int offb = mNEnt-numb*mBlk;
113memcpy((mPtr[numb]+offb*mNVar), x, mNVar*sizeof(r_4));
114mNEnt++;
115return;
116}
117
118
119/* --Methode-- */
120//++
121float NTuple::GetVal(int n, int k) const
122//
123// Retourne la valeur de la variable `k' de l'evenement `n'.
124//--
125{
126if (n >= mNEnt) return(BADVAL);
127if ( (k < 0) || (k >= mNVar) ) return(BADVAL);
128int numb = n/mBlk;
129int offb = n-numb*mBlk;
130return(*(mPtr[numb]+offb*mNVar+k));
131}
132
133
134/* --Methode-- */
135//++
136int NTuple::IndexNom(const char* nom) const
137//
138// Retourne le numero de la variable de nom `nom'.
139//--
140{
141int i;
142for(i=0; i<mNVar; i++)
143 if ( strcmp(nom, mNames+i*LENNAME1) == 0) return(i);
144return(-1);
145}
146
147
148static char nomretour[2*LENNAME1];
149/* --Methode-- */
150//++
151char* NTuple::NomIndex(int k) const
152//
153// Retourne le nom de la variable numero 'k'
154//--
155{
156nomretour[0] = '\0';
157if ((k >= 0) && (k < mNVar)) strcpy(nomretour, mNames+k*LENNAME1);
158return(nomretour);
159}
160
161/* --Methode-- */
162//++
163string NTuple::VarList_C(const char* nomx) const
164//
165// Retourne une chaine de caracteres avec la declaration des noms de
166// variables. si "nomx!=NULL" , des instructions d'affectation
167// a partir d'un tableau "nomx[i]" sont ajoutees.
168//--
169{
170string rets;
171int i;
172for(i=0; i<mNVar; i++) {
173 if ( (i%5 == 0) && (i > 0) ) rets += ";";
174 if (i%5 == 0) rets += "\ndouble ";
175 else rets += ",";
176 rets += mNames+i*LENNAME1;
177 }
178rets += "; \n";
179if (nomx) {
180 char buff[256];
181 for(i=0; i<mNVar; i++) {
182 sprintf(buff,"%s=%s[%d]; ", mNames+i*LENNAME1, nomx, i);
183 rets += buff;
184 if ( (i%3 == 0) && (i > 0) ) rets += "\n";
185 }
186 }
187
188return(rets);
189}
190
191
192/* --Methode-- */
193//++
194r_4* NTuple::GetVec(int n, r_4* ret) const
195//
196// Retourne l'evenement `n' dans le vecteur `ret'.
197//--
198{
199int i;
200if (ret == NULL) ret = mVar;
201if (n >= mNEnt) {
202 for(i=0; i<mNVar; i++) ret[i] = BADVAL;
203 return(ret);
204}
205
206int numb = n/mBlk;
207int offb = n-numb*mBlk;
208memcpy(ret, (mPtr[numb]+offb*mNVar), mNVar*sizeof(r_4));
209return(ret);
210}
211
212/* --Methode-- */
213//++
214void NTuple::GetMinMax(int k, float& min, float& max) const
215//
216// Retourne le minimum et le maximum de la variable `k'.
217//--
218{
219min = 9.e19; max = -9.e19;
220if ( (k < 0) || (k >= mNVar) ) return;
221int jb,ib,i;
222float x;
223i=0;
224for(jb=0; jb< mNBlk; jb++)
225 for(ib=0; ib< mBlk; ib++) {
226 if (i >= mNEnt) break;
227 i++;
228 x = *(mPtr[jb]+ib*mNVar+k);
229 if(i==1) {min = x; max = x;}
230 if (x < min) min = x;
231 if (x > max) max = x;
232 }
233return;
234}
235
236/* --Methode-- */
237//++
238DVList& NTuple::Info()
239//
240// Renvoie une référence sur l'objet DVList Associé
241//--
242{
243if (mInfo == NULL) mInfo = new DVList;
244return(*mInfo);
245}
246
247/* --Methode-- */
248//++
249void NTuple::Print(int num, int nmax) const
250//
251// Imprime `nmax' evenements a partir du numero `num'.
252//--
253{
254int i,j;
255
256printf("Num ");
257for(i=0; i<mNVar; i++) printf("%8s ", mNames+i*LENNAME1);
258putchar('\n');
259
260if (nmax <= 0) nmax = 1;
261if (num < 0) num = 0;
262nmax += num;
263if (nmax > mNEnt) nmax = mNEnt;
264for(i=num; i<nmax; i++) {
265 GetVec(i, NULL);
266 printf("%6d ", i);
267 for(j=0; j<mNVar; j++) printf("%8g ", (float)mVar[j]);
268 putchar('\n');
269}
270return;
271}
272
273/* --Methode-- */
274//++
275void NTuple::Show(ostream& os) const
276//
277// Imprime l'information generale sur le ntuple.
278//--
279{
280os << "NTuple: NVar= " << mNVar << " NEnt=" << mNEnt
281 << " (Blk Sz,Nb= " << mBlk << " ," << mNBlk << ")\n";
282os << " Variables Min Max \n";
283int i;
284float min, max;
285char buff[128];
286for(i=0; i<mNVar; i++) {
287 GetMinMax(i, min, max);
288 sprintf(buff, "%3d %16s %10g %10g \n", i, mNames+i*LENNAME1, min, max);
289 os << (string)buff ;
290 }
291os << endl;
292}
293
294
295/* --Methode-- */
296//++
297void NTuple::WriteSelf(POutPersist& s) const
298//
299// Ecriture ppersist du ntuple.
300//--
301{
302char strg[256];
303if (mInfo) sprintf(strg, "NVar=%6d NEnt=%9d BlkSz=%6d NBlk=%6d HasInfo",
304 (int)mNVar, (int)mNEnt, (int)mBlk, (int)mNBlk);
305else sprintf(strg, "NVar=%6d NEnt=%9d BlkSz=%6d NBlk=%6d ",
306 (int)mNVar, (int)mNEnt, (int)mBlk, (int)mNBlk);
307s.PutLine(strg);
308s.PutI4(mNVar);
309s.PutBytes(mNames, mNVar*LENNAME1);
310s.PutI4(mNEnt);
311s.PutI4(mBlk);
312s.PutI4(mNBlk);
313if (mInfo) s << (*mInfo);
314int jb;
315for(jb=0; jb<mNBlk; jb++)
316 s.PutR4s(mPtr[jb], mNVar*mBlk);
317return;
318}
319
320/* --Methode-- */
321//++
322void NTuple::ReadSelf(PInPersist& s)
323//
324// Lecture ppersist du ntuple.
325//--
326{
327
328Clean();
329
330char strg[256];
331s.GetLine(strg, 255);
332// Pour savoir s'il y avait un DVList Info associe
333bool hadinfo = false;
334if (strncmp(strg+strlen(strg)-7, "HasInfo", 7) == 0) hadinfo = true;
335
336s.GetI4(mNVar);
337mNames = new char[mNVar*LENNAME1];
338mVar = new r_4[mNVar];
339s.GetBytes(mNames, mNVar*LENNAME1);
340s.GetI4(mNEnt);
341s.GetI4(mBlk);
342s.GetI4(mNBlk);
343
344if (hadinfo) { // Lecture eventuelle du DVList Info
345 if (mInfo == NULL) mInfo = new DVList;
346 s >> (*mInfo);
347 }
348
349int jb;
350for(jb=0; jb<mNBlk; jb++) {
351 r_4* pt = new r_4[mNVar*mBlk];
352 mPtr.push_back(pt);
353 s.GetR4s(mPtr[jb], mNVar*mBlk);
354}
355
356}
Note: See TracBrowser for help on using the repository browser.