Skip to content

Commit

Permalink
correct expected result
Browse files Browse the repository at this point in the history
  • Loading branch information
jedwards4b committed Apr 5, 2019
1 parent ee2ed3e commit 44c8819
Show file tree
Hide file tree
Showing 2 changed files with 23 additions and 22 deletions.
9 changes: 5 additions & 4 deletions src/clib/pioc_sc.c
Original file line number Diff line number Diff line change
Expand Up @@ -297,7 +297,7 @@ PIO_Offset GCDblocksize(int arrlen, const PIO_Offset *arr_in)
if (arrlen == 1)
return 1;

/* We can use the array length as the initial value.
/* We can use the array length as the initial value.
* Suppose we have n contiguous blocks with lengths
* b1, b2, ..., bn, then gcd(b1, b2, ..., bn) =
* gcd(b1 + b2 + ... + bn, b1, b2, ..., bn) =
Expand Down Expand Up @@ -330,9 +330,10 @@ PIO_Offset GCDblocksize(int arrlen, const PIO_Offset *arr_in)
blk_len = 1;
}
}

LOG((2, "bsize before last block %ld %d\n",bsize, blk_len));
/* Handle the last block. */
bsize = lgcd(bsize, blk_len);
LOG((2, "bsize after last block %ld %d\n",bsize, blk_len));

return bsize;
}
Expand All @@ -355,7 +356,7 @@ PIO_Offset GCDblocksize(int arrlen, const PIO_Offset *arr_in)
int CalcStartandCount(int pio_type, int ndims, const int *gdims, int num_io_procs,
int myiorank, PIO_Offset *start, PIO_Offset *count, int *num_aiotasks)
{
int minbytes;
int minbytes;
int maxbytes;
int minblocksize; /* Like minbytes, but in data elements. */
int basesize; /* Size in bytes of base data type. */
Expand Down Expand Up @@ -522,6 +523,6 @@ int CalcStartandCount(int pio_type, int ndims, const int *gdims, int num_io_proc

/* Return the number of IO procs used to the caller. */
*num_aiotasks = use_io_procs;

return PIO_NOERR;
}
36 changes: 18 additions & 18 deletions tests/cunit/test_spmd.c
Original file line number Diff line number Diff line change
Expand Up @@ -222,13 +222,13 @@ int test_determine_procs()
#define TWO_COMPONENTS 2
#define THREE_PROCS 3
int ret;

{
int num_io_procs = 1;
int component_count = ONE_COMPONENT;
int num_procs_per_comp[ONE_COMPONENT] = {1};
int *my_proc_list[ONE_COMPONENT];

if ((ret = determine_procs(num_io_procs, component_count, num_procs_per_comp, NULL,
my_proc_list)))
return ret;
Expand All @@ -241,17 +241,17 @@ int test_determine_procs()
free(my_proc_list[c]);
}
}

{
int num_io_procs = 3;
int component_count = TWO_COMPONENTS;
int num_procs_per_comp[TWO_COMPONENTS] = {1, 1};
int *my_proc_list[TWO_COMPONENTS];

if ((ret = determine_procs(num_io_procs, component_count, num_procs_per_comp, NULL,
my_proc_list)))
return ret;

/* Check results and free resources. */
for (int c = 0; c < TWO_COMPONENTS; c++)
{
Expand All @@ -260,17 +260,17 @@ int test_determine_procs()
free(my_proc_list[c]);
}
}

{
int num_io_procs = 3;
int component_count = TWO_COMPONENTS;
int num_procs_per_comp[TWO_COMPONENTS] = {THREE_PROCS, THREE_PROCS};
int *my_proc_list[TWO_COMPONENTS];

if ((ret = determine_procs(num_io_procs, component_count, num_procs_per_comp, NULL,
my_proc_list)))
return ret;

/* Check results and free resources. */
for (int c = 0; c < TWO_COMPONENTS; c++)
{
Expand All @@ -280,7 +280,7 @@ int test_determine_procs()
free(my_proc_list[c]);
}
}

{
int num_io_procs = 3;
int component_count = TWO_COMPONENTS;
Expand All @@ -289,11 +289,11 @@ int test_determine_procs()
int proc_list_2[THREE_PROCS] = {11, 12, 13};
int *proc_list[TWO_COMPONENTS] = {proc_list_1, proc_list_2};
int *my_proc_list[TWO_COMPONENTS];

if ((ret = determine_procs(num_io_procs, component_count, num_procs_per_comp,
(int **)proc_list, my_proc_list)))
return ret;

/* Check results and free resources. */
for (int c = 0; c < TWO_COMPONENTS; c++)
{
Expand All @@ -303,7 +303,7 @@ int test_determine_procs()
free(my_proc_list[c]);
}
}

return PIO_NOERR;
}

Expand Down Expand Up @@ -483,7 +483,7 @@ int test_varlists3()
return ERR_WRONG;
if (get_var_desc(3, &varlist, &var_desc) != PIO_ENOTVAR)
return ERR_WRONG;

return 0;
}

Expand Down Expand Up @@ -664,9 +664,10 @@ int test_CalcStartandCount()
return 0;
}

/* Test the GDCblocksize() function. */
int run_GDCblocksize_tests(MPI_Comm test_comm)
/* Test the GCDblocksize() function. */
int run_GCDblocksize_tests(MPI_Comm test_comm)
{

{
int arrlen = 1;
PIO_Offset arr_in[1] = {0};
Expand Down Expand Up @@ -703,7 +704,7 @@ int run_GDCblocksize_tests(MPI_Comm test_comm)
PIO_Offset blocksize;

blocksize = GCDblocksize(arrlen, arr_in);
if (blocksize != 1)
if (blocksize != 2)
return ERR_WRONG;
}

Expand Down Expand Up @@ -732,7 +733,6 @@ int main(int argc, char **argv)
if ((ret = pio_test_init2(argc, argv, &my_rank, &ntasks, MIN_NTASKS,
TARGET_NTASKS, -1, &test_comm)))
ERR(ERR_INIT);

/* Test code runs on TARGET_NTASKS tasks. The left over tasks do
* nothing. */
if (my_rank < TARGET_NTASKS)
Expand All @@ -746,7 +746,7 @@ int main(int argc, char **argv)
if ((ret = run_sc_tests(test_comm)))
return ret;

if ((ret = run_GDCblocksize_tests(test_comm)))
if ((ret = run_GCDblocksize_tests(test_comm)))
return ret;

if ((ret = run_spmd_tests(test_comm)))
Expand Down

0 comments on commit 44c8819

Please sign in to comment.