source: trunk/src/napi5.c @ 1822

Revision 1822, 66.4 KB checked in by Freddie Akeroyd, 3 months ago (diff)

Fix NXgetdata and NXgetinfo for HDF5 Variable length string datasets.
Refs #329

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
Line 
1/*---------------------------------------------------------------------------
2  NeXus - Neutron & X-ray Common Data Format
3 
4  Application Program Interface (HDF5) Routines
5 
6  Copyright (C) 1997-2011 Mark Koennecke, Przemek Klosowski
7 
8  This library is free software; you can redistribute it and/or
9  modify it under the terms of the GNU Lesser General Public
10  License as published by the Free Software Foundation; either
11  version 2 of the License, or (at your option) any later version.
12 
13  This library is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
16  Lesser General Public License for more details.
17 
18  You should have received a copy of the GNU Lesser General Public
19  License along with this library; if not, write to the Free Software
20  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
21             
22  For further information, see <http://www.nexusformat.org>
23
24----------------------------------------------------------------------------*/
25
26#ifdef HDF5
27
28#include <stdlib.h>
29#include <assert.h>
30#include <string.h>
31#include <time.h>
32
33#include "napi.h"
34#include "napi5.h"
35
36#ifdef H5_VERSION_GE
37#if !H5_VERSION_GE(1,8,0)
38#error HDF5 Version must be 1.8.0 or higher
39#endif
40#endif
41
42#ifdef _MSC_VER
43#define snprintf _snprintf
44#endif /* _MSC_VER */
45
46#define NX_UNKNOWN_GROUP "" /* for when no NX_class attr */
47
48extern  void *NXpData;
49
50  typedef struct __NexusFile5 {
51        struct iStack5 {
52          char irefn[1024];
53          int iVref;
54          hsize_t iCurrentIDX;
55        } iStack5[NXMAXSTACK];
56        struct iStack5 iAtt5;
57        int iFID;
58        int iCurrentG;
59        int iCurrentD;
60        int iCurrentS;
61        int iCurrentT;
62        int iCurrentA;
63        int iNX;
64        int iNXID;
65        int iStackPtr;
66        char *iCurrentLGG;
67        char *iCurrentLD;
68        char name_ref[1024];
69        char name_tmp[1024];
70        char iAccess[2];
71  } NexusFile5, *pNexusFile5;
72
73/*
74   forward declaration of NX5closegroup in order to get rid of a nasty
75   warning
76*/
77NXstatus  NX5closegroup (NXhandle fid);
78/*-------------------------------------------------------------------*/
79
80  static pNexusFile5 NXI5assert(NXhandle fid)
81  {
82    pNexusFile5 pRes;
83 
84    assert(fid != NULL);
85    pRes = (pNexusFile5)fid;
86    assert(pRes->iNXID == NX5SIGNATURE);
87    return pRes;
88  }
89 
90  /*--------------------------------------------------------------------*/
91
92   static void NXI5KillDir (pNexusFile5 self)
93  {
94    self->iStack5[self->iStackPtr].iCurrentIDX = 0;
95  }
96
97static herr_t readStringAttribute(hid_t attr, char** data)
98{
99     int iSize;
100     herr_t iRet;
101     hid_t atype = -1, btype = -1;
102
103        atype = H5Aget_type(attr);
104        *data = NULL;
105        if ( H5Tis_variable_str(atype) )
106        {
107            btype = H5Tget_native_type(atype, H5T_DIR_ASCEND); 
108            iRet = H5Aread(attr, btype, data);
109            H5Tclose(btype);
110        }
111        else
112        {
113            iSize = H5Tget_size(atype);
114            *data = malloc(iSize+1);
115            iRet = H5Aread(attr, atype, *data);
116            (*data)[iSize] = '\0';
117        }
118        H5Tclose(atype);
119    return iRet;
120}
121 
122static herr_t readStringAttributeN(hid_t attr, char* data, int maxlen)
123{
124    herr_t iRet;
125    char* vdat = NULL;
126    iRet = readStringAttribute(attr, &vdat);
127    if (iRet >= 0)
128    {
129        strncpy(data, vdat, maxlen);
130        free(vdat);
131    }
132    data[maxlen-1] = '\0';
133    return iRet;
134}
135  /*--------------------------------------------------------------------*/
136
137  static void NXI5KillAttDir (pNexusFile5 self)
138  {
139    self->iAtt5.iCurrentIDX = 0;
140  }
141/*---------------------------------------------------------------------*/
142static void buildCurrentPath(pNexusFile5 self, char *pathBuffer, 
143                             int pathBufferLen){
144
145  memset(pathBuffer,0,pathBufferLen);
146  if(self->iCurrentG != 0) {
147    strcpy(pathBuffer,"/");
148    if(strlen(self->name_ref) + 1 < pathBufferLen){
149      strcat(pathBuffer, self->name_ref);
150    }
151  }
152  if(self->iCurrentD != 0){
153    strcat(pathBuffer,"/");
154    if(strlen(self->iCurrentLD) + strlen(pathBuffer) < pathBufferLen){
155      strcat(pathBuffer,self->iCurrentLD);
156    }
157  }
158}
159  /* ----------------------------------------------------------------------
160 
161                          Definition of NeXus API
162
163   ---------------------------------------------------------------------*/
164
165  NXstatus  NX5reopen(NXhandle pOrigHandle, NXhandle* pNewHandle) 
166{
167        pNexusFile5 pNew = NULL, pOrig = NULL;
168    *pNewHandle = NULL;
169       pOrig = (pNexusFile5)pOrigHandle;
170    pNew = (pNexusFile5) malloc (sizeof (NexusFile5));
171    if (!pNew) {
172      NXReportError ("ERROR: no memory to create File datastructure");
173      return NX_ERROR;
174    }
175    memset (pNew, 0, sizeof (NexusFile5));
176    pNew->iFID = H5Freopen (pOrig->iFID);
177    if (pNew->iFID <= 0) {
178        NXReportError ("cannot clone file");
179                free (pNew);
180      return NX_ERROR;
181       }
182       strcpy(pNew->iAccess, pOrig->iAccess);
183    pNew->iNXID = NX5SIGNATURE;
184    pNew->iStack5[0].iVref = 0;    /* root! */
185    *pNewHandle = (NXhandle)pNew;
186    return NX_OK;
187}
188
189NXstatus  NX5open(CONSTCHAR *filename, NXaccess am, 
190                                 NXhandle* pHandle)
191  {
192  hid_t attr1,aid1, aid2, iVID;
193  pNexusFile5 pNew = NULL;
194  char pBuffer[512];
195  char *time_buffer = NULL;
196  char version_nr[10];
197  unsigned int vers_major, vers_minor, vers_release, am1 ;
198  hid_t fapl = -1;     
199  int mdc_nelmts;
200  size_t rdcc_nelmts;
201  size_t rdcc_nbytes;
202  double rdcc_w0;
203  unsigned hdf5_majnum, hdf5_minnum, hdf5_relnum;
204
205  *pHandle = NULL;
206
207  if (H5get_libversion(&hdf5_majnum, &hdf5_minnum, &hdf5_relnum) < 0)
208  {
209      NXReportError("ERROR: cannot determine HDF5 library version");
210      return NX_ERROR;
211  }
212  if (hdf5_majnum == 1 && hdf5_minnum < 8)
213  {
214      NXReportError("ERROR: HDF5 library 1.8.0 or higher required");
215      return NX_ERROR;
216  }
217
218  /* mask of any options for now */
219  am = (NXaccess)(am & NXACCMASK_REMOVEFLAGS);
220
221  /* turn off the automatic HDF error handling */ 
222      H5Eset_auto(H5E_DEFAULT, NULL, NULL); 
223#ifdef USE_FTIME
224    struct timeb timeb_struct;
225#endif
226
227
228    pNew = (pNexusFile5) malloc (sizeof (NexusFile5));
229    if (!pNew) {
230      NXReportError("ERROR: not enough memory to create file structure");
231      return NX_ERROR;
232    }
233    memset (pNew, 0, sizeof (NexusFile5));
234
235
236    /* start HDF5 interface */
237    if (am == NXACC_CREATE5) { 
238       fapl = H5Pcreate(H5P_FILE_ACCESS);
239       H5Pget_cache(fapl,&mdc_nelmts,&rdcc_nelmts,&rdcc_nbytes,&rdcc_w0);
240       rdcc_nbytes=(size_t)nx_cacheSize;
241       H5Pset_cache(fapl,mdc_nelmts,rdcc_nelmts,rdcc_nbytes,rdcc_w0);
242       H5Pset_fclose_degree(fapl,H5F_CLOSE_STRONG);
243       am1 = H5F_ACC_TRUNC;
244       pNew->iFID = H5Fcreate (filename, am1, H5P_DEFAULT, fapl);
245    } else {
246       if (am == NXACC_READ) {
247          am1 = H5F_ACC_RDONLY;
248        } else {
249          am1 = H5F_ACC_RDWR;
250        }   
251        fapl = H5Pcreate(H5P_FILE_ACCESS);
252        H5Pset_fclose_degree(fapl,H5F_CLOSE_STRONG);
253        pNew->iFID = H5Fopen (filename, am1, fapl);
254    } 
255    if(fapl != -1) {
256      H5Pclose(fapl);
257    }
258    if (pNew->iFID <= 0) {
259      sprintf (pBuffer, "ERROR: cannot open file: %s", filename);
260      NXReportError( pBuffer);
261      free (pNew);
262      return NX_ERROR;
263    }
264
265/*
266 * need to create global attributes         file_name file_time NeXus_version
267 * at some point for new files
268 */
269    if (am1 != H5F_ACC_RDONLY) 
270    {
271        iVID = H5Gopen(pNew->iFID,"/", H5P_DEFAULT);
272        aid2 = H5Screate(H5S_SCALAR);
273        aid1 = H5Tcopy(H5T_C_S1);
274        H5Tset_size(aid1, strlen(NEXUS_VERSION));
275        if (am1 == H5F_ACC_RDWR)
276        {
277        H5Adelete(iVID, "NeXus_version"); 
278        }
279        attr1= H5Acreate(iVID, "NeXus_version", aid1, aid2, H5P_DEFAULT, H5P_DEFAULT);
280        if (attr1<0)
281        {
282          NXReportError( 
283                         "ERROR: failed to store NeXus_version attribute ");
284          return NX_ERROR;
285        } 
286        if (H5Awrite(attr1, aid1,NEXUS_VERSION)<0)
287        {
288          NXReportError( 
289                         "ERROR: failed to store NeXus_version attribute ");
290          return NX_ERROR;
291        }
292       /* Close attribute dataspace  */
293       H5Tclose(aid1);
294       H5Sclose(aid2); 
295       /* Close attribute */
296       H5Aclose(attr1); 
297       H5Gclose(iVID); 
298    }
299    if (am1 == H5F_ACC_TRUNC) 
300    {
301        iVID = H5Gopen(pNew->iFID,"/", H5P_DEFAULT);
302        aid2 = H5Screate(H5S_SCALAR);
303        aid1 = H5Tcopy(H5T_C_S1);
304        H5Tset_size(aid1, strlen(filename));
305        attr1= H5Acreate(iVID, "file_name", aid1, aid2, H5P_DEFAULT, H5P_DEFAULT);
306        if (attr1 < 0)
307        {
308          NXReportError( 
309                          "ERROR: failed to store file_name attribute ");
310          return NX_ERROR;
311        }
312        if (H5Awrite(attr1, aid1, (char*)filename) < 0)
313        {
314          NXReportError( 
315                          "ERROR: failed to store file_name attribute ");
316          return NX_ERROR;
317        }
318        H5Tclose(aid1);
319        H5Sclose(aid2); 
320        H5Aclose(attr1);
321        /*  ------- library version ------*/
322        H5get_libversion(&vers_major, &vers_minor, &vers_release);
323        sprintf (version_nr, "%d.%d.%d", vers_major,vers_minor,vers_release);
324        aid2=H5Screate(H5S_SCALAR);
325        aid1 = H5Tcopy(H5T_C_S1);
326        H5Tset_size(aid1, strlen(version_nr));
327        attr1= H5Acreate(iVID, "HDF5_Version", aid1, aid2, H5P_DEFAULT, H5P_DEFAULT);
328        if (attr1 < 0)
329        {
330          NXReportError( 
331                          "ERROR: failed to store file_name attribute ");
332          return NX_ERROR;
333        }
334        if (H5Awrite(attr1, aid1, (char*)version_nr) < 0)
335        {
336          NXReportError( 
337                          "ERROR: failed to store file_name attribute ");
338          return NX_ERROR;
339        }
340        H5Tclose(aid1);
341        H5Sclose(aid2); 
342        H5Aclose(attr1);
343        /*----------- file time */
344        time_buffer = NXIformatNeXusTime();
345        if(time_buffer != NULL){
346          aid2=H5Screate(H5S_SCALAR);
347          aid1 = H5Tcopy(H5T_C_S1);
348          H5Tset_size(aid1, strlen(time_buffer));
349          attr1=H5Acreate(iVID, "file_time", aid1, aid2, H5P_DEFAULT, H5P_DEFAULT);
350          if (attr1 < 0)
351          {
352              NXReportError( 
353                          "ERROR: failed to store file_time attribute ");
354              free(time_buffer);
355              return NX_ERROR;
356          }
357          if (H5Awrite(attr1, aid1, time_buffer) < 0)
358          {
359              NXReportError( 
360                          "ERROR: failed to store file_time attribute ");
361              free(time_buffer);
362              return NX_ERROR;
363          }
364          /* Close attribute dataspace */
365          H5Tclose(aid1);
366          H5Sclose(aid2); 
367          /* Close attribute */
368          H5Aclose(attr1); 
369          free(time_buffer);
370        }
371        H5Gclose(iVID); 
372    }
373    /* Set HDFgroup access mode */
374    if (am1 == H5F_ACC_RDONLY) {
375      strcpy(pNew->iAccess,"r");
376    } else {
377      strcpy(pNew->iAccess,"w");
378    }
379    pNew->iNXID = NX5SIGNATURE;
380    pNew->iStack5[0].iVref = 0;    /* root! */
381    *pHandle = (NXhandle)pNew;
382     return NX_OK;
383  }
384 
385  /* ------------------------------------------------------------------------- */
386
387  NXstatus  NX5close (NXhandle* fid)
388  {
389    pNexusFile5 pFile = NULL;
390    herr_t iRet;
391 
392    pFile=NXI5assert(*fid);
393
394    iRet=0;
395    /*
396    printf("HDF5 object count before close: %d\n",
397           H5Fget_obj_count(pFile->iFID,H5F_OBJ_ALL));
398    */
399    iRet = H5Fclose(pFile->iFID);
400
401    /*
402      leave this here: it helps in debugging leakage problems   
403    printf("HDF5 object count after close: %d\n",
404           H5Fget_obj_count(H5F_OBJ_ALL,H5F_OBJ_ALL));
405    printf("HDF5 dataset count after close: %d\n",
406           H5Fget_obj_count(H5F_OBJ_ALL,H5F_OBJ_DATASET));
407    printf("HDF5 group count after close: %d\n",
408           H5Fget_obj_count(H5F_OBJ_ALL,H5F_OBJ_GROUP));
409    printf("HDF5 datatype count after close: %d\n",
410           H5Fget_obj_count(H5F_OBJ_ALL,H5F_OBJ_DATATYPE));
411    printf("HDF5 attribute count after close: %d\n",
412           H5Fget_obj_count(H5F_OBJ_ALL,H5F_OBJ_ATTR));
413    */
414
415    if (iRet < 0) {
416      NXReportError( "ERROR: cannot close HDF file");
417    }
418    /* release memory */
419    NXI5KillDir (pFile);
420    if(pFile->iCurrentLGG != NULL){
421      free(pFile->iCurrentLGG);
422    }
423    if(pFile->iCurrentLD != NULL){
424      free(pFile->iCurrentLD);
425    }
426    free (pFile);
427    *fid = NULL;
428    H5garbage_collect();
429    return NX_OK;
430  }
431
432 /*-----------------------------------------------------------------------*/   
433
434  NXstatus  NX5makegroup (NXhandle fid, CONSTCHAR *name, CONSTCHAR *nxclass) 
435  {
436    pNexusFile5 pFile;
437    herr_t iRet;
438    hid_t iVID;
439    hid_t attr1,aid1, aid2;
440    char pBuffer[1024] = "";
441   
442    pFile = NXI5assert (fid);
443    /* create and configure the group */
444    if (pFile->iCurrentG==0) {
445       snprintf(pBuffer,1023,"/%s",name);
446    } else {
447       snprintf(pBuffer,1023,"/%s/%s",pFile->name_ref,name);
448    }   
449    iRet = H5Gcreate(pFile->iFID,(const char*)pBuffer, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
450    if (iRet < 0) {
451      NXReportError( "ERROR: could not create Group");
452      return NX_ERROR;
453    }
454    iVID = iRet;
455    aid2 = H5Screate(H5S_SCALAR);
456    aid1 = H5Tcopy(H5T_C_S1);
457    H5Tset_size(aid1, strlen(nxclass));
458    attr1= H5Acreate(iVID, "NX_class", aid1, aid2, H5P_DEFAULT, H5P_DEFAULT);
459    if (attr1 < 0) {
460       NXReportError( "ERROR: failed to store class name");
461       return NX_ERROR;
462    }
463    if (H5Awrite(attr1, aid1, (char*)nxclass) < 0) {
464      NXReportError( "ERROR: failed to store class name");
465      return NX_ERROR;
466    }
467    /* close group */
468    iRet=H5Sclose(aid2);
469    iRet=H5Tclose(aid1);
470    iRet=H5Aclose(attr1);
471    iRet=H5Gclose(iVID);
472    return NX_OK;
473  }
474 
475  /*------------------------------------------------------------------------*/
476
477  herr_t attr_check (hid_t loc_id, const char *member_name, const H5A_info_t *unused, void *opdata)
478  {
479    char attr_name[8+1];   /* need to leave space for \0 as well */
480   
481    strcpy(attr_name,"NX_class");
482    return strstr(member_name, attr_name) ? 1 : 0;
483  }
484  /*------------------------------------------------------------------------*/
485  NXstatus  NX5opengroup(NXhandle fid, CONSTCHAR *name, CONSTCHAR *nxclass)
486  {
487
488    pNexusFile5 pFile;
489    hid_t attr1, atype;
490    herr_t iRet;
491    char pBuffer[1024];
492    char data[128];
493         
494    pFile = NXI5assert (fid);
495    if (pFile->iCurrentG == 0) {
496      strcpy(pBuffer,name);
497    } else {
498       sprintf(pBuffer,"%s/%s",pFile->name_tmp,name);
499    }
500    iRet = H5Gopen(pFile->iFID, (const char *)pBuffer, H5P_DEFAULT);
501    if (iRet < 0) {
502      sprintf (pBuffer, "ERROR: group %s does not exist", pFile->name_tmp);
503      NXReportError( pBuffer);
504      return NX_ERROR; 
505    }
506    pFile->iCurrentG = iRet;
507    strcpy(pFile->name_tmp,pBuffer);
508    strcpy(pFile->name_ref,pBuffer);
509
510    if ((nxclass != NULL) && (strcmp(nxclass, NX_UNKNOWN_GROUP) != 0)) {
511        /* check group attribute */
512        iRet=H5Aiterate(pFile->iCurrentG,H5_INDEX_CRT_ORDER,H5_ITER_INC,0,attr_check,NULL);
513        if (iRet < 0) {
514          NXReportError( "ERROR: iterating through attribute list");
515          return NX_ERROR; 
516        } else if (iRet == 1) {
517         /* group attribute was found */
518        } else {
519         /* no group attribute available */
520         NXReportError( "ERROR: no group attribute available");
521         return NX_ERROR;
522        }
523        /* check contents of group attribute */
524        attr1 = H5Aopen_by_name(pFile->iCurrentG, ".", "NX_class", H5P_DEFAULT, H5P_DEFAULT);
525        if (attr1 < 0) {
526              NXReportError( "ERROR: opening NX_class group attribute");
527              return NX_ERROR; 
528        }
529        atype=H5Tcopy(H5T_C_S1);
530        H5Tset_size(atype,sizeof(data)); 
531        iRet = readStringAttributeN(attr1, data, sizeof(data));
532        //iRet = H5Aread(attr1, atype, data);
533        if (strcmp(data, nxclass) == 0) {
534              /* test OK */
535        } else {
536              snprintf(pBuffer, sizeof(pBuffer), "ERROR: group class is not identical: \"%s\" != \"%s\"", data, nxclass);
537              NXReportError(pBuffer);
538              iRet = H5Tclose(atype);
539              iRet = H5Aclose(attr1); 
540              return NX_ERROR; 
541        }         
542        iRet = H5Tclose(atype);
543        iRet = H5Aclose(attr1); 
544    }
545
546    /* maintain stack */
547    pFile->iStackPtr++;
548    pFile->iStack5[pFile->iStackPtr].iVref=pFile->iCurrentG;
549    strcpy(pFile->iStack5[pFile->iStackPtr].irefn,name);
550    pFile->iAtt5.iCurrentIDX=0;
551    pFile->iCurrentD = 0;
552    if(pFile->iCurrentLGG != NULL){
553      free(pFile->iCurrentLGG);
554    }
555    pFile->iCurrentLGG = strdup(name);
556    NXI5KillDir (pFile);
557    return NX_OK;
558  }
559
560  /* ------------------------------------------------------------------- */
561
562  NXstatus  NX5closegroup (NXhandle fid)
563  {
564    pNexusFile5 pFile;
565    int i,ii;
566    char *uname = NULL;
567    char *u1name = NULL;
568   
569    pFile = NXI5assert (fid);
570    /* first catch the trivial case: we are at root and cannot get
571       deeper into a negative directory hierarchy (anti-directory)
572     */
573    if (pFile->iCurrentG == 0) {
574      NXI5KillDir (pFile);
575      return NX_OK;
576    } else {                     
577      /* close the current group and decrement name_ref */
578      H5Gclose (pFile->iCurrentG);
579      i = 0;
580      i = strlen(pFile->iStack5[pFile->iStackPtr].irefn);
581      ii = strlen(pFile->name_ref);
582      if (pFile->iStackPtr>1) {
583         ii=ii-i-1;
584      } else {
585          ii=ii-i;
586      }       
587      if (ii>0) {
588        uname = strdup(pFile->name_ref); 
589        u1name = (char*) malloc((ii+1)*sizeof(char));
590        memset(u1name,0,ii);
591        for (i=0; i<ii; i++) {
592           *(u1name + i) = *(uname + i);
593        }
594        *(u1name + i) = '\0';
595        /*
596        strncpy(u1name, uname, ii);
597        */
598        strcpy(pFile->name_ref,u1name);
599        strcpy(pFile->name_tmp,u1name);
600        free(uname);
601        free(u1name);
602      } else {
603        strcpy(pFile->name_ref,"");
604        strcpy(pFile->name_tmp,"");
605      }
606      NXI5KillDir (pFile);
607      pFile->iStackPtr--;
608      if (pFile->iStackPtr>0) {
609         pFile->iCurrentG=pFile->iStack5[pFile->iStackPtr].iVref;
610      } else {
611         pFile->iCurrentG=0;
612      }
613    }
614    return NX_OK;
615  }
616/*-----------------------------------------------------------------------*/
617static hid_t nxToHDF5Type(int datatype)
618{
619  hid_t type;
620    if (datatype == NX_CHAR) {
621        type=H5T_C_S1;
622    } else if (datatype == NX_INT8) {
623        type=H5T_NATIVE_CHAR;
624    } else if (datatype == NX_UINT8) {
625        type=H5T_NATIVE_UCHAR;
626    } else if (datatype == NX_INT16) {
627        type=H5T_NATIVE_SHORT;
628    } else if (datatype == NX_UINT16) {
629        type=H5T_NATIVE_USHORT;
630    } else if (datatype == NX_INT32) {
631        type=H5T_NATIVE_INT;
632    } else if (datatype == NX_UINT32) {
633        type=H5T_NATIVE_UINT;
634    } else if (datatype == NX_INT64) {
635        type = H5T_NATIVE_INT64;
636    } else if (datatype == NX_UINT64) {
637        type = H5T_NATIVE_UINT64;
638    } else if (datatype == NX_FLOAT32) {
639        type=H5T_NATIVE_FLOAT;
640    } else if (datatype == NX_FLOAT64) {
641        type=H5T_NATIVE_DOUBLE;
642    } else {
643      NXReportError( "ERROR: nxToHDF5Type: unknown type");
644      type = -1;
645    }
646    return type;
647}
648 
649 /* --------------------------------------------------------------------- */
650
651  NXstatus  NX5compmakedata64 (NXhandle fid, CONSTCHAR *name, 
652                                          int datatype, 
653                                          int rank, int64_t dimensions[],
654                                          int compress_type, int64_t chunk_size[])
655  {
656      hid_t datatype1, dataspace, iNew;
657      herr_t iRet;
658      hid_t type, cparms = -1;
659      pNexusFile5 pFile;
660      char pBuffer[256];
661      int i, byte_zahl = 0;
662      hsize_t chunkdims[H5S_MAX_RANK];
663      hsize_t mydim[H5S_MAX_RANK], mydim1[H5S_MAX_RANK]; 
664      hsize_t size[H5S_MAX_RANK];
665      hsize_t maxdims[H5S_MAX_RANK];
666      int compress_level;
667      int unlimiteddim = 0;
668
669      pFile = NXI5assert (fid);
670      if (pFile->iCurrentG <= 0) {
671          sprintf(pBuffer, "ERROR: no group open for makedata on %s", name);
672          NXReportError( pBuffer);
673          return NX_ERROR;
674      }
675
676      if (rank <= 0) {
677          sprintf (pBuffer, "ERROR: invalid rank specified %s", name);
678          NXReportError( pBuffer);
679          return NX_ERROR;
680      }
681
682      type = nxToHDF5Type(datatype);
683
684      /*
685      Check dimensions for consistency. Dimension may be -1
686      thus denoting an unlimited dimension.
687      */
688      for (i = 0; i < rank; i++)
689      {
690          chunkdims[i] = chunk_size[i];
691          mydim[i] = dimensions[i];
692          maxdims[i] = dimensions[i];
693          size[i] = dimensions[i];
694          if (dimensions[i] <= 0) 
695          {
696                mydim[i] = 1;
697                maxdims[i] = H5S_UNLIMITED;
698                size[i] = 1;
699                unlimiteddim = 1;
700          } else {
701                mydim[i] = dimensions[i];
702                maxdims[i] = dimensions[i];
703                size[i] = dimensions[i];
704          }
705      }
706
707      if (datatype == NX_CHAR)
708      {
709          /*
710          *  This assumes string lenght is in the last dimensions and
711          *  the logic must be the same as used in NX5getslab and NX5getinfo
712          *
713          *  search for tests on H5T_STRING
714          */
715          byte_zahl=mydim[rank-1]; 
716          for(i = 0; i < rank; i++)
717          {
718              mydim1[i] = mydim[i];
719                if (dimensions[i] <= 0) 
720                  {
721                      mydim1[0] = 1;
722                       maxdims[0] = H5S_UNLIMITED;
723                 }
724
725          }
726          mydim1[rank-1] = 1;
727          if (mydim[rank-1] > 1)
728          {
729              mydim[rank-1] = maxdims[rank-1] = size[rank-1] = 1;
730          }
731          if (chunkdims[rank-1] > 1) 
732          { 
733              chunkdims[rank-1] = 1; 
734          }
735          dataspace=H5Screate_simple(rank,mydim1,maxdims);
736      } 
737      else 
738      {
739          if (unlimiteddim)
740          {
741              dataspace=H5Screate_simple(rank, mydim, maxdims);
742          } 
743          else 
744          {
745              /* dataset creation */
746              dataspace=H5Screate_simple(rank, mydim, NULL); 
747          }
748      } 
749      datatype1=H5Tcopy(type);
750      if (datatype == NX_CHAR)
751      {
752          H5Tset_size(datatype1, byte_zahl);
753          /*       H5Tset_strpad(H5T_STR_SPACEPAD); */
754      }
755      compress_level = 6;
756      if ( (compress_type / 100) ==  NX_COMP_LZW )
757      {
758          compress_level = compress_type % 100;
759          compress_type = NX_COMP_LZW;
760      }
761      if(compress_type == NX_COMP_LZW)
762      {
763          cparms = H5Pcreate(H5P_DATASET_CREATE);
764          iNew = H5Pset_chunk(cparms,rank,chunkdims);
765          if (iNew < 0) 
766          {
767              NXReportError( "ERROR: size of chunks could not be set");
768              return NX_ERROR;
769          }
770          H5Pset_deflate(cparms,compress_level); 
771          iRet = H5Dcreate(pFile->iCurrentG, (char*)name, datatype1, dataspace, H5P_DEFAULT, cparms, H5P_DEFAULT);   
772      } 
773      else if (compress_type == NX_COMP_NONE) 
774      {
775          if (unlimiteddim)
776          {
777              cparms = H5Pcreate(H5P_DATASET_CREATE);
778              iNew = H5Pset_chunk(cparms,rank,chunkdims);
779              if (iNew < 0) 
780              {
781                  NXReportError( "ERROR: size of chunks could not be set");
782                  return NX_ERROR;
783              }
784              iRet = H5Dcreate(pFile->iCurrentG, (char*)name, datatype1, dataspace, H5P_DEFAULT, cparms, H5P_DEFAULT);   
785          } 
786          else 
787          {
788              iRet = H5Dcreate(pFile->iCurrentG, (char*)name, datatype1, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
789          }               
790      } 
791      else if (compress_type == NX_CHUNK) 
792      {
793          cparms = H5Pcreate(H5P_DATASET_CREATE);
794          iNew = H5Pset_chunk(cparms,rank,chunkdims);
795          if (iNew < 0) 
796          {
797              NXReportError("ERROR: size of chunks could not be set");
798              return NX_ERROR;
799          }
800          iRet = H5Dcreate(pFile->iCurrentG, (char*)name, datatype1, dataspace, H5P_DEFAULT, cparms, H5P_DEFAULT);
801
802      } 
803      else 
804      {
805          NXReportError( 
806              "HDF5 doesn't support selected compression method! Dataset was saved without compression");
807          iRet = H5Dcreate (pFile->iCurrentG, (char*)name, datatype1, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); 
808      }
809      if (iRet < 0) 
810      {
811          NXReportError( "ERROR: creating chunked dataset failed");
812          return NX_ERROR;
813      } 
814      else 
815      {
816          pFile->iCurrentD = iRet;
817      }
818      if (unlimiteddim)
819      {     
820          iNew = H5Dset_extent  (pFile->iCurrentD, size);
821          if (iNew < 0) 
822          {
823              sprintf (pBuffer, "ERROR: cannot create dataset %s",
824                  name);
825              NXReportError( pBuffer);
826              return NX_ERROR;
827          }
828      }
829      if (cparms != -1)
830      {
831          iRet = H5Pclose(cparms);
832      }
833      iRet = H5Sclose(dataspace);
834      iRet = H5Tclose(datatype1);
835      iRet = H5Dclose(pFile->iCurrentD);
836      pFile->iCurrentD = 0;
837      if (iRet < 0)
838      {
839          NXReportError( "ERROR: HDF cannot close dataset");
840          return NX_ERROR;
841      }
842      return NX_OK;
843  }
844
845
846  /* --------------------------------------------------------------------- */
847
848  NXstatus  NX5makedata64 (NXhandle fid, CONSTCHAR *name, int datatype, 
849                                  int rank, int64_t dimensions[])
850  {
851  int64_t chunk_size[H5S_MAX_RANK];
852  int i;
853   
854  NXI5assert(fid);
855   
856  memset(chunk_size,0,H5S_MAX_RANK*sizeof(int64_t));
857  memcpy(chunk_size,dimensions,rank*sizeof(int64_t));
858  for(i = 0; i< rank; i++) {
859  if (dimensions[i] == NX_UNLIMITED || dimensions[i] <= 0){
860         chunk_size[i]= 1;
861  }   
862  }
863  return NX5compmakedata64 (fid, name, datatype, rank, dimensions, NX_COMP_NONE, chunk_size);
864  }
865
866 
867  /* --------------------------------------------------------------------- */
868
869  NXstatus  NX5compress (NXhandle fid, int compress_type)
870  {
871    printf(" NXcompress ERROR: NeXus API  based  on  HDF5  doesn't support\n");
872    printf("                   NXcompress  function!  Using  HDF5 library,\n");
873    printf("                   the NXcompmakedata function can be applied\n"); 
874    printf("                   for compression of data!\n");
875    return NX_ERROR;
876  } 
877
878  /* --------------------------------------------------------------------- */
879
880  NXstatus  NX5opendata (NXhandle fid, CONSTCHAR *name)
881  {
882    pNexusFile5 pFile;
883    char pBuffer[256];
884     
885    pFile = NXI5assert (fid);
886    /* clear pending attribute directories first */
887    NXI5KillAttDir (pFile);
888
889
890    /* find the ID number and open the dataset */
891    pFile->iCurrentD = H5Dopen(pFile->iCurrentG, name, H5P_DEFAULT);
892    if (pFile->iCurrentD < 0) {
893      sprintf (pBuffer, "ERROR: dataset \"%s\" not found at this level", name);
894      NXReportError( pBuffer);
895      return NX_ERROR;
896    }
897     /* find the ID number of datatype */
898    pFile->iCurrentT = H5Dget_type(pFile->iCurrentD);
899    if (pFile->iCurrentT < 0) {
900      NXReportError( "ERROR: error opening dataset");
901      pFile->iCurrentT=0;
902      return NX_ERROR;
903    }
904    /* find the ID number of dataspace */
905    pFile->iCurrentS = H5Dget_space(pFile->iCurrentD);
906    if (pFile->iCurrentS < 0) {
907      NXReportError( "ERROR:HDF error opening dataset");
908      pFile->iCurrentS=0;
909      return NX_ERROR;
910    }
911    if(pFile->iCurrentLD != NULL){
912      free(pFile->iCurrentLD);
913    }
914    pFile->iCurrentLD = strdup(name); 
915    return NX_OK;
916  }
917 
918  /* ----------------------------------------------------------------- */
919
920  NXstatus  NX5closedata (NXhandle fid)
921  {
922    pNexusFile5 pFile;
923    herr_t iRet;
924 
925    pFile = NXI5assert (fid);
926    iRet = H5Sclose(pFile->iCurrentS);
927    iRet = H5Tclose(pFile->iCurrentT);
928    iRet = H5Dclose(pFile->iCurrentD);
929    if (iRet < 0) {
930        NXReportError( "ERROR: cannot end access to dataset");
931        return NX_ERROR;
932     }
933    pFile->iCurrentD=0;
934    return NX_OK;
935  }
936 
937  /* ------------------------------------------------------------------- */
938 
939
940
941  NXstatus  NX5putdata (NXhandle fid, const void *data)
942  {
943    pNexusFile5 pFile;
944    herr_t iRet;
945    int64_t myStart[H5S_MAX_RANK];
946    int64_t mySize[H5S_MAX_RANK];
947    hsize_t thedims[H5S_MAX_RANK],maxdims[H5S_MAX_RANK];
948    int i, rank, unlimiteddim = 0;
949   
950    char pError[512] = "";
951     
952    pFile = NXI5assert (fid);
953    rank = H5Sget_simple_extent_ndims(pFile->iCurrentS);   
954    if (rank < 0)
955    {
956        NXReportError("ERROR: Cannot determine dataset rank");
957        return NX_ERROR;
958    }
959    iRet = H5Sget_simple_extent_dims(pFile->iCurrentS, thedims, maxdims);
960    if (iRet < 0)
961    {
962        NXReportError("ERROR: Cannot determine dataset dimensions");
963        return NX_ERROR;
964    }
965    for(i=0; i<rank; ++i)
966    {
967        myStart[i] = 0;
968        mySize[i] = thedims[i];
969        if (maxdims[i] == H5S_UNLIMITED)
970        {
971            unlimiteddim = 1;
972            myStart[i] = thedims[i] + 1;
973            mySize[i] = 1;
974        }
975    }
976    /* If we are using putdata on an unlimied dimesnion dataset, assume we want to append one single new slab */
977    if ( unlimiteddim )
978    {
979        return NX5putslab64(fid, data, myStart, mySize);
980    }
981    else
982    {
983        iRet = H5Dwrite (pFile->iCurrentD, pFile->iCurrentT, H5S_ALL, H5S_ALL, 
984                     H5P_DEFAULT, data);
985        if (iRet < 0) {
986          snprintf (pError,sizeof(pError)-1, "ERROR: failure to write data");
987          NXReportError( pError);
988          return NX_ERROR;
989        }
990    }
991    return NX_OK;
992  }
993/*------------------------------------------------------------------*/
994static int getAttVID(pNexusFile5 pFile){
995  int vid;
996     if(pFile->iCurrentG == 0 && pFile->iCurrentD == 0){
997       /* global attribute */
998       vid = H5Gopen(pFile->iFID,"/", H5P_DEFAULT);
999     } else if(pFile->iCurrentD != 0) {
1000       /* dataset attribute */
1001       vid = pFile->iCurrentD;
1002     } else {
1003       /* group attribute */;
1004       vid = pFile->iCurrentG;
1005     }
1006     return vid;
1007}
1008/*---------------------------------------------------------------*/
1009static void killAttVID(pNexusFile5 pFile, int vid){
1010  if(pFile->iCurrentG == 0 && pFile->iCurrentD == 0){
1011    H5Gclose(vid);
1012  }
1013}
1014  /* ------------------------------------------------------------------- */
1015
1016  NXstatus  NX5putattr (NXhandle fid, CONSTCHAR *name, const void *data, 
1017
1018                        int datalen, int iType)
1019  {
1020    pNexusFile5 pFile;
1021    hid_t  attr1, aid1, aid2;
1022    hid_t type;
1023    herr_t iRet;
1024    int vid; 
1025
1026    pFile = NXI5assert (fid);
1027
1028    type = nxToHDF5Type(iType);
1029   
1030    /* determine vid */
1031    vid = getAttVID(pFile);
1032    aid2=H5Screate(H5S_SCALAR);
1033    aid1=H5Tcopy(type);
1034    if (iType == NX_CHAR){
1035      H5Tset_size(aid1,datalen); 
1036    }         
1037    iRet = H5Aopen_by_name(vid, ".", name, H5P_DEFAULT, H5P_DEFAULT);
1038    if (iRet>0) {
1039      H5Aclose(iRet);
1040      iRet=H5Adelete(vid,name);
1041      if (iRet<0) {
1042        NXReportError( "ERROR: old attribute cannot be removed! ");
1043        killAttVID(pFile,vid);
1044        return NX_ERROR;
1045      }
1046    }   
1047    attr1 = H5Acreate(vid, name, aid1, aid2, H5P_DEFAULT, H5P_DEFAULT);
1048    if (attr1 < 0) {
1049      NXReportError( "ERROR: attribute cannot created! ");
1050      killAttVID(pFile,vid);
1051      return NX_ERROR;
1052    }   
1053    if (H5Awrite(attr1,aid1,data) < 0) {
1054      NXReportError( "ERROR: failed to store attribute ");
1055      killAttVID(pFile,vid);
1056      return NX_ERROR;
1057    }
1058    /* Close attribute dataspace */
1059    iRet=H5Tclose(aid1);
1060    iRet=H5Sclose(aid2); 
1061    /* Close attribute  */
1062    iRet=H5Aclose(attr1); 
1063    killAttVID(pFile,vid);
1064    return NX_OK;
1065  }
1066 
1067  /* ------------------------------------------------------------------- */
1068 
1069  NXstatus  NX5putslab64 (NXhandle fid, const void *data, const int64_t iStart[], const int64_t iSize[])
1070  {
1071    pNexusFile5 pFile;
1072    int iRet, rank, i;
1073    hsize_t myStart[H5S_MAX_RANK];
1074    hsize_t mySize[H5S_MAX_RANK];
1075    hsize_t size[H5S_MAX_RANK],thedims[H5S_MAX_RANK],maxdims[H5S_MAX_RANK];
1076    hid_t   filespace,dataspace; 
1077    int unlimiteddim = 0;
1078 
1079    pFile = NXI5assert (fid);
1080     /* check if there is an Dataset open */
1081    if (pFile->iCurrentD == 0) {
1082      NXReportError( "ERROR: no dataset open");
1083      return NX_ERROR;
1084    }
1085    rank = H5Sget_simple_extent_ndims(pFile->iCurrentS);   
1086    if (rank < 0) {
1087      NXReportError( "ERROR: cannot get rank");
1088      return NX_ERROR;
1089    }
1090    iRet = H5Sget_simple_extent_dims(pFile->iCurrentS, thedims, maxdims);
1091    if (iRet < 0) {
1092      NXReportError( "ERROR: cannot get dimensions");
1093      return NX_ERROR;
1094    }
1095
1096    for(i = 0; i < rank; i++)
1097    {
1098       myStart[i] = iStart[i];
1099       mySize[i]  = iSize[i];
1100       size[i]    =  iStart[i] + iSize[i];
1101       if (maxdims[i] == H5S_UNLIMITED) {
1102        unlimiteddim = 1;
1103       }
1104    }
1105    if (H5Tget_class(pFile->iCurrentT) == H5T_STRING)
1106    {
1107        mySize[rank - 1] = 1;
1108        myStart[rank - 1] = 0;
1109        size[rank - 1] = 1;
1110    }
1111    dataspace = H5Screate_simple(rank, mySize, NULL);
1112    if (unlimiteddim)
1113    {
1114       for(i = 0; i < rank; i++)
1115       {
1116        if (size[i] < thedims[i]) {
1117                size[i] = thedims[i];
1118        }
1119       } 
1120       iRet = H5Dset_extent(pFile->iCurrentD, size);
1121       if (iRet < 0) 
1122       {
1123           NXReportError( "ERROR: extend slab failed");
1124           return NX_ERROR;
1125       }
1126
1127       filespace = H5Dget_space(pFile->iCurrentD);
1128       
1129       /* define slab */
1130       iRet = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, myStart,
1131                               NULL, mySize, NULL);
1132       /* deal with HDF errors */
1133       if (iRet < 0) 
1134       {
1135       NXReportError( "ERROR: selecting slab failed");
1136       return NX_ERROR;
1137       }
1138       /* write slab */ 
1139       iRet = H5Dwrite(pFile->iCurrentD, pFile->iCurrentT, dataspace, 
1140                    filespace, H5P_DEFAULT, data);
1141       if (iRet < 0)
1142       {
1143           NXReportError( "ERROR: writing slab failed");
1144       }
1145       /* update with new size */
1146       iRet = H5Sclose(pFile->iCurrentS);             
1147       pFile->iCurrentS = filespace;
1148   } else {
1149       /* define slab */
1150       iRet = H5Sselect_hyperslab(pFile->iCurrentS, H5S_SELECT_SET, myStart,
1151                               NULL, mySize, NULL);
1152       /* deal with HDF errors */
1153       if (iRet < 0) 
1154       {
1155       NXReportError( "ERROR: selecting slab failed");
1156       return NX_ERROR;
1157       }
1158       /* write slab */ 
1159       iRet = H5Dwrite(pFile->iCurrentD, pFile->iCurrentT, dataspace, 
1160                    pFile->iCurrentS, H5P_DEFAULT, data);
1161       if (iRet < 0)
1162       {
1163           NXReportError( "ERROR: writing slab failed");
1164       }
1165   }
1166   /* deal with HDF errors */
1167   iRet = H5Sclose(dataspace);
1168   if (iRet < 0) 
1169   {
1170      NXReportError( "ERROR: closing slab failed");
1171      return NX_ERROR;
1172   }
1173    return NX_OK;
1174  }
1175 
1176  /* ------------------------------------------------------------------- */
1177
1178  NXstatus  NX5getdataID (NXhandle fid, NXlink* sRes)
1179  {
1180    pNexusFile5 pFile;
1181    int datalen, type = NX_CHAR;
1182 
1183    pFile = NXI5assert (fid);
1184
1185    /*
1186      we cannot return ID's when no datset is open
1187    */
1188    if(pFile->iCurrentD <= 0){
1189      return NX_ERROR;
1190    }
1191
1192    /*
1193      this means: if the item is already linked: use the target attribute else,
1194      the path to the current node
1195    */
1196    NXMDisableErrorReporting();
1197    datalen = 1024;
1198    memset(&sRes->targetPath,0,datalen*sizeof(char));
1199    if(NX5getattr(fid,"target",&sRes->targetPath,&datalen,&type) != NX_OK)
1200    {
1201      buildCurrentPath(pFile, sRes->targetPath, 1024);
1202    }
1203    NXMEnableErrorReporting();
1204    sRes->linkType = 1;
1205    return NX_OK;
1206  }
1207
1208 /* ------------------------------------------------------------------- */
1209 
1210  NXstatus  NX5printlink (NXhandle fid, NXlink* sLink)
1211  {
1212    NXI5assert(fid);
1213    printf("HDF5 link: targetPath = \"%s\", linkType = \"%d\"\n", sLink->targetPath, sLink->linkType);
1214    return NX_OK;
1215  }
1216/*--------------------------------------------------------------------*/
1217static NXstatus NX5settargetattribute(pNexusFile5 pFile, NXlink *sLink)
1218{
1219  hid_t  dataID, aid2, aid1, attID;
1220  herr_t status;
1221  char name[] = "target";
1222
1223    /*
1224      set the target attribute
1225    */
1226    if(sLink->linkType > 0)
1227    {
1228      dataID = H5Dopen(pFile->iFID,sLink->targetPath, H5P_DEFAULT);
1229    } else {
1230      dataID = H5Gopen(pFile->iFID,sLink->targetPath, H5P_DEFAULT);
1231    }
1232    if(dataID < 0)
1233    {
1234      NXReportError("Internal error, path to link does not exist");
1235      return NX_ERROR;
1236    }
1237    status = H5Aopen_by_name(dataID,".", name, H5P_DEFAULT, H5P_DEFAULT);
1238    if(status > 0)
1239    {
1240      H5Aclose(status);
1241      status = H5Adelete(dataID,name);
1242      if(status < 0)
1243      {
1244        return NX_OK;
1245      }
1246    }
1247    aid2 = H5Screate(H5S_SCALAR);
1248    aid1 = H5Tcopy(H5T_C_S1);
1249    H5Tset_size(aid1,strlen(sLink->targetPath));
1250    attID = H5Acreate(dataID,name,aid1,aid2,H5P_DEFAULT, H5P_DEFAULT);
1251    if(attID < 0)
1252    {
1253        return NX_OK;
1254    }
1255    status = H5Awrite(attID,aid1,sLink->targetPath);
1256    H5Tclose(aid1);
1257    H5Sclose(aid2); 
1258    H5Aclose(attID); 
1259    if(sLink->linkType > 0){
1260      H5Dclose(dataID);
1261    } else {
1262      H5Gclose(dataID);
1263    }
1264    return NX_OK;
1265}
1266/*---------------------------------------------------------------------*/
1267NXstatus NX5makenamedlink(NXhandle fid, CONSTCHAR *name, NXlink *sLink)
1268{
1269    pNexusFile5 pFile;
1270    char        linkTarget[1024];
1271
1272    pFile = NXI5assert (fid);
1273    if (pFile->iCurrentG == 0) { /* root level, can not link here */
1274      return NX_ERROR;
1275    }
1276   
1277    /*
1278      build pathname to link from our current group and the name
1279      of the thing to link
1280    */
1281    if(strlen(pFile->name_ref) + strlen(name) + 2 < 1024)
1282    {
1283      strcpy(linkTarget,"/");
1284      strcat(linkTarget,pFile->name_ref);
1285      strcat(linkTarget,"/");
1286      strcat(linkTarget,name);
1287    }
1288    else 
1289    {
1290      NXReportError("ERROR: path string to long");
1291      return NX_ERROR;
1292    } 
1293
1294    //targetid = H5Oopen(pFile->iFID, sLink->targetPath, H5P_DEFAULT);
1295    H5Lcreate_hard(pFile->iFID, sLink->targetPath, H5L_SAME_LOC, linkTarget, H5P_DEFAULT, H5P_DEFAULT);
1296    //H5Oclose(targetid);
1297
1298    return NX5settargetattribute(pFile,sLink);
1299}
1300 /* ------------------------------------------------------------------- */
1301 
1302  NXstatus  NX5makelink (NXhandle fid, NXlink* sLink)
1303  {
1304    pNexusFile5 pFile;
1305    char        linkTarget[1024];
1306    char        *itemName = NULL;
1307
1308    pFile = NXI5assert (fid);
1309    if (pFile->iCurrentG == 0) { /* root level, can not link here */
1310      return NX_ERROR;
1311    }
1312   
1313    /*
1314      locate name of the element to link
1315    */
1316    itemName = strrchr(sLink->targetPath,'/');
1317    if(itemName == NULL){
1318      NXReportError("ERROR: bad link structure");
1319      return NX_ERROR;
1320    }
1321    itemName++;
1322
1323    /*
1324      build pathname to link from our current group and the name
1325      of the thing to link
1326    */
1327    if(strlen(pFile->name_ref) + strlen(itemName) + 2 < 1024)
1328    {
1329      strcpy(linkTarget,"/");
1330      strcat(linkTarget,pFile->name_ref);
1331      strcat(linkTarget,"/");
1332      strcat(linkTarget,itemName);
1333    }
1334    else 
1335    {
1336      NXReportError("ERROR: path string to long");
1337      return NX_ERROR;
1338    } 
1339
1340    //targetid = H5Oopen(pFile->iFID, sLink->targetPath, H5P_DEFAULT);
1341    H5Lcreate_hard(pFile->iFID, sLink->targetPath, H5L_SAME_LOC, linkTarget, H5P_DEFAULT, H5P_DEFAULT);
1342    //H5Oclose(targetid);
1343
1344    return NX5settargetattribute(pFile,sLink);
1345   }
1346 
1347  /*----------------------------------------------------------------------*/
1348
1349  NXstatus  NX5flush(NXhandle *pHandle)
1350  {
1351    pNexusFile5 pFile = NULL;
1352    herr_t      iRet;
1353   
1354    pFile = NXI5assert (*pHandle);   
1355    if (pFile->iCurrentD != 0)
1356    {   
1357    iRet=H5Fflush(pFile->iCurrentD,H5F_SCOPE_LOCAL);
1358    }
1359    else if (pFile->iCurrentG != 0)
1360    {   
1361    iRet=H5Fflush(pFile->iCurrentG,H5F_SCOPE_LOCAL);
1362    }
1363    else
1364    { 
1365    iRet=H5Fflush(pFile->iFID,H5F_SCOPE_LOCAL);
1366    }
1367    if (iRet < 0){
1368      NXReportError( "ERROR: The object cannot be flushed");
1369      return NX_ERROR; 
1370    }
1371    return NX_OK;
1372  }   
1373 
1374  /*-------------------------------------------------------------------------*/
1375 
1376  /* Operator function. */
1377           
1378  herr_t nxgroup_info(hid_t loc_id, const char *name, const H5L_info_t *statbuf, void *op_data)
1379  {
1380    pinfo self = (pinfo) op_data;
1381    H5O_info_t object_info;
1382    H5Oget_info_by_name(loc_id, name, &object_info, H5P_DEFAULT);
1383    switch ((object_info).type) {
1384      case H5O_TYPE_GROUP:
1385         self->iname = strdup(name);
1386         self->type = H5O_TYPE_GROUP;
1387         break;
1388      case H5O_TYPE_DATASET:
1389         self->iname = strdup(name);
1390         self->type = H5O_TYPE_DATASET;
1391         break;
1392      default:
1393         // TODO defaults to group. not what we would want?
1394         self->type=0;
1395         break;     
1396    }
1397    return 1; 
1398  }
1399  /* --------------------------------------------------------------------- */
1400
1401  /* Operator function. */
1402
1403  herr_t group_info1(hid_t loc_id, const char *name, const H5L_info_t *statbuf, void *opdata)
1404  {
1405    int iNX = *((int*)opdata);
1406    H5O_info_t object_info;
1407    H5Oget_info_by_name(loc_id, name, &object_info, H5P_DEFAULT);
1408    switch ((object_info).type) {
1409      case H5O_TYPE_GROUP:
1410        iNX++;
1411        *((int*)opdata)=iNX;
1412        break;
1413      case H5O_TYPE_DATASET:
1414        iNX++;
1415        *((int*)opdata)=iNX;
1416        break;
1417      default:
1418        break;
1419    }
1420    return 0; 
1421  }
1422 
1423  /*-------------------------------------------------------------------------*/
1424
1425  NXstatus  NX5getgroupinfo_recurse(NXhandle fid, int *iN, NXname pName, NXname pClass)
1426  {
1427    pNexusFile5 pFile;
1428    hid_t       atype, attr_id, grp;
1429    char        data[64];
1430       
1431    pFile = NXI5assert (fid);
1432    /* check if there is a group open */
1433    if (pFile->iCurrentG == 0) {
1434       strcpy (pName, "root");
1435       strcpy (pClass, "NXroot");
1436       pFile->iNX=0;
1437       grp = H5Gopen(pFile->iFID,"/",H5P_DEFAULT);
1438       H5Literate(grp, H5_INDEX_NAME, H5_ITER_INC, 0, group_info1, &pFile->iNX);
1439       H5Gclose(grp);
1440       *iN=pFile->iNX;
1441    }
1442    else {
1443      strcpy (pName,pFile->name_ref);
1444      attr_id = H5Aopen_by_name(pFile->iCurrentG,".", "NX_class", H5P_DEFAULT, H5P_DEFAULT);
1445      if (attr_id<0) {
1446         strcpy(pClass, NX_UNKNOWN_GROUP);
1447      } else {
1448        atype=H5Tcopy(H5T_C_S1);
1449        H5Tset_size(atype,sizeof(data)); 
1450        readStringAttributeN(attr_id, data, sizeof(data));
1451        //H5Aread(attr_id, atype, data);
1452        strcpy(pClass,data);
1453        pFile->iNX=0;
1454        grp = H5Gopen(pFile->iFID,pFile->name_ref,H5P_DEFAULT);
1455        H5Literate(grp, H5_INDEX_NAME, H5_ITER_INC, 0, group_info1, &pFile->iNX);
1456        H5Gclose(grp);
1457        *iN=pFile->iNX;
1458        H5Aclose(attr_id);
1459      }
1460    }
1461    return NX_OK;
1462  }
1463/*---------------------------------------------------------------------------*/
1464static int countObjectsInGroup(hid_t loc_id)
1465{
1466  int count = 0;
1467  H5G_info_t numobj;
1468 
1469  herr_t status;
1470
1471  status = H5Gget_info(loc_id, &numobj);
1472  if(status < 0) {
1473    NXReportError("Internal error, failed to retrive no of objects");
1474    return 0;
1475  }
1476
1477  count = numobj.nlinks;
1478  return count;
1479}
1480/*----------------------------------------------------------------------------*/
1481  NXstatus  NX5getgroupinfo (NXhandle fid, int *iN, NXname pName, NXname pClass)
1482  {
1483    pNexusFile5 pFile;
1484    hid_t atype, attr_id, gid;
1485    char data[64];
1486
1487    pFile = NXI5assert (fid);
1488    /* check if there is a group open */
1489    if (pFile->iCurrentG == 0) {
1490       strcpy (pName, "root");
1491       strcpy (pClass, "NXroot");
1492       gid = H5Gopen(pFile->iFID,"/", H5P_DEFAULT);
1493       *iN = countObjectsInGroup(gid);
1494       H5Gclose(gid);
1495    }
1496    else {
1497      strcpy (pName,pFile->name_ref);
1498      attr_id = H5Aopen_by_name(pFile->iCurrentG,".", "NX_class", H5P_DEFAULT, H5P_DEFAULT);
1499      if (attr_id<0) {
1500         strcpy(pClass, NX_UNKNOWN_GROUP);
1501      } else {
1502        atype=H5Tcopy(H5T_C_S1);
1503        H5Tset_size(atype,sizeof(data)); 
1504        readStringAttributeN(attr_id, data, sizeof(data));
1505        //H5Aread(attr_id, atype, data);
1506        strcpy(pClass,data);
1507        pFile->iNX=0;
1508        *iN = countObjectsInGroup(pFile->iCurrentG);
1509        H5Aclose(attr_id);
1510      }
1511    }
1512    return NX_OK;
1513  }
1514
1515
1516/*-------------------------------------------------------------------------
1517 * Function: hdf5ToNXType
1518 *
1519 * Purpose:     Convert a HDF5 class to a NeXus type;  it handles the following HDF5 classes
1520 *  H5T_STRING
1521 *  H5T_INTEGER
1522 *  H5T_FLOAT
1523 *
1524 * Return: the NeXus type
1525 *
1526 *-------------------------------------------------------------------------
1527 */
1528  static int hdf5ToNXType(H5T_class_t tclass, hid_t atype)
1529  {
1530      int        iPtype = -1;
1531      size_t     size;
1532      H5T_sign_t sign;
1533
1534      if (tclass==H5T_STRING)
1535      {
1536          iPtype=NX_CHAR;
1537      }
1538      else if (tclass==H5T_INTEGER)
1539      {
1540          size=H5Tget_size(atype);
1541          sign=H5Tget_sign(atype);
1542          if (size==1)
1543          {
1544              if (sign==H5T_SGN_2)
1545              {
1546                  iPtype=NX_INT8;
1547              } else {
1548                  iPtype=NX_UINT8;
1549              }
1550          } 
1551          else if (size==2) 
1552          {
1553              if (sign==H5T_SGN_2)
1554              {
1555                  iPtype=NX_INT16;
1556              } else {
1557                  iPtype=NX_UINT16;
1558              }
1559          }
1560          else if (size==4) 
1561          {
1562              if (sign==H5T_SGN_2)
1563              {
1564                  iPtype=NX_INT32;
1565              } else {
1566                  iPtype=NX_UINT32;
1567              }
1568          } 
1569          else if(size == 8)
1570          {
1571              if (sign==H5T_SGN_2)
1572              {
1573                  iPtype=NX_INT64;
1574              } else {
1575                  iPtype=NX_UINT64;
1576              }
1577          }
1578      } 
1579      else if (tclass==H5T_FLOAT)     
1580      {
1581          size=H5Tget_size(atype);
1582          if (size==4)
1583          {
1584              iPtype=NX_FLOAT32;
1585          } 
1586          else if (size==8) 
1587          {
1588              iPtype=NX_FLOAT64;
1589          }
1590      }
1591      if (iPtype == -1)
1592      {
1593          NXReportError( "ERROR: hdf5ToNXtype: invalid type");
1594      }
1595
1596      return iPtype;
1597  }
1598/*--------------------------------------------------------------------------*/
1599  static hid_t h5MemType(hid_t atype)
1600  {
1601      hid_t       memtype_id = -1;
1602      size_t      size;
1603      H5T_sign_t  sign;
1604      H5T_class_t tclass;
1605
1606      tclass = H5Tget_class(atype);
1607
1608      if (tclass==H5T_INTEGER)
1609      {
1610          size=H5Tget_size(atype);
1611          sign=H5Tget_sign(atype);
1612          if (size==1)
1613          {
1614              if (sign==H5T_SGN_2)
1615              {
1616                  memtype_id = H5T_NATIVE_INT8;
1617              } else {
1618                  memtype_id = H5T_NATIVE_UINT8;
1619              }
1620          } 
1621          else if (size==2) 
1622          {
1623              if (sign==H5T_SGN_2)
1624              {
1625                  memtype_id = H5T_NATIVE_INT16;
1626              } else {
1627                  memtype_id = H5T_NATIVE_UINT16; 
1628              }
1629          }
1630          else if (size==4) 
1631          {
1632              if (sign==H5T_SGN_2)
1633              {
1634                  memtype_id = H5T_NATIVE_INT32;
1635              } else {
1636                  memtype_id = H5T_NATIVE_UINT32; 
1637              }
1638          }
1639          else if (size==8) 
1640          {
1641              if (sign==H5T_SGN_2)
1642              {
1643                  memtype_id = H5T_NATIVE_INT64;
1644              } else {
1645                  memtype_id = H5T_NATIVE_UINT64;
1646              }
1647          }
1648      } else if (tclass==H5T_FLOAT)     
1649      {
1650          size=H5Tget_size(atype);
1651          if (size==4)
1652          {
1653              memtype_id = H5T_NATIVE_FLOAT; 
1654          } else if (size==8) {
1655              memtype_id = H5T_NATIVE_DOUBLE;
1656          }
1657      }           
1658      if (memtype_id == -1)
1659      {
1660          NXReportError( "ERROR: h5MemType: invalid type");
1661      }
1662      return memtype_id;
1663  }
1664  /*-------------------------------------------------------------------------*/
1665
1666  NXstatus  NX5getnextentry(NXhandle fid, NXname name, NXname nxclass, int *datatype)
1667  {
1668    pNexusFile5 pFile;
1669    hid_t       grp, attr1,type,atype;
1670    herr_t      iRet;
1671    int         iPtype, i;
1672    hsize_t     idx;
1673    H5T_class_t tclass;
1674    char        data[128];
1675    char        ph_name[1024];
1676    info_type   op_data;
1677    herr_t      iRet_iNX=-1;
1678    char        pBuffer[256];
1679     
1680    pFile = NXI5assert (fid);
1681    op_data.iname = NULL;
1682
1683    /*
1684      iterate to next entry in group list
1685    */
1686    idx=pFile->iStack5[pFile->iStackPtr].iCurrentIDX;
1687    if (strlen(pFile->name_ref) == 0) {
1688       /* root group */
1689       strcpy(pFile->name_ref,"/");
1690    } 
1691    grp = H5Gopen(pFile->iFID, pFile->name_ref, H5P_DEFAULT);
1692    // index can be wrong here
1693    iRet=H5Literate(grp, H5_INDEX_NAME, H5_ITER_INC, &idx, nxgroup_info, &op_data);
1694    H5Gclose(grp);
1695    strcpy(nxclass, NX_UNKNOWN_GROUP);
1696   
1697    /*
1698      figure out the number of items in the current group. We need this in order to
1699      find out if we are at the end of the search.
1700    */
1701    if (pFile->iCurrentG == 0) {
1702        // if pFile->iCurrentG == 0 would not pFile->name_ref be "/" already, so we could skip that if statement ?
1703       pFile->iNX=0;
1704       grp = H5Gopen(pFile->iFID, "/", H5P_DEFAULT);
1705       iRet_iNX=H5Literate(grp, H5_INDEX_NAME, H5_ITER_INC, 0, group_info1, &pFile->iNX);
1706       H5Gclose(grp);
1707    } else {
1708      pFile->iNX=0;
1709      grp = H5Gopen(pFile->iFID, pFile->name_ref, H5P_DEFAULT);
1710        // index can be wrong here
1711      iRet_iNX=H5Literate(grp, H5_INDEX_NAME, H5_ITER_INC, 0, group_info1, &pFile->iNX);
1712      H5Gclose(grp);
1713    }
1714    if (idx == pFile->iNX) {
1715        // why 2?
1716       iRet_iNX = 2;   
1717    } 
1718         
1719    if (iRet > 0)
1720      {
1721        pFile->iStack5[pFile->iStackPtr].iCurrentIDX++;
1722        if (op_data.iname != NULL) {
1723          strcpy(name,op_data.iname);
1724           free(op_data.iname);
1725        } else {
1726          pFile->iStack5[pFile->iStackPtr].iCurrentIDX = 0;   
1727          return NX_EOD;
1728        }         
1729        if (op_data.type == H5O_TYPE_GROUP) {
1730          /*
1731            open group and find class name attribute
1732          */
1733           strcpy(ph_name,"");
1734           for(i = 1; i < (pFile->iStackPtr + 1); i++)
1735           {
1736              strcat(ph_name,pFile->iStack5[i].irefn);
1737              strcat(ph_name,"/");
1738           }
1739           strcat(ph_name,name);
1740           grp = H5Gopen(pFile->iFID,ph_name, H5P_DEFAULT);
1741           if (grp < 0) {
1742              sprintf(pBuffer, "ERROR: group %s does not exist", ph_name);
1743              NXReportError(pBuffer);
1744              return NX_ERROR; 
1745           }
1746           attr1 = H5Aopen_by_name(grp,".",  "NX_class", H5P_DEFAULT, H5P_DEFAULT);
1747           if (attr1 < 0) {
1748              strcpy(nxclass, NX_UNKNOWN_GROUP);
1749           } else {
1750               type=H5T_C_S1;
1751               atype=H5Tcopy(type);
1752               H5Tset_size(atype,sizeof(data)); 
1753               iRet = readStringAttributeN(attr1, data, sizeof(data));
1754               //iRet = H5Aread(attr1, atype, data);
1755               strcpy(nxclass,data);
1756               H5Tclose(atype);
1757               H5Aclose(attr1);
1758           }
1759           H5Gclose(grp);
1760        } else if (op_data.type==H5O_TYPE_DATASET) {
1761          /*
1762            open dataset and find type
1763          */
1764           grp=H5Dopen(pFile->iCurrentG,name, H5P_DEFAULT);
1765           type=H5Dget_type(grp);
1766           atype=H5Tcopy(type);
1767           tclass = H5Tget_class(atype);
1768           iPtype = hdf5ToNXType(tclass, atype);
1769           *datatype=iPtype;
1770           strcpy(nxclass, "SDS");
1771           H5Tclose(atype);
1772            H5Tclose(type);
1773            H5Dclose(grp);
1774         }
1775          return NX_OK;
1776       } else { 
1777         /*
1778           we are at the end of the search: clear the data structure and reset
1779           iCurrentIDX to 0
1780         */
1781          if (iRet_iNX == 2) {
1782            if (op_data.iname != NULL) {
1783               free(op_data.iname);
1784            }
1785            pFile->iStack5[pFile->iStackPtr].iCurrentIDX = 0;   
1786            return NX_EOD;
1787         } 
1788         if (op_data.iname != NULL) {
1789            free(op_data.iname);
1790         } 
1791         NXReportError("ERROR: iterating through group not successful");
1792         return NX_ERROR;             
1793       }
1794   }
1795
1796   /*-------------------------------------------------------------------------*/
1797
1798   NXstatus  NX5getdata (NXhandle fid, void *data)
1799   {
1800     pNexusFile5 pFile;
1801     int i, iStart[H5S_MAX_RANK], status;
1802     hid_t memtype_id; 
1803     H5T_class_t tclass;
1804     hsize_t ndims,dims[H5S_MAX_RANK],len;   
1805     char** vstrdata = NULL;
1806
1807     pFile = NXI5assert (fid);
1808     /* check if there is an Dataset open */
1809     if (pFile->iCurrentD == 0) 
1810     {
1811        NXReportError( "ERROR: no dataset open");
1812        return NX_ERROR;
1813     }
1814     ndims = H5Sget_simple_extent_dims(pFile->iCurrentS, dims, NULL); 
1815     if (ndims <= 0)
1816     {
1817        NXReportError( "ERROR: unable to read dims");
1818        return NX_ERROR;
1819     }
1820     memset (iStart, 0, H5S_MAX_RANK * sizeof(int));
1821     /* map datatypes of other plateforms */
1822     tclass = H5Tget_class(pFile->iCurrentT);
1823     if ( H5Tis_variable_str(pFile->iCurrentT) )
1824     {
1825        vstrdata = (char **) malloc (dims[0] * sizeof (char *));
1826        memtype_id = H5Tcopy(H5T_C_S1);
1827        H5Tset_size(memtype_id, H5T_VARIABLE);
1828        status = H5Dread (pFile->iCurrentD, memtype_id, 
1829                       H5S_ALL, H5S_ALL,H5P_DEFAULT, vstrdata);
1830        ((char*)data)[0] = '\0';
1831        if (status >= 0)
1832        {
1833            for(i=0; i<dims[0]; ++i)
1834            {
1835                if (vstrdata[i] != NULL)
1836                {
1837                    strcat((char*)data, vstrdata[i]);
1838                }
1839            }
1840        }
1841        H5Dvlen_reclaim(memtype_id, pFile->iCurrentS, H5P_DEFAULT, vstrdata);
1842        free(vstrdata);
1843        H5Tclose(memtype_id);
1844     }
1845     else if (tclass==H5T_STRING)
1846     {
1847        len = H5Tget_size(pFile->iCurrentT);
1848        memtype_id = H5Tcopy(H5T_C_S1);
1849        H5Tset_size(memtype_id, len);
1850        status = H5Dread (pFile->iCurrentD, memtype_id, 
1851                       H5S_ALL, H5S_ALL,H5P_DEFAULT, data);
1852        H5Tclose(memtype_id);
1853     }
1854     else 
1855     {
1856       memtype_id = h5MemType(pFile->iCurrentT);
1857       status = H5Dread (pFile->iCurrentD, memtype_id, 
1858                       H5S_ALL, H5S_ALL,H5P_DEFAULT, data);
1859     }
1860     if(status < 0)
1861     {
1862        NXReportError( "ERROR: failed to transfer dataset");
1863        return NX_ERROR;
1864
1865     }
1866     return NX_OK;
1867   }
1868
1869   /*-------------------------------------------------------------------------*/
1870
1871   NXstatus  NX5getinfo64 (NXhandle fid, int *rank, int64_t dimension[], int *iType)
1872   {
1873     pNexusFile5 pFile;
1874     int i, iRank, mType;
1875     hsize_t myDim[H5S_MAX_RANK], vlen_bytes = 0, total_dims_size = 1; 
1876     H5T_class_t tclass;
1877
1878     pFile = NXI5assert (fid);
1879     /* check if there is an Dataset open */
1880     if (pFile->iCurrentD == 0) {
1881       NXReportError( "ERROR: no dataset open");
1882       return NX_ERROR;
1883     }
1884
1885     /* read information */
1886     tclass = H5Tget_class(pFile->iCurrentT);
1887     mType = hdf5ToNXType(tclass,pFile->iCurrentT);
1888     iRank = H5Sget_simple_extent_ndims(pFile->iCurrentS);
1889     H5Sget_simple_extent_dims(pFile->iCurrentS, myDim, NULL);   
1890     for(i=0; i<iRank; ++i)
1891     {
1892        total_dims_size *= myDim[i];
1893     }
1894     /* conversion to proper ints for the platform */ 
1895     *iType = (int)mType;
1896     if (tclass==H5T_STRING && myDim[iRank-1] == 1) {
1897        if ( H5Tis_variable_str(pFile->iCurrentT) )
1898        {
1899            /* this would not work for ragged arrays */
1900            H5Dvlen_get_buf_size(pFile->iCurrentD, pFile->iCurrentT, pFile->iCurrentS, &vlen_bytes);
1901            myDim[iRank-1] = vlen_bytes / total_dims_size;
1902        }
1903        else
1904        {
1905            myDim[iRank-1] = H5Tget_size(pFile->iCurrentT);
1906        }
1907     } 
1908     *rank = (int)iRank;
1909     for (i = 0; i < iRank; i++)
1910     {
1911          dimension[i] = (int)myDim[i];
1912     }
1913     return NX_OK;
1914   }
1915
1916   /*-------------------------------------------------------------------------*/
1917
1918   NXstatus  NX5getslab64 (NXhandle fid, void *data, const int64_t iStart[], const int64_t iSize[])
1919   {
1920     pNexusFile5 pFile;
1921     hsize_t myStart[H5S_MAX_RANK];
1922     hsize_t mySize[H5S_MAX_RANK];
1923     hsize_t mStart[H5S_MAX_RANK];
1924     hid_t   memspace, iRet;
1925     H5T_class_t tclass;
1926     hid_t   memtype_id;
1927     char *tmp_data = NULL;
1928     char *data1;
1929     int i, dims, iRank, mtype = 0;
1930
1931     pFile = NXI5assert (fid);
1932     /* check if there is an Dataset open */
1933     if (pFile->iCurrentD == 0) 
1934       {
1935         NXReportError( "ERROR: no dataset open");
1936         return NX_ERROR;
1937       }
1938     iRank = H5Sget_simple_extent_ndims(pFile->iCurrentS); 
1939     for (i = 0; i < iRank; i++)
1940        {
1941          myStart[i] = (hssize_t)iStart[i];
1942          mySize[i]  = (hsize_t)iSize[i];
1943          mStart[i] = (hsize_t)0;
1944        }
1945      tclass = H5Tget_class(pFile->iCurrentT);
1946      if (tclass == H5T_STRING) {
1947/*
1948 * FAA 24/1/2007: I don't think this will work for multidimensional
1949 * string arrays.
1950 * MK 23/7/2007: You are right Freddie.
1951*/
1952         mtype = NX_CHAR;
1953         if (mySize[0] == 1) {
1954             mySize[0]  = H5Tget_size(pFile->iCurrentT);
1955         }   
1956         tmp_data = (char*) malloc(mySize[0]);
1957         memset(tmp_data,0,sizeof(mySize[0]));
1958         iRet = H5Sselect_hyperslab(pFile->iCurrentS, H5S_SELECT_SET, mStart,
1959                                NULL, mySize, NULL);
1960      } else {
1961         iRet = H5Sselect_hyperslab(pFile->iCurrentS, H5S_SELECT_SET, myStart,
1962                                NULL, mySize, NULL);
1963      }
1964      /* define slab */
1965      /* deal with HDF errors */
1966      if (iRet < 0) 
1967        {
1968          NXReportError( "ERROR: selecting slab failed");
1969          return NX_ERROR;
1970        }
1971
1972      memspace=H5Screate_simple(iRank, mySize, NULL);
1973      iRet = H5Sselect_hyperslab(memspace, H5S_SELECT_SET, mStart,
1974                                NULL, mySize, NULL);
1975      if (iRet < 0) 
1976        {
1977          NXReportError( "ERROR: selecting memspace failed");
1978          return NX_ERROR;
1979        }
1980       /* map datatypes of other plateforms */
1981       if (tclass == H5T_STRING)
1982       {
1983         dims = H5Tget_size(pFile->iCurrentT);
1984         memtype_id = H5Tcopy(H5T_C_S1);
1985         H5Tset_size(memtype_id, dims); 
1986       }
1987       else 
1988       {
1989         memtype_id = h5MemType(pFile->iCurrentT);
1990       }
1991
1992      /* read slab */ 
1993      if (mtype == NX_CHAR) {
1994         iRet = H5Dread(pFile->iCurrentD, memtype_id, H5S_ALL, 
1995                     H5S_ALL, H5P_DEFAULT,tmp_data);
1996          data1 = tmp_data + myStart[0];
1997          strncpy((char*)data,data1,(hsize_t)iSize[0]);
1998          free(tmp_data);           
1999      } else {   
2000         iRet = H5Dread(pFile->iCurrentD, memtype_id, memspace, 
2001                     pFile->iCurrentS, H5P_DEFAULT,data);
2002      }   
2003      /* cleanup */
2004      if (tclass == H5T_STRING) { /* we used H5Tcopy */
2005         H5Tclose(memtype_id);
2006      }
2007      H5Sclose(memspace);
2008
2009      if (iRet < 0) 
2010
2011        {
2012          NXReportError( "ERROR: reading slab failed");
2013          return NX_ERROR;
2014        }
2015      return NX_OK;
2016   }
2017
2018   /*-------------------------------------------------------------------------*/
2019
2020   /* Operator function. */
2021
2022   herr_t attr_info(hid_t loc_id, const char *name, const H5A_info_t *unused, void *opdata)
2023   {
2024     *((char**)opdata)=strdup(name);
2025     return 1;
2026   }
2027
2028   NXstatus  NX5getnextattr (NXhandle fileid, NXname pName, int *iLength, int *iType)
2029   {
2030     pNexusFile5 pFile;
2031     hid_t attr_id;
2032     hid_t atype, aspace;
2033     herr_t iRet;
2034     int iPType,rank;
2035     char *iname = NULL, *vlen_str = NULL; 
2036     hsize_t idx, intern_idx=-1;
2037     int vid;
2038     H5O_info_t oinfo;
2039
2040     pFile = NXI5assert (fileid);
2041
2042     vid = getAttVID(pFile);
2043
2044     pName[0] = '\0';
2045     idx=pFile->iAtt5.iCurrentIDX;
2046     iRet=0;
2047 
2048     H5Oget_info(vid, &oinfo);
2049     intern_idx=oinfo.num_attrs;
2050     if(intern_idx == idx) {
2051       killAttVID(pFile,vid);
2052       return NX_EOD;
2053     }
2054
2055     if (intern_idx > idx) {
2056       iRet=H5Aiterate(vid,H5_INDEX_CRT_ORDER,H5_ITER_INC,&idx,attr_info,&iname);
2057     } else {
2058       iRet=0;
2059     } 
2060     intern_idx=-1;
2061     if (iRet < 0) {
2062           NXReportError( "ERROR: iterating through attribute list");
2063           killAttVID(pFile,vid);
2064           return NX_ERROR; 
2065     } 
2066     pFile->iAtt5.iCurrentIDX++;
2067     if (iname != NULL) {
2068       if(strcmp(iname, "NX_class") == 0 && pFile->iCurrentG != 0) {
2069         /*
2070           skip NXclass attribute which is internal
2071         */
2072         killAttVID(pFile, vid);
2073         return NX5getnextattr(fileid, pName, iLength, iType);
2074       }
2075       strcpy(pName, iname);
2076       free(iname);
2077       iname = NULL;
2078     } else {
2079       strcpy(pName,"What is this?");
2080     }
2081     pFile->iCurrentA = H5Aopen_by_name(vid, ".", pName, H5P_DEFAULT, H5P_DEFAULT);
2082     atype  = H5Aget_type(pFile->iCurrentA);
2083     aspace = H5Aget_space(pFile->iCurrentA);
2084     rank = H5Sget_simple_extent_ndims(aspace);
2085     attr_id = H5Tget_class(atype);
2086     if (attr_id==H5T_STRING) {
2087       iPType=NX_CHAR;
2088       readStringAttribute(pFile->iCurrentA, &vlen_str);
2089       rank = strlen(vlen_str);
2090       free(vlen_str);
2091     }
2092     if (rank == 0) {
2093       rank++;
2094     } 
2095     iPType = hdf5ToNXType(attr_id,atype);
2096     *iType=iPType;
2097     *iLength=rank;
2098     H5Tclose(atype);
2099     H5Sclose(aspace);
2100     H5Aclose(pFile->iCurrentA);
2101 
2102     H5Oget_info(vid, &oinfo);
2103     intern_idx=oinfo.num_attrs;
2104
2105     killAttVID(pFile,vid);
2106     return NX_OK;
2107   }
2108 /*-------------------------------------------------------------------------*/
2109
2110   NXstatus  NX5getattr (NXhandle fid, char *name, 
2111                         void *data, int* datalen, int* iType)
2112   {
2113     pNexusFile5 pFile;
2114     int iNew, vid;
2115     herr_t iRet;
2116     hid_t type, atype = -1;
2117     char pBuffer[256];
2118
2119     pFile = NXI5assert (fid);
2120
2121     type = nxToHDF5Type(*iType);
2122
2123     vid = getAttVID(pFile);
2124     iNew = H5Aopen_by_name(vid, ".", name, H5P_DEFAULT, H5P_DEFAULT);
2125     if (iNew < 0) {
2126       sprintf (pBuffer, "ERROR: attribute \"%s\" not found", name);
2127       killAttVID(pFile,vid);
2128       NXReportError( pBuffer);
2129       return NX_ERROR;
2130     }
2131     pFile->iCurrentA = iNew;
2132     /* finally read the data */
2133     if (type==H5T_C_S1)
2134     {
2135        atype = H5Aget_type(pFile->iCurrentA);
2136        iRet = readStringAttributeN(pFile->iCurrentA, data, *datalen);
2137        *datalen = strlen((char*)data);
2138     } else {
2139       iRet = H5Aread(pFile->iCurrentA, type, data);
2140       *datalen=1;
2141     }
2142
2143     if (iRet < 0) {
2144       sprintf (pBuffer, "ERROR: could not read attribute data for \"%s\"", name);
2145       NXReportError( pBuffer);
2146       killAttVID(pFile,vid);
2147       return NX_ERROR;
2148     }
2149
2150     H5Aclose(pFile->iCurrentA);
2151
2152     killAttVID(pFile,vid);
2153     if (type==H5T_C_S1)
2154     {
2155       H5Tclose(atype);
2156     }
2157     return NX_OK;
2158   } 
2159
2160   /*-------------------------------------------------------------------------*/
2161
2162   NXstatus  NX5getattrinfo (NXhandle fid, int *iN)
2163   {
2164     pNexusFile5 pFile;
2165     hid_t idx;
2166     int vid;
2167     H5O_info_t oinfo;
2168   
2169     pFile = NXI5assert (fid);
2170     idx=0;
2171     *iN = idx;
2172
2173     vid = getAttVID(pFile);
2174
2175     H5Oget_info(vid, &oinfo);
2176     idx=oinfo.num_attrs;
2177     if (idx > 0) {
2178       if(pFile->iCurrentG > 0 && pFile->iCurrentD == 0){
2179         *iN = idx -1; 
2180       } else {
2181         *iN = idx;
2182       }
2183     } else {
2184       *iN = 0;   
2185     }
2186     killAttVID(pFile,vid);
2187     return NX_OK;
2188   }
2189
2190
2191   /*-------------------------------------------------------------------------*/
2192   NXstatus  NX5getgroupID (NXhandle fileid, NXlink* sRes)
2193  {
2194    pNexusFile5 pFile;
2195    int datalen, type = NX_CHAR;
2196 
2197    pFile = NXI5assert (fileid);
2198    if (pFile->iCurrentG == 0) {
2199      return NX_ERROR;
2200    } 
2201    else {
2202      /*
2203        this means: if the item is already linked: use the target attribute, else
2204        the path to the current node
2205      */
2206      NXMDisableErrorReporting();
2207      datalen = 1024;
2208      memset(sRes->targetPath,0,datalen*sizeof(char));
2209      if(NX5getattr(fileid,"target",sRes->targetPath,&datalen,&type) != NX_OK){
2210        buildCurrentPath(pFile,sRes->targetPath,1024);
2211      }
2212      NXMEnableErrorReporting();
2213      sRes->linkType = 0;
2214      return NX_OK;
2215    }
2216    /* not reached */
2217    return NX_ERROR;
2218  } 
2219 
2220  /* ------------------------------------------------------------------- */
2221
2222   NXstatus  NX5nativeexternallink(NXhandle fileid, const char* name, const char* externalfile, const char* remotetarget)
2223  {
2224     herr_t iRet;
2225     pNexusFile5 pFile;
2226     hid_t openwhere;
2227
2228     pFile = NXI5assert(fileid);
2229
2230      if (pFile->iCurrentG <= 0) {
2231          openwhere = pFile->iFID;
2232      } else {
2233          openwhere = pFile->iCurrentG;
2234      }
2235
2236     iRet = H5Lcreate_external(externalfile, remotetarget, openwhere, name, H5P_DEFAULT, H5P_DEFAULT);
2237     if (iRet < 0) {
2238       NXReportError("ERROR: making external link failed");
2239       return NX_ERROR;
2240     }
2241     return NX_OK;
2242  }
2243  /* ------------------------------------------------------------------- */
2244
2245   NXstatus  NX5nativeinquirefile(NXhandle fileid, char* externalfile, const int filenamelen)
2246  {
2247     pNexusFile5 pFile;
2248     ssize_t name_size;
2249     hid_t openthing;
2250
2251     pFile = NXI5assert(fileid);
2252      if (pFile->iCurrentD > 0) {
2253          openthing = pFile->iCurrentD;
2254      } else if (pFile->iCurrentG > 0) {
2255          openthing = pFile->iCurrentG;
2256      } else {
2257          openthing = pFile->iFID;
2258      }
2259
2260     name_size = H5Fget_name(openthing, externalfile, filenamelen);
2261
2262     // Check for failure again
2263     if( name_size < 0 ) {
2264         NXReportError("ERROR: retrieving file name");
2265         return NX_ERROR;
2266     }
2267     return NX_OK;
2268  }
2269  /* ------------------------------------------------------------------- */
2270
2271   NXstatus  NX5nativeisexternallink(NXhandle fileid, const char* name, char* url, const int urllen)
2272  {
2273     pNexusFile5 pFile;
2274     herr_t ret;
2275     H5L_info_t link_buff;
2276     char linkval_buff[1024];
2277     const char *filepath = NULL, *objpath = NULL;
2278     size_t val_size;
2279
2280     pFile = NXI5assert(fileid);
2281     memset(url, 0, urllen);
2282
2283     ret = H5Lget_info(pFile->iFID, name, &link_buff, H5P_DEFAULT);
2284     if (ret < 0 || link_buff.type != H5L_TYPE_EXTERNAL) {
2285       return NX_ERROR;
2286     }
2287
2288     val_size = link_buff.u.val_size;
2289     if (val_size > sizeof(linkval_buff)) {
2290       NXReportError("ERROR: linkval_buff too small");
2291       return NX_ERROR;
2292     }
2293
2294     ret = H5Lget_val(pFile->iFID, name, linkval_buff, val_size, H5P_DEFAULT);
2295     if (ret < 0) {
2296       NXReportError("ERROR: H5Lget_val failed");
2297       return NX_ERROR;
2298     }
2299
2300     ret = H5Lunpack_elink_val(linkval_buff, val_size, NULL, &filepath, &objpath);
2301     if (ret < 0) {
2302       NXReportError("ERROR: H5Lunpack_elink_val failed");
2303       return NX_ERROR;
2304     }
2305
2306    snprintf(url, urllen-1, "nxfile://%s#%s", filepath, objpath);
2307    return NX_OK;
2308     
2309  }
2310  /* ------------------------------------------------------------------- */
2311
2312  NXstatus  NX5sameID (NXhandle fileid, NXlink* pFirstID, NXlink* pSecondID)
2313  {
2314    NXI5assert(fileid);
2315    if ((strcmp(pFirstID->targetPath,pSecondID->targetPath) == 0)){
2316       return NX_OK;
2317    } else {
2318       return NX_ERROR;
2319    }
2320  }
2321 
2322 /*-------------------------------------------------------------------------*/
2323 
2324  NXstatus  NX5initattrdir (NXhandle fid)
2325  {
2326    pNexusFile5 pFile;
2327       
2328    pFile = NXI5assert (fid);
2329    NXI5KillAttDir (pFile);
2330    return NX_OK;
2331  }
2332
2333  /*-------------------------------------------------------------------------*/
2334 
2335  NXstatus  NX5initgroupdir (NXhandle fid)
2336  {
2337    pNexusFile5 pFile;
2338       
2339    pFile = NXI5assert (fid);
2340    NXI5KillDir (pFile);
2341    return NX_OK;
2342  }
2343/*------------------------------------------------------------------------*/
2344void NX5assignFunctions(pNexusFunction fHandle)
2345{
2346               
2347      fHandle->nxclose=NX5close;
2348      fHandle->nxreopen=NX5reopen;
2349      fHandle->nxflush=NX5flush;
2350      fHandle->nxmakegroup=NX5makegroup;
2351      fHandle->nxopengroup=NX5opengroup;
2352      fHandle->nxclosegroup=NX5closegroup;
2353      fHandle->nxmakedata64=NX5makedata64;
2354      fHandle->nxcompmakedata64=NX5compmakedata64;
2355      fHandle->nxcompress=NX5compress;
2356      fHandle->nxopendata=NX5opendata;
2357      fHandle->nxclosedata=NX5closedata;
2358      fHandle->nxputdata=NX5putdata;
2359      fHandle->nxputattr=NX5putattr;
2360      fHandle->nxputslab64=NX5putslab64;   
2361      fHandle->nxgetdataID=NX5getdataID;
2362      fHandle->nxmakelink=NX5makelink;
2363      fHandle->nxmakenamedlink=NX5makenamedlink;
2364      fHandle->nxgetdata=NX5getdata;
2365      fHandle->nxgetinfo64=NX5getinfo64;
2366      fHandle->nxgetnextentry=NX5getnextentry;
2367      fHandle->nxgetslab64=NX5getslab64;
2368      fHandle->nxgetnextattr=NX5getnextattr;
2369      fHandle->nxgetattr=NX5getattr;
2370      fHandle->nxgetattrinfo=NX5getattrinfo;
2371      fHandle->nxgetgroupID=NX5getgroupID;
2372      fHandle->nxgetgroupinfo=NX5getgroupinfo;
2373      fHandle->nxsameID=NX5sameID;
2374      fHandle->nxinitgroupdir=NX5initgroupdir;
2375      fHandle->nxinitattrdir=NX5initattrdir;
2376      fHandle->nxprintlink=NX5printlink;
2377      fHandle->nxnativeexternallink=NX5nativeexternallink;
2378      fHandle->nxnativeinquirefile=NX5nativeinquirefile;
2379      fHandle->nxnativeisexternallink=NX5nativeisexternallink;
2380}
2381
2382#endif /* HDF5 */
Note: See TracBrowser for help on using the repository browser.