1998-08-17 23:13:45 +08:00
|
|
|
|
/*
|
2001-04-06 01:29:14 +08:00
|
|
|
|
* Copyright (C) 1998-2001 NCSA
|
|
|
|
|
* All rights reserved.
|
1998-08-17 23:13:45 +08:00
|
|
|
|
*
|
|
|
|
|
* Programmer: rky 980813
|
|
|
|
|
*
|
|
|
|
|
* Purpose: Functions to read/write directly between app buffer and file.
|
|
|
|
|
*
|
|
|
|
|
* Beware of the ifdef'ed print statements.
|
|
|
|
|
* I didn't make them portable.
|
|
|
|
|
*/
|
|
|
|
|
|
2000-10-11 04:55:54 +08:00
|
|
|
|
#define H5F_PACKAGE /*suppress error about including H5Fpkg */
|
2000-10-10 15:43:38 +08:00
|
|
|
|
#define H5S_PACKAGE /*suppress error about including H5Spkg */
|
|
|
|
|
|
2001-04-06 01:29:14 +08:00
|
|
|
|
#include "H5private.h"
|
|
|
|
|
#include "H5Eprivate.h"
|
|
|
|
|
#include "H5Fpkg.h" /* Ugly, but necessary for the MPIO I/O accesses */
|
2001-07-20 23:01:18 +08:00
|
|
|
|
#include "H5FDprivate.h" /* Necessary for the H5FD_write & H5FD_read prototypes.. */
|
2002-06-05 05:29:05 +08:00
|
|
|
|
#include "H5Iprivate.h" /* Object IDs */
|
2001-08-11 06:30:01 +08:00
|
|
|
|
#include "H5Pprivate.h" /* Property Lists */
|
2001-04-06 01:29:14 +08:00
|
|
|
|
#include "H5Spkg.h"
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
2001-04-06 01:29:14 +08:00
|
|
|
|
#include "H5FDmpio.h" /*the MPIO file driver */
|
1999-08-11 04:21:32 +08:00
|
|
|
|
|
1999-12-17 22:37:22 +08:00
|
|
|
|
#ifndef H5_HAVE_PARALLEL
|
1998-08-19 22:07:14 +08:00
|
|
|
|
/*
|
|
|
|
|
* The H5S_mpio_xxxx functions are for parallel I/O only and are
|
1999-12-17 22:37:22 +08:00
|
|
|
|
* valid only when H5_HAVE_PARALLEL is #defined. This empty #ifndef
|
1998-08-19 22:07:14 +08:00
|
|
|
|
* body is used to allow this source file be included in the serial
|
|
|
|
|
* distribution.
|
|
|
|
|
* Some compilers/linkers may complain about "empty" object file.
|
|
|
|
|
* If that happens, uncomment the following statement to pacify
|
|
|
|
|
* them.
|
|
|
|
|
*/
|
|
|
|
|
/* const hbool_t H5S_mpio_avail = FALSE; */
|
1999-12-17 22:37:22 +08:00
|
|
|
|
#else /* H5_HAVE_PARALLEL */
|
1998-08-17 23:13:45 +08:00
|
|
|
|
/* Interface initialization */
|
1999-04-15 05:48:05 +08:00
|
|
|
|
#define PABLO_MASK H5Sall_mask
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#define INTERFACE_INIT NULL
|
2001-08-15 06:09:56 +08:00
|
|
|
|
static int interface_initialize_g = 0;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
1998-09-01 05:30:31 +08:00
|
|
|
|
static herr_t
|
|
|
|
|
H5S_mpio_all_type( const H5S_t *space, const size_t elmt_size,
|
|
|
|
|
/* out: */
|
|
|
|
|
MPI_Datatype *new_type,
|
2001-07-20 23:01:18 +08:00
|
|
|
|
size_t *count,
|
1998-09-01 05:30:31 +08:00
|
|
|
|
hbool_t *is_derived_type );
|
|
|
|
|
static herr_t
|
|
|
|
|
H5S_mpio_hyper_type( const H5S_t *space, const size_t elmt_size,
|
|
|
|
|
/* out: */
|
|
|
|
|
MPI_Datatype *new_type,
|
2001-07-20 23:01:18 +08:00
|
|
|
|
size_t *count,
|
1998-09-01 05:30:31 +08:00
|
|
|
|
hbool_t *is_derived_type );
|
|
|
|
|
static herr_t
|
|
|
|
|
H5S_mpio_space_type( const H5S_t *space, const size_t elmt_size,
|
|
|
|
|
/* out: */
|
|
|
|
|
MPI_Datatype *new_type,
|
2001-07-20 23:01:18 +08:00
|
|
|
|
size_t *count,
|
1998-09-01 05:30:31 +08:00
|
|
|
|
hbool_t *is_derived_type );
|
|
|
|
|
static herr_t
|
2002-04-26 01:56:56 +08:00
|
|
|
|
H5S_mpio_spaces_xfer(H5F_t *f, const H5O_layout_t *layout,
|
|
|
|
|
H5P_genplist_t *dc_plist, size_t elmt_size,
|
1998-09-01 05:30:31 +08:00
|
|
|
|
const H5S_t *file_space, const H5S_t *mem_space,
|
2002-04-04 01:07:14 +08:00
|
|
|
|
hid_t dxpl_id, void *buf/*out*/, const hbool_t do_write);
|
1998-09-01 05:30:31 +08:00
|
|
|
|
|
1998-08-17 23:13:45 +08:00
|
|
|
|
/*-------------------------------------------------------------------------
|
|
|
|
|
* Function: H5S_mpio_all_type
|
|
|
|
|
*
|
|
|
|
|
* Purpose: Translate an HDF5 "all" selection into an MPI type.
|
|
|
|
|
*
|
|
|
|
|
* Return: non-negative on success, negative on failure.
|
|
|
|
|
*
|
|
|
|
|
* Outputs: *new_type the MPI type corresponding to the selection
|
|
|
|
|
* *count how many objects of the new_type in selection
|
|
|
|
|
* (useful if this is the buffer type for xfer)
|
|
|
|
|
* *is_derived_type 0 if MPI primitive type, 1 if derived
|
|
|
|
|
*
|
|
|
|
|
* Programmer: rky 980813
|
|
|
|
|
*
|
|
|
|
|
* Modifications:
|
|
|
|
|
*
|
|
|
|
|
*-------------------------------------------------------------------------
|
|
|
|
|
*/
|
2002-05-29 02:17:12 +08:00
|
|
|
|
static herr_t
|
1998-08-28 11:59:23 +08:00
|
|
|
|
H5S_mpio_all_type( const H5S_t *space, const size_t elmt_size,
|
1998-08-17 23:13:45 +08:00
|
|
|
|
/* out: */
|
|
|
|
|
MPI_Datatype *new_type,
|
2001-07-20 23:01:18 +08:00
|
|
|
|
size_t *count,
|
1998-08-17 23:13:45 +08:00
|
|
|
|
hbool_t *is_derived_type )
|
|
|
|
|
{
|
|
|
|
|
hsize_t total_bytes;
|
2001-07-20 23:01:18 +08:00
|
|
|
|
unsigned u;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
2002-05-29 02:17:12 +08:00
|
|
|
|
FUNC_ENTER_NOINIT(H5S_mpio_all_type);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* Check args */
|
|
|
|
|
assert (space);
|
|
|
|
|
|
|
|
|
|
/* Just treat the entire extent as a block of bytes */
|
1998-08-28 11:59:23 +08:00
|
|
|
|
total_bytes = (hsize_t)elmt_size;
|
2001-07-20 23:01:18 +08:00
|
|
|
|
for (u=0; u<space->extent.u.simple.rank; ++u)
|
|
|
|
|
total_bytes *= space->extent.u.simple.size[u];
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* fill in the return values */
|
|
|
|
|
*new_type = MPI_BYTE;
|
2001-07-20 23:01:18 +08:00
|
|
|
|
H5_CHECK_OVERFLOW(total_bytes, hsize_t, size_t);
|
2001-11-28 00:29:13 +08:00
|
|
|
|
*count = (size_t)total_bytes;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
*is_derived_type = 0;
|
|
|
|
|
|
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-06-20 05:49:01 +08:00
|
|
|
|
HDfprintf(stdout, "Leave %s total_bytes=%Hu\n", FUNC, total_bytes );
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#endif
|
|
|
|
|
FUNC_LEAVE (SUCCEED);
|
1999-08-11 04:21:32 +08:00
|
|
|
|
}
|
|
|
|
|
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/*-------------------------------------------------------------------------
|
|
|
|
|
* Function: H5S_mpio_hyper_type
|
|
|
|
|
*
|
|
|
|
|
* Purpose: Translate an HDF5 hyperslab selection into an MPI type.
|
|
|
|
|
*
|
|
|
|
|
* Return: non-negative on success, negative on failure.
|
|
|
|
|
*
|
|
|
|
|
* Outputs: *new_type the MPI type corresponding to the selection
|
|
|
|
|
* *count how many objects of the new_type in selection
|
|
|
|
|
* (useful if this is the buffer type for xfer)
|
|
|
|
|
* *is_derived_type 0 if MPI primitive type, 1 if derived
|
|
|
|
|
*
|
|
|
|
|
* Programmer: rky 980813
|
|
|
|
|
*
|
1999-07-01 00:30:56 +08:00
|
|
|
|
* Modifications: ppw 990401
|
[svn-r2646] Purpose:
Bug fix (done by Kim Yates)
Description:
The optimized mpio code was broken and when read was done, it hanged.
Solution:
H5FDmpio.c:
In H5FD_mpio_write, moved the 16-line block of code in which
all procs other than p0 skip the actual write
to be just before the call to MPI_File_write_at.
Previously, the values of the local vars that controlled
"allsame" were not always set correctly when the moved block
was reached.
H5S.c:
Changed default value of H5_mpi_opt_types_g to TRUE, so that
the MPI-IO hyperslab code is executed by default in parallel HDF5,
rather than executing the serial hyperslab code.
H5Smpio.c:
In function H5S_mpio_hyper_type, added a call to free
an intermediate type. Cures a small memory leak.
Added code for cases of empty hyperslab
Changed displacements to be MPI_Aint
Platforms tested:
modi4 -64: worked fine with mpich 1.2.0 but failed with the messages
saying it ran out of entries for MPI_Types during the collective_read
test. After tracing the code all the way to the collective read, all
MPI Types have been freed properly. It aborted with the above message
when it executed the line
if (MPI_SUCCESS!= MPI_File_read_at_all(file->f, mpi_off, buf, size_i, buf_type, &mpi_stat ))
Could not see any problem with this line. It could be a bug in the
SGI version of MPI.
2000-10-10 10:32:45 +08:00
|
|
|
|
* rky, ppw 2000-09-26 Freed old type after creating struct type.
|
|
|
|
|
* rky 2000-10-05 Changed displacements to be MPI_Aint.
|
|
|
|
|
* rky 2000-10-06 Added code for cases of empty hyperslab.
|
2000-11-17 12:48:58 +08:00
|
|
|
|
* akc, rky 2000-11-16 Replaced hard coded dimension size with
|
|
|
|
|
* H5S_MAX_RANK.
|
1998-08-17 23:13:45 +08:00
|
|
|
|
*
|
|
|
|
|
*-------------------------------------------------------------------------
|
|
|
|
|
*/
|
2002-05-29 02:17:12 +08:00
|
|
|
|
static herr_t
|
1998-08-28 11:59:23 +08:00
|
|
|
|
H5S_mpio_hyper_type( const H5S_t *space, const size_t elmt_size,
|
1998-08-17 23:13:45 +08:00
|
|
|
|
/* out: */
|
|
|
|
|
MPI_Datatype *new_type,
|
2001-07-20 23:01:18 +08:00
|
|
|
|
size_t *count,
|
1998-08-17 23:13:45 +08:00
|
|
|
|
hbool_t *is_derived_type )
|
|
|
|
|
{
|
|
|
|
|
struct dim { /* less hassle than malloc/free & ilk */
|
2001-07-20 23:01:18 +08:00
|
|
|
|
hssize_t start;
|
|
|
|
|
hsize_t strid;
|
|
|
|
|
hsize_t block;
|
|
|
|
|
hsize_t xtent;
|
|
|
|
|
hsize_t count;
|
2000-11-17 12:48:58 +08:00
|
|
|
|
} d[H5S_MAX_RANK];
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
1999-08-11 04:21:32 +08:00
|
|
|
|
int i, err, new_rank, num_to_collapse;
|
2001-12-13 22:04:47 +08:00
|
|
|
|
herr_t ret_value = SUCCEED;
|
2000-11-17 12:48:58 +08:00
|
|
|
|
int offset[H5S_MAX_RANK];
|
|
|
|
|
int max_xtent[H5S_MAX_RANK];
|
1998-08-17 23:13:45 +08:00
|
|
|
|
H5S_hyper_dim_t *diminfo; /* [rank] */
|
2001-08-15 06:09:56 +08:00
|
|
|
|
int rank;
|
2000-11-17 12:48:58 +08:00
|
|
|
|
int block_length[2];
|
|
|
|
|
MPI_Datatype inner_type, outer_type, old_type[2];
|
[svn-r2646] Purpose:
Bug fix (done by Kim Yates)
Description:
The optimized mpio code was broken and when read was done, it hanged.
Solution:
H5FDmpio.c:
In H5FD_mpio_write, moved the 16-line block of code in which
all procs other than p0 skip the actual write
to be just before the call to MPI_File_write_at.
Previously, the values of the local vars that controlled
"allsame" were not always set correctly when the moved block
was reached.
H5S.c:
Changed default value of H5_mpi_opt_types_g to TRUE, so that
the MPI-IO hyperslab code is executed by default in parallel HDF5,
rather than executing the serial hyperslab code.
H5Smpio.c:
In function H5S_mpio_hyper_type, added a call to free
an intermediate type. Cures a small memory leak.
Added code for cases of empty hyperslab
Changed displacements to be MPI_Aint
Platforms tested:
modi4 -64: worked fine with mpich 1.2.0 but failed with the messages
saying it ran out of entries for MPI_Types during the collective_read
test. After tracing the code all the way to the collective read, all
MPI Types have been freed properly. It aborted with the above message
when it executed the line
if (MPI_SUCCESS!= MPI_File_read_at_all(file->f, mpi_off, buf, size_i, buf_type, &mpi_stat ))
Could not see any problem with this line. It could be a bug in the
SGI version of MPI.
2000-10-10 10:32:45 +08:00
|
|
|
|
MPI_Aint extent_len, displacement[2];
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
2002-05-29 02:17:12 +08:00
|
|
|
|
FUNC_ENTER_NOINIT(H5S_mpio_hyper_type);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* Check and abbreviate args */
|
|
|
|
|
assert (space);
|
[svn-r2646] Purpose:
Bug fix (done by Kim Yates)
Description:
The optimized mpio code was broken and when read was done, it hanged.
Solution:
H5FDmpio.c:
In H5FD_mpio_write, moved the 16-line block of code in which
all procs other than p0 skip the actual write
to be just before the call to MPI_File_write_at.
Previously, the values of the local vars that controlled
"allsame" were not always set correctly when the moved block
was reached.
H5S.c:
Changed default value of H5_mpi_opt_types_g to TRUE, so that
the MPI-IO hyperslab code is executed by default in parallel HDF5,
rather than executing the serial hyperslab code.
H5Smpio.c:
In function H5S_mpio_hyper_type, added a call to free
an intermediate type. Cures a small memory leak.
Added code for cases of empty hyperslab
Changed displacements to be MPI_Aint
Platforms tested:
modi4 -64: worked fine with mpich 1.2.0 but failed with the messages
saying it ran out of entries for MPI_Types during the collective_read
test. After tracing the code all the way to the collective read, all
MPI Types have been freed properly. It aborted with the above message
when it executed the line
if (MPI_SUCCESS!= MPI_File_read_at_all(file->f, mpi_off, buf, size_i, buf_type, &mpi_stat ))
Could not see any problem with this line. It could be a bug in the
SGI version of MPI.
2000-10-10 10:32:45 +08:00
|
|
|
|
assert(sizeof(MPI_Aint) >= sizeof(elmt_size));
|
1998-11-11 02:04:07 +08:00
|
|
|
|
diminfo = space->select.sel_info.hslab.diminfo;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
assert (diminfo);
|
|
|
|
|
rank = space->extent.u.simple.rank;
|
|
|
|
|
assert (rank >= 0);
|
2001-07-20 23:01:18 +08:00
|
|
|
|
if (0==rank)
|
|
|
|
|
goto empty;
|
|
|
|
|
if (0==elmt_size)
|
|
|
|
|
goto empty;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* make a local copy of the dimension info so we can transform them */
|
2000-11-17 21:42:39 +08:00
|
|
|
|
assert(rank<=H5S_MAX_RANK); /* within array bounds */
|
1998-08-17 23:13:45 +08:00
|
|
|
|
for ( i=0; i<rank; ++i) {
|
2001-07-20 23:01:18 +08:00
|
|
|
|
d[i].start = diminfo[i].start;
|
|
|
|
|
d[i].strid = diminfo[i].stride;
|
|
|
|
|
d[i].block = diminfo[i].block;
|
|
|
|
|
d[i].count = diminfo[i].count;
|
|
|
|
|
d[i].xtent = space->extent.u.simple.size[i];
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-07-20 23:01:18 +08:00
|
|
|
|
HDfprintf(stdout, "hyper_type: start=%Hd stride=%Hu count=%Hu "
|
|
|
|
|
"block=%Hu xtent=%Hu",
|
|
|
|
|
d[i].start, d[i].strid, d[i].count, d[i].block, d[i].xtent );
|
|
|
|
|
if (i==0)
|
|
|
|
|
HDfprintf(stdout, " rank=%d\n", rank );
|
|
|
|
|
else
|
|
|
|
|
HDfprintf(stdout, "\n" );
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#endif
|
2001-07-20 23:01:18 +08:00
|
|
|
|
if (0==d[i].block)
|
|
|
|
|
goto empty;
|
|
|
|
|
if (0==d[i].count)
|
|
|
|
|
goto empty;
|
|
|
|
|
if (0==d[i].xtent)
|
|
|
|
|
goto empty;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
}
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
|
|
|
|
/**********************************************************************
|
|
|
|
|
Compute array "offset[rank]" which gives the offsets for a multi-
|
|
|
|
|
dimensional array with dimensions "d[i].xtent" (i=0,1,...,rank-1).
|
|
|
|
|
**********************************************************************/
|
|
|
|
|
offset[rank-1] = 1;
|
|
|
|
|
max_xtent[rank-1] = d[rank-1].xtent;
|
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-07-20 23:01:18 +08:00
|
|
|
|
i=rank-1;
|
|
|
|
|
HDfprintf(stdout, " offset[%2d]=%d; max_xtent[%2d]=%d\n",
|
1999-07-01 00:30:56 +08:00
|
|
|
|
i, offset[i], i, max_xtent[i]);
|
|
|
|
|
#endif
|
|
|
|
|
for (i=rank-2; i>=0; --i) {
|
|
|
|
|
offset[i] = offset[i+1]*d[i+1].xtent;
|
|
|
|
|
max_xtent[i] = max_xtent[i+1]*d[i].xtent;
|
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-07-20 23:01:18 +08:00
|
|
|
|
HDfprintf(stdout, " offset[%2d]=%d; max_xtent[%2d]=%d\n",
|
1999-07-01 00:30:56 +08:00
|
|
|
|
i, offset[i], i, max_xtent[i]);
|
|
|
|
|
#endif
|
|
|
|
|
}
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* Create a type covering the selected hyperslab.
|
|
|
|
|
* Multidimensional dataspaces are stored in row-major order.
|
|
|
|
|
* The type is built from the inside out, going from the
|
|
|
|
|
* fastest-changing (i.e., inner) dimension * to the slowest (outer). */
|
|
|
|
|
|
|
|
|
|
/* Optimization: check for contiguous inner dimensions.
|
|
|
|
|
* Supposing the dimensions were numbered from 1 to rank, we find that
|
|
|
|
|
*
|
|
|
|
|
* dim d=rank is contiguous if:
|
|
|
|
|
* stride[d] = block[d]
|
|
|
|
|
* and count[d] * block[d] = extent.u.simple.size[d]
|
|
|
|
|
*
|
|
|
|
|
* (i.e., there's no overlap or gaps and the entire extent is filled.)
|
|
|
|
|
*
|
|
|
|
|
* dim d (1<=d<rank) is contiguous if:
|
|
|
|
|
* dim d+1 is contiguous
|
|
|
|
|
* and stride[d] = block[d]
|
|
|
|
|
* and count[d] * block[d] = extent.u.simple.size[d]
|
|
|
|
|
*
|
|
|
|
|
* There is also a weak sense in which the first noncollapsible dim
|
|
|
|
|
* is contiguous if it consists of a single unbroken range,
|
|
|
|
|
* and we also take advantage of that.
|
|
|
|
|
*/
|
|
|
|
|
|
|
|
|
|
/* figure out how many dimensions we can eliminate */
|
|
|
|
|
/* This loop examines contiguity from the inside out. */
|
|
|
|
|
for ( i=0; i<rank; ++i) {
|
1998-08-28 11:59:23 +08:00
|
|
|
|
if ((d[rank-i-1].strid != d[rank-i-1].block)
|
2001-07-20 23:01:18 +08:00
|
|
|
|
|| (d[rank-i-1].count*d[rank-i-1].block) != space->extent.u.simple.size[rank-i-1])
|
|
|
|
|
break;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
} /* end for */
|
2001-07-20 23:01:18 +08:00
|
|
|
|
|
1998-08-28 11:59:23 +08:00
|
|
|
|
num_to_collapse = i;
|
2001-07-20 23:01:18 +08:00
|
|
|
|
if (num_to_collapse == rank)
|
|
|
|
|
num_to_collapse--;
|
1998-08-28 11:59:23 +08:00
|
|
|
|
|
1998-08-17 23:13:45 +08:00
|
|
|
|
assert(0<=num_to_collapse && num_to_collapse<rank);
|
|
|
|
|
new_rank = rank - num_to_collapse;
|
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-06-20 05:49:01 +08:00
|
|
|
|
HDfprintf(stdout, "num_to_collapse=%d\n", num_to_collapse);
|
|
|
|
|
HDfprintf(stdout, "hyper_type: new_rank=%d\n", new_rank );
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#endif
|
|
|
|
|
|
1998-08-28 11:59:23 +08:00
|
|
|
|
/* To collapse dims, just transform dimension info (from inner to outer) */
|
1999-07-01 00:30:56 +08:00
|
|
|
|
/* must multiply r.h.s. by "count" --- ppw 03/11/99 */
|
1998-08-28 11:59:23 +08:00
|
|
|
|
for (i=rank-1; i>=new_rank; --i) {
|
2001-07-20 23:01:18 +08:00
|
|
|
|
d[i-1].block *= d[i].count * d[i].block;
|
|
|
|
|
d[i-1].strid *= d[i].count * d[i].block;
|
|
|
|
|
d[i-1].xtent *= d[i].count * d[i].block;
|
|
|
|
|
assert( d[i].start == 0 );
|
|
|
|
|
/* d[i-1].start stays unchanged */
|
|
|
|
|
/* d[i-1].count stays unchanged */
|
1998-08-17 23:13:45 +08:00
|
|
|
|
}
|
1998-08-28 11:59:23 +08:00
|
|
|
|
|
1998-08-27 01:51:19 +08:00
|
|
|
|
/* check for possibility to coalesce blocks of the uncoalesced dimensions */
|
|
|
|
|
for (i=0; i<new_rank; ++i) {
|
|
|
|
|
if (d[i].strid == d[i].block) {
|
2001-07-20 23:01:18 +08:00
|
|
|
|
/* transform smaller blocks to 1 larger block of combined size */
|
|
|
|
|
d[i].strid = d[i].block *= d[i].count;
|
|
|
|
|
d[i].count = 1;
|
|
|
|
|
}
|
1998-08-17 23:13:45 +08:00
|
|
|
|
}
|
|
|
|
|
|
1999-07-01 00:30:56 +08:00
|
|
|
|
/*******************************************************
|
|
|
|
|
* Construct contig type for inner contig dims:
|
|
|
|
|
*******************************************************/
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-06-20 05:49:01 +08:00
|
|
|
|
HDfprintf(stdout, "hyper_type: Making contig type %d MPI_BYTEs\n", elmt_size );
|
2001-07-20 23:01:18 +08:00
|
|
|
|
for (i=new_rank-1; i>=0; --i)
|
|
|
|
|
HDfprintf(stdout, "d[%d].xtent=%Hu \n", i, d[i].xtent);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#endif
|
1998-08-28 11:59:23 +08:00
|
|
|
|
err = MPI_Type_contiguous( (int)elmt_size, MPI_BYTE, &inner_type );
|
2001-07-20 23:01:18 +08:00
|
|
|
|
if (err)
|
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_MPI, FAIL,"couldn't create MPI contiguous type");
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
1999-07-01 00:30:56 +08:00
|
|
|
|
/*******************************************************
|
|
|
|
|
* Construct the type by walking the hyperslab dims
|
|
|
|
|
* from the inside out:
|
|
|
|
|
*******************************************************/
|
1998-08-17 23:13:45 +08:00
|
|
|
|
for ( i=new_rank-1; i>=0; --i) {
|
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-07-20 23:01:18 +08:00
|
|
|
|
HDfprintf(stdout, "hyper_type: Dimension i=%d \n"
|
|
|
|
|
"count=%Hu block=%Hu stride=%Hu\n",
|
|
|
|
|
i, d[i].count, d[i].block, d[i].strid );
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#endif
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-07-20 23:01:18 +08:00
|
|
|
|
HDfprintf(stdout, "hyper_type: i=%d Making vector-type \n", i);
|
1999-07-01 00:30:56 +08:00
|
|
|
|
#endif
|
|
|
|
|
/****************************************
|
|
|
|
|
* Build vector in current dimension:
|
|
|
|
|
****************************************/
|
2001-07-20 23:01:18 +08:00
|
|
|
|
err = MPI_Type_vector ( (int)(d[i].count), /* count */
|
|
|
|
|
(int)(d[i].block), /* blocklength */
|
2000-11-17 12:48:58 +08:00
|
|
|
|
(int)(d[i].strid), /* stride */
|
1999-07-01 00:30:56 +08:00
|
|
|
|
inner_type, /* old type */
|
|
|
|
|
&outer_type ); /* new type */
|
|
|
|
|
|
2001-07-20 23:01:18 +08:00
|
|
|
|
MPI_Type_free( &inner_type );
|
|
|
|
|
if (err)
|
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_MPI, FAIL,"couldn't create MPI vector type");
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
2001-07-20 23:01:18 +08:00
|
|
|
|
displacement[1] = (MPI_Aint)elmt_size * max_xtent[i];
|
|
|
|
|
err = MPI_Type_extent(outer_type, &extent_len);
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
|
|
|
|
/*************************************************
|
|
|
|
|
* Restructure this datatype ("outer_type")
|
|
|
|
|
* so that it still starts at 0, but its extent
|
|
|
|
|
* is the full extent in this dimension.
|
|
|
|
|
*************************************************/
|
2001-07-20 23:01:18 +08:00
|
|
|
|
if ((int)extent_len < displacement[1]) {
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-07-20 23:01:18 +08:00
|
|
|
|
HDfprintf(stdout, "hyper_type: i=%d Extending struct type\n"
|
|
|
|
|
"***displacements: 0, %d\n", i, displacement[1]);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#endif
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
2000-11-11 23:58:12 +08:00
|
|
|
|
#ifdef H5_HAVE_MPI2 /* have MPI-2 */
|
1999-07-01 00:30:56 +08:00
|
|
|
|
err = MPI_Type_create_resized
|
|
|
|
|
( outer_type, /* old type */
|
|
|
|
|
0, /* blocklengths */
|
|
|
|
|
displacement[1], /* displacements */
|
|
|
|
|
&inner_type); /* new type */
|
|
|
|
|
#else /* do not have MPI-2 */
|
|
|
|
|
block_length[0] = 1;
|
|
|
|
|
block_length[1] = 1;
|
|
|
|
|
|
|
|
|
|
displacement[0] = 0;
|
|
|
|
|
|
|
|
|
|
old_type[0] = outer_type;
|
|
|
|
|
old_type[1] = MPI_UB;
|
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-07-20 23:01:18 +08:00
|
|
|
|
HDfprintf(stdout, "hyper_type: i=%d Extending struct type\n"
|
|
|
|
|
"***displacements: 0, %d\n", i, displacement[1]);
|
1999-07-01 00:30:56 +08:00
|
|
|
|
#endif
|
|
|
|
|
err = MPI_Type_struct ( 2, /* count */
|
|
|
|
|
block_length, /* blocklengths */
|
|
|
|
|
displacement, /* displacements */
|
|
|
|
|
old_type, /* old types */
|
|
|
|
|
&inner_type); /* new type */
|
|
|
|
|
#endif
|
|
|
|
|
|
|
|
|
|
MPI_Type_free (&outer_type);
|
2001-07-20 23:01:18 +08:00
|
|
|
|
if (err)
|
1999-07-01 00:30:56 +08:00
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_MPI, FAIL,"couldn't resize MPI vector type");
|
2001-07-20 23:01:18 +08:00
|
|
|
|
}
|
|
|
|
|
else {
|
1999-07-01 00:30:56 +08:00
|
|
|
|
inner_type = outer_type;
|
2001-07-20 23:01:18 +08:00
|
|
|
|
}
|
1998-08-17 23:13:45 +08:00
|
|
|
|
} /* end for */
|
1999-07-01 00:30:56 +08:00
|
|
|
|
/***************************
|
|
|
|
|
* End of loop, walking
|
|
|
|
|
* thru dimensions.
|
|
|
|
|
***************************/
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
/* At this point inner_type is actually the outermost type, even for 0-trip loop */
|
|
|
|
|
|
|
|
|
|
/***************************************************************
|
|
|
|
|
* Final task: create a struct which is a "clone" of the
|
|
|
|
|
* current struct, but displaced according to the d[i].start
|
|
|
|
|
* values given in the hyperslab description:
|
|
|
|
|
***************************************************************/
|
|
|
|
|
displacement[0] = 0;
|
|
|
|
|
for (i=new_rank-1; i>=0; i--)
|
2001-07-20 23:01:18 +08:00
|
|
|
|
displacement[0] += d[i].start * offset[i];
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
|
|
|
|
if (displacement[0] > 0) {
|
2001-07-20 23:01:18 +08:00
|
|
|
|
displacement[0] *= elmt_size;
|
|
|
|
|
block_length[0] = 1;
|
|
|
|
|
old_type[0] = inner_type;
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
|
|
|
|
#ifdef H5Smpi_DEBUG
|
2001-07-20 23:01:18 +08:00
|
|
|
|
HDfprintf(stdout, "hyper_type: Making final struct\n***count=1:\n");
|
1999-07-01 00:30:56 +08:00
|
|
|
|
HDfprintf(stdout, "\tblocklength[0]=%d; displacement[0]=%d\n",
|
|
|
|
|
block_length[0], displacement[0]);
|
|
|
|
|
#endif
|
|
|
|
|
|
2001-07-20 23:01:18 +08:00
|
|
|
|
err = MPI_Type_struct( 1, /* count */
|
|
|
|
|
block_length, /* blocklengths */
|
1999-07-01 00:30:56 +08:00
|
|
|
|
displacement, /* displacements */
|
|
|
|
|
old_type, /* old type */
|
|
|
|
|
new_type ); /* new type */
|
|
|
|
|
|
2001-07-20 23:01:18 +08:00
|
|
|
|
if (err)
|
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_MPI, FAIL,"couldn't create MPI struct type");
|
|
|
|
|
|
|
|
|
|
err=MPI_Type_free (&old_type[0]);
|
|
|
|
|
if (err)
|
[svn-r2646] Purpose:
Bug fix (done by Kim Yates)
Description:
The optimized mpio code was broken and when read was done, it hanged.
Solution:
H5FDmpio.c:
In H5FD_mpio_write, moved the 16-line block of code in which
all procs other than p0 skip the actual write
to be just before the call to MPI_File_write_at.
Previously, the values of the local vars that controlled
"allsame" were not always set correctly when the moved block
was reached.
H5S.c:
Changed default value of H5_mpi_opt_types_g to TRUE, so that
the MPI-IO hyperslab code is executed by default in parallel HDF5,
rather than executing the serial hyperslab code.
H5Smpio.c:
In function H5S_mpio_hyper_type, added a call to free
an intermediate type. Cures a small memory leak.
Added code for cases of empty hyperslab
Changed displacements to be MPI_Aint
Platforms tested:
modi4 -64: worked fine with mpich 1.2.0 but failed with the messages
saying it ran out of entries for MPI_Types during the collective_read
test. After tracing the code all the way to the collective read, all
MPI Types have been freed properly. It aborted with the above message
when it executed the line
if (MPI_SUCCESS!= MPI_File_read_at_all(file->f, mpi_off, buf, size_i, buf_type, &mpi_stat ))
Could not see any problem with this line. It could be a bug in the
SGI version of MPI.
2000-10-10 10:32:45 +08:00
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_MPI, FAIL,"couldn't resize MPI vector type");
|
1999-07-01 00:30:56 +08:00
|
|
|
|
}
|
|
|
|
|
else {
|
2001-07-20 23:01:18 +08:00
|
|
|
|
*new_type = inner_type;
|
1999-07-01 00:30:56 +08:00
|
|
|
|
}
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
err = MPI_Type_commit( new_type );
|
2001-07-20 23:01:18 +08:00
|
|
|
|
if (err)
|
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_MPI, FAIL,"couldn't commit MPI vector type");
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* fill in the remaining return values */
|
|
|
|
|
*count = 1; /* only have to move one of these suckers! */
|
|
|
|
|
*is_derived_type = 1;
|
2001-12-13 02:40:09 +08:00
|
|
|
|
HGOTO_DONE(SUCCEED);
|
[svn-r2646] Purpose:
Bug fix (done by Kim Yates)
Description:
The optimized mpio code was broken and when read was done, it hanged.
Solution:
H5FDmpio.c:
In H5FD_mpio_write, moved the 16-line block of code in which
all procs other than p0 skip the actual write
to be just before the call to MPI_File_write_at.
Previously, the values of the local vars that controlled
"allsame" were not always set correctly when the moved block
was reached.
H5S.c:
Changed default value of H5_mpi_opt_types_g to TRUE, so that
the MPI-IO hyperslab code is executed by default in parallel HDF5,
rather than executing the serial hyperslab code.
H5Smpio.c:
In function H5S_mpio_hyper_type, added a call to free
an intermediate type. Cures a small memory leak.
Added code for cases of empty hyperslab
Changed displacements to be MPI_Aint
Platforms tested:
modi4 -64: worked fine with mpich 1.2.0 but failed with the messages
saying it ran out of entries for MPI_Types during the collective_read
test. After tracing the code all the way to the collective read, all
MPI Types have been freed properly. It aborted with the above message
when it executed the line
if (MPI_SUCCESS!= MPI_File_read_at_all(file->f, mpi_off, buf, size_i, buf_type, &mpi_stat ))
Could not see any problem with this line. It could be a bug in the
SGI version of MPI.
2000-10-10 10:32:45 +08:00
|
|
|
|
|
|
|
|
|
empty:
|
|
|
|
|
/* special case: empty hyperslab */
|
|
|
|
|
*new_type = MPI_BYTE;
|
|
|
|
|
*count = 0;
|
|
|
|
|
*is_derived_type = 0;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
[svn-r2646] Purpose:
Bug fix (done by Kim Yates)
Description:
The optimized mpio code was broken and when read was done, it hanged.
Solution:
H5FDmpio.c:
In H5FD_mpio_write, moved the 16-line block of code in which
all procs other than p0 skip the actual write
to be just before the call to MPI_File_write_at.
Previously, the values of the local vars that controlled
"allsame" were not always set correctly when the moved block
was reached.
H5S.c:
Changed default value of H5_mpi_opt_types_g to TRUE, so that
the MPI-IO hyperslab code is executed by default in parallel HDF5,
rather than executing the serial hyperslab code.
H5Smpio.c:
In function H5S_mpio_hyper_type, added a call to free
an intermediate type. Cures a small memory leak.
Added code for cases of empty hyperslab
Changed displacements to be MPI_Aint
Platforms tested:
modi4 -64: worked fine with mpich 1.2.0 but failed with the messages
saying it ran out of entries for MPI_Types during the collective_read
test. After tracing the code all the way to the collective read, all
MPI Types have been freed properly. It aborted with the above message
when it executed the line
if (MPI_SUCCESS!= MPI_File_read_at_all(file->f, mpi_off, buf, size_i, buf_type, &mpi_stat ))
Could not see any problem with this line. It could be a bug in the
SGI version of MPI.
2000-10-10 10:32:45 +08:00
|
|
|
|
done:
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#ifdef H5Smpi_DEBUG
|
[svn-r2646] Purpose:
Bug fix (done by Kim Yates)
Description:
The optimized mpio code was broken and when read was done, it hanged.
Solution:
H5FDmpio.c:
In H5FD_mpio_write, moved the 16-line block of code in which
all procs other than p0 skip the actual write
to be just before the call to MPI_File_write_at.
Previously, the values of the local vars that controlled
"allsame" were not always set correctly when the moved block
was reached.
H5S.c:
Changed default value of H5_mpi_opt_types_g to TRUE, so that
the MPI-IO hyperslab code is executed by default in parallel HDF5,
rather than executing the serial hyperslab code.
H5Smpio.c:
In function H5S_mpio_hyper_type, added a call to free
an intermediate type. Cures a small memory leak.
Added code for cases of empty hyperslab
Changed displacements to be MPI_Aint
Platforms tested:
modi4 -64: worked fine with mpich 1.2.0 but failed with the messages
saying it ran out of entries for MPI_Types during the collective_read
test. After tracing the code all the way to the collective read, all
MPI Types have been freed properly. It aborted with the above message
when it executed the line
if (MPI_SUCCESS!= MPI_File_read_at_all(file->f, mpi_off, buf, size_i, buf_type, &mpi_stat ))
Could not see any problem with this line. It could be a bug in the
SGI version of MPI.
2000-10-10 10:32:45 +08:00
|
|
|
|
HDfprintf(stdout, "Leave %s, count=%Hu is_derived_type=%d\n",
|
|
|
|
|
FUNC, *count, *is_derived_type );
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#endif
|
2001-12-13 22:04:47 +08:00
|
|
|
|
FUNC_LEAVE (ret_value);
|
1999-08-11 04:21:32 +08:00
|
|
|
|
}
|
1999-07-01 00:30:56 +08:00
|
|
|
|
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/*-------------------------------------------------------------------------
|
|
|
|
|
* Function: H5S_mpio_space_type
|
|
|
|
|
*
|
|
|
|
|
* Purpose: Translate an HDF5 dataspace selection into an MPI type.
|
|
|
|
|
* Currently handle only hyperslab and "all" selections.
|
|
|
|
|
*
|
|
|
|
|
* Return: non-negative on success, negative on failure.
|
|
|
|
|
*
|
|
|
|
|
* Outputs: *new_type the MPI type corresponding to the selection
|
|
|
|
|
* *count how many objects of the new_type in selection
|
|
|
|
|
* (useful if this is the buffer type for xfer)
|
|
|
|
|
* *is_derived_type 0 if MPI primitive type, 1 if derived
|
|
|
|
|
*
|
|
|
|
|
* Programmer: rky 980813
|
|
|
|
|
*
|
|
|
|
|
* Modifications:
|
|
|
|
|
*
|
|
|
|
|
*-------------------------------------------------------------------------
|
|
|
|
|
*/
|
2002-05-29 02:17:12 +08:00
|
|
|
|
static herr_t
|
1998-08-17 23:13:45 +08:00
|
|
|
|
H5S_mpio_space_type( const H5S_t *space, const size_t elmt_size,
|
|
|
|
|
/* out: */
|
|
|
|
|
MPI_Datatype *new_type,
|
2001-07-20 23:01:18 +08:00
|
|
|
|
size_t *count,
|
1998-08-17 23:13:45 +08:00
|
|
|
|
hbool_t *is_derived_type )
|
|
|
|
|
{
|
|
|
|
|
int err;
|
|
|
|
|
herr_t ret_value = SUCCEED;
|
|
|
|
|
|
2002-05-29 02:17:12 +08:00
|
|
|
|
FUNC_ENTER_NOINIT(H5S_mpio_space_type);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* Check args */
|
|
|
|
|
assert (space);
|
|
|
|
|
|
|
|
|
|
/* Creat MPI type based on the kind of selection */
|
|
|
|
|
switch (space->extent.type) {
|
|
|
|
|
case H5S_SCALAR:
|
|
|
|
|
/* not yet implemented */
|
2001-07-20 23:01:18 +08:00
|
|
|
|
ret_value = FAIL;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
break;
|
|
|
|
|
|
|
|
|
|
case H5S_SIMPLE:
|
|
|
|
|
switch(space->select.type) {
|
|
|
|
|
case H5S_SEL_NONE:
|
|
|
|
|
case H5S_SEL_ALL:
|
2001-07-20 23:01:18 +08:00
|
|
|
|
err = H5S_mpio_all_type( space, elmt_size,
|
|
|
|
|
/* out: */ new_type, count, is_derived_type );
|
|
|
|
|
if (err<0)
|
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_BADTYPE, FAIL,"couldn't convert \"all\" selection to MPI type");
|
1998-08-17 23:13:45 +08:00
|
|
|
|
break;
|
|
|
|
|
|
|
|
|
|
case H5S_SEL_POINTS:
|
2001-07-20 23:01:18 +08:00
|
|
|
|
/* not yet implemented */
|
|
|
|
|
ret_value = FAIL;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
break;
|
|
|
|
|
|
|
|
|
|
case H5S_SEL_HYPERSLABS:
|
2001-07-20 23:01:18 +08:00
|
|
|
|
err = H5S_mpio_hyper_type( space, elmt_size,
|
|
|
|
|
/* out: */ new_type, count, is_derived_type );
|
|
|
|
|
if (err)
|
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_BADTYPE, FAIL,"couldn't convert \"all\" selection to MPI type");
|
1998-08-17 23:13:45 +08:00
|
|
|
|
break;
|
|
|
|
|
|
|
|
|
|
default:
|
|
|
|
|
assert("unknown selection type" && 0);
|
|
|
|
|
break;
|
|
|
|
|
} /* end switch */
|
|
|
|
|
break;
|
|
|
|
|
|
|
|
|
|
case H5S_COMPLEX:
|
|
|
|
|
/* not yet implemented */
|
2001-07-20 23:01:18 +08:00
|
|
|
|
ret_value = FAIL;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
break;
|
|
|
|
|
|
|
|
|
|
default:
|
|
|
|
|
assert("unknown data space type" && 0);
|
|
|
|
|
break;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
FUNC_LEAVE (ret_value);
|
1999-08-11 04:21:32 +08:00
|
|
|
|
}
|
|
|
|
|
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/*-------------------------------------------------------------------------
|
|
|
|
|
* Function: H5S_mpio_spaces_xfer
|
|
|
|
|
*
|
|
|
|
|
* Purpose: Use MPI-IO to transfer data efficiently
|
|
|
|
|
* directly between app buffer and file.
|
|
|
|
|
*
|
|
|
|
|
* Return: non-negative on success, negative on failure.
|
|
|
|
|
*
|
|
|
|
|
* Programmer: rky 980813
|
|
|
|
|
*
|
|
|
|
|
* Modifications:
|
2000-11-26 14:56:29 +08:00
|
|
|
|
* rky 980918
|
|
|
|
|
* Added must_convert parameter to let caller know we can't optimize
|
|
|
|
|
* the xfer.
|
2000-11-28 02:25:18 +08:00
|
|
|
|
*
|
2000-11-26 14:56:29 +08:00
|
|
|
|
* Albert Cheng, 001123
|
|
|
|
|
* Include the MPI_type freeing as part of cleanup code.
|
1998-09-19 07:08:09 +08:00
|
|
|
|
*
|
2002-04-04 01:07:14 +08:00
|
|
|
|
* QAK - 2002/04/02
|
|
|
|
|
* Removed the must_convert parameter and move preconditions to
|
|
|
|
|
* H5S_mpio_opt_possible() routine
|
|
|
|
|
*
|
1998-08-17 23:13:45 +08:00
|
|
|
|
*-------------------------------------------------------------------------
|
|
|
|
|
*/
|
2002-04-26 01:56:56 +08:00
|
|
|
|
static herr_t
|
|
|
|
|
H5S_mpio_spaces_xfer(H5F_t *f, const H5O_layout_t *layout,
|
|
|
|
|
H5P_genplist_t UNUSED *dc_plist, size_t elmt_size,
|
1998-08-17 23:13:45 +08:00
|
|
|
|
const H5S_t *file_space, const H5S_t *mem_space,
|
1999-08-11 04:21:32 +08:00
|
|
|
|
hid_t dxpl_id, void *buf /*out*/,
|
1998-08-17 23:13:45 +08:00
|
|
|
|
const hbool_t do_write )
|
|
|
|
|
{
|
1998-08-27 01:51:19 +08:00
|
|
|
|
herr_t ret_value = SUCCEED;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
int err;
|
|
|
|
|
haddr_t disp, addr;
|
2001-07-20 04:31:37 +08:00
|
|
|
|
size_t mpi_count;
|
2001-07-20 23:01:18 +08:00
|
|
|
|
size_t mpi_buf_count, mpi_unused_count;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
MPI_Datatype mpi_buf_type, mpi_file_type;
|
2000-11-26 14:56:29 +08:00
|
|
|
|
hbool_t mbt_is_derived=0,
|
|
|
|
|
mft_is_derived=0;
|
2002-06-05 05:29:05 +08:00
|
|
|
|
#if 0
|
2001-11-21 03:07:22 +08:00
|
|
|
|
H5P_genplist_t *plist; /* Property list pointer */
|
2002-06-05 05:29:05 +08:00
|
|
|
|
#endif /* 0 */
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
2002-05-29 02:17:12 +08:00
|
|
|
|
FUNC_ENTER_NOINIT(H5S_mpio_spaces_xfer);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* Check args */
|
|
|
|
|
assert (f);
|
|
|
|
|
assert (layout);
|
|
|
|
|
assert (file_space);
|
|
|
|
|
assert (mem_space);
|
|
|
|
|
assert (buf);
|
2000-03-24 10:12:44 +08:00
|
|
|
|
assert (IS_H5FD_MPIO(f));
|
2001-08-11 04:47:05 +08:00
|
|
|
|
/* Make certain we have the correct type of property list */
|
|
|
|
|
assert(H5I_GENPROP_LST==H5I_get_type(dxpl_id));
|
2001-11-21 03:07:22 +08:00
|
|
|
|
assert(TRUE==H5P_isa_class(dxpl_id,H5P_DATASET_XFER));
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
2000-11-28 02:25:18 +08:00
|
|
|
|
/*
|
|
|
|
|
* For collective data transfer only since this would eventually
|
|
|
|
|
* call H5FD_mpio_setup to do setup to eveually call MPI_File_set_view
|
|
|
|
|
* in H5FD_mpio_read or H5FD_mpio_write. MPI_File_set_view is a
|
|
|
|
|
* collective call. Letting independent data transfer use this
|
|
|
|
|
* route would result in hanging.
|
|
|
|
|
*/
|
|
|
|
|
#if 0
|
|
|
|
|
/* For now, the checking is being done in
|
|
|
|
|
* H5D_write and H5D_read before it is called because
|
|
|
|
|
* the following block of code, though with the right idea, is not
|
|
|
|
|
* correct yet.
|
|
|
|
|
*/
|
2001-08-11 04:47:05 +08:00
|
|
|
|
{
|
|
|
|
|
/* Get the transfer mode */
|
|
|
|
|
H5FD_mpio_dxpl_t *dx;
|
|
|
|
|
hid_t driver_id; /* VFL driver ID */
|
|
|
|
|
|
2001-11-21 03:07:22 +08:00
|
|
|
|
if(NULL == (plist = H5I_object(dxpl_id)))
|
|
|
|
|
HRETURN_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not a file access property list");
|
|
|
|
|
|
2001-08-11 04:47:05 +08:00
|
|
|
|
/* Get the driver ID */
|
2001-11-21 03:07:22 +08:00
|
|
|
|
if(H5P_get(plist, H5D_XFER_VFL_ID_NAME, &driver_id)<0)
|
2001-08-11 04:47:05 +08:00
|
|
|
|
HGOTO_ERROR (H5E_PLIST, H5E_CANTGET, FAIL, "Can't retrieve VFL driver ID");
|
|
|
|
|
|
|
|
|
|
/* Get the driver information */
|
2001-11-21 03:07:22 +08:00
|
|
|
|
if(H5P_get(plist, H5D_XFER_VFL_INFO_NAME, &dx)<0)
|
2001-08-11 04:47:05 +08:00
|
|
|
|
HGOTO_ERROR (H5E_PLIST, H5E_CANTGET, FAIL, "Can't retrieve VFL driver info");
|
2000-11-28 02:25:18 +08:00
|
|
|
|
}
|
|
|
|
|
#endif
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* create the MPI buffer type */
|
|
|
|
|
err = H5S_mpio_space_type( mem_space, elmt_size,
|
|
|
|
|
/* out: */
|
|
|
|
|
&mpi_buf_type,
|
|
|
|
|
&mpi_buf_count,
|
|
|
|
|
&mbt_is_derived );
|
|
|
|
|
if (MPI_SUCCESS != err)
|
2000-11-26 14:56:29 +08:00
|
|
|
|
HGOTO_ERROR(H5E_DATASPACE, H5E_BADTYPE, FAIL,"couldn't create MPI buf type");
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* create the MPI file type */
|
|
|
|
|
err = H5S_mpio_space_type( file_space, elmt_size,
|
|
|
|
|
/* out: */
|
|
|
|
|
&mpi_file_type,
|
|
|
|
|
&mpi_unused_count,
|
|
|
|
|
&mft_is_derived );
|
|
|
|
|
if (MPI_SUCCESS != err)
|
2000-11-26 14:56:29 +08:00
|
|
|
|
HGOTO_ERROR(H5E_DATASPACE, H5E_BADTYPE, FAIL,"couldn't create MPI file type");
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* calculate the absolute base addr (i.e., the file view disp) */
|
1999-07-29 03:37:35 +08:00
|
|
|
|
disp = f->shared->base_addr + layout->addr;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#ifdef H5Smpi_DEBUG
|
[svn-r2646] Purpose:
Bug fix (done by Kim Yates)
Description:
The optimized mpio code was broken and when read was done, it hanged.
Solution:
H5FDmpio.c:
In H5FD_mpio_write, moved the 16-line block of code in which
all procs other than p0 skip the actual write
to be just before the call to MPI_File_write_at.
Previously, the values of the local vars that controlled
"allsame" were not always set correctly when the moved block
was reached.
H5S.c:
Changed default value of H5_mpi_opt_types_g to TRUE, so that
the MPI-IO hyperslab code is executed by default in parallel HDF5,
rather than executing the serial hyperslab code.
H5Smpio.c:
In function H5S_mpio_hyper_type, added a call to free
an intermediate type. Cures a small memory leak.
Added code for cases of empty hyperslab
Changed displacements to be MPI_Aint
Platforms tested:
modi4 -64: worked fine with mpich 1.2.0 but failed with the messages
saying it ran out of entries for MPI_Types during the collective_read
test. After tracing the code all the way to the collective read, all
MPI Types have been freed properly. It aborted with the above message
when it executed the line
if (MPI_SUCCESS!= MPI_File_read_at_all(file->f, mpi_off, buf, size_i, buf_type, &mpi_stat ))
Could not see any problem with this line. It could be a bug in the
SGI version of MPI.
2000-10-10 10:32:45 +08:00
|
|
|
|
HDfprintf(stdout, "spaces_xfer: disp=%Hu\n", disp );
|
1998-08-17 23:13:45 +08:00
|
|
|
|
#endif
|
|
|
|
|
|
|
|
|
|
/* Effective address determined by base addr and the MPI file type */
|
1999-07-29 03:37:35 +08:00
|
|
|
|
addr = 0;
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
1999-08-11 04:21:32 +08:00
|
|
|
|
/*
|
|
|
|
|
* Pass buf type, file type, and absolute base address (i.e., the file
|
|
|
|
|
* view disp) to the file driver. Request a dataspace transfer (instead
|
|
|
|
|
* of an elementary byteblock transfer).
|
|
|
|
|
*/
|
|
|
|
|
H5FD_mpio_setup(f->shared->lf, mpi_buf_type, mpi_file_type, disp, 1);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/* transfer the data */
|
2001-07-20 04:31:37 +08:00
|
|
|
|
H5_CHECK_OVERFLOW(mpi_buf_count, hsize_t, size_t);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
mpi_count = (size_t)mpi_buf_count;
|
|
|
|
|
if (do_write) {
|
2000-09-06 11:40:21 +08:00
|
|
|
|
err = H5FD_write(f->shared->lf, H5FD_MEM_DRAW, dxpl_id, addr, mpi_count, buf);
|
1998-08-28 11:59:23 +08:00
|
|
|
|
if (err) {
|
|
|
|
|
HRETURN_ERROR(H5E_IO, H5E_WRITEERROR, FAIL,"MPI write failed");
|
|
|
|
|
}
|
1998-08-17 23:13:45 +08:00
|
|
|
|
} else {
|
2000-10-25 13:54:05 +08:00
|
|
|
|
err = H5FD_read (f->shared->lf, H5FD_MEM_DRAW, dxpl_id, addr, mpi_count, buf);
|
1998-08-28 11:59:23 +08:00
|
|
|
|
if (err) {
|
|
|
|
|
HRETURN_ERROR(H5E_IO, H5E_READERROR, FAIL,"MPI read failed");
|
|
|
|
|
}
|
1998-08-17 23:13:45 +08:00
|
|
|
|
}
|
|
|
|
|
|
2000-11-26 14:56:29 +08:00
|
|
|
|
done:
|
1998-08-17 23:13:45 +08:00
|
|
|
|
/* free the MPI buf and file types */
|
|
|
|
|
if (mbt_is_derived) {
|
|
|
|
|
err = MPI_Type_free( &mpi_buf_type );
|
|
|
|
|
if (MPI_SUCCESS != err) {
|
1999-08-11 04:21:32 +08:00
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_MPI, FAIL,
|
|
|
|
|
"unable to free MPI file type");
|
1998-08-17 23:13:45 +08:00
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
if (mft_is_derived) {
|
|
|
|
|
err = MPI_Type_free( &mpi_file_type );
|
|
|
|
|
if (MPI_SUCCESS != err) {
|
1999-08-11 04:21:32 +08:00
|
|
|
|
HRETURN_ERROR(H5E_DATASPACE, H5E_MPI, FAIL,
|
|
|
|
|
"unable to free MPI file type");
|
1998-08-17 23:13:45 +08:00
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
1998-08-27 01:51:19 +08:00
|
|
|
|
FUNC_LEAVE (ret_value);
|
1999-08-11 04:21:32 +08:00
|
|
|
|
}
|
|
|
|
|
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/*-------------------------------------------------------------------------
|
|
|
|
|
* Function: H5S_mpio_spaces_read
|
|
|
|
|
*
|
|
|
|
|
* Purpose: MPI-IO function to read directly from app buffer to file.
|
|
|
|
|
*
|
|
|
|
|
* Return: non-negative on success, negative on failure.
|
|
|
|
|
*
|
|
|
|
|
* Programmer: rky 980813
|
|
|
|
|
*
|
|
|
|
|
* Modifications:
|
|
|
|
|
*
|
1998-09-19 07:08:09 +08:00
|
|
|
|
* rky 980918
|
|
|
|
|
* Added must_convert parameter to let caller know we can't optimize the xfer.
|
|
|
|
|
*
|
2002-04-04 01:07:14 +08:00
|
|
|
|
* QAK - 2002/04/02
|
|
|
|
|
* Removed the must_convert parameter and move preconditions to
|
|
|
|
|
* H5S_mpio_opt_possible() routine
|
|
|
|
|
*
|
1998-08-17 23:13:45 +08:00
|
|
|
|
*-------------------------------------------------------------------------
|
|
|
|
|
*/
|
|
|
|
|
herr_t
|
2002-04-26 01:56:56 +08:00
|
|
|
|
H5S_mpio_spaces_read(H5F_t *f, const H5O_layout_t *layout,
|
|
|
|
|
H5P_genplist_t *dc_plist, size_t elmt_size,
|
1998-08-17 23:13:45 +08:00
|
|
|
|
const H5S_t *file_space, const H5S_t *mem_space,
|
2002-04-04 01:07:14 +08:00
|
|
|
|
hid_t dxpl_id, void *buf/*out*/)
|
1998-08-17 23:13:45 +08:00
|
|
|
|
{
|
|
|
|
|
herr_t ret_value = FAIL;
|
|
|
|
|
|
2002-05-29 23:07:55 +08:00
|
|
|
|
FUNC_ENTER_NOAPI(H5S_mpio_spaces_read, FAIL);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
2002-04-26 01:56:56 +08:00
|
|
|
|
ret_value = H5S_mpio_spaces_xfer(f, layout, dc_plist, elmt_size,
|
1999-08-11 04:21:32 +08:00
|
|
|
|
file_space, mem_space, dxpl_id,
|
2002-04-04 01:07:14 +08:00
|
|
|
|
buf, 0/*read*/);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
FUNC_LEAVE (ret_value);
|
1999-08-11 04:21:32 +08:00
|
|
|
|
}
|
|
|
|
|
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
/*-------------------------------------------------------------------------
|
|
|
|
|
* Function: H5S_mpio_spaces_write
|
|
|
|
|
*
|
|
|
|
|
* Purpose: MPI-IO function to write directly from app buffer to file.
|
|
|
|
|
*
|
|
|
|
|
* Return: non-negative on success, negative on failure.
|
|
|
|
|
*
|
|
|
|
|
* Programmer: rky 980813
|
|
|
|
|
*
|
|
|
|
|
* Modifications:
|
|
|
|
|
*
|
1998-09-19 07:08:09 +08:00
|
|
|
|
* rky 980918
|
|
|
|
|
* Added must_convert parameter to let caller know we can't optimize the xfer.
|
|
|
|
|
*
|
2002-04-04 01:07:14 +08:00
|
|
|
|
* QAK - 2002/04/02
|
|
|
|
|
* Removed the must_convert parameter and move preconditions to
|
|
|
|
|
* H5S_mpio_opt_possible() routine
|
|
|
|
|
*
|
1998-08-17 23:13:45 +08:00
|
|
|
|
*-------------------------------------------------------------------------
|
|
|
|
|
*/
|
|
|
|
|
herr_t
|
2002-04-26 01:56:56 +08:00
|
|
|
|
H5S_mpio_spaces_write(H5F_t *f, const H5O_layout_t *layout,
|
|
|
|
|
H5P_genplist_t *dc_plist, size_t elmt_size,
|
2002-02-14 23:57:48 +08:00
|
|
|
|
const H5S_t *file_space, const H5S_t *mem_space,
|
2002-04-04 01:07:14 +08:00
|
|
|
|
hid_t dxpl_id, const void *buf)
|
1998-08-17 23:13:45 +08:00
|
|
|
|
{
|
|
|
|
|
herr_t ret_value = FAIL;
|
|
|
|
|
|
2002-05-29 23:07:55 +08:00
|
|
|
|
FUNC_ENTER_NOAPI(H5S_mpio_spaces_write, FAIL);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
2002-04-26 01:56:56 +08:00
|
|
|
|
ret_value = H5S_mpio_spaces_xfer(f, layout, dc_plist, elmt_size,
|
1999-08-11 04:21:32 +08:00
|
|
|
|
file_space, mem_space, dxpl_id,
|
2002-04-04 01:07:14 +08:00
|
|
|
|
(void*)buf, 1/*write*/);
|
1998-08-17 23:13:45 +08:00
|
|
|
|
|
|
|
|
|
FUNC_LEAVE (ret_value);
|
1999-08-11 04:21:32 +08:00
|
|
|
|
}
|
1998-08-19 22:07:14 +08:00
|
|
|
|
|
2002-04-04 01:07:14 +08:00
|
|
|
|
|
|
|
|
|
/*-------------------------------------------------------------------------
|
|
|
|
|
* Function: H5S_mpio_opt_possible
|
|
|
|
|
*
|
|
|
|
|
* Purpose: Checks if an direct I/O transfer is possible between memory and
|
|
|
|
|
* the file.
|
|
|
|
|
*
|
|
|
|
|
* Return: Success: Non-negative: TRUE or FALSE
|
|
|
|
|
* Failure: Negative
|
|
|
|
|
*
|
|
|
|
|
* Programmer: Quincey Koziol
|
|
|
|
|
* Wednesday, April 3, 2002
|
|
|
|
|
*
|
|
|
|
|
* Modifications:
|
|
|
|
|
*
|
|
|
|
|
*-------------------------------------------------------------------------
|
|
|
|
|
*/
|
|
|
|
|
htri_t
|
|
|
|
|
H5S_mpio_opt_possible( const H5S_t *mem_space, const H5S_t *file_space, const unsigned flags)
|
|
|
|
|
{
|
|
|
|
|
htri_t c1,c2; /* Flags whether a selection is optimizable */
|
|
|
|
|
htri_t ret_value=TRUE;
|
|
|
|
|
|
2002-05-29 23:07:55 +08:00
|
|
|
|
FUNC_ENTER_NOAPI(H5S_all_opt_possible, FAIL);
|
2002-04-04 01:07:14 +08:00
|
|
|
|
|
|
|
|
|
/* Check args */
|
|
|
|
|
assert(mem_space);
|
|
|
|
|
assert(file_space);
|
|
|
|
|
|
|
|
|
|
/* Check whether these are both simple dataspaces */
|
|
|
|
|
if (H5S_SIMPLE!=mem_space->extent.type || H5S_SIMPLE!=file_space->extent.type)
|
|
|
|
|
HGOTO_DONE(FALSE);
|
|
|
|
|
|
|
|
|
|
/* Check whether both selections are "regular" */
|
|
|
|
|
c1=H5S_select_regular(file_space);
|
|
|
|
|
c2=H5S_select_regular(mem_space);
|
|
|
|
|
if(c1==FAIL || c2==FAIL)
|
|
|
|
|
HGOTO_ERROR(H5E_DATASPACE, H5E_BADRANGE, FAIL, "invalid check for single selection blocks");
|
|
|
|
|
if(c1==FALSE || c2==FALSE)
|
|
|
|
|
HGOTO_DONE(FALSE);
|
|
|
|
|
|
|
|
|
|
/* Can't currently handle point selections */
|
|
|
|
|
if (H5S_SEL_POINTS==mem_space->select.type || H5S_SEL_POINTS==file_space->select.type)
|
|
|
|
|
HGOTO_DONE(FALSE);
|
|
|
|
|
|
|
|
|
|
/* Dataset storage must be contiguous currently */
|
2002-04-04 01:37:02 +08:00
|
|
|
|
if ((flags&H5S_CONV_STORAGE_MASK)!=H5S_CONV_STORAGE_CONTIGUOUS)
|
2002-04-04 01:07:14 +08:00
|
|
|
|
HGOTO_DONE(FALSE);
|
|
|
|
|
|
|
|
|
|
/* Parallel I/O conversion flag must be set */
|
|
|
|
|
if(!(flags&H5S_CONV_PAR_IO_POSSIBLE))
|
|
|
|
|
HGOTO_DONE(FALSE);
|
|
|
|
|
|
|
|
|
|
done:
|
|
|
|
|
FUNC_LEAVE(ret_value);
|
|
|
|
|
} /* H5S_mpio_opt_possible() */
|
|
|
|
|
|
1999-12-17 22:37:22 +08:00
|
|
|
|
#endif /* H5_HAVE_PARALLEL */
|