H5Part.c 54.8 KB
Newer Older
gsell's avatar
gsell committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168
/*
	H5Part C API	
*/

/*
  H5PartSetStep()
  H5PartGetNumParticles()
  H5PartSetView()

  H5PartHasView()
*/

#include <stdio.h>
#include <stdlib.h>
#include <stdarg.h>	/* va_arg - System dependent ?! */
#include <string.h>
#include <unistd.h>
#include <errno.h>
#include <fcntl.h>
#include <hdf5.h>

#include "H5PartTypes.h"
#include "H5Part.h"
#include "H5PartPrivate.h"
#include "H5PartErrors.h"

/********* Private Variable Declarations *************/

static unsigned			_debug = 0;
static h5part_int64_t		_errno = H5PART_SUCCESS;
static h5part_error_handler	_err_handler = H5PartDefaultErrorHandler;
static char *__funcname = "NONE";

/********** Declaration of private functions ******/

static h5part_int64_t
_init(
	void
	);

static h5part_int64_t
_file_is_valid (
	const H5PartFile *f
	);

static herr_t
_h5_error_handler (
	void *
	);

/*========== File Opening/Closing ===============*/

/*!
  Opens file with specified filename. 

  If you open with flag \c H5PART_WRITE, it will truncate any
  file with the specified filename and start writing to it. If 
  you open with \c H5PART_APPEND, then you can append new timesteps.
  If you open with \c H5PART_READ, then it will open the file
  readonly.

  The typical extension for these files is \c .h5.
  
  H5PartFile should be treated as an essentially opaque
  datastructure.  It acts as the file handle, but internally
  it maintains several key state variables associated with 
  the file.

  \return	File handle or \c NULL
 */

H5PartFile*
H5PartOpenFileParallel (
	const char *filename,	/*!< [in] The name of the data file to open. */
	unsigned flags		/*!< [in] The access mode for the file. */
#ifdef PARALLEL_IO
	,MPI_Comm comm		/*!< [in] MPI communicator */
#endif
) {
	SET_FNAME ( "H5PartOpenFileParallel" );

	if ( _init() < 0 ) {
		HANDLE_H5PART_INIT_ERR;
		return NULL;
	}
	_errno = H5PART_SUCCESS;
	H5PartFile *f = NULL;

	f = (H5PartFile*) malloc( sizeof (H5PartFile) );
	if( f == NULL ) {
		HANDLE_H5PART_NOMEM_ERR;
		goto error_cleanup;
	}
	memset (f, 0, sizeof (H5PartFile));
	f->xfer_prop = f->create_prop = f->access_prop = H5P_DEFAULT;

#ifdef PARALLEL_IO
	/* for the SP2... perhaps different for linux */
	MPI_Info info = MPI_INFO_NULL;

	if (MPI_Comm_size (comm, &f->nprocs) != MPI_SUCCESS) {
		HANDLE_MPI_COMM_SIZE_ERR;
		goto error_cleanup;
	}
	if (MPI_Comm_rank (comm, &f->myproc) != MPI_SUCCESS) {
		HANDLE_MPI_COMM_RANK_ERR;
		goto error_cleanup;
	}

	f->pnparticles = malloc (f->nprocs * sizeof (h5part_int64_t));
	if (f->pnparticles == NULL) {
		HANDLE_H5PART_NOMEM_ERR;
		goto error_cleanup;
	}
		
	f->access_prop = H5Pcreate (H5P_FILE_ACCESS);
	if (f->access_prop < 0) {
		HANDLE_H5P_CREATE_ERR;
		goto error_cleanup;
	}

	if (H5Pset_fapl_mpio (f->access_prop, comm, info) < 0) {
		HANDLE_H5P_SET_FAPL_MPIO_ERR;
		goto error_cleanup;
	}
		
	/* create_prop: tunable parameters like blocksize and btree sizes */
	/* f->create_prop = H5Pcreate(H5P_FILE_CREATE); */
	f->create_prop = H5P_DEFAULT;

	/* currently create_prop is empty */
	/* xfer_prop:  also used for parallel I/O, during actual writes
	   rather than the access_prop which is for file creation. */
	f->xfer_prop = H5Pcreate (H5P_DATASET_XFER);
	if (f->xfer_prop < 0) {
		HANDLE_H5P_CREATE_ERR;
		goto error_cleanup;
	}

	if (H5Pset_dxpl_mpio (f->xfer_prop, H5FD_MPIO_COLLECTIVE) < 0) {
		HANDLE_H5P_SET_DXPL_MPIO_ERR;
		goto error_cleanup;
	}

	f->comm = comm;
#endif
	if ( flags == H5PART_READ ) {
		f->file = H5Fopen (filename, H5F_ACC_RDONLY, f->access_prop);
	}
	else if ( flags == H5PART_WRITE ){
		f->file = H5Fcreate (filename, H5F_ACC_TRUNC, f->create_prop,
				     f->access_prop);
	}
	else if ( flags == H5PART_APPEND ) {
		int fd = open (filename, O_RDONLY, 0);
		if ( (fd == -1) && (errno == ENOENT) ) {
			f->file = H5Fcreate(filename, H5F_ACC_TRUNC,
					    f->create_prop, f->access_prop);
		}
		else if (fd != -1) {
			close (fd);
			f->file = H5Fopen (filename, H5F_ACC_RDWR,
					   f->access_prop);
			/*
			  The following function call returns an error,
			  if f->file < 0. But we can safely ignore this.
			*/
			f->timestep = _H5Part_get_num_objects_matching_pattern(
gsell's avatar
gsell committed
169
				f->file, "/", H5G_GROUP, H5PART_PARTICLES_GROUP );
gsell's avatar
gsell committed
170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189
			if ( f->timestep < 0 ) goto error_cleanup;
		}
	}
	else {
		HANDLE_H5PART_FILE_ACCESS_TYPE_ERR ( flags );
		goto error_cleanup;
	}

	if (f->file < 0) {
		HANDLE_H5F_OPEN_ERR ( filename, flags );
		goto error_cleanup;
	}
	f->mode = flags;
	f->timegroup = 0;
	f->shape = 0;
	f->diskshape = H5S_ALL;
	f->memshape = H5S_ALL;
	f->viewstart = -1;
	f->viewend = -1;

gsell's avatar
gsell committed
190 191 192 193 194 195
	_H5Part_print_debug (
		"Proc[%d]: Opened file \"%s\" val=%lld",
		f->myproc,
		filename,
		(long long)(size_t)f );

gsell's avatar
gsell committed
196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277
	return f;

 error_cleanup:
	if (f != NULL ) {
		if (f->pnparticles != NULL) {
			free (f->pnparticles);
		}
		free (f);
	}
	return NULL;
}

/*!
  Opens file with specified filename. 

  If you open with flag \c H5PART_WRITE, it will truncate any
  file with the specified filename and start writing to it. If 
  you open with \c H5PART_APPEND, then you can append new timesteps.
  If you open with \c H5PART_READ, then it will open the file
  readonly.

  The typical extension for these files is \c .h5.
  
  H5PartFile should be treated as an essentially opaque
  datastructure.  It acts as the file handle, but internally
  it maintains several key state variables associated with 
  the file.

  \return	File handle or \c NULL
 */

H5PartFile*
H5PartOpenFile (
	const char *filename,	/*!< [in] The name of the data file to open. */
	unsigned flags		/*!< [in] The access mode for the file. */
	) {

	SET_FNAME ( "H5PartOpenFile" );
	if ( _init() < 0 ) {
		HANDLE_H5PART_INIT_ERR;
		return NULL;
	}

	_errno = H5PART_SUCCESS;
	H5PartFile *f = NULL;

	f = (H5PartFile*) malloc( sizeof (H5PartFile) );
	if( f == NULL ) {
		HANDLE_H5PART_NOMEM_ERR;
		goto error_cleanup;
	}
	memset (f, 0, sizeof (H5PartFile));
	f->xfer_prop = f->create_prop = f->access_prop = H5P_DEFAULT;

#ifdef PARALLEL_IO
	f->pnparticles = 0;
	f->comm = MPI_COMM_WORLD;
	f->nprocs = 1;
	f->myproc = 0;
#endif
	if (flags == H5PART_READ) {
		f->file = H5Fopen (filename, H5F_ACC_RDONLY, f->access_prop);
	}
	else if (flags == H5PART_WRITE){
		f->file = H5Fcreate (filename, H5F_ACC_TRUNC, f->create_prop,
				     f->access_prop);
	}
	else if (flags == H5PART_APPEND) {
		int fd = open (filename, O_RDONLY, 0);
		if ( (fd == -1) && (errno == ENOENT) ) {
			f->file = H5Fcreate(filename, H5F_ACC_TRUNC,
					    f->create_prop, f->access_prop);
		}
		else if (fd != -1) {
			close (fd);
			f->file = H5Fopen (filename, H5F_ACC_RDWR,
					   f->access_prop);
			/*
			  The following function call returns an error,
			  if f->file < 0. But we can safely ignore it
			*/
			f->timestep = _H5Part_get_num_objects_matching_pattern(
gsell's avatar
gsell committed
278
				f->file, "/", H5G_GROUP, H5PART_PARTICLES_GROUP );
gsell's avatar
gsell committed
279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299
			if ( f->timestep < 0 ) goto error_cleanup;
		}
	}
	else {
		HANDLE_H5PART_FILE_ACCESS_TYPE_ERR ( flags );
		goto error_cleanup;
	}

	if (f->file < 0) {
		HANDLE_H5F_OPEN_ERR ( filename, flags );
		goto error_cleanup;
	}

	f->mode = flags;
	f->timegroup = 0;
	f->shape = 0;
	f->diskshape = H5S_ALL;
	f->memshape = H5S_ALL;
	f->viewstart = -1;
	f->viewend = -1;

gsell's avatar
gsell committed
300 301 302 303 304
	_H5Part_print_debug (
		"Opened file \"%s\" val=%lld",
		filename,
		(long long)(size_t)f );

gsell's avatar
gsell committed
305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491
	return f;

 error_cleanup:
	if (f != NULL ) {
		if (f->pnparticles != NULL) {
			free (f->pnparticles);
		}
		free (f);
	}
	return NULL;
}

/*!
  Checks if a file was successfully opened.

  \return	\c H5PART_SUCCESS or error code
 */
static h5part_int64_t
_file_is_valid (
	const H5PartFile *f	/*!< filehandle  to check validity of */
	) {

	if( f == NULL )
		return H5PART_ERR_BADFD;
	else if(f->file > 0)
		return H5PART_SUCCESS;
	else
		return H5PART_ERR_BADFD;
}

/*!
  Closes an open file.

  \return	\c H5PART_SUCCESS or error code
*/
h5part_int64_t
H5PartCloseFile (
	H5PartFile *f		/*!< [in] filehandle of the file to close */
	) {

	SET_FNAME ( "H5PartCloseFile" );
	herr_t r = 0;
	_errno = H5PART_SUCCESS;

	CHECK_FILEHANDLE ( f );

	if ( f->block && f->close_block ) {
		(*f->close_block) ( f );
		f->block = NULL;
		f->close_block = NULL;
	}

	if( f->shape > 0 ) {
		r = H5Sclose( f->shape );
		if ( r < 0 ) HANDLE_H5S_CLOSE_ERR;
		f->shape = 0;
	}
	if( f->timegroup > 0 ) {
		r = H5Gclose( f->timegroup );
		if ( r < 0 ) HANDLE_H5G_CLOSE_ERR;
		f->timegroup = 0;
	}
	if( f->diskshape != H5S_ALL ) {
		r = H5Sclose( f->diskshape );
		if ( r < 0 ) HANDLE_H5S_CLOSE_ERR;
		f->diskshape = 0;
	}
	if( f->xfer_prop != H5P_DEFAULT ) {
		r = H5Pclose( f->xfer_prop );
		if ( r < 0 ) HANDLE_H5P_CLOSE_ERR ( "f->xfer_prop" );
		f->xfer_prop = H5P_DEFAULT;
	}
	if( f->access_prop != H5P_DEFAULT ) {
		r = H5Pclose( f->access_prop );
		if ( r < 0 ) HANDLE_H5P_CLOSE_ERR ( "f->access_prop" );
		f->access_prop = H5P_DEFAULT;
	}  
	if( f->create_prop != H5P_DEFAULT ) {
		r = H5Pclose( f->create_prop );
		if ( r < 0 ) HANDLE_H5P_CLOSE_ERR ( "f->create_prop" );
		f->create_prop = H5P_DEFAULT;
	}
	if ( f->file ) {
		r = H5Fclose( f->file );
		if ( r < 0 ) HANDLE_H5F_CLOSE_ERR;
		f->file = 0;
	}
	if( f->pnparticles ) {
		free( f->pnparticles );
	}
	free( f );

	return _errno;
}

/*============== File Writing Functions ==================== */
/*!
  Set number of particles for current time-step.

  This function's sole purpose is to prevent 
  needless creation of new HDF5 DataSpace handles if the number of 
  particles is invariant throughout the simulation. That's its only reason 
  for existence. After you call this subroutine, all subsequent 
  operations will assume this number of particles will be written.


  \return	\c H5PART_SUCCESS or error code
 */
h5part_int64_t
H5PartSetNumParticles (
	H5PartFile *f,			/*!< [in] Handle to open file */
	h5part_int64_t nparticles	/*!< [in] Number of particles */
	) {

	SET_FNAME ( "H5PartSetNumParticles" );
	int r;
#ifdef PARALLEL_IO
#ifdef HDF5V160
	hssize_t start[1];
#else
	hsize_t start[1];
#endif

	hsize_t stride[1];
	hsize_t count[1];
	hsize_t total;
	hsize_t dmax = H5S_UNLIMITED;
	register int i;
#endif

	CHECK_FILEHANDLE( f );

#ifndef PARALLEL_IO
	/*
	  if we are not using parallel-IO, there is enough information
	   to know that we can short circuit this routine.  However,
	   for parallel IO, this is going to cause problems because
	   we don't know if things have changed globally
	*/
	if ( f->nparticles == nparticles ) {
		return H5PART_SUCCESS;
	}
#endif
	if ( f->diskshape != H5S_ALL ) {
		r = H5Sclose( f->diskshape );
		if ( r < 0 ) return HANDLE_H5S_CLOSE_ERR;
		f->diskshape = H5S_ALL;
	}
	if(f->memshape != H5S_ALL) {
		r = H5Sclose( f->memshape );
		if ( r < 0 ) return HANDLE_H5S_CLOSE_ERR;
		f->memshape = H5S_ALL;
	}
	if( f->shape ) {
		r = H5Sclose(f->shape);
		if ( r < 0 ) return HANDLE_H5S_CLOSE_ERR;
	}
	f->nparticles =(hsize_t) nparticles;
#ifndef PARALLEL_IO
	f->shape = H5Screate_simple (1,
				     &(f->nparticles),
				     NULL);
	if ( f->shape < 0 ) HANDLE_H5S_CREATE_SIMPLE_ERR ( f->nparticles );

#else /* PARALLEL_IO */
	/*
	  The Gameplan here is to declare the overall size of the on-disk
	  data structure the same way we do for the serial case.  But
	  then we must have additional "DataSpace" structures to define
	  our in-memory layout of our domain-decomposed portion of the particle
	  list as well as a "selection" of a subset of the on-disk 
	  data layout that will be written in parallel to mutually exclusive
	  regions by all of the processors during a parallel I/O operation.
	  These are f->shape, f->memshape and f->diskshape respectively.
	*/

	/*
	  acquire the number of particles to be written from each MPI process
	*/
	r = MPI_Allgather (
		&nparticles, 1, MPI_LONG_LONG,
		f->pnparticles, 1, MPI_LONG_LONG,
		f->comm);
	if ( r != MPI_SUCCESS) {
		return HANDLE_MPI_ALLGATHER_ERR;
	}
	if ( f->myproc == 0 ) {
gsell's avatar
gsell committed
492
		_H5Part_print_debug ( "Particle offsets:" );
gsell's avatar
gsell committed
493
		for(i=0;i<f->nprocs;i++) 
gsell's avatar
gsell committed
494 495
			_H5Part_print_debug ( "\tnp=%lld",
					      (long long) f->pnparticles[i] );
gsell's avatar
gsell committed
496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555
	}
	/* should I create a selection here? */

	/* compute start offsets */
	stride[0] = 1;
	start[0] = 0;
	for (i=0; i<f->myproc; i++) {
		start[0] += f->pnparticles[i];
	}
	


        /* compute total nparticles */
	total = 0;
	for (i=0; i < f->nprocs; i++) {
		total += f->pnparticles[i];
	}

	/* declare overall datasize */
	f->shape = H5Screate_simple (1, &total, &total);
	if (f->shape < 0) return HANDLE_H5S_CREATE_SIMPLE_ERR ( total );


	/* declare overall data size  but then will select a subset */
	f->diskshape = H5Screate_simple (1, &total, &total);
	if (f->diskshape < 0) return HANDLE_H5S_CREATE_SIMPLE_ERR ( total );

	/* declare local memory datasize */
	f->memshape = H5Screate_simple (1, &(f->nparticles), &dmax);
	if (f->memshape < 0)
		return HANDLE_H5S_CREATE_SIMPLE_ERR ( f->nparticles );

	count[0] = nparticles;
	r = H5Sselect_hyperslab (
		f->diskshape,
		H5S_SELECT_SET,
		start,
		stride,
		count, NULL );
	if ( r < 0 ) return HANDLE_H5S_SELECT_HYPERSLAB_ERR;

	if ( f->timegroup < 0 ) {
		r = H5PartSetStep ( f, 0 );
		if ( r < 0 ) return HANDLE_H5PART_SETSTEP_ERR( r, 0 );
		
	}
#endif
	return H5PART_SUCCESS;
}

static h5part_int64_t
_write_data (
	H5PartFile *f,		/*!< IN: Handle to open file */
	char *name,		/*!< IN: Name to associate array with */
	void *array,		/*!< IN: Array to commit to disk */
	hid_t type		/*!< IN: Type of data */
	) {
	herr_t herr;
	hid_t dataset_id;

gsell's avatar
gsell committed
556
	_H5Part_print_debug ( "Create a dataset[%s] mounted on the "
gsell's avatar
gsell committed
557
			      "timestep %lld",
gsell's avatar
gsell committed
558
			      name, (long long)f->timestep );
gsell's avatar
gsell committed
559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699

	dataset_id = H5Dcreate ( 
		f->timegroup,
		name,
		type,
		f->shape,
		H5P_DEFAULT );
	if ( dataset_id < 0 )
		return HANDLE_H5D_CREATE_ERR ( name, f->timestep );

	herr = H5Dwrite (
		dataset_id,
		type,
		f->memshape,
		f->diskshape,
		H5P_DEFAULT,
		array );
	if ( herr < 0 ) return HANDLE_H5D_WRITE_ERR ( name, f->timestep );

	herr = H5Dclose ( dataset_id );
	if ( herr < 0 ) return HANDLE_H5D_CLOSE_ERR;

	return H5PART_SUCCESS;
}

/*!
  Write array of 64 bit floating point data to file.

  After setting the number of particles with \c H5PartSetNumParticles() and
  the current timestep using \c H5PartSetStep(), you can start writing datasets
  into the file. Each dataset has a name associated with it (chosen by the
  user) in order to facilitate later retrieval. The name of the dataset is
  specified in the parameter \c name, which must be a null-terminated string.

  There are no restrictions on naming of datasets, but it is useful to arrive
  at some common naming convention when sharing data with other groups.

  The writing routines also implicitly store the datatype of the array so that
  the array can be reconstructed properly on other systems with incompatible
  type representations.

  All data that is written after setting the timestep is associated with that
  timestep. While the number of particles can change for each timestep, you
  cannot change the number of particles in the middle of a given timestep.

  The data is committed to disk before the routine returns.

  \return	\c H5PART_SUCCESS or error code
 */
h5part_int64_t
H5PartWriteDataFloat64 (
	H5PartFile *f,		/*!< [in] Handle to open file */
	char *name,		/*!< [in] Name to associate array with */
	h5part_float64_t *array	/*!< [in] Array to commit to disk */
	) {

	SET_FNAME ( "H5PartWriteDataFloat64" );
	herr_t herr;

	CHECK_FILEHANDLE ( f );
	CHECK_WRITABLE_MODE( f );
	CHECK_TIMEGROUP( f );

	herr = _write_data ( f, name, (void*)array, H5T_NATIVE_DOUBLE );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/*!
  Write array of 64 bit integer data to file.

  After setting the number of particles with \c H5PartSetNumParticles() and
  the current timestep using \c H5PartSetStep(), you can start writing datasets
  into the file. Each dataset has a name associated with it (chosen by the
  user) in order to facilitate later retrieval. The name of the dataset is
  specified in the parameter \c name, which must be a null-terminated string.

  There are no restrictions on naming of datasets, but it is useful to arrive
  at some common naming convention when sharing data with other groups.

  The writing routines also implicitly store the datatype of the array so that
  the array can be reconstructed properly on other systems with incompatible
  type representations.

  All data that is written after setting the timestep is associated with that
  timestep. While the number of particles can change for each timestep, you
  cannot change the number of particles in the middle of a given timestep.

  The data is committed to disk before the routine returns.

  \return	\c H5PART_SUCCESS or error code
 */
h5part_int64_t
H5PartWriteDataInt64 (
	H5PartFile *f,		/*!< [in] Handle to open file */
	char *name,		/*!< [in] Name to associate array with */
	h5part_int64_t *array	/*!< [in] Array to commit to disk */
	) {

	SET_FNAME ( "H5PartOpenWriteDataInt64" );

	herr_t herr;

	CHECK_FILEHANDLE ( f );
	CHECK_WRITABLE_MODE( f );
	CHECK_TIMEGROUP( f );

	herr = _write_data ( f, name, (void*)array, H5T_NATIVE_INT64 );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/********************** reading and writing attribute ************************/

/********************** private functions to handle attributes ***************/
/*!
   Normalize HDF5 type
*/
hid_t
_H5Part_normalize_h5_type (
	hid_t type
	) {
	H5T_class_t tclass = H5Tget_class ( type );
	int size = H5Tget_size ( type );

	switch ( tclass ){
	case H5T_INTEGER:
		if ( size==8 ) {
			return H5T_NATIVE_INT64;
		}
		else if ( size==1 ) {
			return H5T_NATIVE_CHAR;
		}
		break;
	case H5T_FLOAT:
		return H5T_NATIVE_DOUBLE;
	default:
		; /* NOP */
	}
gsell's avatar
gsell committed
700
	_H5Part_print_warn ( "Unknown type %d", (int)type );
gsell's avatar
gsell committed
701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998

	return -1;
}

h5part_int64_t
_H5Part_read_attrib (
	hid_t id,
	const char *attrib_name,
	void *attrib_value
	) {

	herr_t herr;
	hid_t attrib_id;
	hid_t space_id;
	hid_t type_id;
	hid_t mytype;
	hsize_t nelem;

	attrib_id = H5Aopen_name ( id, attrib_name );
	if ( attrib_id <= 0 ) return HANDLE_H5A_OPEN_NAME_ERR( attrib_name );

	mytype = H5Aget_type ( attrib_id );
	if ( mytype < 0 ) return HANDLE_H5A_GET_TYPE_ERR;

	space_id = H5Aget_space ( attrib_id );
	if ( space_id < 0 ) return HANDLE_H5A_GET_SPACE_ERR;

	nelem = H5Sget_simple_extent_npoints ( space_id );
	if ( nelem < 0 ) return HANDLE_H5S_GET_SIMPLE_EXTENT_NPOINTS_ERR;

	type_id = _H5Part_normalize_h5_type ( mytype );

	herr = H5Aread (attrib_id, type_id, attrib_value );
	if ( herr < 0 ) return HANDLE_H5A_READ_ERR;

	herr = H5Sclose ( space_id );
	if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;

	herr = H5Tclose ( mytype );
	if ( herr < 0 ) return HANDLE_H5T_CLOSE_ERR;

	herr = H5Aclose ( attrib_id );
	if ( herr < 0 ) return HANDLE_H5A_CLOSE_ERR;

	return H5PART_SUCCESS;
}

h5part_int64_t
_H5Part_write_attrib (
	hid_t id,
	const char *attrib_name,
	const hid_t attrib_type,
	const void *attrib_value,
	const hsize_t attrib_nelem
	) {

	herr_t herr;
	hid_t space_id;
	hid_t attrib_id;

	space_id = H5Screate_simple (1, &attrib_nelem, NULL);
	if ( space_id < 0 )
		return HANDLE_H5S_CREATE_SIMPLE_ERR ( attrib_nelem );

	attrib_id = H5Acreate ( 
		id,
		attrib_name,
		attrib_type,
		space_id,
		H5P_DEFAULT );
	if ( attrib_id < 0 ) return HANDLE_H5A_CREATE_ERR ( attrib_name );

	herr = H5Awrite ( attrib_id, attrib_type, attrib_value);
	if ( herr < 0 ) return HANDLE_H5A_WRITE_ERR ( attrib_name );

	herr = H5Aclose ( attrib_id );
	if ( herr < 0 ) return HANDLE_H5A_CLOSE_ERR;

	herr = H5Sclose ( space_id );
	if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;

	return H5PART_SUCCESS;
}

h5part_int64_t
_H5Part_get_attrib_info (
	hid_t id,
	const h5part_int64_t attrib_idx,
	char *attrib_name,
	const h5part_int64_t len_attrib_name,
	h5part_int64_t *attrib_type,
	h5part_int64_t *attrib_nelem
	) {

	herr_t herr;
	hid_t attrib_id;
	hid_t mytype;
	hid_t space_id;

	attrib_id = H5Aopen_idx ( id, attrib_idx );
	if ( attrib_id < 0 ) return HANDLE_H5A_OPEN_IDX_ERR ( attrib_idx );

	mytype = H5Aget_type ( attrib_id );
	if ( attrib_id < 0 ) return HANDLE_H5A_GET_TYPE_ERR;

	space_id =  H5Aget_space ( attrib_id );
	if ( attrib_id < 0 ) return HANDLE_H5A_GET_SPACE_ERR;

	if ( attrib_nelem ) {
		*attrib_nelem = H5Sget_simple_extent_npoints ( space_id );
		if ( *attrib_nelem < 0 )
			return HANDLE_H5S_GET_SIMPLE_EXTENT_NPOINTS_ERR;
	}
	if ( attrib_name ) {
		herr = H5Aget_name (
			attrib_id,
			len_attrib_name,
			attrib_name );
		if ( herr < 0 ) return HANDLE_H5A_GET_NAME_ERR;
	}
	if ( attrib_type ) {
		*attrib_type = _H5Part_normalize_h5_type ( mytype );
	}
	herr = H5Sclose ( space_id );
 	if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;

	herr = H5Tclose ( mytype );
	if ( herr < 0 ) return HANDLE_H5T_CLOSE_ERR;

	herr = H5Aclose ( attrib_id);
	if ( herr < 0 ) return HANDLE_H5A_CLOSE_ERR;

	return H5PART_SUCCESS;
}

/********************** attribute API ****************************************/

/*!
  Writes a string attribute bound to a file.

  This function creates a new attribute \c name with the string \c value as
  content. The attribute is bound to the file associated with the file handle 
  \c f.

  If the attribute already exists an error will be returned. There
  is currently no way to change the content of an existing attribute.

  \return	\c H5PART_SUCCESS or error code   
*/

h5part_int64_t
H5PartWriteFileAttribString (
	H5PartFile *f,		/*!< [in] Handle to open file */
	const char *attrib_name,/*!< [in] Name of attribute to create */
	const char *attrib_value/*!< [in] Value of attribute */ 
	) {

	SET_FNAME ( "H5PartWriteFileAttribString" );

	CHECK_FILEHANDLE ( f );
	CHECK_WRITABLE_MODE( f );

	hid_t group_id = H5Gopen(f->file,"/");
	if ( group_id < 0 ) return HANDLE_H5G_OPEN_ERR( "/" );

	h5part_int64_t herr = _H5Part_write_attrib (
		group_id,
		attrib_name,
		H5T_NATIVE_CHAR,
		attrib_value,
		strlen ( attrib_value ) + 1 );
	if ( herr < 0 ) return herr;

	herr = H5Gclose ( group_id );
	if ( herr < 0 ) return HANDLE_H5G_CLOSE_ERR;

	return H5PART_SUCCESS;
}

/*!
  Writes a string attribute bound to the current time-step.

  This function creates a new attribute \c name with the string \c value as
  content. The attribute is bound to the current time step in the file given
  by the file handle \c f.

  If the attribute already exists an error will be returned. There
  is currently no way to change the content of an existing attribute.

  \return	\c H5PART_SUCCESS or error code   
*/

h5part_int64_t
H5PartWriteStepAttribString (
	H5PartFile *f,		/*!< [in] Handle to open file */
	const char *attrib_name,/*!< [in] Name of attribute to create */
	const char *attrib_value/*!< [in] Value of attribute */ 
	) {

	SET_FNAME ( "H5PartWriteStepAttribString" );

	CHECK_FILEHANDLE ( f );
	CHECK_WRITABLE_MODE( f );
	CHECK_TIMEGROUP( f );

	h5part_int64_t herr = _H5Part_write_attrib (
		f->timegroup,
		attrib_name,
		H5T_NATIVE_CHAR,
		attrib_value,
		strlen ( attrib_value ) + 1 );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/*!
  Writes a attribute bound to the current time-step.

  This function creates a new attribute \c name with the string \c value as
  content. The attribute is bound to the current time step in the file given
  by the file handle \c f.

  The value of the attribute is given the parameter \c type, which must be one
  of \c H5T_NATIVE_DOUBLE, \c H5T_NATIVE_INT64 of \c H5T_NATIVE_CHAR, the array
  \c value and the number of elements \c nelem in the array.

  If the attribute already exists an error will be returned. There
  is currently no way to change the content of an existing attribute.

  \return	\c H5PART_SUCCESS or error code   
*/

h5part_int64_t
H5PartWriteStepAttrib (
	H5PartFile *f,			/*!< [in] Handle to open file */
	const char *attrib_name,	/*!< [in] Name of attribute */
	const h5part_int64_t attrib_type,/*!< [in] Type of value. */
	const void *attrib_value,	/*!< [in] Value of attribute */ 
	const h5part_int64_t attrib_nelem/*!< [in] Number of elements */
	){

	SET_FNAME ( "H5PartWriteStepAttrib" );

	herr_t herr;

	CHECK_FILEHANDLE ( f );
	CHECK_WRITABLE_MODE( f );
	CHECK_TIMEGROUP( f );

	herr = _H5Part_write_attrib (
		f->timegroup,
		attrib_name,
		attrib_type,
		attrib_value,
		attrib_nelem );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/*!
  Writes a attribute bound to a file.

  This function creates a new attribute \c name with the string \c value as
  content. The attribute is bound to the file file given by the file handle
  \c f.

  The value of the attribute is given the parameter \c type, which must be one
  of H5T_NATIVE_DOUBLE, H5T_NATIVE_INT64 of H5T_NATIVE_CHAR, the array \c value
  and the number of elements \c nelem in the array.

  If the attribute already exists an error will be returned. There
  is currently no way to change the content of an existing attribute.

  \return	\c H5PART_SUCCESS or error code   
*/

h5part_int64_t
H5PartWriteFileAttrib (
	H5PartFile *f,			/*!< [in] Handle to open file */
	const char *attrib_name,	/*!< [in] Name of attribute */
	const h5part_int64_t attrib_type,/*!< [in] Type of value. */
	const void *attrib_value,	/*!< [in] Value of attribute */ 
	const h5part_int64_t attrib_nelem/*!< [in] Number of elements */
	) {

	SET_FNAME ( "H5PartWriteFileAttrib" );

	herr_t herr;
	hid_t group_id;

	CHECK_FILEHANDLE ( f );
	CHECK_WRITABLE_MODE ( f );

	group_id = H5Gopen(f->file,"/");
	if ( group_id < 0 ) return HANDLE_H5G_OPEN_ERR( "/" );

gsell's avatar
gsell committed
999
	herr = _H5Part_write_attrib (
gsell's avatar
gsell committed
1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248
		group_id,
		attrib_name,
		attrib_type,
		attrib_value,
		attrib_nelem );
	if ( herr < 0 ) return herr;

	herr = H5Gclose ( group_id );
	if ( herr < 0 ) return HANDLE_H5G_CLOSE_ERR;

	return H5PART_SUCCESS;
}


/*!
  Gets the number of attributes bound to the current step.

  \return	Number of attributes bound to current time step or error code.
*/
h5part_int64_t
H5PartGetNumStepAttribs (
	H5PartFile *f			/*!< [in] Handle to open file */
	) {

	SET_FNAME ( "H5PartGetNumStepAttribs" );
	h5part_int64_t nattribs;

	CHECK_FILEHANDLE ( f );

	nattribs = H5Aget_num_attrs(f->timegroup);
	if ( nattribs < 0 ) HANDLE_H5A_GET_NUM_ATTRS_ERR;

	return nattribs;
}

/*!
  Gets the number of attributes bound to the file.

  \return	Number of attributes bound to file \c f or error code.
*/
h5part_int64_t
H5PartGetNumFileAttribs (
	H5PartFile *f			/*!< [in] Handle to open file */
	) {

	SET_FNAME ( "H5PartGetNumFileAttribs" );
	herr_t herr;
	h5part_int64_t nattribs;

	CHECK_FILEHANDLE ( f );

	hid_t group_id = H5Gopen ( f->file, "/" );
	if ( group_id < 0 ) HANDLE_H5G_OPEN_ERR ( "/" );

	nattribs = H5Aget_num_attrs ( group_id );
	if ( nattribs < 0 ) HANDLE_H5A_GET_NUM_ATTRS_ERR;

	herr = H5Gclose ( group_id );
	if ( herr < 0 ) HANDLE_H5G_CLOSE_ERR;
	return nattribs;
}

/*!
  Gets the name, type and number of elements of the step attribute
  specified by its index.

  This function can be used to retrieve all attributes bound to the
  current time-step by looping from \c 0 to the number of attribute
  minus one.  The number of attributes bound to the current
  time-step can be queried by calling the function
  \c H5PartGetNumStepAttribs().

  \return	\c H5PART_SUCCESS or error code 
*/
h5part_int64_t
H5PartGetStepAttribInfo (
	H5PartFile *f,			/*!< [in]  Handle to open file */
	const h5part_int64_t attrib_idx,/*!< [in]  Index of attribute to
					           get infos about */
	char *attrib_name,		/*!< [out] Name of attribute */
	const h5part_int64_t len_of_attrib_name,
					/*!< [in]  length of buffer \c name */
	h5part_int64_t *attrib_type,	/*!< [out] Type of value. */
	h5part_int64_t *attrib_nelem	/*!< [out] Number of elements */
	) {
	
	SET_FNAME ( "H5PartGetStepAttribInfo" );
	hid_t herr;

	CHECK_FILEHANDLE( f );

	herr = _H5Part_get_attrib_info (
		f->timegroup,
		attrib_idx,
		attrib_name,
		len_of_attrib_name,
		attrib_type,
		attrib_nelem );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/*!
  Gets the name, type and number of elements of the file attribute
  specified by its index.

  This function can be used to retrieve all attributes bound to the
  file \c f by looping from \c 0 to the number of attribute minus
  one.  The number of attributes bound to file \c f can be queried
  by calling the function \c H5PartGetNumFileAttribs().

  \return	\c H5PART_SUCCESS or error code 
*/

h5part_int64_t
H5PartGetFileAttribInfo (
	H5PartFile *f,			/*!< [in]  Handle to open file */
	const h5part_int64_t attrib_idx,/*!< [in]  Index of attribute to get
					           infos about */
	char *attrib_name,		/*!< [out] Name of attribute */
	const h5part_int64_t len_of_attrib_name,
					/*!< [in]  length of buffer \c name */
	h5part_int64_t *attrib_type,	/*!< [out] Type of value. */
	h5part_int64_t *attrib_nelem	/*!< [out] Number of elements */
	) {

	SET_FNAME ( "H5PartGetFileAttribInfo" );
	hid_t group_id;
	herr_t herr;

	CHECK_FILEHANDLE( f );

	group_id = H5Gopen(f->file,"/");
	if ( group_id < 0 ) return HANDLE_H5G_OPEN_ERR( "/" );

	herr = _H5Part_get_attrib_info (
		group_id,
		attrib_idx,
		attrib_name,
		len_of_attrib_name,
		attrib_type,
		attrib_nelem );
	if ( herr < 0 ) return herr;

	herr = H5Gclose ( group_id );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/*!
  Reads an attribute bound to current time-step.

  \return \c H5PART_SUCCESS or error code 
*/

h5part_int64_t
H5PartReadStepAttrib (
	H5PartFile *f,			/*!< [in]  Handle to open file */
	char *attrib_name,		/*!< [in] Name of attribute to read */
	void *attrib_value		/*!< [out] Value of attribute */
	) {

	SET_FNAME ( "H5PartReadStepAttrib" );

	hid_t herr;

	CHECK_FILEHANDLE( f );

	herr = _H5Part_read_attrib ( f->timegroup, attrib_name, attrib_value );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/*!
  Reads an attribute bound to file \c f.

  \return \c H5PART_SUCCESS or error code 
*/

h5part_int64_t
H5PartReadFileAttrib ( 
	H5PartFile *f,
	char *attrib_name,
	void *attrib_value
	) {

	SET_FNAME ( "H5PartReadFileAttrib" );

	hid_t group_id;
	hid_t herr;

	CHECK_FILEHANDLE( f );

	group_id = H5Gopen(f->file,"/");
	if ( group_id < 0 ) return HANDLE_H5G_OPEN_ERR( "/" );

	herr = _H5Part_read_attrib ( group_id, attrib_name, attrib_value );
	if ( herr < 0 ) return herr;

	herr = H5Gclose ( group_id );
	if ( herr < 0 ) return HANDLE_H5G_CLOSE_ERR;

	return H5PART_SUCCESS;
}


/*================== File Reading Routines =================*/
/*
  H5PartSetStep:


  So you use this to random-access the file for a particular timestep.
  Failure to explicitly set the timestep on each read will leave you
  stuck on the same timestep for *all* of your reads.  That is to say
  the writes auto-advance the file pointer, but the reads do not
  (they require explicit advancing by selecting a particular timestep).
*/

/*!
  Set the current time-step.

  When writing data to a file the current time step must be set first
  (even if there is only one). In write-mode this function creates a new
  time-step! You are not allowed to step to an already existing time-step.
  This prevents you from overwriting existing data. Another consequence is,
  that you \b must write all data before going to the next time-step.

  In read-mode you can use this function to random-access the file for a
  particular timestep.

  \return \c H5PART_SUCCESS or error code 
*/

h5part_int64_t
H5PartSetStep (
	H5PartFile *f,			/*!< [in]  Handle to open file */
	h5part_int64_t step		/*!< [in]  Time-step to set. */
	) {

	SET_FNAME ( "H5PartSetStep" );
	herr_t r;

	char name[128];

	CHECK_FILEHANDLE ( f );

gsell's avatar
gsell committed
1249
	_H5Part_print_info ( "Proc[%d]: Set step to #%lld for file %d",
gsell's avatar
gsell committed
1250
			     f->myproc, (long long)step, (int) f );
gsell's avatar
gsell committed
1251
	sprintf ( name, "%s#%lld", H5PART_PARTICLES_GROUP, (long long) step );
gsell's avatar
gsell committed
1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268
	r = H5Gget_objinfo( f->file, name, 1, NULL );
	if ( ( (f->mode == H5PART_APPEND) || (f->mode == H5PART_WRITE) )
	     && ( r >= 0 ) ) {
		return HANDLE_H5PART_STEP_EXISTS_ERR ( step );
	}

	if(f->timegroup>0) {
		r = H5Gclose(f->timegroup); if ( r < 0 ) HANDLE_H5G_CLOSE_ERR;
	}
	f->timegroup = -1;
	f->timestep = step;

	if(f->mode==H5PART_READ) {
		f->timegroup = H5Gopen(f->file,name); 
		if ( f->timegroup < 0 ) return HANDLE_H5G_OPEN_ERR( name );
	}
	else {
gsell's avatar
gsell committed
1269 1270
		_H5Part_print_debug ( "Proc[%d]: create step #%lld", 
				      f->myproc, (long long)step );
gsell's avatar
gsell committed
1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405

		f->timegroup = H5Gcreate(f->file,name,0);
		if ( f->timegroup < 0 ) return HANDLE_H5G_CREATE_ERR ( name );
	}

	return H5PART_SUCCESS;
}

/********************** query file structure *********************************/


herr_t
_H5Part_iteration_operator (
	hid_t group_id,		/*!< [in]  group id */
	const char *member_name,/*!< [in]  group name */
	void *operator_data	/*!< [in,out] data passed to the iterator */
	) {

	struct _iter_op_data *data = operator_data;
	herr_t herr;
	H5G_stat_t objinfo;

	herr = H5Gget_objinfo ( group_id, member_name, 1, &objinfo );
	if ( herr < 0 ) return HANDLE_H5G_GET_OBJINFO_ERR ( member_name );

	if ( objinfo.type != data->type )
		return 0;	/* don't count, continue iteration */

	if ( data->name && (data->stop_idx == data->count) ) {
		memset ( data->name, 0, data->len );
		strncpy ( data->name, member_name, data->len-1 );
		
		return 1;	/* stop iteration */
	}
	/*
	  count only if pattern is NULL or member name matches
	*/
	if ( !data->pattern ||
	     (strncmp (member_name, data->pattern, strlen(data->pattern)) == 0)
	      ) {
		data->count++;
	}
	return 0;		/* continue iteration */
}

h5part_int64_t
_H5Part_get_num_objects (
	hid_t group_id,
	const char *group_name,
	const hid_t type
	) {

	return _H5Part_get_num_objects_matching_pattern (
		group_id,
		group_name,
		type,
		NULL );
}

h5part_int64_t
_H5Part_get_num_objects_matching_pattern (
	hid_t group_id,
	const char *group_name,
	const hid_t type,
	char * const pattern
	) {

	h5part_int64_t herr;
	int idx = 0;
	struct _iter_op_data data;

	memset ( &data, 0, sizeof ( data ) );
	data.type = type;
	data.pattern = pattern;

	herr = H5Giterate ( group_id, group_name, &idx,
			    _H5Part_iteration_operator, &data );
	if ( herr < 0 ) return herr;
	
	return data.count;
}

h5part_int64_t
_H5Part_get_object_name (
	hid_t group_id,
	const char *group_name,
	const hid_t type,
	const h5part_int64_t idx,
	char *obj_name,
	const h5part_int64_t len_obj_name
	) {

	herr_t herr;
	struct _iter_op_data data;
	int iterator_idx = 0;

	memset ( &data, 0, sizeof ( data ) );
	data.stop_idx = idx;
	data.type = type;
	data.name = obj_name;
	data.len = len_obj_name;

	herr = H5Giterate ( group_id, group_name, &iterator_idx,
			    _H5Part_iteration_operator,
			    &data );
	if ( herr < 0 ) return (h5part_int64_t)herr;

	if ( herr == 0 ) HANDLE_H5PART_NOENTRY_ERR( group_name,
						    type, idx );

	return H5PART_SUCCESS;
}

/*!
  Get the number of time-steps that are currently stored in the file
  \c f.

  It works for both reading and writing of files, but is probably
  only typically used when you are reading.

  \return	number of time-steps or error code
*/
h5part_int64_t
H5PartGetNumSteps (
	H5PartFile *f			/*!< [in]  Handle to open file */
	) {

	SET_FNAME ( "H5PartGetNumSteps" );

	CHECK_FILEHANDLE( f );

	return _H5Part_get_num_objects_matching_pattern (
		f->file,
		"/",
		H5G_GROUP,
gsell's avatar
gsell committed
1406
		H5PART_PARTICLES_GROUP );
gsell's avatar
gsell committed
1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425
}

/*!
  Get the number of datasets that are stored at the current time-step.

  \return	number of datasets in current timestep or error code
*/

h5part_int64_t
H5PartGetNumDatasets (
	H5PartFile *f			/*!< [in]  Handle to open file */
	) {

	SET_FNAME ( "H5PartGetNumDatasets" );

	char stepname[128];

	CHECK_FILEHANDLE( f );

gsell's avatar
gsell committed
1426 1427
	sprintf ( stepname, "%s#%lld",
		  H5PART_PARTICLES_GROUP, (long long) f->timestep );
gsell's avatar
gsell committed
1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454

	return _H5Part_get_num_objects ( f->file, stepname, H5G_DATASET );
}

/*!
  This reads the name of a dataset specified by it's index in the current
  time-step.

  If the number of datasets is \c n, the range of \c _index is \c 0 to \c n-1.

  \result	\c H5PART_SUCCESS
*/
h5part_int64_t
H5PartGetDatasetName (
	H5PartFile *f,			/*!< [in]  Handle to open file */
	const h5part_int64_t idx,	/*!< [in]  Index of the dataset */
	char *name,			/*!< [out] Name of dataset */
	const h5part_int64_t len_of_name/*!< [in]  Size of buffer \c name */
	) {

	SET_FNAME ( "H5PartGetDatasetName" );

	char stepname[128];

	CHECK_FILEHANDLE ( f );
	CHECK_TIMEGROUP ( f );

gsell's avatar
gsell committed
1455 1456
	sprintf ( stepname, "%s#%lld",
		  H5PART_PARTICLES_GROUP, (long long)f->timestep);
gsell's avatar
gsell committed
1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495

	return _H5Part_get_object_name (
		f->file,
		stepname,
		H5G_DATASET,
		idx,
		name,
		len_of_name );
}

/*!
  Gets the name, type and number of elements of a dataset specified by it's
  index in the current time-step.

  Type is one of \c H5T_NATIVE_DOUBLE or \c H5T_NATIVE_INT64.

  \return	\c H5PART_SUCCESS
*/
h5part_int64_t
H5PartGetDatasetInfo (
	H5PartFile *f,		/*!< [in]  Handle to open file */
	const h5part_int64_t idx,/*!< [in]  Index of the dataset */
	char *dataset_name,	/*!< [out] Name of dataset */
	const h5part_int64_t len_dataset_name,
				/*!< [in]  Size of buffer \c dataset_name */
	h5part_int64_t *type,	/*!< [out] Type of data in dataset */
	h5part_int64_t *nelem	/*!< [out] Number of elements. */
	) {

	SET_FNAME ( "H5PartGetDatasetInfo" );

	herr_t herr;
	hid_t dataset_id;
	hid_t mytype;
	char step_name[128];

	CHECK_FILEHANDLE ( f );
	CHECK_TIMEGROUP ( f );

gsell's avatar
gsell committed
1496 1497
	sprintf ( step_name, "%s#%lld",
		  H5PART_PARTICLES_GROUP, (long long)f->timestep);
gsell's avatar
gsell committed
1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553

	herr = _H5Part_get_object_name (
		f->timegroup,
		step_name,
		H5G_DATASET,
		idx,
		dataset_name,
		len_dataset_name );
	if ( herr < 0 ) return herr;

	*nelem = H5PartGetNumParticles(f);
	if ( *nelem < 0 ) return *nelem;

	dataset_id = H5Dopen ( f->timegroup, dataset_name );
	if ( dataset_id < 0 ) HANDLE_H5D_OPEN_ERR ( dataset_name );

	mytype = H5Dget_type ( dataset_id );
	if ( mytype < 0 ) HANDLE_H5D_GET_TYPE_ERR;

	if(type)
		*type = (h5part_int64_t) _H5Part_normalize_h5_type ( mytype );

	herr = H5Tclose(mytype);
	if ( herr < 0 ) HANDLE_H5T_CLOSE_ERR;

	herr = H5Dclose(dataset_id);
	if ( herr < 0 ) HANDLE_H5D_CLOSE_ERR;

	return H5PART_SUCCESS;
}

/*!
  \intenal 

*/
static hid_t
_get_diskshape_for_reading (
	H5PartFile *f,
	hid_t dataset
	) {

	herr_t r;

	CHECK_FILEHANDLE( f );

	hid_t space = H5Dget_space(dataset);
	if ( space < 0 ) return HANDLE_H5D_GET_SPACE_ERR;

	if ( H5PartHasView(f) ){ 
		hsize_t stride;
		hsize_t count;
#ifdef HDF5V160
		hssize_t start;
#else
		hsize_t start;
#endif
gsell's avatar
gsell committed
1554
		_H5Part_print_debug ( "Selection is available" );
gsell's avatar
gsell committed
1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574

		/* so, is this selection inclusive or exclusive? */
		start = f->viewstart;
		count = f->viewend - f->viewstart; /* to be inclusive */
		stride=1;

		/* now we select a subset */
		if ( f->diskshape > 0 ) {
			r = H5Sselect_hyperslab (
				f->diskshape, H5S_SELECT_SET,
				&start, &stride, &count, NULL);
			if ( r < 0 ) return HANDLE_H5S_SELECT_HYPERSLAB_ERR;
		}
		/* now we select a subset */
		r = H5Sselect_hyperslab (
			space,H5S_SELECT_SET,
			&start, &stride, &count, NULL );
		if ( r < 0 ) return HANDLE_H5S_SELECT_HYPERSLAB_ERR;

		_H5Part_print_debug (
gsell's avatar
gsell committed
1575
			"Selection: range=%d:%d, npoints=%d s=%d",
gsell's avatar
gsell committed
1576 1577 1578 1579
			(int)f->viewstart,(int)f->viewend,
			(int)H5Sget_simple_extent_npoints(space),
			(int)H5Sget_select_npoints(space) );
	} else {
gsell's avatar
gsell committed
1580
		_H5Part_print_debug ( "Selection" );
gsell's avatar
gsell committed
1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638
	}
	return space;
}

static hid_t
_get_memshape_for_reading (
	H5PartFile *f,
	hid_t dataset
	) {

	hid_t r;

	CHECK_FILEHANDLE( f );
 
	if(H5PartHasView(f)) {
		hsize_t dmax=H5S_UNLIMITED;
		hsize_t len = f->viewend - f->viewstart;
		r = H5Screate_simple(1,&len,&dmax);
		if ( r < 0 ) return HANDLE_H5S_CREATE_SIMPLE_ERR ( len );
		return r;
	}
	else {
		return H5S_ALL;
	}
}

/*!
  This gets the number of particles stored in the current timestep. 
  It will arbitrarily select a time-step if you haven't already set
  the timestep with \c H5PartSetStep().

  \return	number of particles in current timestep or an error
		code.
 */
h5part_int64_t
H5PartGetNumParticles (
	H5PartFile *f			/*!< [in]  Handle to open file */
	) {

	SET_FNAME ( "H5PartGetNumParticles" );

	h5part_int64_t herr;
	hid_t space_id;
	hid_t dataset_id;
	char dataset_name[128];
	char step_name[128];
	hsize_t nparticles;

	CHECK_FILEHANDLE( f );

	if(f->timegroup<0) {
		h5part_int64_t step = (f->timestep<0 ? 0 : f->timestep );
		herr = H5PartSetStep ( f, step );
		if ( herr < 0 )
			return HANDLE_H5PART_SETSTEP_ERR ( herr, step );
	}

	/* Get first dataset in current time-step */
gsell's avatar
gsell committed
1639 1640
	sprintf ( step_name, "%s#%lld",
		  H5PART_PARTICLES_GROUP, (long long) f->timestep );
gsell's avatar
gsell committed
1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697
	herr = _H5Part_get_object_name ( f->file, step_name, H5G_DATASET, 0,
					 dataset_name, sizeof (dataset_name) );
	if ( herr < 0 ) return herr;

	dataset_id = H5Dopen ( f->timegroup, dataset_name );
	if ( dataset_id < 0 )
		return HANDLE_H5D_OPEN_ERR ( dataset_name );

	space_id = _get_diskshape_for_reading ( f, dataset_id );
	if ( space_id < 0 ) return (h5part_int64_t)space_id;

	if ( H5PartHasView ( f ) ) {
		nparticles = H5Sget_select_npoints ( space_id );
		if ( nparticles < 0 ) return HANDLE_H5S_GET_SELECT_NPOINTS_ERR;
	}
	else {
		nparticles = H5Sget_simple_extent_npoints ( space_id );
		if ( nparticles < 0 )
			return HANDLE_H5S_GET_SIMPLE_EXTENT_NPOINTS_ERR;
	}
	if ( space_id != H5S_ALL ) {
		herr = H5Sclose ( space_id );
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
	}
	herr = H5Dclose ( dataset_id );
	if ( herr < 0 ) return HANDLE_H5D_CLOSE_ERR;

	return (h5part_int64_t) nparticles;
}


/*!
  For parallel I/O or for subsetting operations on the datafile, the
  \c H5PartSetView() function allows you to define a subset of the total
  particle dataset to read.  The concept of "view" works for both serial
  and for parallel I/O.  The "view" will remain in effect until a new view
  is set, or the number of particles in a dataset changes, or the view is
  "unset" by calling \c H5PartSetView(file,-1,-1);

  Before you set a view, the \c H5PartGetNumParticles() will return the
  total number of particles in the current time-step (even for the parallel
  reads).  However, after you set a view, it will return the number of
  particles contained in the view.

  The range is inclusive (the start and the end index).

  \return	\c H5PART_SUCCESS or error code
*/
h5part_int64_t
H5PartSetView (
	H5PartFile *f,			/*!< [in]  Handle to open file */
	h5part_int64_t start,		/*!< [in]  Start particle */
	h5part_int64_t end		/*!< [in]  End particle */
	) {

	SET_FNAME ( "H5PartSetView" );

gsell's avatar
gsell committed
1698
	h5part_int64_t herr = 0;
gsell's avatar
gsell committed
1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713
	hsize_t total;
	hsize_t stride;
	hsize_t dmax=H5S_UNLIMITED;

#ifdef HDF5V160
	hssize_t range[2];
#else
	hsize_t range[2];
#endif

	CHECK_FILEHANDLE( f );

	if(f->mode==H5PART_WRITE || f->mode==H5PART_APPEND)
		return HANDLE_H5PART_FILE_ACCESS_TYPE_ERR ( f->mode );

gsell's avatar
gsell committed
1714 1715
	_H5Part_print_debug ( "Set view (%lld,%lld).",
			      (long long)start,(long long)end);
gsell's avatar
gsell committed
1716 1717 1718 1719 1720

	/* if there is already a view selected, lets destroy it */ 
	f->viewstart = -1;
	f->viewend = -1;
	if ( f->shape != 0 ){
gsell's avatar
gsell committed
1721 1722
		herr = H5Sclose(f->shape);
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
gsell's avatar
gsell committed
1723 1724 1725
		f->shape=0;
	}
	if(f->diskshape!=0 && f->diskshape!=H5S_ALL){
gsell's avatar
gsell committed
1726 1727
		herr = H5Sclose(f->diskshape);
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
gsell's avatar
gsell committed
1728 1729 1730 1731
		f->diskshape=H5S_ALL;
	}
	f->diskshape = H5S_ALL;
	if(f->memshape!=0 && f->memshape!=H5S_ALL){
gsell's avatar
gsell committed
1732 1733
		herr = H5Sclose ( f->memshape );
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
gsell's avatar
gsell committed
1734 1735 1736
		f->memshape=H5S_ALL;
	}
	if(start==-1 && end==-1) {
gsell's avatar
gsell committed
1737
		_H5Part_print_debug( "Early Termination: Unsetting View" );
gsell's avatar
gsell committed
1738 1739 1740 1741 1742 1743 1744
		return H5PART_SUCCESS; /* all done */
	}
	/* for now, we interpret start=-1 to mean 0 and 
	   end==-1 to mean end of file */
	total = (hsize_t) H5PartGetNumParticles(f);
	if ( total < 0 ) return HANDLE_H5PART_GET_NUM_PARTICLES_ERR ( total );

gsell's avatar
gsell committed
1745
	_H5Part_print_debug ( "Total nparticles=%lld", (long long)total );
gsell's avatar
gsell committed
1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756
	if(start==-1) start=0;
	if(end==-1) end=total; /* can we trust nparticles (no)? 
				  fortunately, view has been reset
				  so H5PartGetNumParticles will tell
				  us the total number of particles. */

	/* so, is this selection inclusive or exclusive? 
	   it appears to be inclusive for both ends of the range.
	*/
	if(end<start) {
		_H5Part_print_warn (
gsell's avatar
gsell committed
1757 1758 1759
			"Nonfatal error. "
			"End of view (%lld) is less than start (%lld).",
			(long long)end, (long long)start );
gsell's avatar
gsell committed
1760 1761 1762 1763 1764 1765 1766 1767
		end = start; /* ensure that we don't have a range error */
	}
	range[0]=start;
	range[1]=end;
	/* setting up the new view */
	f->viewstart=range[0]; /* inclusive start */
	f->viewend=range[1]; /* inclusive end */
	f->nparticles=range[1]-range[0];
gsell's avatar
gsell committed
1768 1769
	_H5Part_print_debug ( "Range is now %lld:%lld",
			      (long long)range[0], (long long)range[1]);
gsell's avatar
gsell committed
1770 1771 1772 1773 1774 1775 1776 1777 1778 1779 1780 1781 1782 1783 1784 1785 1786 1787
	/* OK, now we must create a selection from this */
	
	/* declare overall datasize */
	f->shape = H5Screate_simple ( 1, &total, &total );
	if ( f->shape < 0 )
		return HANDLE_H5S_CREATE_SIMPLE_ERR ( total );

	/* declare overall data size  but then will select a subset */
	f->diskshape= H5Screate_simple ( 1, &total, &total );
	if ( f->diskshape < 0 )
		return HANDLE_H5S_CREATE_SIMPLE_ERR ( total );

	/* declare local memory datasize */
	f->memshape = H5Screate_simple(1,&(f->nparticles),&dmax);
	if ( f->memshape < 0 )
		return HANDLE_H5S_CREATE_SIMPLE_ERR ( f->nparticles );;

	stride=1;
gsell's avatar
gsell committed
1788
	herr = H5Sselect_hyperslab ( 
gsell's avatar
gsell committed
1789 1790 1791 1792 1793 1794
		f->diskshape,
		H5S_SELECT_SET,
		range,
		&stride,
		&total,
		NULL );
gsell's avatar
gsell committed
1795
	if ( herr < 0 ) return HANDLE_H5S_SELECT_HYPERSLAB_ERR;
gsell's avatar
gsell committed
1796 1797 1798 1799 1800 1801 1802 1803 1804 1805 1806 1807 1808 1809 1810 1811 1812 1813 1814 1815 1816 1817 1818 1819 1820 1821 1822 1823 1824 1825 1826 1827 1828 1829 1830 1831 1832 1833 1834 1835 1836 1837 1838 1839 1840 1841 1842 1843 1844 1845 1846 1847 1848 1849 1850 1851 1852 1853 1854 1855 1856 1857 1858 1859 1860 1861 1862 1863 1864 1865 1866 1867 1868 1869 1870 1871 1872 1873 1874 1875 1876 1877 1878 1879 1880 1881 1882 1883 1884 1885 1886 1887 1888 1889 1890 1891 1892 1893 1894 1895 1896 1897 1898 1899 1900 1901 1902 1903 1904 1905 1906 1907 1908 1909 1910 1911 1912 1913 1914 1915 1916 1917 1918 1919 1920 1921 1922 1923 1924 1925 1926 1927 1928 1929 1930 1931 1932 1933 1934 1935 1936 1937 1938 1939 1940 1941 1942 1943 1944 1945 1946 1947 1948 1949 1950 1951 1952 1953 1954 1955 1956 1957 1958 1959 1960 1961 1962 1963 1964 1965 1966 1967 1968 1969 1970 1971 1972 1973 1974 1975 1976 1977 1978 1979 1980 1981 1982 1983 1984 1985 1986 1987 1988 1989 1990 1991 1992 1993 1994 1995 1996 1997 1998 1999 2000 2001 2002 2003 2004 2005 2006 2007 2008 2009 2010 2011 2012 2013 2014 2015 2016 2017 2018 2019 2020 2021 2022 2023 2024 2025 2026 2027 2028 2029 2030 2031 2032 2033 2034 2035 2036 2037 2038 2039 2040 2041 2042 2043 2044 2045 2046 2047 2048 2049 2050 2051 2052 2053 2054 2055 2056 2057 2058 2059 2060 2061 2062 2063 2064 2065 2066 2067 2068 2069 2070 2071 2072 2073 2074 2075 2076 2077 2078 2079 2080 2081 2082 2083 2084 2085 2086 2087 2088 2089 2090 2091 2092 2093 2094


	return H5PART_SUCCESS;
}

/*!
   Allows you to query the current view. Start and End
   will be \c -1 if there is no current view established.
   Use \c H5PartHasView() to see if the view is smaller than the
   total dataset.

   \return       the number of elements in the view 
*/
h5part_int64_t
H5PartGetView (
	H5PartFile *f,			/*!< [in]  Handle to open file */
	h5part_int64_t *start,		/*!< [out]  Start particle */
	h5part_int64_t *end		/*!< [out]  End particle */
	) {

	SET_FNAME ( "H5PartGetView" );

	h5part_int64_t range[2];
	h5part_int64_t viewend;

	CHECK_FILEHANDLE( f );

	range[0] = (f->viewstart>=0) ? f->viewstart : 0;

	if ( f->viewend >= 0 ) {
		viewend = f->viewend;
	}
	else {
		viewend = H5PartGetNumParticles(f);
		if ( viewend < 0 )
			return HANDLE_H5PART_GET_NUM_PARTICLES_ERR ( viewend );
	}

	range[1] = viewend;
	if(start) {
		*start=range[0];
	}
	if(end) {
		*end=range[1];
	}

	return range[1]-range[0];
}

/*!
  If it is too tedious to manually set the start and end coordinates
  for a view, the \c H5SetCanonicalView() will automatically select an
  appropriate domain decomposition of the data arrays for the degree
  of parallelism and set the "view" accordingly.

  \return		H5PART_SUCCESS or error code
*/
h5part_int64_t
H5PartSetCanonicalView (
	H5PartFile *f			/*!< [in]  Handle to open file */
	) {

	SET_FNAME ( "H5PartSetCanonicalView" );

	h5part_int64_t r;

	CHECK_FILEHANDLE( f );

	if(f->mode==H5PART_WRITE || f->mode==H5PART_APPEND)
		return HANDLE_H5PART_FILE_ACCESS_TYPE_ERR ( f->mode );

	/* if a read_only file, search for on-disk canonical view */
	/* if this view does not exist, then if MPI, subdivide by numprocs */
	/* else, "unset" any existing View */

	/* unset the view */
	r = H5PartSetView(f,-1,-1);
	if ( r < 0 ) return HANDLE_H5PART_SET_VIEW_ERR( r, -1, -1 );
#ifdef PARALLEL_IO
	{
		h5part_int64_t start = 0;
		h5part_int64_t end = 0;
		h5part_int64_t n = 0;
		int i = 0;
		
		if ( f->timegroup < 0 ) {
			/* set to first step in file */
			r = H5PartSetStep(f,0);
			if ( r < 0 ) return HANDLE_H5PART_SETSTEP_ERR ( r, 0 );
		}
		n = H5PartGetNumParticles(f);
		if ( n < 0 ) return HANDLE_H5PART_GET_NUM_PARTICLES_ERR ( n );
		/* 
		   now lets query the attributes for this group to see if there
		   is a 'pnparticles' group that contains the offsets for the
		   processors.
		*/
		/* try to read pnparticles right off of the disk */
		if ( _H5Part_read_attrib (
			     f->timegroup,
			     "pnparticles", f->pnparticles ) < 0) {
			/*
			  automatically subdivide the view into NP mostly
			  equal pieces
			*/

			n /= f->nprocs;
			for ( i=0; i<f->nprocs; i++ ) {
				f->pnparticles[i] = n;
			}
		}

		/* now we set the view for this processor */
		for ( i = 0; i < f->myproc; i++ ){
			start += f->pnparticles[i];
		}
		end = start + f->pnparticles[f->myproc] - 1;
		r = H5PartSetView ( f, start, end );
		if ( r < 0 )
			return HANDLE_H5PART_SET_VIEW_ERR ( r, start, end );
	}
#endif
	/* the canonical view is to see everything if this is serial
	   so there is nothing left to do */
	return H5PART_SUCCESS;
}

static h5part_int64_t
_read_data (
	H5PartFile *f,		/*!< [in] Handle to open file */
	char *name,		/*!< [in] Name to associate dataset with */
	void *array,		/*!< [out] Array of data */
	hid_t type
	) {

	herr_t herr;
	hid_t dataset_id;
	hid_t space_id;
	hid_t memspace_id;

	if ( ! f->timegroup ) {
		herr = H5PartSetStep ( f, f->timestep );
		if ( herr < 0 )
			return HANDLE_H5PART_SETSTEP_ERR ( herr, f->timestep );
	}
	dataset_id = H5Dopen ( f->timegroup, name );
	if ( dataset_id < 0 ) return HANDLE_H5D_OPEN_ERR ( name );

	space_id = _get_diskshape_for_reading ( f, dataset_id );
	if ( space_id < 0 ) return (h5part_int64_t)space_id;

	memspace_id = _get_memshape_for_reading ( f, dataset_id );
	if ( memspace_id < 0 ) return (h5part_int64_t)memspace_id;

	herr = H5Dread (
		dataset_id,
		type,
		memspace_id,		/* shape/size of data in memory (the
					   complement to disk hyperslab) */
		space_id,		/* shape/size of data on disk 
					   (get hyperslab if needed) */
		H5P_DEFAULT,		/* ignore... its for parallel reads */
		array );
	if ( herr < 0 ) return HANDLE_H5D_READ_ERR ( name, f->timestep );

	if ( space_id != H5S_ALL ) {
		herr = H5Sclose (space_id );
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
	}

	if ( memspace_id != H5S_ALL )
		herr = H5Sclose ( memspace_id );
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;

	herr = H5Dclose ( dataset_id );
	if ( herr < 0 ) return HANDLE_H5D_CLOSE_ERR;
	
	return H5PART_SUCCESS;
}

/*!
  Read array of 64 bit floating point data from file.

  When retrieving datasets from disk, you ask for them
  by name. There are no restrictions on naming of arrays,
  but it is useful to arrive at some common naming
  convention when sharing data with other groups.

  \return	\c H5PART_SUCCESS or error code
*/
h5part_int64_t
H5PartReadDataFloat64 (
	H5PartFile *f,		/*!< [in] Handle to open file */
	char *name,		/*!< [in] Name to associate dataset with */
	h5part_float64_t *array	/*!< [out] Array of data */
	) {

	SET_FNAME ( "H5PartReadDataFloat64" );

	h5part_int64_t herr;

	CHECK_FILEHANDLE( f );

	herr = _read_data ( f, name, array, H5T_NATIVE_DOUBLE );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/**
  Read array of 64 bit floating point data from file.

  When retrieving datasets from disk, you ask for them
  by name. There are no restrictions on naming of arrays,
  but it is useful to arrive at some common naming
  convention when sharing data with other groups.

  \return	\c H5PART_SUCCESS or error code
*/
h5part_int64_t
H5PartReadDataInt64 (
	H5PartFile *f,		/*!< [in] Handle to open file */
	char *name,		/*!< [in] Name to associate dataset with */
	h5part_int64_t *array	/*!< [out] Array of data */
	) {

	SET_FNAME ( "H5PartReadDataInt64" );

	h5part_int64_t herr;

	CHECK_FILEHANDLE( f );

	herr = _read_data ( f, name, array, H5T_NATIVE_INT64 );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/*!
  This is the mongo read function that pulls in all of the data for a
  given timestep in one shot. It also takes the timestep as an argument
  and will call \c H5PartSetStep() internally so that you don't have to 
  make that call separately.

  \note
  See also \c H5PartReadDataInt64() and \c H5PartReadDataFloat64() if you want
  to just read in one of the many datasets.

  \return	\c H5PART_SUCCESS or error code
*/
h5part_int64_t
H5PartReadParticleStep (
	H5PartFile *f,		/*!< [in]  Handle to open file */
	h5part_int64_t step,	/*!< [in]  Step to read */
	h5part_float64_t *x,	/*!< [out] Buffer for dataset named "x" */
	h5part_float64_t *y,	/*!< [out] Buffer for dataset named "y" */
	h5part_float64_t *z,	/*!< [out] Buffer for dataset named "z" */
	h5part_float64_t *px,	/*!< [out] Buffer for dataset named "px" */
	h5part_float64_t *py,	/*!< [out] Buffer for dataset named "py" */
	h5part_float64_t *pz,	/*!< [out] Buffer for dataset named "pz" */
	h5part_int64_t *id	/*!< [out] Buffer for dataset named "id" */
	) {

	SET_FNAME ( "H5PartReadParticleStep" );
	h5part_int64_t herr;

	CHECK_FILEHANDLE( f );

	herr = H5PartSetStep(f,step);
	if ( herr < 0 ) return HANDLE_H5PART_SETSTEP_ERR ( herr, step );

	herr = _read_data ( f, "x", (void*)x, H5T_NATIVE_DOUBLE );
	if ( herr < 0 ) return herr;

	herr = _read_data ( f, "y", (void*)y, H5T_NATIVE_DOUBLE );
	if ( herr < 0 ) return herr;

	herr = _read_data ( f, "z", (void*)z, H5T_NATIVE_DOUBLE );
	if ( herr < 0 ) return herr;

	herr = _read_data ( f, "px", (void*)px, H5T_NATIVE_DOUBLE );
	if ( herr < 0 ) return herr;

	herr = _read_data ( f, "py", (void*)py, H5T_NATIVE_DOUBLE );
	if ( herr < 0 ) return herr;

	herr = _read_data ( f, "pz", (void*)pz, H5T_NATIVE_DOUBLE );
	if ( herr < 0 ) return herr;

	herr = _read_data ( f, "id", (void*)id, H5T_NATIVE_INT64 );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/****************** error handling ******************/

h5part_int64_t
H5PartSetVerbosityLevel (
gsell's avatar
gsell committed
2095
	h5part_int64_t level
gsell's avatar
gsell committed
2096 2097 2098 2099 2100 2101 2102 2103 2104 2105 2106 2107 2108 2109 2110 2111 2112 2113 2114 2115 2116 2117 2118 2119 2120 2121 2122 2123 2124 2125 2126 2127 2128 2129 2130 2131 2132 2133 2134 2135
	) {

	_debug = level;
	return H5PART_SUCCESS;
}

h5part_int64_t
H5PartSetErrorHandler (
	h5part_error_handler handler
	) {
	_err_handler = handler;
	return H5PART_SUCCESS;
}

h5part_error_handler
H5PartGetErrorHandler (
	void
	) {
	return _err_handler;
}

h5part_int64_t
H5PartGetErrno (
	void
	) {
	return _errno;
}

h5part_int64_t
H5PartDefaultErrorHandler (
	const char *funcname,
	const h5part_int64_t eno,
	const char *fmt,
	...
	) {

	_errno = eno;
	if ( _debug > 0 ) {
		va_list ap;
		va_start ( ap, fmt );
gsell's avatar
gsell committed
2136 2137
		_H5Part_vprint_error ( fmt, ap );
		va_end ( ap );
gsell's avatar
gsell committed
2138 2139 2140 2141 2142 2143 2144 2145 2146 2147 2148 2149 2150 2151 2152 2153 2154 2155 2156 2157 2158 2159 2160 2161 2162 2163 2164 2165 2166 2167 2168 2169 2170 2171 2172 2173 2174 2175 2176 2177 2178 2179 2180 2181 2182 2183
	}
	return _errno;
}

h5part_int64_t
H5PartAbortErrorHandler (
	const char *funcname,
	const h5part_int64_t eno,
	const char *fmt,
	...
	) {

	_errno = eno;
	if ( _debug > 0 ) {
		va_list ap;
		va_start ( ap, fmt );
		fprintf ( stderr, "%s: ", funcname );
		vfprintf ( stderr, fmt, ap );
		fprintf ( stderr, "\n" );
	}
	exit (-_errno);
}


static h5part_int64_t
_init ( void ) {
	static int __init = 0;

	herr_t r5;
	if ( ! __init ) {
		r5 = H5Eset_auto ( _h5_error_handler, NULL );
		if ( r5 < 0 ) return H5PART_ERR_INIT;
	}
	__init = 1;
	return H5PART_SUCCESS;
}

static herr_t
_h5_error_handler ( void* unused ) {
	
	if ( _debug >= 5 ) {
		H5Eprint (stderr);
	}
	return 0;
}

gsell's avatar
gsell committed
2184 2185 2186 2187 2188 2189 2190 2191 2192 2193 2194 2195 2196 2197
static void
_vprint (
	FILE* f,
	const char *prefix,
	const char *fmt,
	va_list ap
	) {
	char *fmt2 = malloc( strlen ( prefix ) +strlen ( fmt ) + strlen ( __funcname ) + 16 );
	if ( fmt2 == NULL ) return;
	sprintf ( fmt2, "%s: %s: %s\n", prefix, __funcname, fmt ); 
	vfprintf ( stderr, fmt2, ap );
	free ( fmt2 );
}

gsell's avatar
gsell committed
2198 2199 2200 2201 2202 2203 2204
void
_H5Part_vprint_error (
	const char *fmt,
	va_list ap
	) {

	if ( _debug < 1 ) return;
gsell's avatar
gsell committed
2205
	_vprint ( stderr, "E", fmt, ap );
gsell's avatar
gsell committed
2206 2207 2208 2209 2210 2211 2212 2213 2214 2215 2216 2217 2218 2219 2220 2221 2222 2223 2224 2225 2226
}

void
_H5Part_print_error (
	const char *fmt,
	...
	) {

	va_list ap;
	va_start ( ap, fmt );
	_H5Part_vprint_error ( fmt, ap );
	va_end ( ap );
}

void
_H5Part_vprint_warn (
	const char *fmt,
	va_list ap
	) {

	if ( _debug < 2 ) return;
gsell's avatar
gsell committed
2227
	_vprint ( stderr, "W", fmt, ap );
gsell's avatar
gsell committed
2228 2229 2230 2231 2232 2233 2234 2235 2236 2237 2238 2239 2240 2241 2242 2243 2244 2245 2246 2247 2248
}

void
_H5Part_print_warn (
	const char *fmt,
	...
	) {

	va_list ap;
	va_start ( ap, fmt );
	_H5Part_vprint_warn ( fmt, ap );
	va_end ( ap );
}

void
_H5Part_vprint_info (
	const char *fmt,
	va_list ap
	) {

	if ( _debug < 3 ) return;
gsell's avatar
gsell committed
2249
	_vprint ( stdout, "I", fmt, ap );
gsell's avatar
gsell committed
2250 2251 2252 2253 2254 2255 2256 2257 2258 2259 2260 2261 2262 2263 2264 2265 2266 2267 2268 2269 2270
}

void
_H5Part_print_info (
	const char *fmt,
	...
	) {

	va_list ap;
	va_start ( ap, fmt );
	_H5Part_vprint_info ( fmt, ap );
	va_end ( ap );
}

void
_H5Part_vprint_debug (
	const char *fmt,
	va_list ap
	) {

	if ( _debug < 4 ) return;
gsell's avatar
gsell committed
2271
	_vprint ( stdout, "D", fmt, ap );
gsell's avatar
gsell committed
2272 2273 2274 2275 2276 2277 2278 2279 2280 2281 2282 2283 2284 2285 2286 2287 2288 2289 2290 2291 2292 2293 2294 2295 2296 2297 2298
}

void
_H5Part_print_debug (
	const char *fmt,
	...
	) {

	va_list ap;
	va_start ( ap, fmt );
	_H5Part_vprint_debug ( fmt, ap );
	va_end ( ap );
}

void
_H5Part_set_funcname (
	char  * const fname
	) {
	__funcname = fname;
}

const char *
_H5Part_get_funcname (
	void
	) {
	return __funcname;
}