libxmp/libxmpf in Omni Compiler  1.3.4
xmp_io.c File Reference
#include <stdio.h>
#include <string.h>
#include <stdlib.h>
#include <mpi.h>
#include "xmp.h"
#include "xmp_constant.h"
#include "xmp_data_struct.h"
#include "xmp_io_sys.h"
#include "xmp_internal.h"
Include dependency graph for xmp_io.c:

Macros

#define MPI_PORTABLE_PLATFORM_H
 
#define MPI_TYPE_CREATE_RESIZED1   MPI_Type_create_resized
 
#define func_m(p, q)   ((q) >= 0 ? -(q)/(p) : ((p) >= 0 ? (-(q)+(p)-1)/(p) : (-(q)-(p)-1)/(p) ))
 
#define RP_DIMS   (rp_dims)
 
#define RP_LB(i)   (rp_lb_addr[(i)])
 
#define RP_UB(i)   (rp_ub_addr[(i)])
 
#define RP_STEP(i)   (rp_step_addr[(i)])
 
#define RP_DIMS   (rp_dims)
 
#define RP_LB(i)   (rp_lb_addr[(i)])
 
#define RP_UB(i)   (rp_ub_addr[(i)])
 
#define RP_STEP(i)   (rp_step_addr[(i)])
 
#define RP_DIMS   (rp_dims)
 
#define RP_LB(i)   (rp_lb_addr[(i)])
 
#define RP_UB(i)   (rp_ub_addr[(i)])
 
#define RP_STEP(i)   (rp_step_addr[(i)])
 
#define RP_DIMS   (rp_dims)
 
#define RP_LB(i)   (rp_lb_addr[(i)])
 
#define RP_UB(i)   (rp_ub_addr[(i)])
 
#define RP_STEP(i)   (rp_step_addr[(i)])
 
#define RP_DIMS   (rp_dims)
 
#define RP_LB(i)   (rp_lb_addr[(i)])
 
#define RP_UB(i)   (rp_ub_addr[(i)])
 
#define RP_STEP(i)   (rp_step_addr[(i)])
 

Functions

void _XMP_fatal (char *msg)
 
xmp_range_txmp_allocate_range (int n_dim)
 
void xmp_set_range (xmp_range_t *rp, int i_dim, int lb, int length, int step)
 
void xmp_free_range (xmp_range_t *rp)
 
xmp_file_txmp_fopen_all (const char *fname, const char *amode)
 
int xmp_fclose_all (xmp_file_t *pstXmp_file)
 
int xmp_fseek (xmp_file_t *pstXmp_file, long long offset, int whence)
 
int xmp_fseek_shared (xmp_file_t *pstXmp_file, long long offset, int whence)
 
long long xmp_ftell (xmp_file_t *pstXmp_file)
 
long long xmp_ftell_shared (xmp_file_t *pstXmp_file)
 
long long xmp_file_sync_all (xmp_file_t *pstXmp_file)
 
ssize_t xmp_fread_all (xmp_file_t *pstXmp_file, void *buffer, size_t size, size_t count)
 
ssize_t xmp_fwrite_all (xmp_file_t *pstXmp_file, void *buffer, size_t size, size_t count)
 
int xmp_fread_darray_unpack (xmp_file_t *fp, xmp_desc_t apd, xmp_range_t *rp)
 
ssize_t xmp_fread_darray_all (xmp_file_t *pstXmp_file, xmp_desc_t apd, xmp_range_t *rp)
 
int xmp_fwrite_darray_pack (xmp_file_t *fp, xmp_desc_t apd, xmp_range_t *rp)
 
ssize_t xmp_fwrite_darray_all (xmp_file_t *pstXmp_file, xmp_desc_t apd, xmp_range_t *rp)
 
ssize_t xmp_fread_shared (xmp_file_t *pstXmp_file, void *buffer, size_t size, size_t count)
 
ssize_t xmp_fwrite_shared (xmp_file_t *pstXmp_file, void *buffer, size_t size, size_t count)
 
ssize_t xmp_fread (xmp_file_t *pstXmp_file, void *buffer, size_t size, size_t count)
 
ssize_t xmp_fwrite (xmp_file_t *pstXmp_file, void *buffer, size_t size, size_t count)
 
int xmp_file_set_view_all (xmp_file_t *pstXmp_file, long long disp, xmp_desc_t apd, xmp_range_t *rp)
 
int xmp_file_clear_view_all (xmp_file_t *pstXmp_file, long long disp)
 

Macro Definition Documentation

◆ func_m

#define func_m (   p,
 
)    ((q) >= 0 ? -(q)/(p) : ((p) >= 0 ? (-(q)+(p)-1)/(p) : (-(q)-(p)-1)/(p) ))

◆ MPI_PORTABLE_PLATFORM_H

#define MPI_PORTABLE_PLATFORM_H

◆ MPI_TYPE_CREATE_RESIZED1

#define MPI_TYPE_CREATE_RESIZED1   MPI_Type_create_resized

◆ RP_DIMS [1/5]

#define RP_DIMS   (rp_dims)

◆ RP_DIMS [2/5]

#define RP_DIMS   (rp_dims)

◆ RP_DIMS [3/5]

#define RP_DIMS   (rp_dims)

◆ RP_DIMS [4/5]

#define RP_DIMS   (rp_dims)

◆ RP_DIMS [5/5]

#define RP_DIMS   (rp_dims)

◆ RP_LB [1/5]

#define RP_LB (   i)    (rp_lb_addr[(i)])

◆ RP_LB [2/5]

#define RP_LB (   i)    (rp_lb_addr[(i)])

◆ RP_LB [3/5]

#define RP_LB (   i)    (rp_lb_addr[(i)])

◆ RP_LB [4/5]

#define RP_LB (   i)    (rp_lb_addr[(i)])

◆ RP_LB [5/5]

#define RP_LB (   i)    (rp_lb_addr[(i)])

◆ RP_STEP [1/5]

#define RP_STEP (   i)    (rp_step_addr[(i)])

◆ RP_STEP [2/5]

#define RP_STEP (   i)    (rp_step_addr[(i)])

◆ RP_STEP [3/5]

#define RP_STEP (   i)    (rp_step_addr[(i)])

◆ RP_STEP [4/5]

#define RP_STEP (   i)    (rp_step_addr[(i)])

◆ RP_STEP [5/5]

#define RP_STEP (   i)    (rp_step_addr[(i)])

◆ RP_UB [1/5]

#define RP_UB (   i)    (rp_ub_addr[(i)])

◆ RP_UB [2/5]

#define RP_UB (   i)    (rp_ub_addr[(i)])

◆ RP_UB [3/5]

#define RP_UB (   i)    (rp_ub_addr[(i)])

◆ RP_UB [4/5]

#define RP_UB (   i)    (rp_ub_addr[(i)])

◆ RP_UB [5/5]

#define RP_UB (   i)    (rp_ub_addr[(i)])

Function Documentation

◆ _XMP_fatal()

void _XMP_fatal ( char *  msg)
43 {
44  fprintf(stderr, "[RANK:%d] XcalableMP runtime error: %s\n", _XMP_world_rank, msg);
45  MPI_Abort(MPI_COMM_WORLD, 1);
46 }
Here is the caller graph for this function:

◆ xmp_allocate_range()

xmp_range_t* xmp_allocate_range ( int  n_dim)
1268 {
1269 #ifdef CHECK_POINT
1270  fprintf(stderr, "IO:START(xmp_allocate_range)\n");
1271 #endif /* CHECK_POINT */
1272  xmp_range_t *rp = NULL;
1273  if (n_dim <= 0){ return rp; }
1274  rp = (xmp_range_t *)malloc(sizeof(xmp_range_t));
1275  rp->dims = n_dim;
1276  rp->lb = (int*)malloc(sizeof(int)*rp->dims);
1277  rp->ub = (int*)malloc(sizeof(int)*rp->dims);
1278  rp->step = (int*)malloc(sizeof(int)*rp->dims);
1279  if(!rp->lb || !rp->ub || !rp->step){ return rp; }
1280 #ifdef CHECK_POINT
1281  fprintf(stderr, "IO:END (xmp_allocate_range)\n");
1282 #endif /* CHECK_POINT */
1283  return rp;
1284 }

◆ xmp_fclose_all()

int xmp_fclose_all ( xmp_file_t pstXmp_file)
1503 {
1504 #ifdef CHECK_POINT
1505  fprintf(stderr, "IO:START(xmp_fclose_all)\n");
1506 #endif /* CHECK_POINT */
1507  // check argument
1508  if (pstXmp_file == NULL) { return 1; }
1509 
1510  // file close
1511  if (MPI_File_close(&(pstXmp_file->fh)) != MPI_SUCCESS)
1512  {
1513  free(pstXmp_file);
1514 #ifdef CHECK_POINT
1515  fprintf(stderr, "IO:END (xmp_fclose_all)\n");
1516 #endif /* CHECK_POINT */
1517  return 2;
1518  }
1519  free(pstXmp_file);
1520 #ifdef CHECK_POINT
1521  fprintf(stderr, "IO:END (xmp_fclose_all)\n");
1522 #endif /* CHECK_POINT */
1523  return 0;
1524 }

◆ xmp_file_clear_view_all()

int xmp_file_clear_view_all ( xmp_file_t pstXmp_file,
long long  disp 
)
3616 {
3617  // check argument
3618  if (pstXmp_file == NULL) { return 1; }
3619  if (disp < 0) { return 1; }
3620 
3621  // initialize view
3622  if (MPI_File_set_view(pstXmp_file->fh,
3623  disp,
3624  MPI_BYTE,
3625  MPI_BYTE,
3626  "native",
3627  MPI_INFO_NULL) != MPI_SUCCESS)
3628  {
3629  return 1;
3630  }
3631 
3632  return 0;
3633 }

◆ xmp_file_set_view_all()

int xmp_file_set_view_all ( xmp_file_t pstXmp_file,
long long  disp,
xmp_desc_t  apd,
xmp_range_t rp 
)
3229 {
3230  int i = 0;
3231  int mpiRet; // return value of MPI functions
3232  int lower; // lower bound accessed by this node
3233  int upper; // upper bound accessed by this node
3234  long contiguous_size; // contiguous size
3235  MPI_Datatype dataType[2];
3236  MPI_Aint tmp1, type_size;
3237  xmp_desc_t tempd;
3238  int rp_dims;
3239  int *rp_lb_addr = NULL;
3240  int *rp_ub_addr = NULL;
3241  int *rp_step_addr = NULL;
3242  int array_ndims;
3243  size_t array_type_size;
3244  //int ierr;
3245 
3246  int rank, nproc;
3247  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
3248  MPI_Comm_size(MPI_COMM_WORLD, &nproc);
3249 
3250 #ifdef CHECK_POINT
3251  fprintf(stderr, "IO:START(xmp_file_set_view_all): rank=%d\n", rank);
3252 #endif /* CHECK_POINT */
3253 
3254  // check argument
3255  if (pstXmp_file == NULL) { return 1001; }
3256  if (apd == NULL) { return 1002; }
3257  if (rp == NULL) { return 1004; }
3258  if (disp < 0) { return 1005; }
3259 
3260  /*ierr =*/ xmp_align_template(apd, &tempd);
3261  if (tempd == NULL){ return 1006; }
3262  /*ierr =*/ xmp_array_ndims(apd, &array_ndims);
3263  array_type_size = xmp_array_type_size(apd);
3264 
3265  rp_dims = _xmp_range_get_dims(rp);
3266  rp_lb_addr = _xmp_range_get_lb_addr(rp);
3267  rp_ub_addr = _xmp_range_get_ub_addr(rp);
3268  rp_step_addr = _xmp_range_get_step_addr(rp);
3269  if (!rp_lb_addr || !rp_ub_addr || !rp_step_addr){ return 1007; }
3270 #define RP_DIMS (rp_dims)
3271 #define RP_LB(i) (rp_lb_addr[(i)])
3272 #define RP_UB(i) (rp_ub_addr[(i)])
3273 #define RP_STEP(i) (rp_step_addr[(i)])
3274 
3275  // check number of dimensions
3276  if (array_ndims != RP_DIMS) { return 1008; }
3277 
3278 #ifdef DEBUG
3279 fprintf(stderr, "VIEW(%d/%d) dims=%d\n", rank, nproc, RP_DIMS);
3280 #endif
3281 
3282  // create basic data type
3283  MPI_Type_contiguous(array_type_size, MPI_BYTE, &dataType[0]);
3284 
3285  // loop for each dimension
3286  for (i = RP_DIMS - 1; i >= 0; i--)
3287  {
3288  int par_lower_i = xmp_array_gcllbound(apd, i+1);
3289  int par_upper_i = xmp_array_gclubound_tmp(apd, i+1);
3290  int align_manner_i = xmp_align_format(apd, i+1);
3291 #ifdef DEBUG
3292  fprintf(stderr, "xmp_file_set_view_all: myrank=%d: i=%d: "
3293  "align_manner_i=%d bw_i=%d par_lower_i=%d par_upper_i=%d\n",
3294  rank, i,
3295  xmp_align_format(apd, i+1),
3296  xmp_align_size(apd, i+1),
3297  xmp_array_gcllbound(apd, i+1),
3298  xmp_array_gclubound_tmp(apd, i+1));
3299 #endif /* DEBUG */
3300 
3301  // get extent of data type
3302  mpiRet =MPI_Type_get_extent(dataType[0], &tmp1, &type_size);
3303  if (mpiRet != MPI_SUCCESS) { return -1009; }
3304 
3305  int byte_dataType0; MPI_Type_size(dataType[0], &byte_dataType0);
3306 #ifdef DEBUG
3307  fprintf(stderr, "xmp_file_set_view_all: rank=%d: i=%d align_manner_i=%d type_size=%ld byte_dataType0=%d\n",
3308  rank, i, align_manner_i, (long)type_size, byte_dataType0);
3309 #endif /* DEBUG */
3310 
3311 #ifdef DEBUG
3312 fprintf(stderr, "VIEW(%d/%d) (lb,ub,step)=(%d,%d,%d)\n",
3313  rank, nproc, RP_LB(i), RP_UB(i), RP_STEP(i));
3314 fprintf(stderr, "VIEW(%d/%d) (par_lower,par_upper)=(%d,%d)\n",
3315  rank, nproc, par_lower_i, par_upper_i);
3316 #endif
3317  // no distribution
3318  if (align_manner_i == _XMP_N_ALIGN_NOT_ALIGNED ||
3319  align_manner_i == _XMP_N_ALIGN_DUPLICATION)
3320  {
3321  // contiguous size
3322  contiguous_size = (RP_UB(i) - RP_LB(i)) / RP_STEP(i) + 1;
3323 
3324  // create basic data type
3325  mpiRet = MPI_Type_contiguous(contiguous_size, dataType[0], &dataType[1]);
3326 
3327  // free MPI_Datatype out of use
3328  MPI_Type_free(&dataType[0]);
3329  dataType[0] = dataType[1];
3330 
3331  // on error in MPI_Type_contiguous
3332  if (mpiRet != MPI_SUCCESS) { return 1010; }
3333 
3334 #ifdef DEBUG
3335 fprintf(stderr, "VIEW(%d/%d) NOT_ALIGNED\n", rank, nproc);
3336 fprintf(stderr, "VIEW(%d/%d) contiguous_size=%d\n", rank, nproc, contiguous_size);
3337 #endif
3338  }
3339  // block distribution
3340  else if (align_manner_i == _XMP_N_ALIGN_BLOCK)
3341  {
3342  long space_size;
3343  long total_size;
3344 
3345  // increment is positive
3346  if (RP_STEP(i) >= 0)
3347  {
3348  // lower > upper
3349  if (RP_LB(i) > RP_UB(i))
3350  {
3351  return 1011; /* return 1; *//* MODIFIED */
3352  }
3353  // upper after distribution < lower
3354  else if (par_upper_i < RP_LB(i))
3355  {
3356  contiguous_size = space_size = 0;
3357  }
3358  // lower after distribution > upper
3359  else if (par_lower_i > RP_UB(i))
3360  {
3361  contiguous_size = space_size = 0;
3362  }
3363  // other
3364  else
3365  {
3366  // lower in this node
3367  lower
3368  = (par_lower_i > RP_LB(i)) ?
3369  RP_LB(i) + ((par_lower_i - 1 - RP_LB(i)) / RP_STEP(i) + 1) * RP_STEP(i)
3370  : RP_LB(i);
3371 
3372  // upper in this node
3373  upper
3374  = (par_upper_i < RP_UB(i)) ?
3375  par_upper_i : RP_UB(i);
3376 
3377  // contiguous size
3378  contiguous_size = (upper - lower) / RP_STEP(i) + 1;
3379 
3380  // space size
3381  space_size
3382  = ((lower - RP_LB(i)) / RP_STEP(i)) * type_size;
3383 
3384 /* fprintf(stderr, "set_view_all: rank = %d: lower = %d upper = %d contiguous_size = %d space_size = %d\n", */
3385 /* rank,lower, upper, contiguous_size, space_size); */
3386 
3387  }
3388 
3389  // total size
3390  total_size
3391  = ((RP_UB(i) - RP_LB(i)) / RP_STEP(i) + 1) * type_size;
3392 
3393  // create basic data type
3394  mpiRet = MPI_Type_contiguous(contiguous_size, dataType[0], &dataType[1]);
3395 
3396  // free MPI_Datatype out of use
3397  MPI_Type_free(&dataType[0]);
3398 
3399  // on error in MPI_Type_contiguous
3400  if (mpiRet != MPI_SUCCESS) { return 1012; }
3401 
3402  {
3403  int byte_datatype1; MPI_Aint lb_datatype1, extent_datatype1;
3404  MPI_Type_size(dataType[1], &byte_datatype1);
3405  MPI_Type_get_extent(dataType[1], &lb_datatype1, &extent_datatype1);
3406  if (extent_datatype1 + space_size > total_size){
3407  _XMP_fatal("xmp_file_set_view_all (block): data type is incorrect");
3408  }
3409  }
3410 
3411  // create new file type
3412  mpiRet = MPI_TYPE_CREATE_RESIZED1(dataType[1],
3413  space_size,
3414  total_size,
3415  &dataType[0]);
3416 
3417  // on error in MPI_Type_create_resized1
3418  if (mpiRet != MPI_SUCCESS) { return 1013; }
3419 
3420  // free MPI_Datatype out of use
3421  MPI_Type_free(&dataType[1]);
3422 
3423  int byte_dataType0; MPI_Aint lb_dataType0, extent_dataType0;
3424  MPI_Type_size(dataType[0], &byte_dataType0);
3425  MPI_Type_get_extent(dataType[0], &lb_dataType0, &extent_dataType0);
3426 #ifdef DEBUG
3427  fprintf(stderr, "set_view_all: after block: myrank=%d: byte_dataType0=%d lb=%ld extent=%ld ; space_size=%ld total_size=%ld\n",
3428  rank, byte_dataType0, (long)lb_dataType0, (long)extent_dataType0,
3429  space_size, total_size);
3430 #endif /* DEBUG */
3431 #ifdef DEBUG
3432 fprintf(stderr, "VIEW(%d/%d) ALIGN_BLOCK\n", rank, nproc );
3433 fprintf(stderr, "VIEW(%d/%d) type_size=%ld\n", rank, nproc , (long)type_size);
3434 fprintf(stderr, "VIEW(%d/%d) contiguous_size=%ld\n", rank, nproc , contiguous_size);
3435 fprintf(stderr, "VIEW(%d/%d) space_size=%ld\n", rank, nproc , space_size);
3436 fprintf(stderr, "VIEW(%d/%d) total_size=%ld\n", rank, nproc , total_size);
3437 fprintf(stderr, "VIEW(%d/%d) (lower,upper)=(%d,%d)\n", rank, nproc , lower, upper);
3438 fprintf(stderr, "\n");
3439 #endif
3440  }
3441  // incremnet is negative
3442  else if (RP_STEP(i) < 0)
3443  {
3444  // lower < upper
3445  if (RP_LB(i) < RP_UB(i))
3446  {
3447  return 1014;
3448  }
3449  // lower after distribution < upper
3450  else if (par_lower_i < RP_UB(i))
3451  {
3452  contiguous_size = space_size = 0;
3453  }
3454  // upper after distribution > lower
3455  else if (par_upper_i > RP_LB(i))
3456  {
3457  contiguous_size = space_size = 0;
3458  }
3459  // other
3460  else
3461  {
3462  // lower in this node
3463  lower
3464  = (par_upper_i < RP_LB(i)) ?
3465  RP_LB(i) - (( RP_LB(i) - par_upper_i - 1) / RP_STEP(i) - 1) * RP_STEP(i)
3466  : RP_LB(i);
3467 
3468  // upper in this node
3469  upper
3470  = (par_lower_i > RP_UB(i)) ?
3471  par_lower_i : RP_UB(i);
3472 
3473  // contiguous size
3474  contiguous_size = (upper - lower) / RP_STEP(i) + 1;
3475 
3476  // space size
3477 /* space_size */
3478 /* = ((lower - RP_LB(i)) / RP_STEP(i)) * type_size; */
3479  space_size
3480  = ( - (upper - RP_UB(i)) / RP_STEP(i)) * type_size;
3481  }
3482 
3483  // create basic data type
3484  mpiRet = MPI_Type_contiguous(contiguous_size, dataType[0], &dataType[1]);
3485 
3486  // total size
3487  total_size
3488  = ((RP_UB(i) - RP_LB(i)) / RP_STEP(i) + 1) * type_size;
3489 
3490  // free MPI_Datatype out of use
3491  MPI_Type_free(&dataType[0]);
3492 
3493  // on error in MPI_Type_contiguous
3494  if (mpiRet != MPI_SUCCESS) { return 1015; }
3495 
3496  {
3497  int byte_datatype1; MPI_Aint lb_datatype1, extent_datatype1;
3498  MPI_Type_size(dataType[1], &byte_datatype1);
3499  MPI_Type_get_extent(dataType[1], &lb_datatype1, &extent_datatype1);
3500  if (extent_datatype1 + space_size > total_size){
3501  _XMP_fatal("xmp_file_set_view_all (block): data type is incorrect");
3502  }
3503  }
3504 
3505  // create new file type
3506  mpiRet = MPI_TYPE_CREATE_RESIZED1(dataType[1],
3507  space_size,
3508  total_size,
3509  &dataType[0]);
3510 
3511  // on error in MPI_Type_create_resized1
3512  if (mpiRet != MPI_SUCCESS) { return 1016; }
3513 
3514  // free MPI_Datatype out of use
3515  MPI_Type_free(&dataType[1]);
3516 
3517 #ifdef DEBUG
3518 fprintf(stderr, "VIEW(%d/%d) ALIGN_BLOCK\n", rank, nproc);
3519 fprintf(stderr, "VIEW(%d/%d) contiguous_size=%ld\n", rank, nproc, contiguous_size);
3520 fprintf(stderr, "VIEW(%d/%d) space_size=%ld\n", rank, nproc, space_size);
3521 fprintf(stderr, "VIEW(%d/%d) total_size=%ld\n", rank, nproc, total_size);
3522 fprintf(stderr, "VIEW(%d/%d) (lower,upper)=(%d,%d)\n", rank, nproc, lower, upper);
3523 #endif
3524  }
3525  }
3526  // cyclic or block-cyclic distribution
3527  else if (align_manner_i == _XMP_N_ALIGN_CYCLIC ||
3528  align_manner_i == _XMP_N_ALIGN_BLOCK_CYCLIC)
3529  {
3530  int bw_i = xmp_align_size(apd, i+1);
3531  if (bw_i <= 0){
3532  _XMP_fatal("xmp_file_set_view_all: invalid block width");
3533  return 1021;
3534  }else if(align_manner_i == _XMP_N_ALIGN_CYCLIC && bw_i != 1){
3535  _XMP_fatal("xmp_file_set_view_all: invalid block width for cyclic distribution");
3536  return 1021;
3537  }
3538  int cycle_i = xmp_dist_stride(tempd, i+1);
3539  int ierr = _xmp_io_set_view_block_cyclic(par_lower_i /* in */, par_upper_i /* in */, bw_i /* in */, cycle_i /* in */,
3540  RP_LB(i) /* in */, RP_UB(i) /* in */, RP_STEP(i) /* in */,
3541  dataType[0] /* in */,
3542  &dataType[1] /* out */);
3543  if (ierr != MPI_SUCCESS) { return -1017; }
3544  MPI_Type_free(&dataType[0]);
3545  dataType[0] = dataType[1];
3546  }
3547  // other
3548  else
3549  {
3550  _XMP_fatal("xmp_file_set_view_all: invalid align manner");
3551  return 1018;
3552  } /* align_manner_i */
3553  }
3554 
3555  // commit
3556  mpiRet = MPI_Type_commit(&dataType[0]);
3557 
3558  // on erro in commit
3559  if (mpiRet != MPI_SUCCESS) { return 1019; }
3560 
3561  // set view
3562  {
3563  int byte_dataType0;
3564  MPI_Type_size(dataType[0], &byte_dataType0);
3565 #ifdef DEBUG
3566  fprintf(stderr, "set_view_all: myrank=%d: byte_dataType0=%d\n", rank, byte_dataType0);
3567 #endif /* DEBUG */
3568  if (byte_dataType0 > 0){
3569  mpiRet = MPI_File_set_view(pstXmp_file->fh,
3570  (MPI_Offset)disp,
3571  MPI_BYTE,
3572  dataType[0],
3573  "native",
3574  MPI_INFO_NULL);
3575  }else{
3576  mpiRet = MPI_File_set_view(pstXmp_file->fh,
3577  (MPI_Offset)disp,
3578  MPI_BYTE,
3579  MPI_BYTE, /* dummy */
3580  "native",
3581  MPI_INFO_NULL);
3582  }
3583  }
3584  // free MPI_Datatype out of use
3585  MPI_Type_free(&dataType[0]);
3586 
3587  // on erro in set view
3588  if (mpiRet != MPI_SUCCESS) { return 1020; }
3589 
3590 #ifdef CHECK_POINT
3591  fprintf(stderr, "IO:END (xmp_file_set_view_all): rank=%d\n", rank);
3592 #endif /* CHECK_POINT */
3593 
3594  return 0;
3595 #undef RP_DIMS
3596 #undef RP_LB
3597 #undef RP_UB
3598 #undef RP_STEP
3599 }
Here is the call graph for this function:

◆ xmp_file_sync_all()

long long xmp_file_sync_all ( xmp_file_t pstXmp_file)
1689 {
1690  // check argument
1691  if (pstXmp_file == NULL) { return 1; }
1692 
1693  // sync
1694  if (MPI_File_sync(pstXmp_file->fh) != MPI_SUCCESS)
1695  {
1696  return 1;
1697  }
1698 
1699  // barrier
1700  MPI_Barrier(MPI_COMM_WORLD);
1701 
1702  return 0;
1703 }

◆ xmp_fopen_all()

xmp_file_t* xmp_fopen_all ( const char *  fname,
const char *  amode 
)

mode analysis

1385 {
1386  xmp_file_t *pstXmp_file = NULL;
1387  int iMode = 0;
1388  size_t modelen = 0;
1389 #ifdef CHECK_POINT
1390  fprintf(stderr, "IO:START(xmp_fopen_all)\n");
1391 #endif /* CHECK_POINT */
1392  // allocate
1393  pstXmp_file = malloc(sizeof(xmp_file_t));
1394  if (pstXmp_file == NULL) { return NULL; }
1395  memset(pstXmp_file, 0x00, sizeof(xmp_file_t));
1396 
1400  modelen = strlen(amode);
1401  // mode has single character
1402  if (modelen == 1)
1403  {
1404  if (strncmp(amode, "r", modelen) == 0)
1405  {
1406  iMode = MPI_MODE_RDONLY;
1407  }
1408  else if (strncmp(amode, "w", modelen) == 0)
1409  {
1410  iMode = (MPI_MODE_WRONLY | MPI_MODE_CREATE);
1411  }
1412  else if (strncmp(amode, "a", modelen) == 0)
1413  {
1414  iMode = (MPI_MODE_RDWR | MPI_MODE_CREATE | MPI_MODE_APPEND);
1415  pstXmp_file->is_append = 0x01;
1416  }
1417  else
1418  {
1419  goto ErrorExit;
1420  }
1421  }
1422  // mode has two characters
1423  else if (modelen == 2)
1424  {
1425  if (strncmp(amode, "r+", modelen) == 0)
1426  {
1427  iMode = MPI_MODE_RDWR;
1428  }
1429  else if (strncmp(amode, "w+", modelen) == 0)
1430  {
1431  iMode = (MPI_MODE_RDWR | MPI_MODE_CREATE);
1432  }
1433  else if (strncmp(amode, "a+", modelen) == 0 ||
1434  strncmp(amode, "ra", modelen) == 0 ||
1435  strncmp(amode, "ar", modelen) == 0)
1436  {
1437  iMode = (MPI_MODE_RDWR | MPI_MODE_CREATE);
1438  pstXmp_file->is_append = 0x01;
1439  }
1440  else if (strncmp(amode, "rw", modelen) == 0 ||
1441  strncmp(amode, "wr", modelen) == 0)
1442  {
1443  goto ErrorExit;
1444  }
1445  else
1446  {
1447  goto ErrorExit;
1448  }
1449  }
1450  // mode has more than two characters
1451  else
1452  {
1453  goto ErrorExit;
1454  }
1455 
1456  // file open
1457  if (MPI_File_open(MPI_COMM_WORLD,
1458  (char*)fname,
1459  iMode,
1460  MPI_INFO_NULL,
1461  &(pstXmp_file->fh)) != MPI_SUCCESS)
1462  {
1463  goto ErrorExit;
1464  }
1465 
1466  // if "W" or "W+", then set file size to zero
1467  if ((iMode == (MPI_MODE_WRONLY | MPI_MODE_CREATE) ||
1468  iMode == (MPI_MODE_RDWR | MPI_MODE_CREATE)) &&
1469  pstXmp_file->is_append == 0x00)
1470  {
1471  if (MPI_File_set_size(pstXmp_file->fh, 0) != MPI_SUCCESS)
1472  {
1473  goto ErrorExit;
1474  }
1475  }
1476 
1477  // normal return
1478  return pstXmp_file;
1479 
1480 // on error
1481 ErrorExit:
1482  if (pstXmp_file != NULL)
1483  {
1484  free(pstXmp_file);
1485  }
1486 #ifdef CHECK_POINT
1487  fprintf(stderr, "IO:END (xmp_fopen_all)\n");
1488 #endif /* CHECK_POINT */
1489  return NULL;
1490 }

◆ xmp_fread()

ssize_t xmp_fread ( xmp_file_t pstXmp_file,
void *  buffer,
size_t  size,
size_t  count 
)
3128 {
3129  MPI_Status status;
3130  int readCount;
3131 
3132  // check argument
3133  if (pstXmp_file == NULL) { return -1; }
3134  if (buffer == NULL) { return -1; }
3135  if (size < 1) { return -1; }
3136  if (count < 1) { return -1; }
3137 
3138  // read
3139  if (MPI_File_read(pstXmp_file->fh, buffer, size * count, MPI_BYTE, &status) != MPI_SUCCESS)
3140  {
3141  return -1;
3142  }
3143 
3144  // number of bytes read
3145  if (MPI_Get_count(&status, MPI_BYTE, &readCount) != MPI_SUCCESS)
3146  {
3147  return -1;
3148  }
3149 
3150  return readCount;
3151 }

◆ xmp_fread_all()

ssize_t xmp_fread_all ( xmp_file_t pstXmp_file,
void *  buffer,
size_t  size,
size_t  count 
)
1720 {
1721  MPI_Status status;
1722  int readCount;
1723 #ifdef CHECK_POINT
1724  fprintf(stderr, "IO:START(xmp_fread_all)\n");
1725 #endif /* CHECK_POINT */
1726  // check argument
1727  if (pstXmp_file == NULL) { return -1; }
1728  if (buffer == NULL) { return -1; }
1729  if (size < 1) { return -1; }
1730  if (count < 1) { return -1; }
1731 
1732  // read
1733  if (MPI_File_read_all(pstXmp_file->fh,
1734  buffer, size * count,
1735  MPI_BYTE,
1736  &status) != MPI_SUCCESS)
1737  {
1738 #ifdef CHECK_POINT
1739  fprintf(stderr, "IO:END (xmp_fread_all)\n");
1740 #endif /* CHECK_POINT */
1741  return -1;
1742  }
1743 
1744  // number of bytes read
1745  if (MPI_Get_count(&status, MPI_BYTE, &readCount) != MPI_SUCCESS)
1746  {
1747 #ifdef CHECK_POINT
1748  fprintf(stderr, "IO:END (xmp_fread_all)\n");
1749 #endif /* CHECK_POINT */
1750  return -1;
1751  }
1752 
1753 #ifdef CHECK_POINT
1754  fprintf(stderr, "IO:END (xmp_fread_all)\n");
1755 #endif /* CHECK_POINT */
1756  return readCount;
1757 }

◆ xmp_fread_darray_all()

ssize_t xmp_fread_darray_all ( xmp_file_t pstXmp_file,
xmp_desc_t  apd,
xmp_range_t rp 
)
2094 {
2095  MPI_Status status; // MPI status
2096  int readCount; // read bytes
2097  int mpiRet; // return value of MPI functions
2098  long contiguous_size; // contiguous size
2099  long space_size; // space size
2100  long total_size; // total size
2101  MPI_Aint tmp1, type_size;
2102  MPI_Datatype dataType[2];
2103  int i = 0;
2104  xmp_desc_t tempd;
2105  int rp_dims;
2106  int *rp_lb_addr = NULL;
2107  int *rp_ub_addr = NULL;
2108  int *rp_step_addr = NULL;
2109  int array_ndims;
2110  size_t array_type_size;
2111  int typesize_int;
2112  //int ierr;
2113 
2114  int rank, nproc;
2115  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
2116  MPI_Comm_size(MPI_COMM_WORLD, &nproc);
2117 
2118 #ifdef CHECK_POINT
2119  fprintf(stderr, "IO:START(xmp_fread_darray_all): rank=%d\n", rank);
2120 #endif /* CHECK_POINT */
2121 
2122  // check argument
2123  if (pstXmp_file == NULL) { return -1; }
2124  if (apd == NULL) { return -1; }
2125  if (rp == NULL) { return -1; }
2126 
2127  /*ierr =*/ xmp_align_template(apd, &tempd);
2128  if (tempd == NULL){ return -1; }
2129  /*ierr =*/ xmp_array_ndims(apd, &array_ndims);
2130  array_type_size = xmp_array_type_size(apd);
2131 
2132  rp_dims = _xmp_range_get_dims(rp);
2133  rp_lb_addr = _xmp_range_get_lb_addr(rp);
2134  rp_ub_addr = _xmp_range_get_ub_addr(rp);
2135  rp_step_addr = _xmp_range_get_step_addr(rp);
2136  if (!rp_lb_addr || !rp_ub_addr || !rp_step_addr){ return -1; }
2137 #define RP_DIMS (rp_dims)
2138 #define RP_LB(i) (rp_lb_addr[(i)])
2139 #define RP_UB(i) (rp_ub_addr[(i)])
2140 #define RP_STEP(i) (rp_step_addr[(i)])
2141 
2142  // check number of dimensions
2143  if (array_ndims != RP_DIMS) { return -1; }
2144 
2145  /* case unpack is required */
2146  for (i = RP_DIMS - 1; i >= 0; i--){
2147  if(RP_STEP(i) < 0){
2148  int ret = xmp_fread_darray_unpack(pstXmp_file, apd, rp);
2149  return ret;
2150  }
2151  }
2152 
2153 #ifdef DEBUG
2154 fprintf(stderr, "READ(%d/%d) dims=%d\n", rank, nproc, RP_DIMS);
2155 #endif
2156 
2157  // create basic data type
2158  MPI_Type_contiguous(array_type_size, MPI_BYTE, &dataType[0]);
2159 
2160  // loop for each dimension
2161  for (i = RP_DIMS - 1; i >= 0; i--)
2162  {
2163  int par_lower_i = xmp_array_gcllbound(apd, i+1);
2164  int par_upper_i = xmp_array_gclubound_tmp(apd, i+1);
2165  int align_manner_i = xmp_align_format(apd, i+1);
2166  int local_lower_i = xmp_array_lcllbound(apd, i+1);
2167  int alloc_size_i;
2168 
2169  /*ierr =*/ xmp_array_lsize(apd, i+1, &alloc_size_i);
2170 #ifdef DEBUG
2171 fprintf(stderr, "READ(%d/%d) (lb,ub,step)=(%d,%d,%d)\n",
2172  rank, nproc, RP_LB(i), RP_UB(i), RP_STEP(i));
2173 fprintf(stderr, "READ(%d/%d) (par_lower,par_upper)=(%d,%d)\n",
2174  rank, nproc, par_lower_i, par_upper_i);
2175 #endif
2176  // no distribution
2177  if (align_manner_i == _XMP_N_ALIGN_NOT_ALIGNED ||
2178  align_manner_i == _XMP_N_ALIGN_DUPLICATION)
2179  {
2180  // upper after distribution < lower
2181  if (par_upper_i < RP_LB(i)) { return -1; }
2182  // lower after distribution > upper
2183  if (par_lower_i > RP_UB(i)) { return -1; }
2184 
2185  // incremnet is negative
2186  if ( RP_STEP(i) < 0)
2187  {
2188  }
2189  // incremnet is positive
2190  else
2191  {
2192  // contiguous size
2193  contiguous_size = (RP_UB(i) - RP_LB(i)) / RP_STEP(i) + 1;
2194 
2195  // get extent of data type
2196  mpiRet =MPI_Type_get_extent(dataType[0], &tmp1, &type_size);
2197  if (mpiRet != MPI_SUCCESS) { return -1; }
2198 
2199  // create basic data type
2200  mpiRet = MPI_Type_create_hvector(contiguous_size,
2201  1,
2202  type_size * RP_STEP(i),
2203  dataType[0],
2204  &dataType[1]);
2205 
2206  // free MPI_Datatype out of use
2207  MPI_Type_free(&dataType[0]);
2208 
2209  // on error in MPI_Type_create_hvector
2210  if (mpiRet != MPI_SUCCESS) { return -1; }
2211 
2212  // total size
2213  total_size
2214  = (par_upper_i
2215  - par_lower_i + 1)
2216  * type_size;
2217 
2218  // space size
2219  space_size
2220  = (RP_LB(i) - par_lower_i)
2221  * type_size;
2222 
2223  // create new file type
2224  mpiRet = MPI_TYPE_CREATE_RESIZED1(dataType[1],
2225  (MPI_Aint)space_size,
2226  (MPI_Aint)total_size,
2227  &dataType[0]);
2228 
2229  // on error in MPI_Type_create_resized1
2230  if (mpiRet != MPI_SUCCESS) { return -1; }
2231 
2232  // free MPI_Datatype out of use
2233  MPI_Type_free(&dataType[1]);
2234 
2235 #ifdef DEBUG
2236 fprintf(stderr, "READ(%d/%d) NOT_ALIGNED\n", rank, nproc);
2237 fprintf(stderr, "READ(%d/%d) contiguous_size=%d\n", rank, nproc, contiguous_size);
2238 fprintf(stderr, "READ(%d/%d) space_size=%d\n", rank, nproc, space_size);
2239 fprintf(stderr, "READ(%d/%d) total_size=%d\n", rank, nproc, total_size);
2240 #endif
2241  }
2242  }
2243  // block distribution
2244  else if (align_manner_i == _XMP_N_ALIGN_BLOCK) {
2245  // increment is negative
2246  if ( RP_STEP(i) < 0) { }
2247  // increment is positive
2248  else {
2249  int lower, upper;
2250  // get extent of data type
2251  mpiRet =MPI_Type_get_extent(dataType[0], &tmp1, &type_size);
2252  if (mpiRet != MPI_SUCCESS) { return -1; }
2253 
2254  // upper after distribution < lower
2255  if (par_upper_i < RP_LB(i)) {
2256  contiguous_size = space_size = 0;
2257  }
2258  // lower after distribution > upper
2259  else if (par_lower_i > RP_UB(i)) {
2260  contiguous_size = space_size = 0;
2261  }
2262  // other
2263  else {
2264  // lower in this node
2265  lower = (par_lower_i > RP_LB(i)) ?
2266  RP_LB(i) + ((par_lower_i - 1 - RP_LB(i)) / RP_STEP(i) + 1) * RP_STEP(i)
2267  : RP_LB(i);
2268 
2269  // upper in this node
2270  upper = (par_upper_i < RP_UB(i)) ?
2271  par_upper_i
2272  : RP_UB(i);
2273 
2274  // contiguous size
2275  contiguous_size = (upper - lower + RP_STEP(i)) / RP_STEP(i);
2276 
2277  // space size
2278  space_size = (local_lower_i + (lower - par_lower_i)) * type_size;
2279  }
2280 
2281  // create basic data type
2282  mpiRet = MPI_Type_create_hvector(contiguous_size,
2283  1,
2284  type_size * RP_STEP(i),
2285  dataType[0],
2286  &dataType[1]);
2287 
2288  // free MPI_Datatype out of use
2289  MPI_Type_free(&dataType[0]);
2290 
2291  // on error in MPI_Type_create_hvector
2292  if (mpiRet != MPI_SUCCESS) { return -1; }
2293 
2294  // total size
2295  total_size = (alloc_size_i)* type_size;
2296 
2297  // create new file type
2298  mpiRet = MPI_TYPE_CREATE_RESIZED1(dataType[1],
2299  (MPI_Aint)space_size,
2300  (MPI_Aint)total_size,
2301  &dataType[0]);
2302 
2303  // on error in MPI_Type_create_resized1
2304  if (mpiRet != MPI_SUCCESS) { return -1; }
2305 
2306  // free MPI_Datatype out of use
2307  MPI_Type_free(&dataType[1]);
2308 
2309 #ifdef DEBUG
2310  fprintf(stderr, "fread_darray_all: rank = %d: space_size = %d total_size = %d\n",
2311  rank,space_size,total_size);
2312 fprintf(stderr, "READ(%d/%d) ALIGN_BLOCK\n", rank, nproc);
2313 fprintf(stderr, "READ(%d/%d) contiguous_size=%d\n", rank, nproc, contiguous_size);
2314 fprintf(stderr, "READ(%d/%d) space_size=%d\n", rank, nproc, space_size);
2315 fprintf(stderr, "READ(%d/%d) total_size=%d\n", rank, nproc, total_size);
2316 fprintf(stderr, "READ(%d/%d) (lower,upper)=(%d,%d)\n", rank, nproc, lower, upper);
2317 #endif
2318  }
2319  }
2320  // cyclic or block-cyclic distribution
2321  else if (align_manner_i == _XMP_N_ALIGN_CYCLIC ||
2322  align_manner_i == _XMP_N_ALIGN_BLOCK_CYCLIC) {
2323  int bw_i = xmp_align_size(apd, i+1);
2324  if (bw_i <= 0) {
2325  _XMP_fatal("xmp_fread_darray_all: invalid block width");
2326  return -1;
2327  } else if(align_manner_i == _XMP_N_ALIGN_CYCLIC && bw_i != 1) {
2328  _XMP_fatal("xmp_fread_darray_all: invalid block width for cyclic distribution");
2329  return -1;
2330  }
2331  int cycle_i = xmp_dist_stride(tempd, i+1);
2332  int ierr = _xmp_io_write_read_block_cyclic(par_lower_i /* in */,
2333  par_upper_i /* in */,
2334  bw_i /* in */,
2335  cycle_i /* in */,
2336  RP_LB(i) /* in */,
2337  RP_UB(i) /* in */,
2338  RP_STEP(i) /* in */,
2339  local_lower_i /* in */,
2340  alloc_size_i /* in */,
2341  dataType[0] /* in */,
2342  &dataType[1] /* out */);
2343  if (ierr != MPI_SUCCESS) { return -1; }
2344  MPI_Type_free(&dataType[0]);
2345  dataType[0] = dataType[1];
2346  }
2347  // other
2348  else {
2349  _XMP_fatal("xmp_fread_darray_all: invalid align manner");
2350  return -1;
2351  } /* align_manner_i */
2352  }
2353 
2354  // commit
2355  mpiRet = MPI_Type_commit(&dataType[0]);
2356 
2357  // on erro in commit
2358  if (mpiRet != MPI_SUCCESS) { return 1; }
2359 
2360  char *array_addr;
2361  xmp_array_laddr(apd, (void **)&array_addr);
2362 
2363  // read
2364  MPI_Type_size(dataType[0], &typesize_int);
2365 #ifdef DEBUG
2366  fprintf(stderr, "fread_darray_all: rank=%d: typesize_int = %d\n",rank,typesize_int);
2367 #endif /* DEBUG */
2368 
2369  if(typesize_int > 0){
2370  if (MPI_File_read_all(pstXmp_file->fh,
2371  array_addr,
2372  1,
2373  dataType[0],
2374  &status)
2375  != MPI_SUCCESS){ return -1; }
2376  }else{
2377  if (MPI_File_read_all(pstXmp_file->fh,
2378  array_addr,
2379  0, /* dummy */
2380  MPI_BYTE, /* dummy */
2381  &status)
2382  != MPI_SUCCESS){ return -1; }
2383  }
2384 
2385 #ifdef DEBUG
2386  fprintf(stderr, "CP(aft fread_darray_all) [%d/%d]\n", rank, nproc);
2387 #endif
2388  // free MPI_Datatype out of use
2389  MPI_Type_free(&dataType[0]);
2390 
2391  // number of bytes read
2392  if (MPI_Get_count(&status, MPI_BYTE, &readCount) != MPI_SUCCESS)
2393  {
2394 #ifdef CHECK_POINT
2395  fprintf(stderr, "IO:END (xmp_fread_darray_all): rank=%d\n", rank);
2396 #endif /* CHECK_POINT */
2397  return -1;
2398  }
2399 #ifdef DEBUG
2400  fprintf(stderr, "CP(finish: xmp_fread_darray_all) [%d/%d]\n", rank, nproc);
2401 #endif
2402 #ifdef CHECK_POINT
2403  fprintf(stderr, "IO:END (xmp_fread_darray_all): rank=%d\n", rank);
2404 #endif /* CHECK_POINT */
2405  return readCount;
2406 #undef RP_DIMS
2407 #undef RP_LB
2408 #undef RP_UB
2409 #undef RP_STEP
2410 }
Here is the call graph for this function:

◆ xmp_fread_darray_unpack()

int xmp_fread_darray_unpack ( xmp_file_t fp,
xmp_desc_t  apd,
xmp_range_t rp 
)
1841 {
1842  MPI_Status status;
1843  char *array_addr=NULL;
1844  char *buf=NULL;
1845  char *cp;
1846  int *lb=NULL;
1847  int *ub=NULL;
1848  int *step=NULL;
1849  int *cnt=NULL;
1850  long buf_size, j;
1851  int ret=0;
1852  long disp;
1853  long size;
1854  long array_size;
1855  int i;
1856  xmp_desc_t tempd = NULL;
1857  int **bc2_result = NULL;
1858  size_t array_type_size;
1859  int rp_dims = 0;
1860  int *rp_lb_addr = NULL;
1861  int *rp_ub_addr = NULL;
1862  int *rp_step_addr = NULL;
1863  int array_ndims;
1864  //int ierr;
1865 
1866 #ifdef CHECK_POINT
1867  fprintf(stderr, "IO:START(xmp_fread_darray_unpack)\n");
1868 #endif /* CHECK_POINT */
1869 
1870  // check argument
1871  if (fp == NULL){ ret = -1; goto FunctionExit; }
1872  if (apd == NULL){ ret = -1; goto FunctionExit; }
1873  if (rp == NULL){ ret = -1; goto FunctionExit; }
1874 
1875  /*ierr = */xmp_align_template(apd, &tempd);
1876  if (tempd == NULL){ ret = -1; goto FunctionExit; }
1877  /*ierr =*/ xmp_array_ndims(apd, &array_ndims);
1878 
1879  rp_dims = _xmp_range_get_dims(rp);
1880  rp_lb_addr = _xmp_range_get_lb_addr(rp);
1881  rp_ub_addr = _xmp_range_get_ub_addr(rp);
1882  rp_step_addr = _xmp_range_get_step_addr(rp);
1883  if (!rp_lb_addr || !rp_ub_addr || !rp_step_addr){ ret = -1; goto FunctionExit; }
1884 #define RP_DIMS (rp_dims)
1885 #define RP_LB(i) (rp_lb_addr[(i)])
1886 #define RP_UB(i) (rp_ub_addr[(i)])
1887 #define RP_STEP(i) (rp_step_addr[(i)])
1888 
1889  // check number of dimensions
1890  if (array_ndims != RP_DIMS){ ret = -1; goto FunctionExit; }
1891 
1892  /* allocate arrays for the number of rotations */
1893  lb = (int*)malloc(sizeof(int)*RP_DIMS);
1894  ub = (int*)malloc(sizeof(int)*RP_DIMS);
1895  step = (int*)malloc(sizeof(int)*RP_DIMS);
1896  cnt = (int*)malloc(sizeof(int)*RP_DIMS);
1897  if(!lb || !ub || !step || !cnt){
1898  ret = -1;
1899  goto FunctionExit;
1900  }
1901  bc2_result = (int**)malloc(sizeof(int*)*RP_DIMS);
1902  if(!bc2_result){
1903  ret = -1;
1904  goto FunctionExit;
1905  }
1906  for(i=0; i<RP_DIMS; i++){ bc2_result[i]=NULL; }
1907 
1908  /* calculate the number of rotations */
1909  buf_size = 1;
1910  for(i=0; i<RP_DIMS; i++){
1911  int par_lower_i = xmp_array_gcllbound(apd, i+1);
1912  int par_upper_i = xmp_array_gclubound_tmp(apd, i+1);
1913  int align_manner_i = xmp_align_format(apd, i+1);
1914  /* error check */
1915  if(RP_STEP(i) > 0 && RP_LB(i) > RP_UB(i)){
1916  ret = -1;
1917  goto FunctionExit;
1918  }
1919  if(RP_STEP(i) < 0 && RP_LB(i) < RP_UB(i)){
1920  ret = -1;
1921  goto FunctionExit;
1922  }
1923  if (align_manner_i == _XMP_N_ALIGN_NOT_ALIGNED ||
1924  align_manner_i == _XMP_N_ALIGN_DUPLICATION) {
1925  lb[i] = RP_LB(i);
1926  ub[i] = RP_UB(i);
1927  step[i] = RP_STEP(i);
1928  cnt[i] = (ub[i]-lb[i]+step[i])/step[i];
1929  } else if(align_manner_i == _XMP_N_ALIGN_BLOCK){
1930  if(RP_STEP(i) > 0){
1931  if(par_upper_i < RP_LB(i) ||
1932  par_lower_i > RP_UB(i)){
1933  lb[i] = 1;
1934  ub[i] = 0;
1935  step[i] = 1;
1936  } else {
1937  lb[i] = (par_lower_i > RP_LB(i))?
1938  RP_LB(i)+((par_lower_i-1-RP_LB(i))/RP_STEP(i)+1)*RP_STEP(i):
1939  RP_LB(i);
1940  ub[i] = (par_upper_i < RP_UB(i)) ?
1941  par_upper_i:
1942  RP_UB(i);
1943  step[i] = RP_STEP(i);
1944  }
1945  } else {
1946  if(par_upper_i < RP_UB(i) ||
1947  par_lower_i > RP_LB(i)){
1948  lb[i] = 1;
1949  ub[i] = 0;
1950  step[i] = 1;
1951  } else {
1952  lb[i] = (par_upper_i < RP_LB(i))?
1953  RP_LB(i)-((RP_LB(i)-par_upper_i-1)/RP_STEP(i)-1)*RP_STEP(i):
1954  RP_LB(i);
1955  ub[i] = (par_lower_i > RP_UB(i))?
1956  par_lower_i:
1957  RP_UB(i);
1958  step[i] = RP_STEP(i);
1959  }
1960  }
1961  cnt[i] = (ub[i]-lb[i]+step[i])/step[i];
1962 
1963  } else if(align_manner_i == _XMP_N_ALIGN_CYCLIC||
1964  align_manner_i == _XMP_N_ALIGN_BLOCK_CYCLIC){
1965  int bw_i = xmp_align_size(apd, i+1);
1966  if (bw_i <= 0){
1967  _XMP_fatal("xmp_fread_darray_unpack: invalid block width");
1968  ret = -1; goto FunctionExit;
1969  }else if(align_manner_i == _XMP_N_ALIGN_CYCLIC && bw_i != 1){
1970  _XMP_fatal("xmp_fread_darray_unpack: invalid block width for cyclic distribution");
1971  ret = -1; goto FunctionExit;
1972  }
1973  int cycle_i = xmp_dist_stride(tempd, i+1);
1974  int zzcnt; int *zzptr;
1975  int ierr;
1976  ierr = _xmp_io_pack_unpack_block_cyclic_aux1(par_lower_i /* in */, par_upper_i /* in */, bw_i /* in */, cycle_i /* in */,
1977  RP_LB(i) /* in */, RP_UB(i) /* in */, RP_STEP(i) /* in */,
1978  &zzcnt /* out */, &zzptr /* out */);
1979  if (ierr != MPI_SUCCESS){ ret = -1; goto FunctionExit; }
1980  cnt[i] = zzcnt;
1981  bc2_result[i] = zzptr;
1982 
1983  } else {
1984  ret = -1;
1985  goto FunctionExit;
1986  }
1987  cnt[i] = (cnt[i]>0)? cnt[i]: 0;
1988  buf_size *= cnt[i];
1989  }
1990 
1991  array_type_size = xmp_array_type_size(apd);
1992 
1993  /* allocate buffer */
1994  if(buf_size == 0){
1995  buf = (char*)malloc(array_type_size);
1996  } else {
1997  buf = (char*)malloc(buf_size * array_type_size);
1998  }
1999  if(!buf){
2000  ret = -1;
2001  goto FunctionExit;
2002  }
2003 
2004  // read
2005  if(buf_size > 0){
2006  if (MPI_File_read_all(fp->fh, buf, buf_size * array_type_size, MPI_BYTE, &status) != MPI_SUCCESS){
2007  ret = -1;
2008  goto FunctionExit;
2009  }
2010  }else{
2011  if (MPI_File_read_all(fp->fh, buf, 0, MPI_BYTE, &status) != MPI_SUCCESS){
2012  ret = -1;
2013  goto FunctionExit;
2014  }
2015  }
2016  // number of bytes written
2017  if (MPI_Get_count(&status, MPI_BYTE, &ret) != MPI_SUCCESS){
2018  ret = -1;
2019  goto FunctionExit;
2020  }
2021 
2022  /* unpack data */
2023  cp = buf;
2024  for(j=0; j<buf_size; j++){
2025  disp = 0;
2026  size = 1;
2027  array_size = 1;
2028  for(i=RP_DIMS-1; i>=0; i--){
2029  int par_lower_i = xmp_array_gcllbound(apd, i+1);
2030  int align_manner_i = xmp_align_format(apd, i+1);
2031  int ser_size_i = xmp_array_gsize(apd, i+1);
2032  int local_lower_i = xmp_array_lcllbound(apd, i+1);
2033  ub[i] = (j/size)%cnt[i];
2034  if (align_manner_i == _XMP_N_ALIGN_NOT_ALIGNED ||
2035  align_manner_i == _XMP_N_ALIGN_DUPLICATION) {
2036  disp += (lb[i]+ub[i]*step[i])*array_size;
2037  array_size *= ser_size_i;
2038 
2039  } else if(align_manner_i == _XMP_N_ALIGN_BLOCK){
2040  disp += (lb[i] + ub[i]*step[i] + local_lower_i - par_lower_i)*array_size;
2041 
2042  } else if(align_manner_i == _XMP_N_ALIGN_CYCLIC ||
2043  align_manner_i == _XMP_N_ALIGN_BLOCK_CYCLIC){
2044  int local_index;
2045  int ierr = _xmp_io_pack_unpack_block_cyclic_aux2(ub[i], bc2_result[i],
2046  &local_index);
2047  if (ierr != MPI_SUCCESS){ ret = -1; goto FunctionExit; }
2048  disp += (local_index + local_lower_i) * array_size;
2049  } /* align_manner_i */
2050  size *= cnt[i];
2051  } /* i */
2052  disp *= array_type_size;
2053  memcpy(array_addr+disp, cp, array_type_size);
2054  cp += array_type_size;
2055  } /* j */
2056 
2057  FunctionExit:
2058  if(buf) free(buf);
2059  if(lb) free(lb);
2060  if(ub) free(ub);
2061  if(step) free(step);
2062  if(cnt) free(cnt);
2063  if(bc2_result){
2064  for(i=0; i<RP_DIMS; i++){ if(bc2_result[i]){ free(bc2_result[i]); } }
2065  }
2066 
2067 #ifdef CHECK_POINT
2068  fprintf(stderr, "IO:END (xmp_fread_darray_unpack)\n");
2069 #endif /* CHECK_POINT */
2070  return ret;
2071 #undef RP_DIMS
2072 #undef RP_LB
2073 #undef RP_UB
2074 #undef RP_STEP
2075 }
Here is the call graph for this function:

◆ xmp_fread_shared()

ssize_t xmp_fread_shared ( xmp_file_t pstXmp_file,
void *  buffer,
size_t  size,
size_t  count 
)
3031 {
3032  MPI_Status status;
3033  int readCount;
3034 
3035  // check argument
3036  if (pstXmp_file == NULL) { return -1; }
3037  if (buffer == NULL) { return -1; }
3038  if (size < 1) { return -1; }
3039  if (count < 1) { return -1; }
3040 
3041  // read
3042  if (MPI_File_read_shared(pstXmp_file->fh, buffer, size * count, MPI_BYTE, &status) != MPI_SUCCESS)
3043  {
3044  return -1;
3045  }
3046 
3047  // number of bytes read
3048  if (MPI_Get_count(&status, MPI_BYTE, &readCount) != MPI_SUCCESS)
3049  {
3050  return -1;
3051  }
3052 
3053  return readCount;
3054 }

◆ xmp_free_range()

void xmp_free_range ( xmp_range_t rp)
1331 {
1332 #ifdef CHECK_POINT
1333  fprintf(stderr, "IO:START(xmp_free_range)\n");
1334 #endif /* CHECK_POINT */
1335  if (rp == NULL){
1336  return;
1337  }else{
1338  if (rp->lb){ free(rp->lb); }
1339  if (rp->ub){ free(rp->ub); }
1340  if (rp->step){ free(rp->step); }
1341  free(rp);
1342  }
1343 #ifdef CHECK_POINT
1344  fprintf(stderr, "IO:END (xmp_free_range)\n");
1345 #endif /* CHECK_POINT */
1346 }

◆ xmp_fseek()

int xmp_fseek ( xmp_file_t pstXmp_file,
long long  offset,
int  whence 
)
1539 {
1540  int iMpiWhence;
1541 
1542  // checkk argument
1543  if (pstXmp_file == NULL) { return 1; }
1544 
1545  // convert offset to MPI_Offset
1546  switch (whence)
1547  {
1548  case SEEK_SET:
1549  iMpiWhence = MPI_SEEK_SET;
1550  break;
1551  case SEEK_CUR:
1552  iMpiWhence = MPI_SEEK_CUR;
1553  break;
1554  case SEEK_END:
1555  iMpiWhence = MPI_SEEK_END;
1556  break;
1557  default:
1558  return 1;
1559  }
1560 
1561  // file seek
1562  if (MPI_File_seek(pstXmp_file->fh, (MPI_Offset)offset, iMpiWhence) != MPI_SUCCESS)
1563  {
1564  return 1;
1565  }
1566 
1567  return 0;
1568 }

◆ xmp_fseek_shared()

int xmp_fseek_shared ( xmp_file_t pstXmp_file,
long long  offset,
int  whence 
)
1583 {
1584  int iMpiWhence;
1585 
1586  // check argument
1587  if (pstXmp_file == NULL) { return 1; }
1588 
1589  // convert offset to MPI_Offset
1590  switch (whence)
1591  {
1592  case SEEK_SET:
1593  iMpiWhence = MPI_SEEK_SET;
1594  break;
1595  case SEEK_CUR:
1596  iMpiWhence = MPI_SEEK_CUR;
1597  break;
1598  case SEEK_END:
1599  iMpiWhence = MPI_SEEK_END;
1600  break;
1601  default:
1602  return 1;
1603  }
1604 
1605  // file seek
1606  if (MPI_File_seek_shared(pstXmp_file->fh, (MPI_Offset)offset, iMpiWhence) != MPI_SUCCESS)
1607  {
1608  return 1;
1609  }
1610 
1611  return 0;
1612 }

◆ xmp_ftell()

long long xmp_ftell ( xmp_file_t pstXmp_file)
1626 {
1627  MPI_Offset offset;
1628  MPI_Offset disp;
1629 
1630  // check argument
1631  if (pstXmp_file == NULL) { return -1; }
1632 
1633  if (MPI_File_get_position(pstXmp_file->fh, &offset) != MPI_SUCCESS)
1634  {
1635  return -1;
1636  }
1637 
1638  if (MPI_File_get_byte_offset(pstXmp_file->fh, offset, &disp) != MPI_SUCCESS)
1639  {
1640  return -1;
1641  }
1642 
1643  return (long long)disp;
1644 }

◆ xmp_ftell_shared()

long long xmp_ftell_shared ( xmp_file_t pstXmp_file)
1658 {
1659  MPI_Offset offset;
1660  MPI_Offset disp;
1661 
1662  // check argument
1663  if (pstXmp_file == NULL) { return -1; }
1664 
1665  if (MPI_File_get_position_shared(pstXmp_file->fh, &offset) != MPI_SUCCESS)
1666  {
1667  return -1;
1668  }
1669 
1670  if (MPI_File_get_byte_offset(pstXmp_file->fh, offset, &disp) != MPI_SUCCESS)
1671  {
1672  return -1;
1673  }
1674 
1675  return (long long)disp;
1676 }

◆ xmp_fwrite()

ssize_t xmp_fwrite ( xmp_file_t pstXmp_file,
void *  buffer,
size_t  size,
size_t  count 
)
3168 {
3169  MPI_Status status;
3170  int writeCount;
3171 
3172  // check argument
3173  if (pstXmp_file == NULL) { return -1; }
3174  if (buffer == NULL) { return -1; }
3175  if (size < 1) { return -1; }
3176  if (count < 1) { return -1; }
3177 
3178  // if file open is "r+", then move pointer to end
3179  if (pstXmp_file->is_append)
3180  {
3181  if (MPI_File_seek(pstXmp_file->fh,
3182  (MPI_Offset)0,
3183  MPI_SEEK_END) != MPI_SUCCESS)
3184  {
3185  return -1;
3186  }
3187 
3188  pstXmp_file->is_append = 0x00;
3189  }
3190 
3191  // write
3192  if (MPI_File_write(pstXmp_file->fh,
3193  buffer,
3194  size * count,
3195  MPI_BYTE,
3196  &status) != MPI_SUCCESS)
3197  {
3198  return -1;
3199  }
3200 
3201  // number of bytes written
3202  if (MPI_Get_count(&status, MPI_BYTE, &writeCount) != MPI_SUCCESS)
3203  {
3204  return -1;
3205  }
3206 
3207  return writeCount;
3208 }

◆ xmp_fwrite_all()

ssize_t xmp_fwrite_all ( xmp_file_t pstXmp_file,
void *  buffer,
size_t  size,
size_t  count 
)
1775 {
1776  MPI_Status status;
1777  int writeCount;
1778 #ifdef CHECK_POINT
1779  fprintf(stderr, "IO:START(xmp_fwrite_all)\n");
1780 #endif /* CHECK_POINT */
1781  // check argument
1782  if (pstXmp_file == NULL) { return -1; }
1783  if (buffer == NULL) { return -1; }
1784  if (size < 1) { return -1; }
1785  if (count < 1) { return -1; }
1786 
1787  // if file open is "r+", then move pointer to end
1788  if (pstXmp_file->is_append)
1789  {
1790  if (MPI_File_seek(pstXmp_file->fh,
1791  (MPI_Offset)0,
1792  MPI_SEEK_END) != MPI_SUCCESS)
1793  {
1794  return -1;
1795  }
1796 
1797  pstXmp_file->is_append = 0x00;
1798  }
1799 
1800  // write
1801  if (MPI_File_write_all(pstXmp_file->fh, buffer, size * count, MPI_BYTE, &status) != MPI_SUCCESS)
1802  {
1803 #ifdef CHECK_POINT
1804  fprintf(stderr, "IO:END (xmp_fwrite_all)\n");
1805 #endif /* CHECK_POINT */
1806  return -1;
1807  }
1808 
1809  // number of bytes written
1810  if (MPI_Get_count(&status, MPI_BYTE, &writeCount) != MPI_SUCCESS)
1811  {
1812 #ifdef CHECK_POINT
1813  fprintf(stderr, "IO:END (xmp_fwrite_all)\n");
1814 #endif /* CHECK_POINT */
1815  return -1;
1816  }
1817 #ifdef CHECK_POINT
1818  fprintf(stderr, "IO:END (xmp_fwrite_all)\n");
1819 #endif /* CHECK_POINT */
1820  return writeCount;
1821 }

◆ xmp_fwrite_darray_all()

ssize_t xmp_fwrite_darray_all ( xmp_file_t pstXmp_file,
xmp_desc_t  apd,
xmp_range_t rp 
)
2688 {
2689  MPI_Status status; // MPI status
2690  int writeCount; // write btye
2691  int mpiRet; // return value of MPI functions
2692  long contiguous_size; // contiguous size
2693  long space_size; // space size
2694  long total_size; // total size
2695  MPI_Aint tmp1, type_size;
2696  MPI_Datatype dataType[2];
2697  int i = 0;
2698  xmp_desc_t tempd;
2699  int rp_dims;
2700  int *rp_lb_addr = NULL;
2701  int *rp_ub_addr = NULL;
2702  int *rp_step_addr = NULL;
2703  int array_ndims;
2704  size_t array_type_size;
2705  int typesize_int;
2706  //int ierr;
2707 
2708  int rank, nproc;
2709  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
2710  MPI_Comm_size(MPI_COMM_WORLD, &nproc);
2711 
2712 #ifdef CHECK_POINT
2713  fprintf(stderr, "IO:START(xmp_fwrite_darray_all): rank=%d\n", rank);
2714 #endif /* CHECK_POINT */
2715 
2716  // check argument
2717  if (pstXmp_file == NULL) { return -1101; }
2718  if (apd == NULL) { return -1103; }
2719  if (rp == NULL) { return -1104; }
2720 
2721  /*ierr =*/ xmp_align_template(apd, &tempd);
2722  if (tempd == NULL){ return -1105; }
2723  /*ierr =*/ xmp_array_ndims(apd, &array_ndims);
2724  array_type_size = xmp_array_type_size(apd);
2725 
2726  rp_dims = _xmp_range_get_dims(rp);
2727  rp_lb_addr = _xmp_range_get_lb_addr(rp);
2728  rp_ub_addr = _xmp_range_get_ub_addr(rp);
2729  rp_step_addr = _xmp_range_get_step_addr(rp);
2730  if (!rp_lb_addr || !rp_ub_addr || !rp_step_addr){ return -1106; }
2731 #define RP_DIMS (rp_dims)
2732 #define RP_LB(i) (rp_lb_addr[(i)])
2733 #define RP_UB(i) (rp_ub_addr[(i)])
2734 #define RP_STEP(i) (rp_step_addr[(i)])
2735 
2736  // check number of dimensions
2737  if (array_ndims != RP_DIMS) { return -1107; }
2738 
2739 #ifdef DEBUG
2740 fprintf(stderr, "WRITE(%d/%d) dims=%d\n",rank, nproc, RP_DIMS);
2741 #endif
2742 
2743  /* case pack is required */
2744  for (i = RP_DIMS - 1; i >= 0; i--){
2745  if(RP_STEP(i) < 0){
2746  int ret = xmp_fwrite_darray_pack(pstXmp_file, apd, rp);
2747  return ret;
2748  }
2749  }
2750 
2751  // create basic data type
2752  MPI_Type_contiguous(array_type_size, MPI_BYTE, &dataType[0]);
2753 
2754  // loop for each dimension
2755  for (i = RP_DIMS - 1; i >= 0; i--)
2756  {
2757  int par_lower_i = xmp_array_gcllbound(apd, i+1);
2758  int par_upper_i = xmp_array_gclubound_tmp(apd, i+1);
2759  int align_manner_i = xmp_align_format(apd, i+1);
2760  int ser_lower_i = xmp_array_gcglbound(apd, i+1);
2761  int ser_upper_i = xmp_array_gcgubound(apd, i+1);
2762  int local_lower_i = xmp_array_lcllbound(apd, i+1);
2763  int alloc_size_i;
2764  //int ierr;
2765  /*ierr =*/ xmp_array_lsize(apd, i+1, &alloc_size_i);
2766 /* int local_upper_i = xmp_array_lclubound(apd, i+1); */
2767 /* int shadow_size_lo_i = xmp_array_lshadow(apd, i+1); */
2768 /* int shadow_size_hi_i = xmp_array_ushadow(apd, i+1); */
2769 #ifdef DEBUG
2770 fprintf(stderr, "WRITE(%d/%d) (lb,ub,step)=(%d,%d,%d)\n",
2771  rank, nproc, RP_LB(i), RP_UB(i), RP_STEP(i));
2772 fprintf(stderr, "WRITE(%d/%d) (par_lower,par_upper)=(%d,%d)\n",
2773  rank, nproc, par_lower_i, par_upper_i);
2774 /* fprintf(stderr, "WRITE(%d/%d) (local_lower,local_upper,alloc_size)=(%d,%d,%d)\n", */
2775 /* rank, nproc, local_lower_i, local_upper_i, alloc_size_i); */
2776 /* fprintf(stderr, "WRITE(%d/%d) (shadow_size_lo,shadow_size_hi)=(%d,%d)\n", */
2777 /* rank, nproc, shadow_size_lo_i, shadow_size_hi_i); */
2778 #endif
2779 
2780  // no distribution
2781  if (align_manner_i == _XMP_N_ALIGN_NOT_ALIGNED ||
2782  align_manner_i == _XMP_N_ALIGN_DUPLICATION)
2783  {
2784  // upper after distribution < lower
2785  if (par_upper_i < RP_LB(i)) { return -1108; }
2786  // lower after distribution > upper
2787  if (par_lower_i > RP_UB(i)) { return -1109; }
2788 
2789  // incremnet is negative
2790  if ( RP_STEP(i) < 0)
2791  {
2792  }
2793  // incremnet is positive
2794  else
2795  {
2796  // contiguous size
2797  contiguous_size = (RP_UB(i) - RP_LB(i)) / RP_STEP(i) + 1;
2798 
2799  // get extent of data type
2800  mpiRet =MPI_Type_get_extent(dataType[0], &tmp1, &type_size);
2801  if (mpiRet != MPI_SUCCESS) { return -1110; }
2802 
2803  // create basic data type
2804  mpiRet = MPI_Type_create_hvector(contiguous_size,
2805  1,
2806  type_size * RP_STEP(i),
2807  dataType[0],
2808  &dataType[1]);
2809 
2810  // free MPI_Datatype out of use
2811  MPI_Type_free(&dataType[0]);
2812 
2813  // on error in MPI_Type_contiguous
2814  if (mpiRet != MPI_SUCCESS) { return -1111; }
2815 
2816  // total size
2817  total_size
2818  = (ser_upper_i
2819  - ser_lower_i + 1)
2820  * type_size;
2821 
2822  // space size
2823  space_size
2824  = (RP_LB(i) - par_lower_i)
2825  * type_size;
2826 
2827  // create new file type
2828  mpiRet = MPI_TYPE_CREATE_RESIZED1(dataType[1],
2829  (MPI_Aint)space_size,
2830  (MPI_Aint)total_size,
2831  &dataType[0]);
2832 
2833  // on error in MPI_Type_create_resized1
2834  if (mpiRet != MPI_SUCCESS) { return -1112; }
2835 
2836  // free MPI_Datatype out of use
2837  MPI_Type_free(&dataType[1]);
2838 
2839 #ifdef DEBUG
2840 fprintf(stderr, "WRITE(%d/%d) NOT_ALIGNED\n",rank, nproc);
2841 fprintf(stderr, "WRITE(%d/%d) type_size=%ld\n",rank, nproc, (long)type_size);
2842 fprintf(stderr, "WRITE(%d/%d) contiguous_size=%ld\n",rank, nproc, contiguous_size);
2843 fprintf(stderr, "WRITE(%d/%d) space_size=%ld\n",rank, nproc, space_size);
2844 fprintf(stderr, "WRITE(%d/%d) total_size=%ld\n",rank, nproc, total_size);
2845 #endif
2846  }
2847  }
2848  // block distribution
2849  else if (align_manner_i == _XMP_N_ALIGN_BLOCK)
2850  {
2851  // increment is negative
2852  if ( RP_STEP(i) < 0)
2853  {
2854  }
2855  // increment is positive
2856  else
2857  {
2858  int lower, upper;
2859  // get extent of data type
2860  mpiRet =MPI_Type_get_extent(dataType[0], &tmp1, &type_size);
2861  if (mpiRet != MPI_SUCCESS) { return -1113; }
2862 
2863  // upper after distribution < lower
2864  if (par_upper_i < RP_LB(i))
2865  {
2866  contiguous_size = space_size = 0;
2867  }
2868  // lower after distribution > upper
2869  else if (par_lower_i > RP_UB(i)) {
2870  contiguous_size = space_size = 0;
2871  }
2872  // other
2873  else {
2874  // lower in this node
2875  lower = (par_lower_i > RP_LB(i)) ?
2876  RP_LB(i) + ((par_lower_i - 1 - RP_LB(i)) / RP_STEP(i) + 1) * RP_STEP(i)
2877  : RP_LB(i);
2878 
2879  // upper in this node
2880  upper = (par_upper_i < RP_UB(i)) ?
2881  par_upper_i
2882  : RP_UB(i);
2883 
2884  // contiguous size
2885  contiguous_size = (upper - lower + RP_STEP(i)) / RP_STEP(i);
2886 
2887  if(lower > upper){ type_size = 0; }
2888  // space size
2889  space_size = (local_lower_i + (lower - par_lower_i)) * type_size;
2890  }
2891 
2892  // create basic data type
2893  mpiRet = MPI_Type_create_hvector(contiguous_size,
2894  1,
2895  type_size * RP_STEP(i),
2896  dataType[0],
2897  &dataType[1]);
2898 
2899  // free MPI_Datatype out of use
2900  MPI_Type_free(&dataType[0]);
2901 
2902  // on error in MPI_Type_create_hvector
2903  if (mpiRet != MPI_SUCCESS) { return -1114; }
2904 
2905  // total size
2906  total_size = (alloc_size_i)* type_size;
2907 
2908  // create new file type
2909  mpiRet = MPI_TYPE_CREATE_RESIZED1(dataType[1],
2910  (MPI_Aint)space_size,
2911  (MPI_Aint)total_size,
2912  &dataType[0]);
2913 
2914  // on error in MPI_Type_create_resized1
2915  if (mpiRet != MPI_SUCCESS) { return -1115; }
2916 
2917  // free MPI_Datatype out of use
2918  MPI_Type_free(&dataType[1]);
2919 
2920 #ifdef DEBUG
2921 fprintf(stderr, "WRITE(%d/%d) ALIGN_BLOCK\n",rank, nproc);
2922 fprintf(stderr, "WRITE(%d/%d) type_size=%ld\n",rank, nproc, (long)type_size);
2923 fprintf(stderr, "WRITE(%d/%d) contiguous_size=%ld\n",rank, nproc, contiguous_size);
2924 fprintf(stderr, "WRITE(%d/%d) space_size=%ld\n",rank, nproc, space_size);
2925 fprintf(stderr, "WRITE(%d/%d) total_size=%ld\n",rank, nproc, total_size);
2926 fprintf(stderr, "WRITE(%d/%d) (lower,upper)=(%d,%d)\n",rank, nproc, lower, upper);
2927 #endif
2928  }
2929  }
2930  // cyclic or block-cyclic distribution
2931  else if (align_manner_i == _XMP_N_ALIGN_CYCLIC ||
2932  align_manner_i == _XMP_N_ALIGN_BLOCK_CYCLIC)
2933  {
2934  int bw_i = xmp_align_size(apd, i+1);
2935  if (bw_i <= 0){
2936  _XMP_fatal("xmp_fwrite_darray_all: invalid block width");
2937  return -1122;
2938  }else if(align_manner_i == _XMP_N_ALIGN_CYCLIC && bw_i != 1){
2939  _XMP_fatal("xmp_fwrite_darray_all: invalid block width for cyclic distribution");
2940  return -1122;
2941  }
2942  int cycle_i = xmp_dist_stride(tempd, i+1);
2943  int ierr = _xmp_io_write_read_block_cyclic(par_lower_i /* in */, par_upper_i /* in */, bw_i /* in */, cycle_i /* in */,
2944  RP_LB(i) /* in */, RP_UB(i) /* in */, RP_STEP(i) /* in */,
2945  local_lower_i /* in */,
2946  alloc_size_i /* in */,
2947  dataType[0] /* in */,
2948  &dataType[1] /* out */);
2949  if (ierr != MPI_SUCCESS) { return -1117; }
2950  MPI_Type_free(&dataType[0]);
2951  dataType[0] = dataType[1];
2952  }
2953  // other
2954  else
2955  {
2956  _XMP_fatal("xmp_fwrite_darray_all: invalid align manner");
2957  return -1118;
2958  } /* align_manner_i */
2959  }
2960 
2961  // commit
2962  mpiRet = MPI_Type_commit(&dataType[0]);
2963 
2964  // on error in commit
2965  if (mpiRet != MPI_SUCCESS) { return 1119; }
2966 
2967  char *array_addr;
2968  xmp_array_laddr(apd, (void **)&array_addr);
2969 
2970  // write
2971  MPI_Type_size(dataType[0], &typesize_int);
2972 #ifdef DEBUG
2973  fprintf(stderr, "fwrite_darray_all: rank=%d: typesize_int = %d\n",rank,typesize_int);
2974 #endif /* DEBUG */
2975  {
2976  if(typesize_int > 0){
2977  if ((MPI_File_write_all(pstXmp_file->fh,
2978  array_addr,
2979  1,
2980  dataType[0],
2981  &status))
2982  != MPI_SUCCESS){ return -1120; }
2983  }else{
2984  if ((MPI_File_write_all(pstXmp_file->fh,
2985  array_addr,
2986  0, /* dummy */
2987  MPI_BYTE, /* dummy */
2988  &status))
2989  != MPI_SUCCESS){ return -1120; }
2990  }
2991  }
2992  // free MPI_Datatype out of use
2993  MPI_Type_free(&dataType[0]);
2994 
2995  // number of btyes written
2996  if (MPI_Get_count(&status, MPI_BYTE, &writeCount) != MPI_SUCCESS)
2997  {
2998 #ifdef CHECK_POINT
2999  fprintf(stderr, "IO:END (xmp_fwrite_darray_all): rank=%d\n", rank);
3000 #endif /* CHECK_POINT */
3001  return -1121;
3002  }
3003 #ifdef DEBUG
3004  if(rank==0){fprintf(stderr, "-------------------- fwrite_darray_all: NORMAL END\n");}
3005 #endif /* DEBUG */
3006 #ifdef CHECK_POINT
3007  fprintf(stderr, "IO:END (xmp_fwrite_darray_all): rank=%d\n", rank);
3008 #endif /* CHECK_POINT */
3009  return writeCount;
3010 #undef RP_DIMS
3011 #undef RP_LB
3012 #undef RP_UB
3013 #undef RP_STEP
3014 }
Here is the call graph for this function:

◆ xmp_fwrite_darray_pack()

int xmp_fwrite_darray_pack ( xmp_file_t fp,
xmp_desc_t  apd,
xmp_range_t rp 
)
2430 {
2431  MPI_Status status;
2432  char *array_addr;
2433  char *buf=NULL;
2434  char *cp;
2435  int *lb=NULL;
2436  int *ub=NULL;
2437  int *step=NULL;
2438  int *cnt=NULL;
2439  long buf_size, j;
2440  int ret=0;
2441  long disp;
2442  long size;
2443  long array_size;
2444  int i;
2445  size_t array_type_size;
2446  xmp_desc_t tempd = NULL;
2447  int **bc2_result = NULL;
2448  int rp_dims = 0;
2449  int *rp_lb_addr = NULL;
2450  int *rp_ub_addr = NULL;
2451  int *rp_step_addr = NULL;
2452  int array_ndims;
2453  //int ierr;
2454 
2455  int myrank, nprocs;
2456  MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
2457  MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
2458 
2459 #ifdef CHECK_POINT
2460  fprintf(stderr, "IO:START(xmp_fwrite_darray_pack): rank=%d\n", myrank);
2461 #endif /* CHECK_POINT */
2462 
2463  /*ierr =*/ xmp_align_template(apd, &tempd);
2464  if (tempd == NULL){ ret = -1; goto FunctionExit; }
2465  array_type_size = xmp_array_type_size(apd);
2466  /*ierr =*/ xmp_array_ndims(apd, &array_ndims);
2467 
2468  rp_dims = _xmp_range_get_dims(rp);
2469  rp_lb_addr = _xmp_range_get_lb_addr(rp);
2470  rp_ub_addr = _xmp_range_get_ub_addr(rp);
2471  rp_step_addr = _xmp_range_get_step_addr(rp);
2472  if (!rp_lb_addr || !rp_ub_addr || !rp_step_addr){ ret = -1; goto FunctionExit; }
2473 #define RP_DIMS (rp_dims)
2474 #define RP_LB(i) (rp_lb_addr[(i)])
2475 #define RP_UB(i) (rp_ub_addr[(i)])
2476 #define RP_STEP(i) (rp_step_addr[(i)])
2477 
2478  // check number of dimensions
2479  if (array_ndims != RP_DIMS){ ret = -1; goto FunctionExit; }
2480 
2481  /* allocate arrays for the number of rotations */
2482  lb = (int*)malloc(sizeof(int)*RP_DIMS);
2483  ub = (int*)malloc(sizeof(int)*RP_DIMS);
2484  step = (int*)malloc(sizeof(int)*RP_DIMS);
2485  cnt = (int*)malloc(sizeof(int)*RP_DIMS);
2486  if(!lb || !ub || !step || !cnt){
2487  ret = -1;
2488  goto FunctionExit;
2489  }
2490  bc2_result = (int**)malloc(sizeof(int*)*RP_DIMS);
2491  if(!bc2_result){
2492  ret = -1;
2493  goto FunctionExit;
2494  }
2495  for(i=0; i<RP_DIMS; i++){ bc2_result[i]=NULL; }
2496 
2497  /* calculate the number of rotaions */
2498  buf_size = 1;
2499  for(i=0; i<RP_DIMS; i++){
2500  int par_lower_i = xmp_array_gcllbound(apd, i+1);
2501  int par_upper_i = xmp_array_gclubound_tmp(apd, i+1);
2502  int align_manner_i = xmp_align_format(apd, i+1);
2503 
2504  /* error check */
2505  if(RP_STEP(i) > 0 && RP_LB(i) > RP_UB(i)){
2506  ret = -1;
2507  goto FunctionExit;
2508  }
2509  if(RP_STEP(i) < 0 && RP_LB(i) < RP_UB(i)){
2510  ret = -1;
2511  goto FunctionExit;
2512  }
2513  if (align_manner_i == _XMP_N_ALIGN_NOT_ALIGNED ||
2514  align_manner_i == _XMP_N_ALIGN_DUPLICATION) {
2515  lb[i] = RP_LB(i);
2516  ub[i] = RP_UB(i);
2517  step[i] = RP_STEP(i);
2518  cnt[i] = (ub[i]-lb[i]+step[i])/step[i];
2519 
2520  } else if(align_manner_i == _XMP_N_ALIGN_BLOCK){
2521  if(RP_STEP(i) > 0){
2522  if(par_upper_i < RP_LB(i) ||
2523  par_lower_i > RP_UB(i)){
2524  lb[i] = 1;
2525  ub[i] = 0;
2526  step[i] = 1;
2527  } else {
2528  lb[i] = (par_lower_i > RP_LB(i))?
2529  RP_LB(i)+((par_lower_i-1-RP_LB(i))/RP_STEP(i)+1)*RP_STEP(i):
2530  RP_LB(i);
2531  ub[i] = (par_upper_i < RP_UB(i)) ?
2532  par_upper_i:
2533  RP_UB(i);
2534  step[i] = RP_STEP(i);
2535  }
2536  } else {
2537  if(par_upper_i < RP_UB(i) ||
2538  par_lower_i > RP_LB(i)){
2539  lb[i] = 1;
2540  ub[i] = 0;
2541  step[i] = 1;
2542  } else {
2543  lb[i] = (par_upper_i < RP_LB(i))?
2544  RP_LB(i)-((RP_LB(i)-par_upper_i-1)/RP_STEP(i)-1)*RP_STEP(i):
2545  RP_LB(i);
2546  ub[i] = (par_lower_i > RP_UB(i))?
2547  par_lower_i:
2548  RP_UB(i);
2549  step[i] = RP_STEP(i);
2550  }
2551  }
2552  cnt[i] = (ub[i]-lb[i]+step[i])/step[i];
2553 
2554  } else if(align_manner_i == _XMP_N_ALIGN_CYCLIC ||
2555  align_manner_i == _XMP_N_ALIGN_BLOCK_CYCLIC){
2556  int bw_i = xmp_align_size(apd, i+1);
2557  if (bw_i <= 0){
2558  _XMP_fatal("xmp_fwrite_darray_pack: invalid block width");
2559  ret = -1; goto FunctionExit;
2560  }else if(align_manner_i == _XMP_N_ALIGN_CYCLIC && bw_i != 1){
2561  _XMP_fatal("xmp_fwrite_darray_pack: invalid block width for cyclic distribution");
2562  ret = -1; goto FunctionExit;
2563  }
2564  int cycle_i = xmp_dist_stride(tempd, i+1);
2565  int zzcnt; int *zzptr;
2566  int ierr;
2567  ierr = _xmp_io_pack_unpack_block_cyclic_aux1(par_lower_i /* in */, par_upper_i /* in */, bw_i /* in */, cycle_i /* in */,
2568  RP_LB(i) /* in */, RP_UB(i) /* in */, RP_STEP(i) /* in */,
2569  &zzcnt /* out */, &zzptr /* out */);
2570  if (ierr != MPI_SUCCESS){ ret = -1; goto FunctionExit; }
2571  cnt[i] = zzcnt;
2572  bc2_result[i] = zzptr;
2573 
2574  } else {
2575  ret = -1;
2576  goto FunctionExit;
2577  }
2578  cnt[i] = (cnt[i]>0)? cnt[i]: 0;
2579  buf_size *= cnt[i];
2580  }
2581 
2582  /* allocate buffer */
2583  if(buf_size == 0){
2584  buf = (char*)malloc(array_type_size);
2585  } else {
2586  buf = (char*)malloc(buf_size * array_type_size);
2587  }
2588  if(!buf){
2589  ret = -1;
2590  goto FunctionExit;
2591  }
2592 
2593  /* pack data */
2594  cp = buf;
2595  xmp_array_laddr(apd, (void **)&array_addr);
2596  for(j=0; j<buf_size; j++){
2597  disp = 0;
2598  size = 1;
2599  array_size = 1;
2600  for(i=RP_DIMS-1; i>=0; i--){
2601  int par_lower_i = xmp_array_gcllbound(apd, i+1);
2602  int align_manner_i = xmp_align_format(apd, i+1);
2603  int local_lower_i = xmp_array_lcllbound(apd, i+1);
2604  int ser_size_i = xmp_array_gsize(apd, i+1);
2605  int alloc_size_i;
2606  //int ierr;
2607  /*ierr =*/ xmp_array_lsize(apd, i+1, &alloc_size_i);
2608  ub[i] = (j/size)%cnt[i];
2609  if (align_manner_i == _XMP_N_ALIGN_NOT_ALIGNED ||
2610  align_manner_i == _XMP_N_ALIGN_DUPLICATION) {
2611  disp += (lb[i]+ub[i]*step[i])*array_size;
2612  array_size *= ser_size_i;
2613 
2614  } else if(align_manner_i == _XMP_N_ALIGN_BLOCK){
2615  disp += (lb[i]+ub[i]*step[i] + local_lower_i - par_lower_i)*array_size;
2616  array_size *= alloc_size_i;
2617 
2618  } else if(align_manner_i == _XMP_N_ALIGN_CYCLIC ||
2619  align_manner_i == _XMP_N_ALIGN_BLOCK_CYCLIC){
2620  int local_index;
2621  int ierr = _xmp_io_pack_unpack_block_cyclic_aux2(ub[i] /* in */, bc2_result[i] /* in */,
2622  &local_index /* out */);
2623  if (ierr != MPI_SUCCESS){ ret = -1; goto FunctionExit; }
2624  disp += (local_index + local_lower_i) * array_size;
2625  array_size *= alloc_size_i;
2626  } /* align_manner_i */
2627  size *= cnt[i];
2628  } /* i */
2629  disp *= array_type_size;
2630  memcpy(cp, array_addr+disp, array_type_size);
2631  cp += array_type_size;
2632  } /* j */
2633 
2634  // write
2635  if(buf_size > 0){
2636  if (MPI_File_write_all(fp->fh, buf, buf_size * array_type_size, MPI_BYTE, &status) != MPI_SUCCESS){
2637  ret = -1;
2638  goto FunctionExit;
2639  }
2640  }else{
2641  if (MPI_File_write_all(fp->fh, buf, 0, MPI_BYTE, &status) != MPI_SUCCESS){
2642  ret = -1;
2643  goto FunctionExit;
2644  }
2645  }
2646  // number of bytes written
2647  if (MPI_Get_count(&status, MPI_BYTE, &ret) != MPI_SUCCESS) {
2648  ret = -1;
2649  goto FunctionExit;
2650  }
2651 
2652  FunctionExit:
2653  if(buf) free(buf);
2654  if(lb) free(lb);
2655  if(ub) free(ub);
2656  if(step) free(step);
2657  if(cnt) free(cnt);
2658  if(bc2_result){
2659  for(i=0; i<RP_DIMS; i++){ if(bc2_result[i]){ free(bc2_result[i]); } }
2660  }
2661 #ifdef CHECK_POINT
2662  fprintf(stderr, "IO:END (xmp_fwrite_darray_pack): rank=%d\n", myrank);
2663 #endif /* CHECK_POINT */
2664  return ret;
2665 #undef RP_DIMS
2666 #undef RP_LB
2667 #undef RP_UB
2668 #undef RP_STEP
2669 }

◆ xmp_fwrite_shared()

ssize_t xmp_fwrite_shared ( xmp_file_t pstXmp_file,
void *  buffer,
size_t  size,
size_t  count 
)
3071 {
3072  MPI_Status status;
3073  int writeCount;
3074 
3075  // check argument
3076  if (pstXmp_file == NULL) { return -1; }
3077  if (buffer == NULL) { return -1; }
3078  if (size < 1) { return -1; }
3079  if (count < 1) { return -1; }
3080 
3081  // if file open is "r+", then move pointer to end
3082  if (pstXmp_file->is_append)
3083  {
3084  if (MPI_File_seek_shared(pstXmp_file->fh,
3085  (MPI_Offset)0,
3086  MPI_SEEK_END) != MPI_SUCCESS)
3087  {
3088  return -1;
3089  }
3090 
3091  pstXmp_file->is_append = 0x00;
3092  }
3093 
3094  // write
3095  if (MPI_File_write_shared(pstXmp_file->fh,
3096  buffer,
3097  size * count,
3098  MPI_BYTE,
3099  &status) != MPI_SUCCESS)
3100  {
3101  return -1;
3102  }
3103 
3104  // number of bytes written
3105  if (MPI_Get_count(&status, MPI_BYTE, &writeCount) != MPI_SUCCESS)
3106  {
3107  return -1;
3108  }
3109 
3110  return writeCount;
3111 }

◆ xmp_set_range()

void xmp_set_range ( xmp_range_t rp,
int  i_dim,
int  lb,
int  length,
int  step 
)
1299 {
1300 #ifdef CHECK_POINT
1301  fprintf(stderr, "IO:START(xmp_set_range)\n");
1302 #endif /* CHECK_POINT */
1303  if (rp == NULL){ _XMP_fatal("xmp_set_range: descriptor is NULL"); }
1304  if (step == 0){ _XMP_fatal("xmp_set_range: step must be non-zero"); }
1305  if (i_dim-1 < 0 || i_dim-1 >= rp->dims){ _XMP_fatal("xmp_set_range: i_dim is out of range"); }
1306  if(!rp->lb || !rp->ub || !rp->step){ _XMP_fatal("xmp_set_range: null pointer"); }
1307  if (step != 0 && length == 0){ _XMP_fatal("xmp_set_range: invalid combination of length and step\n"); }
1308  if (length < 0){ _XMP_fatal("xmp_set_range: length must be >= 0\n"); }
1309  rp->lb[i_dim-1] = lb;
1310  rp->step[i_dim-1] = step;
1311  if(step > 0){
1312  rp->ub[i_dim-1] = lb + length - 1;
1313  }else if (step < 0){
1314  rp->ub[i_dim-1] = lb - length + 1;
1315  }else{
1316  _XMP_fatal("xmp_set_range: step must be non-zero");
1317  }
1318 #ifdef CHECK_POINT
1319  fprintf(stderr, "IO:END (xmp_set_range)\n");
1320 #endif /* CHECK_POINT */
1321 }
Here is the call graph for this function:
MPI_TYPE_CREATE_RESIZED1
#define MPI_TYPE_CREATE_RESIZED1
Definition: xmp_io.c:12
xmp_range_t::ub
int * ub
Definition: xmp_io.h:19
xmp_range_t::step
int * step
Definition: xmp_io.h:20
xmp_array_type_size
size_t xmp_array_type_size(xmp_desc_t d)
Definition: xmp_lib.c:117
xmp_align_size
int xmp_align_size(xmp_desc_t d, int dim)
Definition: xmp_lib.c:262
xmp_range_t::lb
int * lb
Definition: xmp_io.h:18
xmp_file_t::is_append
char is_append
Definition: xmp_io.h:12
RP_STEP
#define RP_STEP(i)
xmp_desc_t
void * xmp_desc_t
Definition: xmp.h:29
RP_DIMS
#define RP_DIMS
xmp_array_laddr
int xmp_array_laddr(xmp_desc_t d, void **laddr)
Definition: xmp_lib.c:173
xmp_array_ndims
int xmp_array_ndims(xmp_desc_t d, int *ndims)
Definition: xmp_lib.c:96
_XMP_world_rank
int _XMP_world_rank
Definition: xmp_world.c:9
xmp_range_t
Definition: xmp_io.h:15
xmp_array_gsize
int xmp_array_gsize(xmp_desc_t d, int dim)
Definition: xmp_lib.c:123
xmp_align_format
int xmp_align_format(xmp_desc_t d, int dim)
Definition: xmp_lib.c:256
xmp_array_gcglbound
int xmp_array_gcglbound(xmp_desc_t d, int dim)
Definition: xmp_lib.c:161
xmp_array_lcllbound
int xmp_array_lcllbound(xmp_desc_t d, int dim)
Definition: xmp_lib.c:149
xmp_array_gcgubound
int xmp_array_gcgubound(xmp_desc_t d, int dim)
Definition: xmp_lib.c:167
xmp_file_t::fh
MPI_File fh
Definition: xmp_io.h:10
xmp_file_t
Definition: xmp_io.h:9
_XMP_N_ALIGN_CYCLIC
#define _XMP_N_ALIGN_CYCLIC
Definition: xmp_constant.h:38
_XMP_N_ALIGN_BLOCK
#define _XMP_N_ALIGN_BLOCK
Definition: xmp_constant.h:37
_XMP_N_ALIGN_DUPLICATION
#define _XMP_N_ALIGN_DUPLICATION
Definition: xmp_constant.h:36
_XMP_N_ALIGN_NOT_ALIGNED
#define _XMP_N_ALIGN_NOT_ALIGNED
Definition: xmp_constant.h:35
RP_LB
#define RP_LB(i)
_XMP_N_ALIGN_BLOCK_CYCLIC
#define _XMP_N_ALIGN_BLOCK_CYCLIC
Definition: xmp_constant.h:39
xmp_array_lsize
int xmp_array_lsize(xmp_desc_t d, int dim, int *lsize)
Definition: xmp_lib.c:129
xmp_align_template
int xmp_align_template(xmp_desc_t d, xmp_desc_t *dt)
Definition: xmp_lib.c:323
xmp_fwrite_darray_pack
int xmp_fwrite_darray_pack(xmp_file_t *fp, xmp_desc_t apd, xmp_range_t *rp)
Definition: xmp_io.c:2426
xmp_dist_stride
int xmp_dist_stride(xmp_desc_t d, int dim)
Definition: xmp_lib.c:409
xmp_range_t::dims
int dims
Definition: xmp_io.h:17
_XMP_fatal
void _XMP_fatal(char *msg)
Definition: xmp_util.c:42
RP_UB
#define RP_UB(i)
xmp_array_gcllbound
int xmp_array_gcllbound(xmp_desc_t d, int dim)
Definition: xmp_lib.c:137
xmp_fread_darray_unpack
int xmp_fread_darray_unpack(xmp_file_t *fp, xmp_desc_t apd, xmp_range_t *rp)
Definition: xmp_io.c:1837