H5Part.c 88.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
<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>
22
<LI> Mark Howison (NERSC/LBNL)</LI>
Kurt Stockinger's avatar
Kurt Stockinger committed
23
</UL>
gsell's avatar
gsell committed
24 25 26 27


Papers: 

Kurt Stockinger's avatar
Kurt Stockinger committed
28 29 30 31 32 33 34
<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>


35
<LI> A. Adelmann, R.D. Ryne, J. Shalf, C. Siegerist, "H5Part: A Portable High Performance Parallel Data Interface for Particle Simulations," <i>
Kurt Stockinger's avatar
Kurt Stockinger committed
36 37 38 39 40 41 42
<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>

gsell's avatar
gsell committed
43 44
*/

gsell's avatar
gsell committed
45 46

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

*/
/*!
gsell's avatar
gsell committed
51
  \ingroup h5part_c_api
52
  \defgroup h5part_open 	File Opening and Closing
gsell's avatar
gsell committed
53 54
*/
/*!
gsell's avatar
gsell committed
55
  \ingroup h5part_c_api
56
  \defgroup h5part_model	Setting up the Data Model
gsell's avatar
gsell committed
57 58
*/  
/*!
gsell's avatar
gsell committed
59
  \ingroup h5part_c_api
60
  \defgroup h5part_data		Readind and Writing Datasets
gsell's avatar
gsell committed
61 62
*/  
/*!
gsell's avatar
gsell committed
63
  \ingroup h5part_c_api
gsell's avatar
gsell committed
64 65 66
  \defgroup h5part_attrib	Reading and Writing Attributes
*/
/*!
gsell's avatar
gsell committed
67
  \ingroup h5part_c_api
gsell's avatar
gsell committed
68 69 70 71 72 73 74 75
  \defgroup h5part_errhandle	Error Handling
*/
/*!
  \internal
  \defgroup h5partkernel H5Part private functions 
*/


gsell's avatar
gsell committed
76 77 78 79 80 81 82 83
#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
84 85 86 87 88 89 90 91
#ifndef WIN32
#include <unistd.h>
#else /* WIN32 */
#include <io.h>
#define open  _open
#define close _close
#endif /* WIN32 */

gsell's avatar
gsell committed
92 93 94 95
#include "H5Part.h"
#include "H5PartPrivate.h"
#include "H5PartErrors.h"

96 97 98
/********* Global Variable Declarations *************/
h5part_error_handler	_err_handler = H5PartReportErrorHandler;

gsell's avatar
gsell committed
99 100
/********* Private Variable Declarations *************/

101
static unsigned			_debug = H5PART_VERB_ERROR;
gsell's avatar
gsell committed
102
static h5part_int64_t		_h5part_errno = H5PART_SUCCESS;
103
static char *__funcname;
gsell's avatar
gsell committed
104 105 106 107 108 109 110 111 112

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

static h5part_int64_t
_init(
	void
	);

static h5part_int64_t
113 114
_reset_view (
	H5PartFile *f
gsell's avatar
gsell committed
115 116
	);

gsell's avatar
gsell committed
117 118 119
/*
  error handler for hdf5
*/
gsell's avatar
gsell committed
120 121
static herr_t
_h5_error_handler (
122
#ifndef H5_USE_16_API
123 124
	hid_t,
#endif
gsell's avatar
gsell committed
125 126 127 128 129
	void *
	);

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

gsell's avatar
gsell committed
130 131
static H5PartFile*
_H5Part_open_file (
gsell's avatar
gsell committed
132
	const char *filename,	/*!< [in] The name of the data file to open. */
133
	const char flags,	/*!< [in] The access mode for the file. */
gsell's avatar
gsell committed
134
	MPI_Comm comm,		/*!< [in] MPI communicator */
135 136 137 138
	int f_parallel,		/*!< [in] 0 for serial io otherwise parallel */
	h5part_int64_t align	/*!< [in] Number of bytes for setting alignment,
					  metadata block size, etc.
					  Set to 0 to disable. */
gsell's avatar
gsell committed
139
	) {
140

gsell's avatar
gsell committed
141
	_h5part_errno = H5PART_SUCCESS;
gsell's avatar
gsell committed
142 143 144 145 146 147 148 149
	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
150

151 152
	f->flags = flags;

153 154
	/* set default step name */
	strncpy ( f->groupname_step, H5PART_GROUPNAME_STEP, H5PART_STEPNAME_LEN );
gsell's avatar
gsell committed
155 156
	f->stepno_width = 0;

157
	f->dxfer_prop = f->dcreate_prop = f->fcreate_prop = H5P_DEFAULT;
158

159 160
	f->faccess_prop = H5Pcreate (H5P_FILE_ACCESS);
	if (f->faccess_prop < 0) {
161 162 163
		HANDLE_H5P_CREATE_ERR;
		goto error_cleanup;
	}
gsell's avatar
gsell committed
164

gsell's avatar
gsell committed
165
	if ( f_parallel ) {
gsell's avatar
gsell committed
166
#ifdef PARALLEL_IO
gsell's avatar
gsell committed
167
		MPI_Info info = MPI_INFO_NULL;
168

gsell's avatar
gsell committed
169 170 171 172 173 174 175 176
		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
177

gsell's avatar
gsell committed
178 179
		f->pnparticles =
		  (h5part_int64_t*) malloc (f->nprocs * sizeof (h5part_int64_t));
gsell's avatar
gsell committed
180 181 182 183
		if (f->pnparticles == NULL) {
			HANDLE_H5PART_NOMEM_ERR;
			goto error_cleanup;
		}
gsell's avatar
gsell committed
184

185 186 187 188 189
		/* select the HDF5 VFD */
		if (flags & H5PART_VFD_MPIPOSIX) {
			if (f->myproc == 0) {
				_H5Part_print_info ( "Selecting MPI-POSIX VFD" );
			}
190
			if (H5Pset_fapl_mpiposix ( f->faccess_prop, comm, 0 ) < 0) {
191 192 193
				HANDLE_H5P_SET_FAPL_ERR;
				goto error_cleanup;
			}
gsell's avatar
gsell committed
194
		}
195 196 197 198
		else {
			if (f->myproc == 0) {
				_H5Part_print_info ( "Selecting MPI-IO VFD" );
			}
199
			if (H5Pset_fapl_mpio ( f->faccess_prop, comm, info ) < 0) {
200 201 202
				HANDLE_H5P_SET_FAPL_ERR;
				goto error_cleanup;
			}
203
			if (flags & H5PART_VFD_MPIIO_IND) {
204 205 206
				if (f->myproc == 0) {
					_H5Part_print_info ( "Using independent mode" );
				}
207
			} else {
208 209 210
				if (f->myproc == 0) {
					_H5Part_print_info ( "Using collective mode" );
				}
211 212
				f->dxfer_prop = H5Pcreate (H5P_DATASET_XFER);
				if (f->dxfer_prop < 0) {
213 214 215
					HANDLE_H5P_CREATE_ERR;
					goto error_cleanup;
				}
216
				if (H5Pset_dxpl_mpio ( f->dxfer_prop, H5FD_MPIO_COLLECTIVE ) < 0) {
217 218 219
					HANDLE_H5P_SET_DXPL_MPIO_ERR;
					goto error_cleanup;
				}
220
			}
gsell's avatar
gsell committed
221
		}
gsell's avatar
gsell committed
222

223 224
		if ( flags & H5PART_FS_LUSTRE )
		{
225
			/* extend the btree size so that metadata pieces are
226 227 228 229 230 231 232 233 234
			 * close to the alignment value */
			unsigned int btree_ik = (align - 256) / 96;
			unsigned int btree_bytes = 64 + 96*btree_ik;
			if ( btree_bytes > align ) {
				HANDLE_H5PART_INVALID_ERR(
					"btree_ik", btree_ik );
				goto error_cleanup;
			}

235 236
			if (f->myproc == 0) {
				_H5Part_print_info (
237 238
					"Setting HDF5 btree parameter to %u",
					btree_ik );
239
				_H5Part_print_info (
240 241
					"Extending HDF5 btree size to %u bytes at rank 3",
					btree_bytes );
242
			}
243

244 245 246 247 248 249 250
			f->fcreate_prop = H5Pcreate(H5P_FILE_CREATE);
			if ( f->fcreate_prop < 0 ) {
				HANDLE_H5P_CREATE_ERR;
				goto error_cleanup;
			}

			H5Pset_istore_k (f->fcreate_prop, btree_ik);
251

252
#ifndef H5_USE_16_API
253 254 255
			/* defer metadata cache flushing until file close */
			H5AC_cache_config_t cache_config;
			cache_config.version = H5AC__CURR_CACHE_CONFIG_VERSION;
256
			H5Pget_mdc_config (f->faccess_prop, &cache_config);
257 258 259 260 261 262
			cache_config.set_initial_size = 1;
			cache_config.initial_size = 16 * 1024 * 1024;
			cache_config.evictions_enabled = 0;
			cache_config.incr_mode = H5C_incr__off;
			cache_config.flash_incr_mode = H5C_flash_incr__off;
			cache_config.decr_mode = H5C_decr__off;
263
			H5Pset_mdc_config (f->faccess_prop, &cache_config);
264
#else // H5_USE_16_API
265 266
			_H5Part_print_info (
					"Unable to defer metadata write: need HDF5 1.8");
267
#endif // H5_USE_16_API
268
		}
269

270
		f->comm = comm;
271
#endif // PARALLEL_IO
gsell's avatar
gsell committed
272
	} else {
gsell's avatar
gsell committed
273
		f->comm = 0;
gsell's avatar
gsell committed
274 275
		f->nprocs = 1;
		f->myproc = 0;
gsell's avatar
gsell committed
276 277
		f->pnparticles = 
			(h5part_int64_t*) malloc (f->nprocs * sizeof (h5part_int64_t));
gsell's avatar
gsell committed
278
	}
279 280 281

	if ( align != 0 ) {
		if (f->myproc == 0) {
282 283 284
			_H5Part_print_info (
				"Setting HDF5 alignment to %ld bytes",
				align );
285
		}
286
		if (H5Pset_alignment ( f->faccess_prop, 0, align ) < 0) {
287 288 289 290
			HANDLE_H5P_SET_FAPL_ERR;
			goto error_cleanup;
		}
		if (f->myproc == 0) {
291 292 293
			_H5Part_print_info (
				"Setting HDF5 meta block to %ld bytes",
				align );
294
		}
295
		if (H5Pset_meta_block_size ( f->faccess_prop, align ) < 0) {
296 297 298 299 300 301 302 303 304 305 306 307 308
			HANDLE_H5P_SET_FAPL_ERR;
			goto error_cleanup;
		}
		/*if (f->myproc == 0) {
			_H5Part_print_info ( "Setting HDF5 sieve buffer to %ld bytes", align );
		}
		if (H5Pset_sieve_buf_size ( f->access_prop, align ) < 0) {
			HANDLE_H5P_SET_FAPL_ERR;
			goto error_cleanup;
		}*/
	}

	if ( flags & H5PART_READ ) {
309
		f->file = H5Fopen(filename, H5F_ACC_RDONLY, f->faccess_prop);
gsell's avatar
gsell committed
310
	}
311
	else if ( flags & H5PART_WRITE ){
312 313
		f->file = H5Fcreate (filename, H5F_ACC_TRUNC, f->fcreate_prop,
				f->faccess_prop);
gsell's avatar
gsell committed
314
		f->empty = 1;
gsell's avatar
gsell committed
315
	}
316
	else if ( flags & H5PART_APPEND ) {
317
		int fd = open(filename, O_RDONLY, 0);
gsell's avatar
gsell committed
318 319
		if ( (fd == -1) && (errno == ENOENT) ) {
			f->file = H5Fcreate(filename, H5F_ACC_TRUNC,
320
					f->fcreate_prop, f->faccess_prop);
gsell's avatar
gsell committed
321
			f->empty = 1;
gsell's avatar
gsell committed
322 323 324
		}
		else if (fd != -1) {
			close (fd);
325
			f->file = H5Fopen(
326
				filename, H5F_ACC_RDWR, f->faccess_prop);
gsell's avatar
gsell committed
327 328 329 330 331
			/*
			  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
332
				f->file, "/", H5G_GROUP, f->groupname_step );
gsell's avatar
gsell committed
333 334 335 336 337 338 339 340 341 342 343 344
			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;
	}
345
	f->nparticles = 0;
gsell's avatar
gsell committed
346
	f->timegroup = -1;
347
	f->shape = H5S_ALL;
gsell's avatar
gsell committed
348 349 350 351
	f->diskshape = H5S_ALL;
	f->memshape = H5S_ALL;
	f->viewstart = -1;
	f->viewend = -1;
352
	f->viewindexed = 0;
353
	f->throttle = 0;
gsell's avatar
gsell committed
354

gsell's avatar
gsell committed
355 356 357 358 359 360
	_H5Part_print_debug (
		"Proc[%d]: Opened file \"%s\" val=%lld",
		f->myproc,
		filename,
		(long long)(size_t)f );

gsell's avatar
gsell committed
361 362 363 364 365 366 367 368 369 370 371 372
	return f;

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

gsell's avatar
gsell committed
373
#ifdef PARALLEL_IO
gsell's avatar
gsell committed
374
/*!
375
  \ingroup h5part_open
gsell's avatar
gsell committed
376

gsell's avatar
gsell committed
377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394
  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
395
H5PartOpenFileParallel (
gsell's avatar
gsell committed
396
	const char *filename,	/*!< [in] The name of the data file to open. */
397
	const char flags,	/*!< [in] The access mode for the file. */
gsell's avatar
gsell committed
398
	MPI_Comm comm		/*!< [in] MPI communicator */
399
	) {
400
	INIT
401 402 403
	SET_FNAME ( "H5PartOpenFileParallel" );

	int f_parallel = 1;	/* parallel i/o */
404
	h5part_int64_t align = 0; /* no alignment tuning */
405 406 407 408 409

	return _H5Part_open_file ( filename, flags, comm, f_parallel, align );
}

/*!
410
  \ingroup h5part_open
411 412 413 414 415 416 417 418 419

  Opens file with specified filename, and also specifices an alignment
  value used for HDF5 tuning parameters.

  \return	File handle or \c NULL
 */
H5PartFile*
H5PartOpenFileParallelAlign (
	const char *filename,	/*!< [in] The name of the data file to open. */
420
	const char flags,	/*!< [in] The access mode for the file. */
421 422 423
	MPI_Comm comm,		/*!< [in] MPI communicator */
	h5part_int64_t align	/*!< [in] Alignment size in bytes. */
	) {
424
	INIT
425
	SET_FNAME ( "H5PartOpenFileParallelAlign" );
426

gsell's avatar
gsell committed
427
	int f_parallel = 1;	/* parallel i/o */
gsell's avatar
gsell committed
428

429
	return _H5Part_open_file ( filename, flags, comm, f_parallel, align );
gsell's avatar
gsell committed
430 431
}
#endif
gsell's avatar
gsell committed
432

gsell's avatar
gsell committed
433
/*!
434
  \ingroup  h5part_open
gsell's avatar
gsell committed
435

gsell's avatar
gsell committed
436
  Opens file with specified filename. 
gsell's avatar
gsell committed
437

gsell's avatar
gsell committed
438 439 440 441 442
  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
443

gsell's avatar
gsell committed
444 445 446 447 448 449
  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
450

gsell's avatar
gsell committed
451 452 453 454 455
  \return	File handle or \c NULL
 */
H5PartFile*
H5PartOpenFile (
	const char *filename,	/*!< [in] The name of the data file to open. */
456
	const char flags	/*!< [in] The access mode for the file. */
gsell's avatar
gsell committed
457
	) {
458
	INIT
gsell's avatar
gsell committed
459
	SET_FNAME ( "H5PartOpenFile" );
gsell's avatar
gsell committed
460

gsell's avatar
gsell committed
461 462
	MPI_Comm comm = 0;	/* dummy */
	int f_parallel = 0;	/* serial open */
463
	int align = 0;		/* no tuning parameters */
gsell's avatar
gsell committed
464

465 466 467 468
	return _H5Part_open_file ( filename, flags, comm, f_parallel, align );
}

/*!
469
  \ingroup h5part_open
470 471 472 473 474 475 476 477 478

  Opens file with specified filename, and also specifices an alignment
  value used for HDF5 tuning parameters.

  \return	File handle or \c NULL
 */
H5PartFile*
H5PartOpenFileAlign (
	const char *filename,	/*!< [in] The name of the data file to open. */
479
	const char flags,	/*!< [in] The access mode for the file. */
480 481
	h5part_int64_t align	/*!< [in] Alignment size in bytes. */
) {
482 483
	INIT
	SET_FNAME ( "H5PartOpenFileAlign" );
484 485 486 487 488

	MPI_Comm comm = 0;	/* dummy */
	int f_parallel = 0;	/* serial open */

	return _H5Part_open_file ( filename, flags, comm, f_parallel, align );
gsell's avatar
gsell committed
489 490 491 492 493 494 495
}

/*!
  Checks if a file was successfully opened.

  \return	\c H5PART_SUCCESS or error code
 */
496 497
h5part_int64_t
_H5Part_file_is_valid (
gsell's avatar
gsell committed
498 499 500 501 502 503 504 505 506 507 508
	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;
}

509 510 511 512 513 514 515
/*!
  \ingroup h5part_open

  Closes an open file.

  \return	\c H5PART_SUCCESS or error code
*/
gsell's avatar
gsell committed
516
h5part_int64_t
517
H5PartCloseFile (
gsell's avatar
gsell committed
518 519 520
	H5PartFile *f		/*!< [in] filehandle of the file to close */
	) {

521
	SET_FNAME ( "H5PartCloseFile" );
gsell's avatar
gsell committed
522
	herr_t r = 0;
523 524 525
	_h5part_errno = H5PART_SUCCESS;

	CHECK_FILEHANDLE ( f );
gsell's avatar
gsell committed
526 527 528 529 530 531 532

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

533 534 535 536 537 538 539 540
#ifdef PARALLEL_IO
	if ( f->multiblock && f->close_multiblock ) {
		(*f->close_multiblock) ( f );
		f->multiblock = NULL;
		f->close_multiblock = NULL;
	}
#endif

541
	if( f->shape != H5S_ALL ) {
gsell's avatar
gsell committed
542 543 544 545
		r = H5Sclose( f->shape );
		if ( r < 0 ) HANDLE_H5S_CLOSE_ERR;
		f->shape = 0;
	}
gsell's avatar
gsell committed
546
	if( f->timegroup >= 0 ) {
gsell's avatar
gsell committed
547 548
		r = H5Gclose( f->timegroup );
		if ( r < 0 ) HANDLE_H5G_CLOSE_ERR;
gsell's avatar
gsell committed
549
		f->timegroup = -1;
gsell's avatar
gsell committed
550 551 552 553 554 555
	}
	if( f->diskshape != H5S_ALL ) {
		r = H5Sclose( f->diskshape );
		if ( r < 0 ) HANDLE_H5S_CLOSE_ERR;
		f->diskshape = 0;
	}
556 557 558 559 560
	if( f->memshape != H5S_ALL ) {
		r = H5Sclose( f->memshape );
		if ( r < 0 ) HANDLE_H5S_CLOSE_ERR;
		f->memshape = 0;
	}
561 562 563 564
	if( f->dxfer_prop != H5P_DEFAULT ) {
		r = H5Pclose( f->dxfer_prop );
		if ( r < 0 ) HANDLE_H5P_CLOSE_ERR ( "f->dxfer_prop" );
		f->dxfer_prop = H5P_DEFAULT;
gsell's avatar
gsell committed
565
	}
566 567 568 569
	if( f->dcreate_prop != H5P_DEFAULT ) {
		r = H5Pclose( f->dcreate_prop );
		if ( r < 0 ) HANDLE_H5P_CLOSE_ERR ( "f->dcreate_prop" );
		f->dcreate_prop = H5P_DEFAULT;
gsell's avatar
gsell committed
570
	}
571 572


gsell's avatar
gsell committed
573 574 575 576 577
	if ( f->file ) {
		r = H5Fclose( f->file );
		if ( r < 0 ) HANDLE_H5F_CLOSE_ERR;
		f->file = 0;
	}
578 579 580 581 582 583 584 585 586 587
	if( f->faccess_prop != H5P_DEFAULT ) {
		r = H5Pclose( f->faccess_prop );
		if ( r < 0 ) HANDLE_H5P_CLOSE_ERR ( "f->faccess_prop" );
		f->faccess_prop = H5P_DEFAULT;
	}  
	if( f->fcreate_prop != H5P_DEFAULT ) {
		r = H5Pclose( f->fcreate_prop );
		if ( r < 0 ) HANDLE_H5P_CLOSE_ERR ( "f->fcreate_prop" );
		f->fcreate_prop = H5P_DEFAULT;
	}
588 589

	/* free memory from H5PartFile struct */
gsell's avatar
gsell committed
590 591 592 593 594
	if( f->pnparticles ) {
		free( f->pnparticles );
	}
	free( f );

gsell's avatar
gsell committed
595
	return _h5part_errno;
gsell's avatar
gsell committed
596 597
}

598 599 600 601 602 603 604
h5part_int64_t
H5PartFileIsValid (
	H5PartFile *f
	) {
	return _H5Part_file_is_valid(f);
}

gsell's avatar
gsell committed
605
/*============== File Writing Functions ==================== */
gsell's avatar
gsell committed
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
h5part_int64_t
_H5Part_get_step_name(
	H5PartFile *f,
	const h5part_int64_t step,
	char *name
	) {

	/* Work around sprintf bug on older systems */
	if (f->stepno_width == 0 && step == 0) {
		sprintf (
			name,
			"%s#%0*lld",
			f->groupname_step, 1, (long long) step );
	}
	else {
		sprintf (
			name,
			"%s#%0*lld",
			f->groupname_step, f->stepno_width, (long long) step );
	}

	return H5PART_SUCCESS;
}

gsell's avatar
gsell committed
631 632 633 634 635 636
h5part_int64_t
H5PartDefineStepName (
	H5PartFile *f,
	const char *name,
	const h5part_int64_t width
	) {
637 638 639 640 641 642 643 644

	CHECK_FILEHANDLE ( f );

	h5part_int64_t len = H5PART_STEPNAME_LEN - width - 2;
	if ( strlen(name) > len ) {
		_H5Part_print_warn (
			"Step name has been truncated to fit within %d chars.",
			H5PART_STEPNAME_LEN );
gsell's avatar
gsell committed
645
	}
646 647

	strncpy ( f->groupname_step, name, len );
gsell's avatar
gsell committed
648
	f->stepno_width = (int)width;
649 650

	_H5Part_print_debug ( "Step name defined as '%s'", f->groupname_step );
gsell's avatar
gsell committed
651 652 653 654
	
	return H5PART_SUCCESS;
}

655 656 657 658 659
static h5part_int64_t
_set_num_particles (
	H5PartFile *f,				/*!< [in] Handle to open file */
	const h5part_int64_t nparticles,	/*!< [in] Number of particles */
	const h5part_int64_t _stride
gsell's avatar
gsell committed
660 661
	) {

662 663 664 665
	int ret;
	h5part_int64_t herr;

	hsize_t count;
gsell's avatar
gsell committed
666
#ifdef HDF5V160
667
	hssize_t start;
gsell's avatar
gsell committed
668
#else
669
	hsize_t start;
gsell's avatar
gsell committed
670
#endif
671 672
	hsize_t stride;
	hsize_t dmax = H5S_UNLIMITED;
gsell's avatar
gsell committed
673

674
#ifdef PARALLEL_IO
gsell's avatar
gsell committed
675 676 677 678
	hsize_t total;
	register int i;
#endif

679 680 681 682 683 684 685 686 687 688 689 690
	if ( nparticles <= 0 )
		return HANDLE_H5PART_INVALID_ERR ( "nparticles", nparticles );

	/* prevent invalid stride value */	
	if (_stride < 1)
	{
		_H5Part_print_warn (
			"Stride < 1 was specified: changing to 1." );
		stride = 1;
	} else {
		stride = (hsize_t) _stride;
	}
gsell's avatar
gsell committed
691 692 693 694 695 696 697 698

#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
	*/
699 700 701
	if ( f->nparticles == nparticles && stride == 1 ) {
		_H5Part_print_debug (
			"Serial mode: skipping unnecessary view creation" );
gsell's avatar
gsell committed
702 703 704
		return H5PART_SUCCESS;
	}
#endif
705 706 707 708

	herr = _reset_view ( f );
	if ( herr < 0 ) return herr;

709 710 711 712 713 714
	if ( f->shape != H5S_ALL ) {
		herr = H5Sclose ( f->shape );
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
		f->shape = H5S_ALL;
	}

715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736
	f->nparticles = (hsize_t) nparticles;

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

	/* we need a hyperslab selection if there is striding
	 * (otherwise, the default H5S_ALL selection is ok)
	 */
	if ( stride > 1 )
	{
		start = 0;
		count = f->nparticles;
		herr = H5Sselect_hyperslab (
			f->memshape,
			H5S_SELECT_SET,
			&start,
			&stride,
			&count, NULL );
		if ( herr < 0 ) return HANDLE_H5S_SELECT_HYPERSLAB_ERR;
gsell's avatar
gsell committed
737
	}
738

gsell's avatar
gsell committed
739
#ifndef PARALLEL_IO
740 741 742 743 744
	count = f->nparticles;
	f->shape = H5Screate_simple ( 1, &count, NULL );
	if ( f->shape < 0 ) HANDLE_H5S_CREATE_SIMPLE_ERR ( count );
	f->viewstart = 0;
	f->viewend   = nparticles - 1; // view range is *inclusive*
gsell's avatar
gsell committed
745 746 747 748 749 750 751 752 753 754 755 756 757 758 759
#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
	*/
760

761
	ret = MPI_Allgather (
gsell's avatar
gsell committed
762 763
		&nparticles, 1, MPI_LONG_LONG,
		f->pnparticles, 1, MPI_LONG_LONG,
764 765 766
		f->comm );
	if ( ret != MPI_SUCCESS) return HANDLE_MPI_ALLGATHER_ERR;

gsell's avatar
gsell committed
767
	if ( f->myproc == 0 ) {
gsell's avatar
gsell committed
768
		_H5Part_print_debug ( "Particle offsets:" );
769
		for ( i=0; i<f->nprocs; i++ ) 
770
			_H5Part_print_debug ( "\t[%d] np=%lld", i,
771
					(long long) f->pnparticles[i] );
gsell's avatar
gsell committed
772 773 774
	}

	/* compute start offsets */
775
	start = 0;
gsell's avatar
gsell committed
776
	for (i=0; i<f->myproc; i++) {
777
		start += f->pnparticles[i];
gsell's avatar
gsell committed
778
	}
779 780
	f->viewstart = start;
	f->viewend   = start + f->nparticles - 1; // view range is *inclusive*
gsell's avatar
gsell committed
781
	
782
	/* compute total nparticles */
gsell's avatar
gsell committed
783 784 785 786 787
	total = 0;
	for (i=0; i < f->nprocs; i++) {
		total += f->pnparticles[i];
	}

788
	/* declare overall datasize */
789
	count = total;
790 791
	f->shape = H5Screate_simple (1, &count, NULL);
	if ( f->shape < 0 ) return HANDLE_H5S_CREATE_SIMPLE_ERR ( count );
gsell's avatar
gsell committed
792 793

	/* declare overall data size  but then will select a subset */
794 795
	f->diskshape = H5Screate_simple (1, &count, NULL);
	if ( f->diskshape < 0 ) return HANDLE_H5S_CREATE_SIMPLE_ERR ( count );
gsell's avatar
gsell committed
796

797 798 799
	count = nparticles;
	stride = 1;
	herr = H5Sselect_hyperslab (
gsell's avatar
gsell committed
800 801
		f->diskshape,
		H5S_SELECT_SET,
802 803 804 805
		&start,
		&stride,
		&count, NULL );
	if ( herr < 0 ) return HANDLE_H5S_SELECT_HYPERSLAB_ERR;
gsell's avatar
gsell committed
806
#endif
807

gsell's avatar
gsell committed
808 809 810
	return H5PART_SUCCESS;
}

811 812 813 814 815 816 817
/*!
  \ingroup h5part_model

  Set the number of particles for the current time-step.
  After you call this subroutine, all subsequent 
  operations will assume this number of particles will be written.

818
  For the parallel library, the \c nparticles value is the number of
819 820 821 822 823 824 825 826 827 828
  particles that the \e individual task will write. You can use
  a different value on different tasks.
  This function uses an \c MPI_Allgather
  call to aggregate each tasks number of particles and determine
  the appropiate offsets. Because of the use of this MPI collective,
  it is advisable to call this function as
  few times as possible when running at large concurrency.

  This function assumes that your particles' data fields are in stored in
  contiguous 1D arrays.
829
  For instance, the fields \e x and \e y for your particles are stored
830 831 832
  in separate arrays \c x[] and \c y[].
  
  If instead you store your particles as tuples, so that the values
833
  are arranged \f$ x_1,y_1,x_2,y_2\f$... than you need to setup striding
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
  (in this case with value 2) using \ref H5PartSetNumParticlesStrided.

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

	SET_FNAME ( "H5PartSetNumParticles" );
	CHECK_FILEHANDLE( f );

	h5part_int64_t herr;
	h5part_int64_t stride = 1;

	herr = _set_num_particles ( f, nparticles, stride );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

/*!
  \ingroup h5part_model

  Set the number of particles for the current time-step.
  After you call this subroutine, all subsequent 
  operations will assume this number of particles will be written.

863
  For the parallel library, the \c nparticles value is the number of
864 865 866 867 868 869 870 871 872
  particles that the \e individual task will write. You can use
  a different value on different tasks.
  This function uses an \c MPI_Allgather
  call to aggregate each tasks number of particles and determine
  the appropiate offsets. Because of the use of this MPI collective,
  it is advisable to call this function as
  few times as possible when running at large concurrency.

  This function assumes that your particles' data fields are
873 874
  stored tuples. For instance, the fields \e x and \e y of your
  particles are arranged \f$x_1,y_1,x_2,y_2\f$... in a single data
875 876 877 878 879 880 881
  array. In this example, the stride value would be 2.
  
  If you instead have a separate array for each fields,
  such as \c x[] and \c y[],
  use \ref H5PartSetNumParticles.

  \return	\c H5PART_SUCCESS or error code
882
*/
883 884 885 886
h5part_int64_t
H5PartSetNumParticlesStrided (
	H5PartFile *f,				/*!< [in] Handle to open file */
	const h5part_int64_t nparticles,	/*!< [in] Number of particles */
887
	const h5part_int64_t stride		/*!< [in] Stride value (e.g. number of fields in the particle array) */
888 889 890 891 892 893 894 895 896 897 898 899 900
	) {

	SET_FNAME ( "H5PartSetNumParticlesStrided" );
	CHECK_FILEHANDLE( f );

	h5part_int64_t herr;

	herr = _set_num_particles ( f, nparticles, stride );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

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
/*!
  \ingroup h5part_model

  Define the chunk \c size and enables chunking in the underlying
  HDF5 layer. When combined with the \c align value in the
  \ref H5PartOpenFileAlign or \ref H5PartOpenFileParallelAlign
  function, this causes each group of \c size particles to be
  padded on disk out to the nearest multiple of \c align bytes.

  Note that this policy wastes disk space, but can improve write
  bandwidth on parallel filesystems that are sensitive to alignment
  to stripe boundaries (e.g. lustre).

  \return	\c H5PART_SUCCESS or error code
*/
h5part_int64_t
H5PartSetChunkSize (
	H5PartFile *f,
	const h5part_int64_t size
	) {

	SET_FNAME ( "H5PartSetChunkSize" );
	CHECK_FILEHANDLE( f );

	if ( f->myproc == 0 )
		_H5Part_print_info (
			"Setting chunk size to %lld elements",
			(long long)size );

	if ( f->dcreate_prop == H5P_DEFAULT ) {
		f->dcreate_prop = H5Pcreate (H5P_DATASET_CREATE);
		if ( f->dcreate_prop < 0 ) return HANDLE_H5P_CREATE_ERR;
	}

	hsize_t hsize = (hsize_t)size;

	herr_t herr = H5Pset_chunk ( f->dcreate_prop, 1, &hsize );
	if ( herr < 0 ) return HANDLE_H5P_SET_CHUNK_ERR;

	return H5PART_SUCCESS;
}

943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960
static void
_normalize_dataset_name (
	const char *name,
	char *name2
	) {

	if ( strlen(name) > H5PART_DATANAME_LEN ) {
		strncpy ( name2, name, H5PART_DATANAME_LEN - 1 );
		name2[H5PART_DATANAME_LEN-1] = '\0';
		_H5Part_print_warn (
			"Dataset name '%s' is longer than maximum %d chars. "
			"Truncated to: '%s'",
			name, H5PART_DATANAME_LEN, name2 );
	} else {
		strcpy ( name2, name );
	}
}

gsell's avatar
gsell committed
961 962 963
static h5part_int64_t
_write_data (
	H5PartFile *f,		/*!< IN: Handle to open file */
gsell's avatar
gsell committed
964 965 966
	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
967
	) {
968

gsell's avatar
gsell committed
969 970 971
	herr_t herr;
	hid_t dataset_id;

972 973 974 975 976
	char name2[H5PART_DATANAME_LEN];
	_normalize_dataset_name ( name, name2 );

	_H5Part_print_debug (
			"Create a dataset[%s] mounted on "
977
			"timestep %lld",
978
			name2, (long long)f->timestep );
gsell's avatar
gsell committed
979

980 981 982 983 984 985 986
	if ( f->shape == H5S_ALL ) {
		_H5Part_print_warn (
			"The view is unset or invalid: please "
			"set the view or specify a number of particles." );
		return HANDLE_H5PART_BAD_VIEW_ERR ( f->viewstart, f->viewend );
	}

987 988
	H5E_BEGIN_TRY
	dataset_id = H5Dopen ( f->timegroup, name2
989
#ifndef H5_USE_16_API
990 991 992 993
		, H5P_DEFAULT
#endif
		);
	H5E_END_TRY
994

995 996 997 998 999
	if ( dataset_id > 0 ) {
		_H5Part_print_warn (
			"Dataset[%s] at timestep %lld "
			"already exists", name2, (long long)f->timestep );
	} else {
1000
		dataset_id = H5Dcreate ( 
1001 1002 1003 1004
			f->timegroup,
			name2,
			type,
			f->shape,
1005
#ifndef H5_USE_16_API
1006
			H5P_DEFAULT,
1007 1008
			f->dcreate_prop,
			H5P_DEFAULT
1009
#else
1010
			f->dcreate_prop
1011
#endif
1012
                );
1013 1014 1015
		if ( dataset_id < 0 )
			return HANDLE_H5D_CREATE_ERR ( name2, f->timestep );
	}
gsell's avatar
gsell committed
1016

1017 1018 1019
#ifdef PARALLEL_IO
	herr = _H5Part_start_throttle ( f );
	if ( herr < 0 ) return herr;
1020
#endif
1021

1022 1023 1024 1025 1026
	herr = H5Dwrite (
		dataset_id,
		type,
		f->memshape,
		f->diskshape,
1027
		f->dxfer_prop,
1028 1029
		array );

1030 1031 1032
#ifdef PARALLEL_IO
	herr = _H5Part_end_throttle ( f );
	if ( herr < 0 ) return herr;
1033
#endif
1034

1035
	if ( herr < 0 ) return HANDLE_H5D_WRITE_ERR ( name2, f->timestep );
gsell's avatar
gsell committed
1036 1037 1038 1039

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

gsell's avatar
gsell committed
1040 1041
	f->empty = 0;

gsell's avatar
gsell committed
1042 1043 1044 1045
	return H5PART_SUCCESS;
}

/*!
1046
  \ingroup h5part_data
gsell's avatar
gsell committed
1047

gsell's avatar
gsell committed
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
  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
1074 1075
	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
1076 1077 1078
	) {

	SET_FNAME ( "H5PartWriteDataFloat64" );
gsell's avatar
gsell committed
1079
	h5part_int64_t herr;
gsell's avatar
gsell committed
1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090

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

1091
/*!
1092
  \ingroup h5part_data
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

  Write array of 32 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
H5PartWriteDataFloat32 (
	H5PartFile *f,		/*!< [in] Handle to open file */
	const char *name,	/*!< [in] Name to associate array with */
	const h5part_float32_t *array	/*!< [in] Array to commit to disk */
	) {

	SET_FNAME ( "H5PartWriteDataFloat32" );
	h5part_int64_t herr;

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

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

	return H5PART_SUCCESS;
}

gsell's avatar
gsell committed
1137
/*!
1138
  \ingroup h5part_data
gsell's avatar
gsell committed
1139

gsell's avatar
gsell committed
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
  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
1166 1167
	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
1168 1169
	) {

1170
	SET_FNAME ( "H5PartWriteDataInt64" );
gsell's avatar
gsell committed
1171

gsell's avatar
gsell committed
1172
	h5part_int64_t herr;
gsell's avatar
gsell committed
1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183

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

1184
/*!
1185
  \ingroup h5part_data
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

  Write array of 32 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
H5PartWriteDataInt32 (
	H5PartFile *f,		/*!< [in] Handle to open file */
	const char *name,	/*!< [in] Name to associate array with */
	const h5part_int32_t *array	/*!< [in] Array to commit to disk */
	) {

	SET_FNAME ( "H5PartWriteDataInt32" );

	h5part_int64_t herr;

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

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

	return H5PART_SUCCESS;
}

gsell's avatar
gsell committed
1231 1232 1233
/********************** reading and writing attribute ************************/

/********************** private functions to handle attributes ***************/
gsell's avatar
gsell committed
1234 1235 1236 1237 1238 1239

/*!
  \ingroup h5partkernel
  @{
*/

1240
h5part_int64_t
1241
_H5Part_make_string_type(
1242
	hid_t *stype,
1243 1244
	int size
	) {
1245 1246 1247 1248 1249
	*stype = H5Tcopy ( H5T_C_S1 );
	if ( *stype < 0 ) return HANDLE_H5T_STRING_ERR;
	herr_t herr = H5Tset_size ( *stype, size );
	if ( herr < 0 ) return HANDLE_H5T_STRING_ERR;
	return H5PART_SUCCESS;
1250 1251
}

gsell's avatar
gsell committed
1252 1253 1254
/*!
   Normalize HDF5 type
*/
1255
h5part_int64_t
gsell's avatar
gsell committed
1256 1257 1258
_H5Part_normalize_h5_type (
	hid_t type
	) {
1259

gsell's avatar
gsell committed
1260 1261 1262
	H5T_class_t tclass = H5Tget_class ( type );
	int size = H5Tget_size ( type );

1263
	switch ( tclass ) {
gsell's avatar
gsell committed
1264
	case H5T_INTEGER:
1265 1266 1267 1268 1269 1270 1271
		if ( size==8 ) {
			return H5PART_INT64;
		}
	  	else if ( size==1 ) {
			return H5PART_CHAR;
		}
		break;
gsell's avatar
gsell committed
1272
	case H5T_FLOAT:
1273 1274 1275 1276 1277 1278 1279
		if ( size==8 ) {
			return H5PART_FLOAT64;
		}
		else if ( size==4 ) {
			return H5PART_FLOAT32;
		}
		break;
1280
	case H5T_STRING:
1281
		return H5PART_STRING;
gsell's avatar
gsell committed
1282
	default:
1283
		;/* NOP */
gsell's avatar
gsell committed
1284
	}
1285
	
1286
	return HANDLE_H5PART_TYPE_ERR;
gsell's avatar
gsell committed
1287 1288 1289 1290 1291 1292 1293 1294 1295 1296
}

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

	herr_t herr;
1297
	h5part_int64_t h5err;
gsell's avatar
gsell committed
1298 1299 1300 1301
	hid_t attrib_id;
	hid_t space_id;
	hid_t type_id;

1302
#ifndef H5_USE_16_API
1303
	if (! H5Aexists ( id, attrib_name )) {
1304
		_H5Part_print_warn ( "Attribute '%s' does not exist!", attrib_name );
1305 1306 1307
	}
	attrib_id = H5Aopen ( id, attrib_name, H5P_DEFAULT );
#else
gsell's avatar
gsell committed
1308
	attrib_id = H5Aopen_name ( id, attrib_name );
1309
#endif
gsell's avatar
gsell committed
1310 1311
	if ( attrib_id <= 0 ) return HANDLE_H5A_OPEN_NAME_ERR( attrib_name );

1312
	type_id = H5Aget_type ( attrib_id );
1313
	if ( type_id < 0 ) return HANDLE_H5A_GET_TYPE_ERR;
gsell's avatar
gsell committed
1314 1315 1316 1317

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

1318
	herr = H5Aread ( attrib_id, type_id, attrib_value );
gsell's avatar
gsell committed
1319 1320 1321 1322 1323
	if ( herr < 0 ) return HANDLE_H5A_READ_ERR;

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

1324
	herr = H5Tclose ( type_id );
gsell's avatar
gsell committed
1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341
	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
	) {

1342
	h5part_int64_t herr;
gsell's avatar
gsell committed
1343 1344
	hid_t space_id;
	hid_t attrib_id;
1345 1346
	hid_t type = attrib_type;

1347 1348 1349 1350 1351
	if ( attrib_type == H5PART_STRING )
	{
		herr = _H5Part_make_string_type ( &type, attrib_nelem );
		if ( herr < 0 ) return herr;

1352 1353 1354 1355 1356
		space_id = H5Screate (H5S_SCALAR);
		if ( space_id < 0 )
			return HANDLE_H5S_CREATE_SCALAR_ERR;
	}
	else {
1357
		space_id = H5Screate_simple ( 1, &attrib_nelem, NULL );
1358 1359 1360
		if ( space_id < 0 )
			return HANDLE_H5S_CREATE_SIMPLE_ERR ( attrib_nelem );
	}
gsell's avatar
gsell committed
1361

1362 1363 1364 1365 1366 1367
	attrib_id = H5Acreate(
		id,
		attrib_name,
		type,
		space_id,
		H5P_DEFAULT
1368
#ifndef H5_USE_16_API
1369
		, H5P_DEFAULT
1370
#endif
1371
		);
1372

gsell's avatar
gsell committed
1373 1374
	if ( attrib_id < 0 ) return HANDLE_H5A_CREATE_ERR ( attrib_name );

1375
	herr = H5Awrite ( attrib_id, type, attrib_value);
gsell's avatar
gsell committed
1376 1377 1378 1379 1380 1381 1382 1383
	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;

1384
	if ( attrib_type == H5PART_STRING ) {
1385 1386 1387 1388
		herr = H5Tclose ( type );
		if ( herr < 0 ) return HANDLE_H5T_CLOSE_ERR;
	}

gsell's avatar
gsell committed
1389 1390 1391
	return H5PART_SUCCESS;
}

1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444
h5part_int64_t
_H5Part_write_file_attrib (
	H5PartFile *f,
	const char *name,
	const hid_t type,
	const void *value,
	const hsize_t nelem
	) {

	hid_t group_id = H5Gopen ( f->file, "/"
#ifndef H5_USE_16_API
				  , H5P_DEFAULT
#endif
				  );

	if ( group_id < 0 ) return HANDLE_H5G_OPEN_ERR( "/" );

	h5part_int64_t herr = _H5Part_write_attrib (
		group_id,
		name,
		type,
		value,
		nelem );
	if ( herr < 0 ) return herr;

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

	return H5PART_SUCCESS;
}

h5part_int64_t
_H5Part_write_step_attrib (
	H5PartFile *f,
	const char *name,
	const hid_t type,
	const void *value,
	const hsize_t nelem
	) {

	CHECK_TIMEGROUP( f );

	h5part_int64_t herr = _H5Part_write_attrib (
		f->timegroup,
		name,
		type,
		value,
		nelem );
	if ( herr < 0 ) return herr;

	return H5PART_SUCCESS;
}

gsell's avatar
gsell committed
1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459
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;

gsell's avatar
gsell committed
1460
	attrib_id = H5Aopen_idx ( id, (unsigned int)attrib_idx );
gsell's avatar
gsell committed
1461 1462 1463
	if ( attrib_id < 0 ) return HANDLE_H5A_OPEN_IDX_ERR ( attrib_idx );

	if ( attrib_nelem ) {
gsell's avatar
gsell committed
1464 1465 1466
		space_id =  H5Aget_space ( attrib_id );
		if ( space_id < 0 ) return HANDLE_H5A_GET_SPACE_ERR;

gsell's avatar
gsell committed
1467 1468 1469
		*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
1470 1471 1472

		herr = H5Sclose ( space_id );
		if ( herr < 0 ) return HANDLE_H5S_CLOSE_ERR;
gsell's avatar
gsell committed
1473 1474 1475 1476
	}
	if ( attrib_name ) {
		herr = H5Aget_name (
			attrib_id,
gsell's avatar
gsell committed
1477
			(size_t)len_attrib_name,
gsell's avatar
gsell committed
1478 1479 1480 1481
			attrib_name );
		if ( herr < 0 ) return HANDLE_H5A_GET_NAME_ERR;
	}
	if ( attrib_type ) {
gsell's avatar
gsell committed
1482 1483
		mytype = H5Aget_type ( attrib_id );
		if ( mytype < 0 ) return HANDLE_H5A_GET_TYPE_ERR;
gsell's avatar
gsell committed
1484

gsell's avatar
gsell committed
1485
		*attrib_type = _H5Part_normalize_h5_type ( mytype );
Marc Howison's avatar