H5Part.c 59.9 KB
Newer Older
gsell's avatar
gsell committed
1 2 3 4 5 6 7 8 9
/*! \mainpage H5Part: A Portable High Performance Parallel Data Interface to HDF5

Particle based simulations of accelerator beam-lines, especially in
six dimensional phase space, generate vast amounts of data. Even
though a subset of statistical information regarding phase space or
analysis needs to be preserved, reading and writing such enormous
restart files on massively parallel supercomputing systems remains
challenging. 

Kurt Stockinger's avatar
Kurt Stockinger committed
10
H5Part consists of Particles and Block structured Fields.
gsell's avatar
gsell committed
11

Kurt Stockinger's avatar
Kurt Stockinger committed
12
Developed by:
gsell's avatar
gsell committed
13

Kurt Stockinger's avatar
Kurt Stockinger committed
14 15 16 17 18 19 20 21 22
<UL>
<LI> Andreas Adelmann (PSI) </LI>
<LI> Achim Gsell (PSI) </LI>
<LI> Benedikt Oswald (PSI) </LI>

<LI> Wes Bethel (NERSC/LBNL)</LI>
<LI> John Shalf (NERSC/LBNL)</LI>
<LI> Cristina Siegerist (NERSC/LBNL)</LI>
</UL>
gsell's avatar
gsell committed
23 24 25 26


Papers: 

Kurt Stockinger's avatar
Kurt Stockinger committed
27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
<UL>
<LI> A. Adelmann, R.D. Ryne, C. Siegerist, J. Shalf,"From Visualization to Data Mining with Large Data Sets," <i>
<a href="http://www.sns.gov/pac05">Particle Accelerator Conference (PAC05)</a></i>, Knoxville TN., May 16-20, 2005. (LBNL-57603)
<a href="http://vis.lbl.gov/Publications/2005/FPAT082.pdf">FPAT082.pdf</a>
</LI>


<LI> A. Adelmann, R.D. Ryne, J. Shalf, C. Siegerist,"H5Part: A Portable High Performance Parallel Data Interface for Particle Simulations," <i>
<a href="http://www.sns.gov/pac05">Particle Accelerator Conference (PAC05)</a></i>, Knoxville TN., May 16-20, 2005.
<a href="http://vis.lbl.gov/Publications/2005/FPAT083.pdf">FPAT083.pdf</a>
</LI>
</UL>

For further information contact: <a href="mailto:h5part@lists.psi.ch">h5part</a>

Last modified on April 19, 2007.

gsell's avatar
gsell committed
44 45
*/

gsell's avatar
gsell committed
46 47

/*!
gsell's avatar
gsell committed
48
  \defgroup h5part_c_api H5Part C API
gsell's avatar
gsell committed
49 50 51

*/
/*!
gsell's avatar
gsell committed
52
  \ingroup h5part_c_api
gsell's avatar
gsell committed
53 54 55
  \defgroup h5part_openclose	File Opening and Closing
*/
/*!
gsell's avatar
gsell committed
56
  \ingroup h5part_c_api
gsell's avatar
gsell committed
57 58 59
  \defgroup h5part_write	File Writing
*/  
/*!
gsell's avatar
gsell committed
60
  \ingroup h5part_c_api
gsell's avatar
gsell committed
61 62 63
  \defgroup h5part_read		File Reading
*/  
/*!
gsell's avatar
gsell committed
64
  \ingroup h5part_c_api
gsell's avatar
gsell committed
65 66 67
  \defgroup h5part_attrib	Reading and Writing Attributes
*/
/*!
gsell's avatar
gsell committed
68
  \ingroup h5part_c_api
gsell's avatar
gsell committed
69 70 71 72 73 74 75 76
  \defgroup h5part_errhandle	Error Handling
*/
/*!
  \internal
  \defgroup h5partkernel H5Part private functions 
*/


gsell's avatar
gsell committed
77 78 79 80 81 82 83 84
#include <stdio.h>
#include <stdlib.h>
#include <stdarg.h>	/* va_arg - System dependent ?! */
#include <string.h>
#include <errno.h>
#include <fcntl.h>
#include <hdf5.h>

gsell's avatar
gsell committed
85 86 87 88 89 90 91 92
#ifndef WIN32
#include <unistd.h>
#else /* WIN32 */
#include <io.h>
#define open  _open
#define close _close
#endif /* WIN32 */

gsell's avatar
gsell committed
93 94 95 96 97 98 99 100
#include "H5PartTypes.h"
#include "H5Part.h"
#include "H5PartPrivate.h"
#include "H5PartErrors.h"

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

static unsigned			_debug = 0;
gsell's avatar
gsell committed
101
static h5part_int64_t		_h5part_errno = H5PART_SUCCESS;
gsell's avatar
gsell committed
102
static h5part_error_handler	_err_handler = H5PartReportErrorHandler;
gsell's avatar
gsell committed
103 104 105 106 107 108 109 110 111 112 113 114 115 116
static char *__funcname = "NONE";

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

static h5part_int64_t
_init(
	void
	);

static h5part_int64_t
_file_is_valid (
	const H5PartFile *f
	);

gsell's avatar
gsell committed
117 118 119
/*
  error handler for hdf5
*/
gsell's avatar
gsell committed
120 121 122 123 124 125 126
static herr_t
_h5_error_handler (
	void *
	);

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

gsell's avatar
gsell committed
127 128
static H5PartFile*
_H5Part_open_file (
gsell's avatar
gsell committed
129
	const char *filename,	/*!< [in] The name of the data file to open. */
gsell's avatar
gsell committed
130 131 132 133
	unsigned flags,		/*!< [in] The access mode for the file. */
	MPI_Comm comm,		/*!< [in] MPI communicator */
	int f_parallel		/*!< [in] 0 for serial io otherwise parallel */
	) {
gsell's avatar
gsell committed
134 135 136 137
	if ( _init() < 0 ) {
		HANDLE_H5PART_INIT_ERR;
		return NULL;
	}
gsell's avatar
gsell committed
138
	_h5part_errno = H5PART_SUCCESS;
gsell's avatar
gsell committed
139 140 141 142 143 144 145 146
	H5PartFile *f = NULL;

	f = (H5PartFile*) malloc( sizeof (H5PartFile) );
	if( f == NULL ) {
		HANDLE_H5PART_NOMEM_ERR;
		goto error_cleanup;
	}
	memset (f, 0, sizeof (H5PartFile));
gsell's avatar
gsell committed
147 148 149 150 151 152 153 154

	f->groupname_step = strdup ( H5PART_GROUPNAME_STEP );
	if( f->groupname_step == NULL ) {
		HANDLE_H5PART_NOMEM_ERR;
		goto error_cleanup;
	}
	f->stepno_width = 0;

gsell's avatar
gsell committed
155 156
	f->xfer_prop = f->create_prop = f->access_prop = H5P_DEFAULT;

gsell's avatar
gsell committed
157
	if ( f_parallel ) {
gsell's avatar
gsell committed
158
#ifdef PARALLEL_IO
gsell's avatar
gsell committed
159 160
		/* for the SP2... perhaps different for linux */
		MPI_Info info = MPI_INFO_NULL;
gsell's avatar
gsell committed
161

162 163 164 165
		/* ks: IBM_large_block_io */
		MPI_Info_create(&info);
		MPI_Info_set(info, "IBM_largeblock_io", "true" );

gsell's avatar
gsell committed
166 167 168 169 170 171 172 173
		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;
		}
gsell's avatar
gsell committed
174

gsell's avatar
gsell committed
175 176
		f->pnparticles =
		  (h5part_int64_t*) malloc (f->nprocs * sizeof (h5part_int64_t));
gsell's avatar
gsell committed
177 178 179 180
		if (f->pnparticles == NULL) {
			HANDLE_H5PART_NOMEM_ERR;
			goto error_cleanup;
		}
gsell's avatar
gsell committed
181
		
gsell's avatar
gsell committed
182 183 184 185 186
		f->access_prop = H5Pcreate (H5P_FILE_ACCESS);
		if (f->access_prop < 0) {
			HANDLE_H5P_CREATE_ERR;
			goto error_cleanup;
		}
gsell's avatar
gsell committed
187

gsell's avatar
gsell committed
188 189 190 191
		if (H5Pset_fapl_mpio (f->access_prop, comm, info) < 0) {
			HANDLE_H5P_SET_FAPL_MPIO_ERR;
			goto error_cleanup;
		}
gsell's avatar
gsell committed
192
		
gsell's avatar
gsell committed
193 194
		/* f->create_prop = H5Pcreate(H5P_FILE_CREATE); */
		f->create_prop = H5P_DEFAULT;
gsell's avatar
gsell committed
195

gsell's avatar
gsell committed
196 197 198 199 200 201 202 203 204 205 206 207 208
		/* 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;
		}
gsell's avatar
gsell committed
209

gsell's avatar
gsell committed
210
		f->comm = comm;
211 212

		MPI_Info_free(&info);
gsell's avatar
gsell committed
213
#endif
gsell's avatar
gsell committed
214
	} else {
gsell's avatar
gsell committed
215
		f->comm = 0;
gsell's avatar
gsell committed
216 217
		f->nprocs = 1;
		f->myproc = 0;
gsell's avatar
gsell committed
218 219
		f->pnparticles = 
			(h5part_int64_t*) malloc (f->nprocs * sizeof (h5part_int64_t));
gsell's avatar
gsell committed
220
	}
gsell's avatar
gsell committed
221 222 223 224 225 226
	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);
gsell's avatar
gsell committed
227
		f->empty = 1;
gsell's avatar
gsell committed
228 229 230 231 232 233
	}
	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);
gsell's avatar
gsell committed
234
			f->empty = 1;
gsell's avatar
gsell committed
235 236 237 238 239 240 241 242 243 244
		}
		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
245
				f->file, "/", H5G_GROUP, f->groupname_step );
gsell's avatar
gsell committed
246 247 248 249 250 251 252 253 254 255 256 257 258
			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;
gsell's avatar
gsell committed
259
	f->timegroup = -1;
gsell's avatar
gsell committed
260 261 262 263 264 265
	f->shape = 0;
	f->diskshape = H5S_ALL;
	f->memshape = H5S_ALL;
	f->viewstart = -1;
	f->viewend = -1;

gsell's avatar
gsell committed
266 267 268 269 270 271
	_H5Part_print_debug (
		"Proc[%d]: Opened file \"%s\" val=%lld",
		f->myproc,
		filename,
		(long long)(size_t)f );

gsell's avatar
gsell committed
272 273 274 275
	return f;

 error_cleanup:
	if (f != NULL ) {
gsell's avatar
gsell committed
276 277 278
		if (f->groupname_step) {
			free (f->groupname_step);
		}
gsell's avatar
gsell committed
279 280 281 282 283 284 285 286
		if (f->pnparticles != NULL) {
			free (f->pnparticles);
		}
		free (f);
	}
	return NULL;
}

gsell's avatar
gsell committed
287
#ifdef PARALLEL_IO
gsell's avatar
gsell committed
288
/*!
gsell's avatar
gsell committed
289
  \ingroup h5part_openclose
gsell's avatar
gsell committed
290

gsell's avatar
gsell committed
291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308
  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*
gsell's avatar
gsell committed
309
H5PartOpenFileParallel (
gsell's avatar
gsell committed
310
	const char *filename,	/*!< [in] The name of the data file to open. */
gsell's avatar
gsell committed
311 312 313 314
	unsigned flags,		/*!< [in] The access mode for the file. */
	MPI_Comm comm		/*!< [in] MPI communicator */
) {
	int f_parallel = 1;	/* parallel i/o */
gsell's avatar
gsell committed
315

gsell's avatar
gsell committed
316 317 318
	return _H5Part_open_file ( filename, flags, comm, f_parallel );
}
#endif
gsell's avatar
gsell committed
319

gsell's avatar
gsell committed
320 321
/*!
  \ingroup  h5part_openclose
gsell's avatar
gsell committed
322

gsell's avatar
gsell committed
323
  Opens file with specified filename. 
gsell's avatar
gsell committed
324

gsell's avatar
gsell committed
325 326 327 328 329
  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.
gsell's avatar
gsell committed
330

gsell's avatar
gsell committed
331 332 333 334 335 336
  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.
gsell's avatar
gsell committed
337

gsell's avatar
gsell committed
338 339
  \return	File handle or \c NULL
 */
gsell's avatar
gsell committed
340

gsell's avatar
gsell committed
341 342 343 344 345
H5PartFile*
H5PartOpenFile (
	const char *filename,	/*!< [in] The name of the data file to open. */
	unsigned flags		/*!< [in] The access mode for the file. */
	) {
gsell's avatar
gsell committed
346

gsell's avatar
gsell committed
347
	SET_FNAME ( "H5PartOpenFile" );
gsell's avatar
gsell committed
348

gsell's avatar
gsell committed
349 350 351 352
	MPI_Comm comm = 0;	/* dummy */
	int f_parallel = 0;	/* serial open */

	return _H5Part_open_file ( filename, flags, comm, f_parallel );
gsell's avatar
gsell committed
353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373
}

/*!
  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;
}

/*!
gsell's avatar
gsell committed
374 375
  \ingroup h5part_openclose

gsell's avatar
gsell committed
376 377 378 379 380 381 382 383 384 385 386
  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;
gsell's avatar
gsell committed
387
	_h5part_errno = H5PART_SUCCESS;
gsell's avatar
gsell committed
388 389 390 391 392 393 394 395 396 397 398 399 400 401

	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;
	}
gsell's avatar
gsell committed
402
	if( f->timegroup >= 0 ) {
gsell's avatar
gsell committed
403 404
		r = H5Gclose( f->timegroup );
		if ( r < 0 ) HANDLE_H5G_CLOSE_ERR;
gsell's avatar
gsell committed
405
		f->timegroup = -1;
gsell's avatar
gsell committed
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
	}
	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;
	}
gsell's avatar
gsell committed
432 433 434
	if (f->groupname_step) {
		free (f->groupname_step);
	}
gsell's avatar
gsell committed
435 436 437 438 439
	if( f->pnparticles ) {
		free( f->pnparticles );
	}
	free( f );

gsell's avatar
gsell committed
440
	return _h5part_errno;
gsell's avatar
gsell committed
441 442 443
}

/*============== File Writing Functions ==================== */
gsell's avatar
gsell committed
444

gsell's avatar
gsell committed
445 446 447 448 449 450 451 452 453 454 455 456 457 458 459
h5part_int64_t
H5PartDefineStepName (
	H5PartFile *f,
	const char *name,
	const h5part_int64_t width
	) {
	f->groupname_step = strdup ( name );
	if( f->groupname_step == NULL ) {
		return HANDLE_H5PART_NOMEM_ERR;
	}
	f->stepno_width = width;
	
	return H5PART_SUCCESS;
}

gsell's avatar
gsell committed
460
/*!
gsell's avatar
gsell committed
461 462
  \ingroup h5part_write

gsell's avatar
gsell committed
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 492 493 494 495 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
  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
	*/
545

gsell's avatar
gsell committed
546 547 548 549 550 551 552 553
	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
554
		_H5Part_print_debug ( "Particle offsets:" );
gsell's avatar
gsell committed
555
		for(i=0;i<f->nprocs;i++) 
gsell's avatar
gsell committed
556 557
			_H5Part_print_debug ( "\tnp=%lld",
					      (long long) f->pnparticles[i] );
gsell's avatar
gsell committed
558 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
	}
	/* 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 ) {
gsell's avatar
gsell committed
598
		r = _H5Part_set_step ( f, 0 );
gsell's avatar
gsell committed
599
		if ( r < 0 ) return r;
gsell's avatar
gsell committed
600 601 602 603 604 605 606 607 608
		
	}
#endif
	return H5PART_SUCCESS;
}

static h5part_int64_t
_write_data (
	H5PartFile *f,		/*!< IN: Handle to open file */
gsell's avatar
gsell committed
609 610 611
	const char *name,	/*!< IN: Name to associate array with */
	const void *array,	/*!< IN: Array to commit to disk */
	const hid_t type	/*!< IN: Type of data */
gsell's avatar
gsell committed
612 613 614 615
	) {
	herr_t herr;
	hid_t dataset_id;

gsell's avatar
gsell committed
616
	_H5Part_print_debug ( "Create a dataset[%s] mounted on the "
gsell's avatar
gsell committed
617
			      "timestep %lld",
gsell's avatar
gsell committed
618
			      name, (long long)f->timestep );
gsell's avatar
gsell committed
619 620 621 622 623 624 625 626 627 628

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

629 630 631 632 633 634 635 636 637
#ifdef COLLECTIVE_IO
	herr = H5Dwrite (
		dataset_id,
		type,
		f->memshape,
		f->diskshape,
		f->xfer_prop,
		array );
#else
gsell's avatar
gsell committed
638 639 640 641 642 643 644
	herr = H5Dwrite (
		dataset_id,
		type,
		f->memshape,
		f->diskshape,
		H5P_DEFAULT,
		array );
645 646
#endif

gsell's avatar
gsell committed
647 648 649 650 651
	if ( herr < 0 ) return HANDLE_H5D_WRITE_ERR ( name, f->timestep );

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

gsell's avatar
gsell committed
652 653
	f->empty = 0;

gsell's avatar
gsell committed
654 655 656 657
	return H5PART_SUCCESS;
}

/*!
gsell's avatar
gsell committed
658 659
  \ingroup h5part_write

gsell's avatar
gsell committed
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
  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 */
gsell's avatar
gsell committed
686 687
	const char *name,	/*!< [in] Name to associate array with */
	const h5part_float64_t *array	/*!< [in] Array to commit to disk */
gsell's avatar
gsell committed
688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703
	) {

	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;
}

/*!
gsell's avatar
gsell committed
704 705
  \ingroup h5part_write

gsell's avatar
gsell committed
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
  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 */
gsell's avatar
gsell committed
732 733
	const char *name,	/*!< [in] Name to associate array with */
	const h5part_int64_t *array	/*!< [in] Array to commit to disk */
gsell's avatar
gsell committed
734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752
	) {

	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 ***************/
gsell's avatar
gsell committed
753 754 755 756 757 758

/*!
  \ingroup h5partkernel
  @{
*/

gsell's avatar
gsell committed
759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782
/*!
   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
783
	_H5Part_print_warn ( "Unknown type %d", (int)type );
gsell's avatar
gsell committed
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

	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 );

	if ( attrib_nelem ) {
gsell's avatar
gsell committed
887 888 889
		space_id =  H5Aget_space ( attrib_id );
		if ( space_id < 0 ) return HANDLE_H5A_GET_SPACE_ERR;

gsell's avatar
gsell committed
890 891 892
		*attrib_nelem = H5Sget_simple_extent_npoints ( space_id );
		if ( *attrib_nelem < 0 )
			return HANDLE_H5S_GET_SIMPLE_EXTENT_NPOINTS_ERR;
gsell's avatar
gsell committed
893 894 895

		herr = H5Sclose ( space_id );
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
gsell's avatar
gsell committed
896 897 898 899 900 901 902 903 904
	}
	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 ) {
gsell's avatar
gsell committed
905 906
		mytype = H5Aget_type ( attrib_id );
		if ( mytype < 0 ) return HANDLE_H5A_GET_TYPE_ERR;
gsell's avatar
gsell committed
907

gsell's avatar
gsell committed
908
		*attrib_type = _H5Part_normalize_h5_type ( mytype );
gsell's avatar
gsell committed
909

gsell's avatar
gsell committed
910 911 912
		herr = H5Tclose ( mytype );
		if ( herr < 0 ) return HANDLE_H5T_CLOSE_ERR;
	}
gsell's avatar
gsell committed
913 914 915 916 917 918 919 920 921
	herr = H5Aclose ( attrib_id);
	if ( herr < 0 ) return HANDLE_H5A_CLOSE_ERR;

	return H5PART_SUCCESS;
}

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

/*!
gsell's avatar
gsell committed
922 923
  \ingroup h5part_attrib

gsell's avatar
gsell committed
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
  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;
}

/*!
gsell's avatar
gsell committed
965 966
  \ingroup h5part_attrib

gsell's avatar
gsell committed
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 999 1000 1001 1002 1003
  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;
}

/*!
gsell's avatar
gsell committed
1004 1005
  \ingroup h5part_attrib

gsell's avatar
gsell committed
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
  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;
}

/*!
gsell's avatar
gsell committed
1051 1052
  \ingroup h5part_attrib

gsell's avatar
gsell committed
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
  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
1089
	herr = _H5Part_write_attrib (
gsell's avatar
gsell committed
1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103
		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;
}

/*!
gsell's avatar
gsell committed
1104 1105
  \ingroup h5part_attrib

gsell's avatar
gsell committed
1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126
  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;
}

/*!
gsell's avatar
gsell committed
1127 1128
  \ingroup h5part_attrib

gsell's avatar
gsell committed
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
  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;
}

/*!
gsell's avatar
gsell committed
1156 1157
  \ingroup h5part_attrib

gsell's avatar
gsell committed
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
  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;
}

/*!
gsell's avatar
gsell committed
1199 1200
  \ingroup h5part_attrib

gsell's avatar
gsell committed
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
  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;
}

/*!
gsell's avatar
gsell committed
1249 1250
  \ingroup h5part_attrib

gsell's avatar
gsell committed
1251 1252 1253 1254 1255 1256 1257
  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 */
gsell's avatar
gsell committed
1258
	const char *attrib_name,	/*!< [in] Name of attribute to read */
gsell's avatar
gsell committed
1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274
	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;
}

/*!
gsell's avatar
gsell committed
1275 1276
  \ingroup h5part_attrib

gsell's avatar
gsell committed
1277 1278 1279 1280 1281 1282 1283
  Reads an attribute bound to file \c f.

  \return \c H5PART_SUCCESS or error code 
*/
h5part_int64_t
H5PartReadFileAttrib ( 
	H5PartFile *f,
gsell's avatar
gsell committed
1284
	const char *attrib_name,
gsell's avatar
gsell committed
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
	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).
*/

gsell's avatar
gsell committed
1320
h5part_int64_t
gsell's avatar
gsell committed
1321
_H5Part_set_step (
gsell's avatar
gsell committed
1322 1323 1324 1325 1326 1327
	H5PartFile *f,			/*!< [in]  Handle to open file */
	const h5part_int64_t step	/*!< [in]  Time-step to set. */
	) {

	char name[128];

gsell's avatar
gsell committed
1328 1329 1330 1331
	sprintf (
		name,
		"%s#%0*lld",
		f->groupname_step, f->stepno_width, (long long) step );
gsell's avatar
gsell committed
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
	herr_t herr = H5Gget_objinfo( f->file, name, 1, NULL );
	if ( (f->mode != H5PART_READ) && ( herr >= 0 ) ) {
		return HANDLE_H5PART_STEP_EXISTS_ERR ( step );
	}

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

	if( f->mode == H5PART_READ ) {
		_H5Part_print_info (
			"Proc[%d]: Set step to #%lld for file %lld",
			f->myproc,
			(long long)step,
			(long long)(size_t) f );

		f->timegroup = H5Gopen ( f->file, name ); 
		if ( f->timegroup < 0 ) return HANDLE_H5G_OPEN_ERR( name );
	}
	else {
		_H5Part_print_debug (
			"Proc[%d]: Create step #%lld for file %lld", 
			f->myproc,
			(long long)step,
			(long long)(size_t) f );

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

	return H5PART_SUCCESS;
}

gsell's avatar
gsell committed
1368
/*!
gsell's avatar
gsell committed
1369 1370
  \ingroup h5part_read

gsell's avatar
gsell committed
1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386
  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 */
gsell's avatar
gsell committed
1387
	const h5part_int64_t step	/*!< [in]  Time-step to set. */
gsell's avatar
gsell committed
1388 1389 1390 1391 1392 1393
	) {

	SET_FNAME ( "H5PartSetStep" );

	CHECK_FILEHANDLE ( f );

gsell's avatar
gsell committed
1394
	return _H5Part_set_step ( f, step );
gsell's avatar
gsell committed
1395 1396 1397 1398
}

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

gsell's avatar
gsell committed
1399 1400
/*!
  \ingroup h5part_kernel
gsell's avatar
gsell committed
1401

gsell's avatar
gsell committed
1402 1403
  Iterator for \c H5Giterate().
*/
gsell's avatar
gsell committed
1404 1405 1406 1407 1408 1409 1410
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 */
	) {

gsell's avatar
gsell committed
1411
	struct _iter_op_data *data = (struct _iter_op_data*)operator_data;
gsell's avatar
gsell committed
1412 1413 1414
	herr_t herr;
	H5G_stat_t objinfo;

gsell's avatar
gsell committed
1415 1416 1417
	if ( data->type != H5G_UNKNOWN ) {
		herr = H5Gget_objinfo ( group_id, member_name, 1, &objinfo );
		if ( herr < 0 ) return HANDLE_H5G_GET_OBJINFO_ERR ( member_name );
gsell's avatar
gsell committed
1418

gsell's avatar
gsell committed
1419 1420 1421
		if ( objinfo.type != data->type )
			return 0;/* don't count, continue iteration */
	}
gsell's avatar
gsell committed
1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439

	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 */
}

gsell's avatar
gsell committed
1440 1441 1442 1443 1444
/*!
  \ingroup h5part_kernel

  Iterator for \c H5Giterate().
*/
gsell's avatar
gsell committed
1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458
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 );
}

gsell's avatar
gsell committed
1459 1460 1461 1462 1463
/*!
  \ingroup h5part_kernel

  Iterator for \c H5Giterate().
*/
gsell's avatar
gsell committed
1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486
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;
}

gsell's avatar
gsell committed
1487 1488 1489 1490 1491
/*!
  \ingroup h5part_kernel

  Iterator for \c H5Giterate().
*/
gsell's avatar
gsell committed
1492 1493 1494 1495 1496 1497 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
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;
}

gsell's avatar
gsell committed
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
/*!
  \ingroup h5part_read

  Query whether a particular step already exists in the file
  \c f.

  It works for both reading and writing of files

  \return      true or false
*/
h5part_int64_t
H5PartHasStep (
	H5PartFile *f,		/*!< [in]  Handle to open file */
	h5part_int64_t step	/*!< [in]  Step number to query */
	) {
  
	SET_FNAME ( "H5PartHasStep" );

	CHECK_FILEHANDLE( f );

	char name[128];
	sprintf ( name, "%s#%0*lld", f->groupname_step, f->stepno_width, (long long) step );
	herr_t herr = H5Gget_objinfo( f->file, name, 1, NULL );

	return ( herr >= 0 );
}


gsell's avatar
gsell committed
1551
/*!
gsell's avatar
gsell committed
1552 1553
  \ingroup h5part_read

gsell's avatar
gsell committed
1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573
  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,
		"/",
gsell's avatar
gsell committed
1574
		H5G_UNKNOWN,
gsell's avatar
gsell committed
1575
		f->groupname_step );
gsell's avatar
gsell committed
1576 1577 1578
}

/*!
gsell's avatar
gsell committed
1579 1580
  \ingroup h5part_read

gsell's avatar
gsell committed
1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596
  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
1597 1598 1599 1600
	sprintf (
		stepname,
		"%s#%0*lld",
		f->groupname_step, f->stepno_width, (long long) f->timestep );
gsell's avatar
gsell committed
1601 1602 1603 1604 1605

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

/*!
gsell's avatar
gsell committed
1606 1607
  \ingroup h5part_read

gsell's avatar
gsell committed
1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629
  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
1630 1631 1632 1633
	sprintf (
		stepname,
		"%s#%0*lld",
		f->groupname_step, f->stepno_width, (long long) f->timestep );
gsell's avatar
gsell committed
1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644

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

/*!
gsell's avatar
gsell committed
1645 1646
  \ingroup h5part_read

gsell's avatar
gsell committed
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
  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
1675 1676 1677 1678
	sprintf (
		step_name,
		"%s#%0*lld",
		f->groupname_step, f->stepno_width, (long long) f->timestep );
gsell's avatar
gsell committed
1679 1680

	herr = _H5Part_get_object_name (
Kurt Stockinger's avatar
Kurt Stockinger committed
1681
		f->file,
gsell's avatar
gsell committed
1682 1683 1684 1685 1686 1687 1688
		step_name,
		H5G_DATASET,
		idx,
		dataset_name,
		len_dataset_name );
	if ( herr < 0 ) return herr;

gsell's avatar
gsell committed
1689
	*nelem = _H5Part_get_num_particles ( f );
gsell's avatar
gsell committed
1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721 1722 1723 1724 1725 1726 1727 1728 1729 1730
	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;
}

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
1731
		_H5Part_print_debug ( "Selection is available" );
gsell's avatar
gsell committed
1732 1733 1734 1735 1736 1737 1738 1739 1740 1741 1742 1743 1744 1745 1746 1747 1748 1749 1750 1751

		/* 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
1752
			"Selection: range=%d:%d, npoints=%d s=%d",
gsell's avatar
gsell committed
1753 1754 1755 1756
			(int)f->viewstart,(int)f->viewend,
			(int)H5Sget_simple_extent_npoints(space),
			(int)H5Sget_select_npoints(space) );
	} else {
gsell's avatar
gsell committed
1757
		_H5Part_print_debug ( "Selection" );
gsell's avatar
gsell committed
1758 1759 1760 1761 1762 1763 1764 1765 1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 1776 1777 1778 1779 1780 1781 1782 1783 1784
	}
	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;
	}
}

h5part_int64_t
gsell's avatar
gsell committed
1785
_H5Part_get_num_particles (
gsell's avatar
gsell committed
1786 1787 1788 1789 1790 1791 1792 1793 1794 1795 1796
	H5PartFile *f			/*!< [in]  Handle to open file */
	) {

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

	/* Get first dataset in current time-step */
gsell's avatar
gsell committed
1797
	sprintf (
gsell's avatar
gsell committed
1798 1799 1800 1801
		step_name,
		"%s#%0*lld",
		f->groupname_step, f->stepno_width, (long long) f->timestep );

gsell's avatar
gsell committed
1802 1803 1804 1805 1806 1807
	herr = _H5Part_get_object_name (
		f->file,
		step_name,
		H5G_DATASET,
		0,
		dataset_name, sizeof (dataset_name) );
gsell's avatar
gsell committed
1808 1809 1810
	if ( herr < 0 ) return herr;

	dataset_id = H5Dopen ( f->timegroup, dataset_name );
gsell's avatar
gsell committed
1811
	if ( dataset_id < 0 ) 
gsell's avatar
gsell committed
1812 1813 1814 1815 1816 1817 1818 1819 1820 1821 1822 1823 1824 1825 1826 1827 1828 1829 1830 1831 1832 1833 1834 1835
		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;
}

gsell's avatar
gsell committed
1836 1837 1838 1839 1840 1841 1842 1843 1844 1845
/*!
  \ingroup h5part_read

  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.
 */
gsell's avatar
gsell committed
1846 1847 1848 1849 1850 1851 1852 1853 1854 1855
h5part_int64_t
H5PartGetNumParticles (
	H5PartFile *f			/*!< [in]  Handle to open file */
	) {

	SET_FNAME ( "H5PartGetNumParticles" );

	CHECK_FILEHANDLE( f );

	if ( f->timegroup < 0 ) {
gsell's avatar
gsell committed
1856
		h5part_int64_t herr = _H5Part_set_step ( f, 0 );
gsell's avatar
gsell committed
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
		if ( herr < 0 ) return herr;
	}

	return _H5Part_get_num_particles ( f );
}

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

	herr_t herr = 0;

	f->viewstart = -1;
	f->viewend = -1;
	if ( f->shape != 0 ){
		herr = H5Sclose(f->shape);
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
		f->shape=0;
	}
	if(f->diskshape!=0 && f->diskshape!=H5S_ALL){
		herr = H5Sclose(f->diskshape);
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
		f->diskshape=H5S_ALL;
	}
	f->diskshape = H5S_ALL;
	if(f->memshape!=0 && f->memshape!=H5S_ALL){
		herr = H5Sclose ( f->memshape );
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
		f->memshape=H5S_ALL;
	}
	return H5PART_SUCCESS;
}

gsell's avatar
gsell committed
1891 1892 1893
/*!
  \ingroup h5part_read
*/
gsell's avatar
gsell committed
1894 1895 1896 1897 1898 1899 1900 1901 1902 1903 1904 1905
h5part_int64_t
H5PartResetView (
 	H5PartFile *f			/*!< [in]  Handle to open file */
	) {
	SET_FNAME ( "H5PartResetView" );

	CHECK_FILEHANDLE( f );
	CHECK_READONLY_MODE ( f );

	return _reset_view ( f );
}

gsell's avatar
gsell committed
1906 1907 1908
/*!
  \ingroup h5part_read
*/
gsell's avatar
gsell committed
1909 1910 1911 1912 1913 1914 1915 1916 1917 1918 1919
h5part_int64_t
H5PartHasView (
 	H5PartFile *f			/*!< [in]  Handle to open file */
	) {
	SET_FNAME ( "H5PartResetView" );

	CHECK_FILEHANDLE( f );
	CHECK_READONLY_MODE ( f );

	return  ( f->viewstart >= 0 ) && ( f->viewend >= 0 );
}
gsell's avatar
gsell committed
1920

gsell's avatar
gsell committed
1921 1922
static h5part_int64_t
_set_view (
gsell's avatar
gsell committed
1923 1924 1925 1926
	H5PartFile *f,			/*!< [in]  Handle to open file */
	h5part_int64_t start,		/*!< [in]  Start particle */
	h5part_int64_t end		/*!< [in]  End particle */
	) {
gsell's avatar
gsell committed
1927
	h5part_int64_t herr = 0;
gsell's avatar
gsell committed
1928
	hsize_t total;
gsell's avatar
gsell committed
1929 1930
	hsize_t stride = 1;
	hsize_t dmax = H5S_UNLIMITED;
gsell's avatar
gsell committed
1931

gsell's avatar
gsell committed
1932 1933 1934
	_H5Part_print_debug (
		"Set view (%lld,%lld).",
		(long long)start,(long long)end);
gsell's avatar
gsell committed
1935

gsell's avatar
gsell committed
1936 1937
	herr = _reset_view ( f );
	if ( herr < 0 ) return herr;
gsell's avatar
gsell committed
1938

gsell's avatar
gsell committed
1939
	if ( start == -1 && end == -1 ) return H5PART_SUCCESS;
gsell's avatar
gsell committed
1940

gsell's avatar
gsell committed
1941 1942 1943
	/*
	  View has been reset so H5PartGetNumParticles will tell
	  us the total number of particles.
gsell's avatar
gsell committed
1944

gsell's avatar
gsell committed
1945 1946 1947 1948
	  For now, we interpret start=-1 to mean 0 and 
	  end==-1 to mean end of file
	*/
	total = (hsize_t) _H5Part_get_num_particles ( f );
gsell's avatar
gsell committed
1949 1950
	if ( total < 0 ) return HANDLE_H5PART_GET_NUM_PARTICLES_ERR ( total );

gsell's avatar
gsell committed
1951 1952 1953
	if ( start == -1 ) start = 0;
	if ( end == -1 )   end = total;

gsell's avatar
gsell committed
1954
	_H5Part_print_debug ( "Total nparticles=%lld", (long long)total );
gsell's avatar
gsell committed
1955 1956 1957 1958

	/* so, is this selection inclusive or exclusive? 
	   it appears to be inclusive for both ends of the range.
	*/
gsell's avatar
gsell committed
1959
	if ( end < start ) {
gsell's avatar
gsell committed
1960
		_H5Part_print_warn (
gsell's avatar
gsell committed
1961 1962 1963
			"Nonfatal error. "
			"End of view (%lld) is less than start (%lld).",
			(long long)end, (long long)start );
gsell's avatar
gsell committed
1964 1965 1966
		end = start; /* ensure that we don't have a range error */
	}
	/* setting up the new view */
gsell's avatar
gsell committed
1967 1968
	f->viewstart =  start;
	f->viewend =    end;
gsell's avatar
gsell committed
1969
	f->nparticles = end - start + 1;
gsell's avatar
gsell committed
1970 1971 1972 1973 1974 1975 1976 1977 1978 1979 1980 1981 1982 1983
	
	/* 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 )
gsell's avatar
gsell committed
1984
		return HANDLE_H5S_CREATE_SIMPLE_ERR ( f->nparticles );
gsell's avatar
gsell committed
1985

gsell's avatar
gsell committed
1986
	herr = H5Sselect_hyperslab ( 
gsell's avatar
gsell committed
1987 1988
		f->diskshape,
		H5S_SELECT_SET,
gsell's avatar
gsell committed
1989
		(hsize_t*)&start,
gsell's avatar
gsell committed
1990 1991 1992
		&stride,
		&total,
		NULL );
gsell's avatar
gsell committed
1993
	if ( herr < 0 ) return HANDLE_H5S_SELECT_HYPERSLAB_ERR;
gsell's avatar
gsell committed
1994 1995 1996 1997

	return H5PART_SUCCESS;
}

gsell's avatar
gsell committed
1998 1999 2000 2001 2002 2003 2004 2005 2006 2007 2008 2009 2010 2011 2012 2013 2014 2015 2016
/*!
  \ingroup h5part_read

  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
*/
gsell's avatar
gsell committed
2017 2018 2019 2020 2021 2022 2023 2024 2025 2026 2027 2028 2029
h5part_int64_t
H5PartSetView (
	H5PartFile *f,			/*!< [in]  Handle to open file */
	const h5part_int64_t start,	/*!< [in]  Start particle */
	const h5part_int64_t end	/*!< [in]  End particle */
	) {

	SET_FNAME ( "H5PartSetView" );

	CHECK_FILEHANDLE( f );
	CHECK_READONLY_MODE ( f );

	if ( f->timegroup < 0 ) {
gsell's avatar
gsell committed
2030
		h5part_int64_t herr = _H5Part_set_step ( f, 0 );
gsell's avatar
gsell committed
2031 2032 2033 2034 2035 2036
		if ( herr < 0 ) return herr;
	}

	return _set_view ( f, start, end );
}

gsell's avatar
gsell committed
2037
/*!
gsell's avatar
gsell committed
2038 2039
  \ingroup h5part_read

gsell's avatar
gsell committed
2040 2041 2042 2043 2044 2045 2046 2047 2048 2049 2050 2051 2052 2053 2054 2055 2056 2057
   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" );

	CHECK_FILEHANDLE( f );

gsell's avatar
gsell committed
2058
	if ( f->timegroup < 0 ) {
gsell's avatar
gsell committed
2059
		h5part_int64_t herr = _H5Part_set_step ( f, 0 );
gsell's avatar
gsell committed
2060 2061 2062 2063 2064 2065 2066 2067
		if ( herr < 0 ) return herr;
	}

	h5part_int64_t viewstart = 0;
	h5part_int64_t viewend = 0;

	if ( f->viewstart >= 0 )
		viewstart = f->viewstart;
gsell's avatar
gsell committed
2068 2069 2070 2071 2072

	if ( f->viewend >= 0 ) {
		viewend = f->viewend;
	}
	else {
gsell's avatar
gsell committed
2073
		viewend = _H5Part_get_num_particles ( f );
gsell's avatar
gsell committed
2074 2075 2076 2077
		if ( viewend < 0 )
			return HANDLE_H5PART_GET_NUM_PARTICLES_ERR ( viewend );
	}

gsell's avatar
gsell committed
2078 2079
	if ( start ) *start = viewstart;
	if ( end )   *end = viewend;
gsell's avatar
gsell committed
2080

gsell's avatar
gsell committed
2081
	return viewend - viewstart;
gsell's avatar
gsell committed
2082 2083 2084
}

/*!
gsell's avatar
gsell committed
2085 2086
  \ingroup h5part_read

gsell's avatar
gsell committed
2087 2088 2089 2090 2091 2092 2093
  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
*/
gsell's avatar
gsell committed
2094 2095 2096 2097 2098 2099 2100
/*
  \note
  There is a bug in this function:
  If (NumParticles % f->nprocs) != 0  then
  the last  (NumParticles % f->nprocs) particles are not handled!
*/

gsell's avatar
gsell committed
2101 2102 2103 2104 2105 2106 2107
h5part_int64_t
H5PartSetCanonicalView (
	H5PartFile *f			/*!< [in]  Handle to open file */
	) {

	SET_FNAME ( "H5PartSetCanonicalView" );

gsell's avatar
gsell committed
2108
	h5part_int64_t herr;
gsell's avatar
gsell committed
2109 2110

	CHECK_FILEHANDLE( f );
gsell's avatar
gsell committed
2111
	CHECK_READONLY_MODE ( f )
gsell's avatar
gsell committed
2112

gsell's avatar
gsell committed
2113
	herr = _reset_view ( f );
gsell's avatar
gsell committed
2114
	if ( herr < 0 ) return HANDLE_H5PART_SET_VIEW_ERR( herr, -1, -1 );
gsell's avatar
gsell committed
2115 2116

#ifdef PARALLEL_IO
gsell's avatar
gsell committed
2117 2118 2119 2120 2121 2122
	h5part_int64_t start = 0;
	h5part_int64_t end = 0;
	h5part_int64_t n = 0;
	int i = 0;
	
	if ( f->timegroup < 0 ) {
gsell's avatar
gsell committed
2123
		herr = _H5Part_set_step ( f, 0 );
gsell's avatar
gsell committed
2124 2125 2126 2127 2128 2129 2130 2131 2132 2133 2134 2135 2136 2137 2138
		if ( herr < 0 ) return herr;
	}
	n = _H5Part_get_num_particles ( 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.
	*/
	if ( _H5Part_read_attrib (
		     f->timegroup,
		     "pnparticles", f->pnparticles ) < 0) {
		/*
		  Attribute "pnparticles" is not available.  So
		  subdivide the view into NP mostly equal pieces
gsell's avatar
gsell committed
2139
		*/
gsell's avatar
gsell committed
2140 2141 2142
		n /= f->nprocs;
		for ( i=0; i<f->nprocs; i++ ) {
			f->pnparticles[i] = n;
gsell's avatar
gsell committed
2143
		}
gsell's avatar
gsell committed
2144
	}
gsell's avatar
gsell committed
2145

gsell's avatar
gsell committed
2146 2147
	for ( i = 0; i < f->myproc; i++ ){
		start += f->pnparticles[i];
gsell's avatar
gsell committed
2148
	}
gsell's avatar
gsell committed
2149
	end = start + f->pnparticles[f->myproc] - 1;
gsell's avatar
gsell committed
2150 2151
	herr = _set_view ( f, start, end );
	if ( herr < 0 ) return HANDLE_H5PART_SET_VIEW_ERR ( herr, start, end );
gsell's avatar
gsell committed
2152

gsell's avatar
gsell committed
2153
#endif
gsell's avatar
gsell committed
2154

gsell's avatar
gsell committed
2155 2156 2157 2158 2159 2160
	return H5PART_SUCCESS;
}

static h5part_int64_t
_read_data (
	H5PartFile *f,		/*!< [in] Handle to open file */