libxmp/libxmpf in Omni Compiler  1.3.4
xmp_gmove.c File Reference
#include <stdarg.h>
#include <stdlib.h>
#include <string.h>
#include "xmp.h"
#include "mpi.h"
#include "xmp_internal.h"
#include "xmp_math_function.h"
Include dependency graph for xmp_gmove.c:

Macros

#define MPI_PORTABLE_PLATFORM_H
 
#define _XMP_SM_GTOL_BLOCK(_i, _m, _w)   (((_i) - (_m)) % (_w))
 
#define _XMP_SM_GTOL_CYCLIC(_i, _m, _P)   (((_i) - (_m)) / (_P))
 
#define _XMP_SM_GTOL_BLOCK_CYCLIC(_b, _i, _m, _P)   (((((_i) - (_m)) / (((_P) * (_b)))) * (_b)) + (((_i) - (_m)) % (_b)))
 
#define XMP_DBG   0
 
#define DBG_RANK   3
 
#define XMP_DBG_OWNER_REGION   0
 

Functions

void * _XMP_get_array_addr (_XMP_array_t *a, int *gidx)
 
void _XMP_gtol_array_ref_triplet (_XMP_array_t *array, int dim_index, int *lower, int *upper, int *stride)
 
void _XMP_calc_gmove_rank_array_SCALAR (_XMP_array_t *array, int *ref_index, int *rank_array)
 
int _XMP_calc_gmove_array_owner_linear_rank_SCALAR (_XMP_array_t *array, int *ref_index)
 
int xmp_calc_gmove_array_owner_linear_rank_scalar_ (_XMP_array_t **a, int *ref_index)
 
unsigned long long _XMP_gmove_bcast_ARRAY (void *dst_addr, int dst_dim, int *dst_l, int *dst_u, int *dst_s, unsigned long long *dst_d, void *src_addr, int src_dim, int *src_l, int *src_u, int *src_s, unsigned long long *src_d, int type, size_t type_size, int root_rank)
 
int _XMP_check_gmove_array_ref_inclusion_SCALAR (_XMP_array_t *array, int array_index, int ref_index)
 
void _XMP_gmove_localcopy_ARRAY (int type, int type_size, void *dst_addr, int dst_dim, int *dst_l, int *dst_u, int *dst_s, unsigned long long *dst_d, void *src_addr, int src_dim, int *src_l, int *src_u, int *src_s, unsigned long long *src_d)
 
int _XMP_calc_global_index_HOMECOPY (_XMP_array_t *dst_array, int dst_dim_index, int *dst_l, int *dst_u, int *dst_s, int *src_l, int *src_u, int *src_s)
 
int _XMP_calc_global_index_BCAST (int dst_dim, int *dst_l, int *dst_u, int *dst_s, _XMP_array_t *src_array, int *src_array_nodes_ref, int *src_l, int *src_u, int *src_s)
 
void _XMP_sendrecv_ARRAY (int type, int type_size, MPI_Datatype *mpi_datatype, _XMP_array_t *dst_array, int *dst_array_nodes_ref, int *dst_lower, int *dst_upper, int *dst_stride, unsigned long long *dst_dim_acc, _XMP_array_t *src_array, int *src_array_nodes_ref, int *src_lower, int *src_upper, int *src_stride, unsigned long long *src_dim_acc)
 
void _XMP_gmove_BCAST_GSCALAR (void *dst_addr, _XMP_array_t *array, int ref_index[])
 
void _XMP_gmove_SENDRECV_GSCALAR (void *dst_addr, void *src_addr, _XMP_array_t *dst_array, _XMP_array_t *src_array, int dst_ref_index[], int src_ref_index[])
 
void _XMP_gmove_calc_unit_size (_XMP_array_t *dst_array, _XMP_array_t *src_array, unsigned long long *alltoall_unit_size, unsigned long long *dst_pack_unit_size, unsigned long long *src_pack_unit_size, unsigned long long *dst_ser_size, unsigned long long *src_ser_size, int dst_block_dim, int src_block_dim)
 
void _XMP_align_local_idx (long long int global_idx, int *local_idx, _XMP_array_t *array, int array_axis, int *rank)
 
void _XMP_gmove_array_array_common (_XMP_gmv_desc_t *gmv_desc_leftp, _XMP_gmv_desc_t *gmv_desc_rightp, int *dst_l, int *dst_u, int *dst_s, unsigned long long *dst_d, int *src_l, int *src_u, int *src_s, unsigned long long *src_d, int mode)
 
unsigned long long _XMP_gtol_calc_offset (_XMP_array_t *a, int g_idx[])
 
void _XMP_copy_scalar_array (char *scalar, _XMP_array_t *a, _XMP_comm_set_t *comm_set[])
 
void _XMP_gmove_INOUT_SCALAR (_XMP_array_t *dst_array, void *scalar,...)
 
void _XMP_gmove_scalar_garray (void *scalar, _XMP_gmv_desc_t *gmv_desc_rightp, int mode)
 
void _XMP_gmove_garray_scalar (_XMP_gmv_desc_t *gmv_desc_leftp, void *scalar, int mode)
 
void _XMP_gmove_garray_garray (_XMP_gmv_desc_t *gmv_desc_leftp, _XMP_gmv_desc_t *gmv_desc_rightp, int mode)
 
void _XMP_gmove_garray_larray (_XMP_gmv_desc_t *gmv_desc_leftp, _XMP_gmv_desc_t *gmv_desc_rightp, int mode)
 
void _XMP_gmove_larray_garray (_XMP_gmv_desc_t *gmv_desc_leftp, _XMP_gmv_desc_t *gmv_desc_rightp, int mode)
 

Variables

void(* _XMP_pack_comm_set )(void *sendbuf, int sendbuf_size, _XMP_array_t *a, _XMP_comm_set_t *comm_set[][_XMP_N_MAX_DIM])
 
void(* _XMP_unpack_comm_set )(void *recvbuf, int recvbuf_size, _XMP_array_t *a, _XMP_comm_set_t *comm_set[][_XMP_N_MAX_DIM])
 
_XMP_nodes_tgmv_nodes
 
int n_gmv_nodes
 
int(* _alloc_size )[_XMP_N_MAX_DIM]
 
int _dim_alloc_size
 

Macro Definition Documentation

◆ _XMP_SM_GTOL_BLOCK

#define _XMP_SM_GTOL_BLOCK (   _i,
  _m,
  _w 
)    (((_i) - (_m)) % (_w))

◆ _XMP_SM_GTOL_BLOCK_CYCLIC

#define _XMP_SM_GTOL_BLOCK_CYCLIC (   _b,
  _i,
  _m,
  _P 
)    (((((_i) - (_m)) / (((_P) * (_b)))) * (_b)) + (((_i) - (_m)) % (_b)))

◆ _XMP_SM_GTOL_CYCLIC

#define _XMP_SM_GTOL_CYCLIC (   _i,
  _m,
  _P 
)    (((_i) - (_m)) / (_P))

◆ DBG_RANK

#define DBG_RANK   3

◆ MPI_PORTABLE_PLATFORM_H

#define MPI_PORTABLE_PLATFORM_H

◆ XMP_DBG

#define XMP_DBG   0

◆ XMP_DBG_OWNER_REGION

#define XMP_DBG_OWNER_REGION   0

Function Documentation

◆ _XMP_align_local_idx()

void _XMP_align_local_idx ( long long int  global_idx,
int *  local_idx,
_XMP_array_t array,
int  array_axis,
int *  rank 
)
1360 {
1361  _XMP_template_t *template = array->align_template;
1362  int template_index = array->info[array_axis].align_template_index;
1363  _XMP_template_chunk_t *chunk = &(template->chunk[template_index]);
1364  _XMP_nodes_info_t *n_info = chunk->onto_nodes_info;
1365  long long base = array->info[array_axis].ser_lower;
1366  long long tbase = template->info[template_index].ser_lower;
1367  int offset = array->info[array_axis].align_subscript + (base - tbase);
1368  int irank, idiv, imod;
1369 
1370  switch(array->info[array_axis].align_manner){
1372  {
1373  *rank=0;
1374  *local_idx = global_idx + offset - base;
1375  }
1376  break;
1377  case _XMP_N_ALIGN_BLOCK:
1378  {
1379  *rank = (global_idx + offset - base) / chunk->par_chunk_width;
1380  *local_idx = (global_idx + offset - base ) - *rank * chunk->par_chunk_width + array->info[array_axis].shadow_size_lo;
1381  idiv = offset / (chunk->par_chunk_width);
1382  if (*rank == idiv){
1383  *local_idx = *local_idx - offset%(chunk->par_chunk_width);
1384  }
1385  }
1386  break;
1387  case _XMP_N_ALIGN_CYCLIC:
1388  {
1389  idiv = offset/n_info->size;
1390  imod = offset%n_info->size;
1391  *rank = (global_idx + offset - base) % n_info->size;
1392  *local_idx = (global_idx + offset - base) / n_info->size;
1393  if (imod > *rank){
1394  *local_idx = *local_idx - (idiv + 1);
1395  } else {
1396  *local_idx = *local_idx - idiv;
1397  }
1398  }
1399  break;
1401  {
1402  int w = chunk->par_width;
1403  idiv = (offset/w)/n_info->size;
1404  int imod1 = (offset/w)%n_info->size;
1405  int imod2 = offset%w;
1406  int off = global_idx + offset - base;
1407  *local_idx = (off / (n_info->size*w)) * w + off%w;
1408 
1409  *rank=(off/w)% (n_info->size);
1410  if (imod1 > 0){
1411  if (imod1 == *rank ){
1412  *local_idx = *local_idx - idiv*w-imod2;
1413  }else if (imod1 > *rank){
1414  *local_idx = *local_idx - (idiv+1)*w;
1415  }
1416  }else if (imod1 == 0){
1417  if (imod1 == *rank ){
1418  *local_idx = *local_idx - idiv*w -imod2;
1419  }else{
1420  *local_idx = *local_idx - idiv*w;
1421  }
1422  }
1423  }
1424  break;
1425  case _XMP_N_ALIGN_GBLOCK:
1426  {
1427  irank=-1;
1428  idiv=0;
1429  for(int i=1;i<(n_info->size+1);i++){
1430  if(global_idx + offset < chunk->mapping_array[i]+ (base - tbase)){
1431  *rank = i - 1;
1432  break;
1433  }
1434  }
1435  *local_idx = global_idx + offset - chunk->mapping_array[*rank]- (base - tbase) + array->info[array_axis].shadow_size_lo;
1436  for(int i=1;i<n_info->size+1;i++){
1437  if(offset < chunk->mapping_array[i]+(base-tbase)){
1438  irank = i - 1;
1439  idiv = offset - (chunk->mapping_array[i-1] + (base - tbase) - base);
1440  break;
1441  }
1442  }
1443  if (*rank == irank){
1444  *local_idx = *local_idx - idiv;
1445  }
1446  }
1447  break;
1449  {
1450  *rank=0;
1451  *local_idx=global_idx - base;
1452  }
1453  break;
1454  default:
1455  _XMP_fatal("_XMP_: unknown chunk dist_manner");
1456  }
1457 
1458 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_calc_global_index_BCAST()

int _XMP_calc_global_index_BCAST ( int  dst_dim,
int *  dst_l,
int *  dst_u,
int *  dst_s,
_XMP_array_t src_array,
int *  src_array_nodes_ref,
int *  src_l,
int *  src_u,
int *  src_s 
)
420  {
421  _XMP_template_t *template = src_array->align_template;
422 
423  int dst_dim_index = 0;
424  int array_dim = src_array->dim;
425  for (int i = 0; i < array_dim; i++) {
426  if (_XMP_M_COUNT_TRIPLETi(src_l[i], src_u[i], src_s[i]) == 1) {
427  continue;
428  }
429 
430  int template_index = src_array->info[i].align_template_index;
431  if (template_index == _XMP_N_NO_ALIGN_TEMPLATE) {
432  do {
433  if (_XMP_M_COUNT_TRIPLETi(dst_l[dst_dim_index], dst_u[dst_dim_index], dst_s[dst_dim_index]) != 1) {
434  dst_dim_index++;
435  break;
436  } else if (dst_dim_index < dst_dim) {
437  dst_dim_index++;
438  } else {
439  _XMP_fatal("wrong assign statement for gmove");
440  }
441  } while (1);
442  } else {
443  int onto_nodes_index = template->chunk[template_index].onto_nodes_index;
444  _XMP_ASSERT(onto_nodes_index != _XMP_N_NO_ONTO_NODES);
445 
446  int array_nodes_index = _XMP_calc_nodes_index_from_inherit_nodes_index(src_array->array_nodes, onto_nodes_index);
447  int rank = src_array_nodes_ref[array_nodes_index];
448 
449  // calc template info
450  int template_lower, template_upper, template_stride;
451  if (!_XMP_calc_template_par_triplet(template, template_index, rank, &template_lower, &template_upper, &template_stride)) {
452  return _XMP_N_INT_FALSE;
453  }
454 
455  do {
456  if (_XMP_M_COUNT_TRIPLETi(dst_l[dst_dim_index], dst_u[dst_dim_index], dst_s[dst_dim_index]) != 1) {
457  if (_XMP_sched_gmove_triplet(template_lower, template_upper, template_stride,
458  src_array, i,
459  &(src_l[i]), &(src_u[i]), &(src_s[i]),
460  &(dst_l[dst_dim_index]), &(dst_u[dst_dim_index]), &(dst_s[dst_dim_index]))) {
461  dst_dim_index++;
462  break;
463  } else {
464  return _XMP_N_INT_FALSE;
465  }
466  } else if (dst_dim_index < dst_dim) {
467  dst_dim_index++;
468  } else {
469  _XMP_fatal("wrong assign statement for gmove");
470  }
471  } while (1);
472  }
473  }
474 
475  return _XMP_N_INT_TRUE;
476 }
Here is the call graph for this function:

◆ _XMP_calc_global_index_HOMECOPY()

int _XMP_calc_global_index_HOMECOPY ( _XMP_array_t dst_array,
int  dst_dim_index,
int *  dst_l,
int *  dst_u,
int *  dst_s,
int *  src_l,
int *  src_u,
int *  src_s 
)
406  {
407  int align_template_index = dst_array->info[dst_dim_index].align_template_index;
408  if (align_template_index != _XMP_N_NO_ALIGN_TEMPLATE) {
409  _XMP_template_chunk_t *tc = &((dst_array->align_template)->chunk[align_template_index]);
410  return _XMP_sched_gmove_triplet(tc->par_lower, tc->par_upper, tc->par_stride,
411  dst_array, dst_dim_index,
412  dst_l, dst_u, dst_s,
413  src_l, src_u, src_s);
414  } else {
415  return _XMP_N_INT_TRUE;
416  }
417 }
Here is the caller graph for this function:

◆ _XMP_calc_gmove_array_owner_linear_rank_SCALAR()

int _XMP_calc_gmove_array_owner_linear_rank_SCALAR ( _XMP_array_t array,
int *  ref_index 
)
216  {
217  _XMP_nodes_t *array_nodes = array->array_nodes;
218  int array_nodes_dim = array_nodes->dim;
219  int rank_array[array_nodes_dim];
220 
221  _XMP_calc_gmove_rank_array_SCALAR(array, ref_index, rank_array);
222 
223  return _XMP_calc_linear_rank_on_target_nodes(array_nodes, rank_array, _XMP_get_execution_nodes());
224 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_calc_gmove_rank_array_SCALAR()

void _XMP_calc_gmove_rank_array_SCALAR ( _XMP_array_t array,
int *  ref_index,
int *  rank_array 
)
197  {
198  _XMP_template_t *template = array->align_template;
199 
200  int array_dim = array->dim;
201  for (int i = 0; i < array_dim; i++) {
202  _XMP_array_info_t *ai = &(array->info[i]);
203  int template_index = ai->align_template_index;
204  if (template_index != _XMP_N_NO_ALIGN_TEMPLATE) {
205  _XMP_template_chunk_t *chunk = &(template->chunk[ai->align_template_index]);
206  int onto_nodes_index = chunk->onto_nodes_index;
207  _XMP_ASSERT(array_nodes_index != _XMP_N_NO_ONTO_NODES);
208  if (onto_nodes_index == -1) continue;
209  int array_nodes_index = _XMP_calc_nodes_index_from_inherit_nodes_index(array->array_nodes, onto_nodes_index);
210  rank_array[array_nodes_index] = _XMP_calc_template_owner_SCALAR(template, template_index,
211  ref_index[i] + ai->align_subscript);
212  }
213  }
214 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_check_gmove_array_ref_inclusion_SCALAR()

int _XMP_check_gmove_array_ref_inclusion_SCALAR ( _XMP_array_t array,
int  array_index,
int  ref_index 
)
309  {
310  _XMP_ASSERT(!(array->align_template)->is_owner);
311 
312  _XMP_array_info_t *ai = &(array->info[array_index]);
314  return _XMP_N_INT_TRUE;
315  } else {
316  int template_ref_index = ref_index + ai->align_subscript;
317  return _XMP_check_template_ref_inclusion(template_ref_index, template_ref_index, 1,
319  }
320 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_copy_scalar_array()

void _XMP_copy_scalar_array ( char *  scalar,
_XMP_array_t a,
_XMP_comm_set_t comm_set[] 
)
4206  {
4207 
4208  int ndims = a->dim;
4209 
4210  char *dst = (char *)a->array_addr_p;
4211 
4212  _XMP_comm_set_t *c[ndims];
4213 
4214  int i[_XMP_N_MAX_DIM];
4215 
4216  switch (ndims){
4217 
4218  case 1:
4219  for (c[0] = comm_set[0]; c[0]; c[0] = c[0]->next){
4220  for (i[0] = c[0]->l; i[0] <= c[0]->u; i[0]++){
4221  memcpy(dst + _XMP_gtol_calc_offset(a, i), scalar, a->type_size);
4222  }}
4223  break;
4224 
4225  case 2:
4226  for (c[1] = comm_set[1]; c[1]; c[1] = c[1]->next){
4227  for (i[1] = c[1]->l; i[1] <= c[1]->u; i[1]++){
4228  for (c[0] = comm_set[0]; c[0]; c[0] = c[0]->next){
4229  for (i[0] = c[0]->l; i[0] <= c[0]->u; i[0]++){
4230  memcpy(dst + _XMP_gtol_calc_offset(a, i), scalar, a->type_size);
4231  }}
4232  }}
4233  break;
4234 
4235  case 3:
4236  for (c[2] = comm_set[2]; c[2]; c[2] = c[2]->next){
4237  for (i[2] = c[2]->l; i[2] <= c[2]->u; i[2]++){
4238  for (c[1] = comm_set[1]; c[1]; c[1] = c[1]->next){
4239  for (i[1] = c[1]->l; i[1] <= c[1]->u; i[1]++){
4240  for (c[0] = comm_set[0]; c[0]; c[0] = c[0]->next){
4241  for (i[0] = c[0]->l; i[0] <= c[0]->u; i[0]++){
4242  memcpy(dst + _XMP_gtol_calc_offset(a, i), scalar, a->type_size);
4243  }}
4244  }}
4245  }}
4246  break;
4247 
4248  case 4:
4249  for (c[3] = comm_set[3]; c[3]; c[3] = c[3]->next){
4250  for (i[3] = c[3]->l; i[3] <= c[3]->u; i[3]++){
4251  for (c[2] = comm_set[2]; c[2]; c[2] = c[2]->next){
4252  for (i[2] = c[2]->l; i[2] <= c[2]->u; i[2]++){
4253  for (c[1] = comm_set[1]; c[1]; c[1] = c[1]->next){
4254  for (i[1] = c[1]->l; i[1] <= c[1]->u; i[1]++){
4255  for (c[0] = comm_set[0]; c[0]; c[0] = c[0]->next){
4256  for (i[0] = c[0]->l; i[0] <= c[0]->u; i[0]++){
4257  memcpy(dst + _XMP_gtol_calc_offset(a, i), scalar, a->type_size);
4258  }}
4259  }}
4260  }}
4261  }}
4262  break;
4263 
4264  case 5:
4265  for (c[4] = comm_set[4]; c[4]; c[4] = c[4]->next){
4266  for (i[4] = c[4]->l; i[4] <= c[4]->u; i[4]++){
4267  for (c[3] = comm_set[3]; c[3]; c[3] = c[3]->next){
4268  for (i[3] = c[3]->l; i[3] <= c[3]->u; i[3]++){
4269  for (c[2] = comm_set[2]; c[2]; c[2] = c[2]->next){
4270  for (i[2] = c[2]->l; i[2] <= c[2]->u; i[2]++){
4271  for (c[1] = comm_set[1]; c[1]; c[1] = c[1]->next){
4272  for (i[1] = c[1]->l; i[1] <= c[1]->u; i[1]++){
4273  for (c[0] = comm_set[0]; c[0]; c[0] = c[0]->next){
4274  for (i[0] = c[0]->l; i[0] <= c[0]->u; i[0]++){
4275  memcpy(dst + _XMP_gtol_calc_offset(a, i), scalar, a->type_size);
4276  }}
4277  }}
4278  }}
4279  }}
4280  }}
4281  break;
4282 
4283  case 6:
4284  for (c[5] = comm_set[5]; c[5]; c[5] = c[5]->next){
4285  for (i[5] = c[5]->l; i[5] <= c[5]->u; i[5]++){
4286  for (c[4] = comm_set[4]; c[4]; c[4] = c[4]->next){
4287  for (i[4] = c[4]->l; i[4] <= c[4]->u; i[4]++){
4288  for (c[3] = comm_set[3]; c[3]; c[3] = c[3]->next){
4289  for (i[3] = c[3]->l; i[3] <= c[3]->u; i[3]++){
4290  for (c[2] = comm_set[2]; c[2]; c[2] = c[2]->next){
4291  for (i[2] = c[2]->l; i[2] <= c[2]->u; i[2]++){
4292  for (c[1] = comm_set[1]; c[1]; c[1] = c[1]->next){
4293  for (i[1] = c[1]->l; i[1] <= c[1]->u; i[1]++){
4294  for (c[0] = comm_set[0]; c[0]; c[0] = c[0]->next){
4295  for (i[0] = c[0]->l; i[0] <= c[0]->u; i[0]++){
4296  memcpy(dst + _XMP_gtol_calc_offset(a, i), scalar, a->type_size);
4297  }}
4298  }}
4299  }}
4300  }}
4301  }}
4302  }}
4303  break;
4304 
4305  case 7:
4306  for (c[6] = comm_set[6]; c[6]; c[6] = c[6]->next){
4307  for (i[6] = c[6]->l; i[6] <= c[6]->u; i[6]++){
4308  for (c[5] = comm_set[5]; c[5]; c[5] = c[5]->next){
4309  for (i[5] = c[5]->l; i[5] <= c[5]->u; i[5]++){
4310  for (c[4] = comm_set[4]; c[4]; c[4] = c[4]->next){
4311  for (i[4] = c[4]->l; i[4] <= c[4]->u; i[4]++){
4312  for (c[3] = comm_set[3]; c[3]; c[3] = c[3]->next){
4313  for (i[3] = c[3]->l; i[3] <= c[3]->u; i[3]++){
4314  for (c[2] = comm_set[2]; c[2]; c[2] = c[2]->next){
4315  for (i[2] = c[2]->l; i[2] <= c[2]->u; i[2]++){
4316  for (c[1] = comm_set[1]; c[1]; c[1] = c[1]->next){
4317  for (i[1] = c[1]->l; i[1] <= c[1]->u; i[1]++){
4318  for (c[0] = comm_set[0]; c[0]; c[0] = c[0]->next){
4319  for (i[0] = c[0]->l; i[0] <= c[0]->u; i[0]++){
4320  memcpy(dst + _XMP_gtol_calc_offset(a, i), scalar, a->type_size);
4321  }}
4322  }}
4323  }}
4324  }}
4325  }}
4326  }}
4327  }}
4328  break;
4329 
4330  default:
4331  _XMP_fatal("wrong array dimension");
4332  }
4333 
4334 }
Here is the call graph for this function:

◆ _XMP_get_array_addr()

void* _XMP_get_array_addr ( _XMP_array_t a,
int *  gidx 
)
48 {
49  int ndims = a->dim;
50  void *ret = a->array_addr_p;
51 
53 
54  for (int i = 0; i < ndims; i++){
55 
56  _XMP_array_info_t *ai = &(a->info[i]);
59  int lidx = 0;
60  int offset;
61  int glb, t_lb, np, w;
62  int l_shadow;
63  size_t type_size = a->type_size;
64 
65  switch (ai->align_manner){
66 
68  lidx = gidx[i] - ai->ser_lower;
69  break;
70 
72  lidx = gidx[i] - ai->ser_lower;
73  break;
74 
75  case _XMP_N_ALIGN_BLOCK:
76  // par_lower is the index of the lower bound of the local section.
77  glb = ai->par_lower;
78  l_shadow = ai->shadow_size_lo;
79  lidx = gidx[i] - glb + l_shadow;
80  break;
81 
83  // assumed that even a cyclic array is distributed equally
84  offset = ai->align_subscript;
85  t_lb = ti->ser_lower;
86  np = ai->par_stride;
87  lidx = (gidx[i] + offset - t_lb) / np;
88  break;
89 
91  // assumed that even a cyclic array is distributed equally
92  offset = ai->align_subscript;
93  t_lb = ti->ser_lower;
94  np = ai->par_stride;
95  w = tc->par_stride;
96  lidx = w * ((gidx[i] + offset - t_lb) / (np * w))
97  + ((gidx[i] + offset - t_lb) % w);
98  break;
99 
100  default:
101  _XMP_fatal("_XMP_get_array_addr: unknown align_manner");
102  }
103 
104  //xmpf_dbg_printf("a->array_addr_p = %p\n", a->array_addr_p);
105  //xmpf_dbg_printf("ret = %p, lidx = %d, ai->dim_acc = %d\n", ret, lidx, ai->dim_acc);
106  ret = (char *)ret + lidx * ai->dim_acc * type_size;
107  }
108 
109  return ret;
110 
111 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_gmove_array_array_common()

void _XMP_gmove_array_array_common ( _XMP_gmv_desc_t gmv_desc_leftp,
_XMP_gmv_desc_t gmv_desc_rightp,
int *  dst_l,
int *  dst_u,
int *  dst_s,
unsigned long long *  dst_d,
int *  src_l,
int *  src_u,
int *  src_s,
unsigned long long *  src_d,
int  mode 
)
2039  {
2040  // NOTE: asynchronous gmove aloways done by _XMP_gmove_1to1
2041  if (!xmp_is_async() && gmv_desc_leftp->is_global && gmv_desc_rightp->is_global && mode == _XMP_N_GMOVE_NORMAL){
2042  if (_XMP_gmove_garray_garray_opt(gmv_desc_leftp, gmv_desc_rightp,
2043  dst_l, dst_u, dst_s, dst_d,
2044  src_l, src_u, src_s, src_d)) return;
2045  // fall through
2046  }
2047 
2048  _XMP_gmove_1to1(gmv_desc_leftp, gmv_desc_rightp, mode);
2049 
2050  return;
2051 
2052 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_gmove_bcast_ARRAY()

unsigned long long _XMP_gmove_bcast_ARRAY ( void *  dst_addr,
int  dst_dim,
int *  dst_l,
int *  dst_u,
int *  dst_s,
unsigned long long *  dst_d,
void *  src_addr,
int  src_dim,
int *  src_l,
int *  src_u,
int *  src_s,
unsigned long long *  src_d,
int  type,
size_t  type_size,
int  root_rank 
)
277  {
278  unsigned long long dst_buffer_elmts = 1;
279  for (int i = 0; i < dst_dim; i++) {
280  dst_buffer_elmts *= _XMP_M_COUNT_TRIPLETi(dst_l[i], dst_u[i], dst_s[i]);
281  }
282 
283  void *buffer = _XMP_alloc(dst_buffer_elmts * type_size);
284 
285  _XMP_nodes_t *exec_nodes = _XMP_get_execution_nodes();
286  _XMP_ASSERT(exec_nodes->is_member);
287 
288  if (root_rank == (exec_nodes->comm_rank)) {
289  unsigned long long src_buffer_elmts = 1;
290  for (int i = 0; i < src_dim; i++) {
291  src_buffer_elmts *= _XMP_M_COUNT_TRIPLETi(src_l[i], src_u[i], src_s[i]);
292  }
293 
294  if (dst_buffer_elmts != src_buffer_elmts) {
295  _XMP_fatal("wrong assign statement for gmove");
296  } else {
297  (*_xmp_pack_array)(buffer, src_addr, type, type_size, src_dim, src_l, src_u, src_s, src_d);
298  }
299  }
300 
301  _XMP_gmove_bcast(buffer, type_size, dst_buffer_elmts, root_rank);
302 
303  (*_xmp_unpack_array)(dst_addr, buffer, type, type_size, dst_dim, dst_l, dst_u, dst_s, dst_d);
304  _XMP_free(buffer);
305 
306  return dst_buffer_elmts;
307 }
Here is the call graph for this function:

◆ _XMP_gmove_BCAST_GSCALAR()

void _XMP_gmove_BCAST_GSCALAR ( void *  dst_addr,
_XMP_array_t array,
int  ref_index[] 
)
687  {
688 
689  _XMP_nodes_t *exec_nodes = _XMP_get_execution_nodes();
690  _XMP_ASSERT(exec_nodes->is_member);
691 
692  void *src_addr = NULL;
693  int type_size = array->type_size;
694 
695  if(_XMP_IS_SINGLE) {
696  src_addr = _XMP_get_array_addr(array, ref_index);
697  memcpy(dst_addr, src_addr, type_size);
698  return;
699  }
700 
701  int root_rank = _XMP_calc_gmove_array_owner_linear_rank_SCALAR(array, ref_index);
702 
703  if (root_rank == exec_nodes->comm_rank){
704  // I am the root.
705  src_addr = _XMP_get_array_addr(array, ref_index);
706  memcpy(dst_addr, src_addr, type_size);
707  }
708 
709  _XMP_gmove_bcast(dst_addr, type_size, 1, root_rank);
710 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_gmove_calc_unit_size()

void _XMP_gmove_calc_unit_size ( _XMP_array_t dst_array,
_XMP_array_t src_array,
unsigned long long *  alltoall_unit_size,
unsigned long long *  dst_pack_unit_size,
unsigned long long *  src_pack_unit_size,
unsigned long long *  dst_ser_size,
unsigned long long *  src_ser_size,
int  dst_block_dim,
int  src_block_dim 
)
1176  {
1177 
1178  int dst_dim=dst_array->dim, src_dim=src_array->dim;
1179  int dst_chunk_size[_XMP_N_MAX_DIM], src_chunk_size[_XMP_N_MAX_DIM];
1180 
1181  *alltoall_unit_size=1;
1182  *dst_pack_unit_size=1;
1183  *src_pack_unit_size=1;
1184  *dst_ser_size=1;
1185  *src_ser_size=1;
1186 
1187  if ((_XMPF_running == 1) && (_XMPC_running == 0)){
1188  for(int i=0; i<dst_dim; i++){
1189  if(i==dst_block_dim){
1190  dst_chunk_size[dst_block_dim]=dst_array->info[dst_block_dim].par_size;
1191  }else if(i==src_block_dim){
1192  dst_chunk_size[src_block_dim]=src_array->info[src_block_dim].par_size;
1193  }else{
1194  dst_chunk_size[i]=dst_array->info[i].ser_size;
1195  }
1196  }
1197  for(int i=0; i<src_dim; i++){
1198  if(i==dst_block_dim){
1199  src_chunk_size[i]=dst_array->info[i].par_size;
1200  }else if(i==src_block_dim){
1201  src_chunk_size[i]=src_array->info[i].par_size;
1202  }else{
1203  src_chunk_size[i]=src_array->info[i].ser_size;
1204  }
1205  *alltoall_unit_size *= src_chunk_size[i];
1206  }
1207  for(int i=0; i<dst_block_dim+1; i++){
1208  *src_pack_unit_size *= src_chunk_size[i];
1209  *src_ser_size *= src_array->info[i].par_size;
1210  }
1211  for(int i=0; i<src_block_dim+1; i++){
1212  *dst_pack_unit_size *= dst_chunk_size[i];
1213  *dst_ser_size *= dst_array->info[i].par_size;
1214  }
1215  }
1216 
1217  if ((_XMPF_running == 0) && (_XMPC_running == 1)){
1218  for(int i=dst_dim-1; i> -1; i--){
1219  if(i==dst_block_dim){
1220  dst_chunk_size[i]=dst_array->info[i].par_size;
1221  }else if(i==src_block_dim){
1222  dst_chunk_size[i]=src_array->info[i].par_size;
1223  }else{
1224  dst_chunk_size[i]=dst_array->info[i].ser_size;
1225  }
1226  }
1227  for(int i=src_dim-1; i>-1; i--){
1228  if(i==dst_block_dim){
1229  src_chunk_size[i]=dst_array->info[i].par_size;
1230  }else if(i==src_block_dim){
1231  src_chunk_size[i]=src_array->info[i].par_size;
1232  }else{
1233  src_chunk_size[i]=src_array->info[i].ser_size;
1234  }
1235  *alltoall_unit_size *= src_chunk_size[i];
1236  }
1237  for(int i=src_array->dim-1; i>dst_block_dim-1; i--){
1238  *src_pack_unit_size *= src_chunk_size[i];
1239  *src_ser_size *= src_array->info[i].par_size;
1240  }
1241  for(int i=src_array->dim-1; i>src_block_dim-1; i--){
1242  *dst_pack_unit_size *= dst_chunk_size[i];
1243  *dst_ser_size *= dst_array->info[i].par_size;
1244  }
1245  }
1246 
1247 }

◆ _XMP_gmove_garray_garray()

void _XMP_gmove_garray_garray ( _XMP_gmv_desc_t gmv_desc_leftp,
_XMP_gmv_desc_t gmv_desc_rightp,
int  mode 
)
4772 {
4773  _XMP_array_t *dst_array = gmv_desc_leftp->a_desc;
4774  _XMP_array_t *src_array = gmv_desc_rightp->a_desc;
4775 
4776  _XMP_ASSERT(src_array->type == type);
4777  _XMP_ASSERT(src_array->type_size == dst_array->type_size);
4778 
4779  //unsigned long long gmove_total_elmts = 0;
4780 
4781  // get dst info
4782  unsigned long long dst_total_elmts = 1;
4783  int dst_dim = dst_array->dim;
4784  int dst_l[dst_dim], dst_u[dst_dim], dst_s[dst_dim];
4785  unsigned long long dst_d[dst_dim];
4786  int dst_scalar_flag = 1;
4787  for (int i = 0; i < dst_dim; i++) {
4788  dst_l[i] = gmv_desc_leftp->lb[i];
4789  dst_u[i] = gmv_desc_leftp->ub[i];
4790  dst_s[i] = gmv_desc_leftp->st[i];
4791  dst_d[i] = dst_array->info[i].dim_acc;
4792  _XMP_normalize_array_section(gmv_desc_leftp, i, &(dst_l[i]), &(dst_u[i]), &(dst_s[i]));
4793  if (dst_s[i] != 0) dst_total_elmts *= _XMP_M_COUNT_TRIPLETi(dst_l[i], dst_u[i], dst_s[i]);
4794  dst_scalar_flag &= (dst_s[i] == 0);
4795  }
4796 
4797  // get src info
4798  unsigned long long src_total_elmts = 1;
4799  int src_dim = src_array->dim;
4800  int src_l[src_dim], src_u[src_dim], src_s[src_dim];
4801  unsigned long long src_d[src_dim];
4802  int src_scalar_flag = 1;
4803  for (int i = 0; i < src_dim; i++) {
4804  src_l[i] = gmv_desc_rightp->lb[i];
4805  src_u[i] = gmv_desc_rightp->ub[i];
4806  src_s[i] = gmv_desc_rightp->st[i];
4807  src_d[i] = src_array->info[i].dim_acc;
4808  _XMP_normalize_array_section(gmv_desc_rightp, i, &(src_l[i]), &(src_u[i]), &(src_s[i]));
4809  if (src_s[i] != 0) src_total_elmts *= _XMP_M_COUNT_TRIPLETi(src_l[i], src_u[i], src_s[i]);
4810  src_scalar_flag &= (src_s[i] == 0);
4811  }
4812 
4813  if (dst_total_elmts != src_total_elmts && !src_scalar_flag){
4814  _XMP_fatal("wrong assign statement for gmove");
4815  }
4816  else {
4817  //gmove_total_elmts = dst_total_elmts;
4818  }
4819 
4820  if (mode == _XMP_N_GMOVE_NORMAL){
4821 
4822  if (dst_scalar_flag && src_scalar_flag){
4823  void *dst_addr = (char *)dst_array->array_addr_p + _XMP_gtol_calc_offset(dst_array, dst_l);
4824  void *src_addr = (char *)src_array->array_addr_p + _XMP_gtol_calc_offset(src_array, src_l);
4825  _XMP_gmove_SENDRECV_GSCALAR(dst_addr, src_addr,
4826  dst_array, src_array,
4827  dst_l, src_l);
4828  return;
4829  }
4830  else if (!dst_scalar_flag && src_scalar_flag){
4831  char *tmp = _XMP_alloc(src_array->type_size);
4832  //char *src_addr = (char *)src_array->array_addr_p + _XMP_gtol_calc_offset(src_array, src_l);
4833  //_XMP_gmove_BCAST_GSCALAR(tmp, src_addr, src_array, src_l);
4834  _XMP_gmove_BCAST_GSCALAR(tmp, src_array, src_l);
4835  _XMP_gmove_gsection_scalar(dst_array, dst_l, dst_u, dst_s, tmp);
4836  _XMP_free(tmp);
4837  return;
4838  }
4839  }
4840 
4841  _XMP_gmove_array_array_common(gmv_desc_leftp, gmv_desc_rightp,
4842  dst_l, dst_u, dst_s, dst_d,
4843  src_l, src_u, src_s, src_d,
4844  mode);
4845 }
Here is the call graph for this function:

◆ _XMP_gmove_garray_larray()

void _XMP_gmove_garray_larray ( _XMP_gmv_desc_t gmv_desc_leftp,
_XMP_gmv_desc_t gmv_desc_rightp,
int  mode 
)
4854 {
4855  _XMP_array_t *dst_array = gmv_desc_leftp->a_desc;
4856 
4857  int type = dst_array->type;
4858  size_t type_size = dst_array->type_size;
4859 
4860  if (!dst_array->is_allocated && mode != _XMP_N_GMOVE_OUT) return;
4861 
4862  // get dst info
4863  unsigned long long dst_total_elmts = 1;
4864  void *dst_addr = dst_array->array_addr_p;
4865  int dst_dim = dst_array->dim;
4866  int dst_l[dst_dim], dst_u[dst_dim], dst_s[dst_dim];
4867  unsigned long long dst_d[dst_dim];
4868  int dst_scalar_flag = 1;
4869  for (int i = 0; i < dst_dim; i++) {
4870  dst_l[i] = gmv_desc_leftp->lb[i];
4871  dst_u[i] = gmv_desc_leftp->ub[i];
4872  dst_s[i] = gmv_desc_leftp->st[i];
4873  dst_d[i] = dst_array->info[i].dim_acc;
4874  _XMP_normalize_array_section(gmv_desc_leftp, i, &(dst_l[i]), &(dst_u[i]), &(dst_s[i]));
4875  if (dst_s[i] != 0) dst_total_elmts *= _XMP_M_COUNT_TRIPLETi(dst_l[i], dst_u[i], dst_s[i]);
4876  dst_scalar_flag &= (dst_s[i] == 0);
4877  }
4878 
4879  // get src info
4880  unsigned long long src_total_elmts = 1;
4881  void *src_addr = gmv_desc_rightp->local_data;
4882  int src_dim = gmv_desc_rightp->ndims;
4883  int src_l[src_dim], src_u[src_dim], src_s[src_dim];
4884  unsigned long long src_d[src_dim];
4885  int src_scalar_flag = 1;
4886  for (int i = 0; i < src_dim; i++) {
4887  src_l[i] = gmv_desc_rightp->lb[i];
4888  src_u[i] = gmv_desc_rightp->ub[i];
4889  src_s[i] = gmv_desc_rightp->st[i];
4890  if (i == 0) src_d[i] = 1;
4891  else src_d[i] = src_d[i-1] * (gmv_desc_rightp->a_ub[i] - gmv_desc_rightp->a_lb[i] + 1);
4892  _XMP_normalize_array_section(gmv_desc_rightp, i, &(src_l[i]), &(src_u[i]), &(src_s[i]));
4893  if (src_s[i] != 0) src_total_elmts *= _XMP_M_COUNT_TRIPLETi(src_l[i], src_u[i], src_s[i]);
4894  src_scalar_flag &= (src_s[i] == 0);
4895  }
4896 
4897  if (dst_total_elmts != src_total_elmts && !src_scalar_flag){
4898  _XMP_fatal("wrong assign statement for gmove");
4899  }
4900 
4901  char *scalar = (char *)src_addr;
4902  if (src_scalar_flag){
4903  for (int i = 0; i < src_dim; i++){
4904  scalar += ((src_l[i] - gmv_desc_rightp->a_lb[i]) * src_d[i] * type_size);
4905  }
4906  }
4907 
4908  if (mode == _XMP_N_GMOVE_OUT){
4909  if (src_scalar_flag){
4910 #ifdef _XMP_MPI3_ONESIDED
4911  _XMP_gmove_inout_scalar(scalar, gmv_desc_leftp, _XMP_N_COARRAY_PUT);
4912 #else
4913  _XMP_fatal("Not supported gmove in/out on non-MPI3 environments");
4914 #endif
4915  }
4916  else {
4917  _XMP_gmove_array_array_common(gmv_desc_leftp, gmv_desc_rightp,
4918  dst_l, dst_u, dst_s, dst_d,
4919  src_l, src_u, src_s, src_d,
4920  mode);
4921  }
4922  return;
4923  }
4924 
4925  if (dst_scalar_flag && src_scalar_flag){
4926  _XMP_gmove_garray_scalar(gmv_desc_leftp, scalar, mode);
4927  return;
4928  }
4929 
4930  if (_XMP_IS_SINGLE) {
4931  for (int i = 0; i < dst_dim; i++) {
4932  _XMP_gtol_array_ref_triplet(dst_array, i, &(dst_l[i]), &(dst_u[i]), &(dst_s[i]));
4933  }
4934 
4935  _XMP_gmove_localcopy_ARRAY(type, type_size,
4936  dst_addr, dst_dim, dst_l, dst_u, dst_s, dst_d,
4937  src_addr, src_dim, src_l, src_u, src_s, src_d);
4938  return;
4939  }
4940 
4941  for (int i = 0; i < src_dim; i++){
4942  src_l[i] -= gmv_desc_rightp->a_lb[i];
4943  src_u[i] -= gmv_desc_rightp->a_lb[i];
4944  }
4945 
4946  // calc index ref
4947  int src_dim_index = 0;
4948  unsigned long long dst_buffer_elmts = 1;
4949  unsigned long long src_buffer_elmts = 1;
4950  for (int i = 0; i < dst_dim; i++) {
4951  int dst_elmts = _XMP_M_COUNT_TRIPLETi(dst_l[i], dst_u[i], dst_s[i]);
4952  if (dst_elmts == 1) {
4953  if(!_XMP_check_gmove_array_ref_inclusion_SCALAR(dst_array, i, dst_l[i])) {
4954  return;
4955  }
4956  } else {
4957  dst_buffer_elmts *= dst_elmts;
4958 
4959  int src_elmts;
4960  do {
4961  src_elmts = _XMP_M_COUNT_TRIPLETi(src_l[src_dim_index], src_u[src_dim_index], src_s[src_dim_index]);
4962  if (src_elmts != 1) {
4963  break;
4964  } else if (src_dim_index < src_dim) {
4965  src_dim_index++;
4966  } else {
4967  _XMP_fatal("wrong assign statement for gmove");
4968  }
4969  } while (1);
4970 
4971  if (_XMP_calc_global_index_HOMECOPY(dst_array, i,
4972  &(dst_l[i]), &(dst_u[i]), &(dst_s[i]),
4973  &(src_l[src_dim_index]), &(src_u[src_dim_index]), &(src_s[src_dim_index]))) {
4974  src_buffer_elmts *= src_elmts;
4975  src_dim_index++;
4976  } else {
4977  return;
4978  }
4979  }
4980 
4981  _XMP_gtol_array_ref_triplet(dst_array, i, &(dst_l[i]), &(dst_u[i]), &(dst_s[i]));
4982  }
4983 
4984  for (int i = src_dim_index; i < src_dim; i++) {
4985  src_buffer_elmts *= _XMP_M_COUNT_TRIPLETi(src_l[i], src_u[i], src_s[i]);
4986  }
4987 
4988  // alloc buffer
4989  if (dst_buffer_elmts != src_buffer_elmts) {
4990  _XMP_fatal("wrong assign statement for gmove");
4991  }
4992 
4993  void *buffer = _XMP_alloc(dst_buffer_elmts * type_size);
4994  (*_xmp_pack_array)(buffer, src_addr, type, type_size, src_dim, src_l, src_u, src_s, src_d);
4995  (*_xmp_unpack_array)(dst_addr, buffer, type, type_size, dst_dim, dst_l, dst_u, dst_s, dst_d);
4996  _XMP_free(buffer);
4997 }
Here is the call graph for this function:

◆ _XMP_gmove_garray_scalar()

void _XMP_gmove_garray_scalar ( _XMP_gmv_desc_t gmv_desc_leftp,
void *  scalar,
int  mode 
)
4733 {
4734  if (mode == _XMP_N_GMOVE_NORMAL){
4735  _XMP_array_t *array = gmv_desc_leftp->a_desc;
4736  int type_size = array->type_size;
4737  void *dst_addr = NULL;
4738  _XMP_nodes_t *exec_nodes = _XMP_get_execution_nodes();
4739 
4740  int ndims = gmv_desc_leftp->ndims;
4741  int lidx[ndims];
4742 
4743  for (int i=0;i<ndims;i++)
4744  lidx[i] = gmv_desc_leftp->lb[i];
4745 
4746  int owner_rank = _XMP_calc_gmove_array_owner_linear_rank_SCALAR(array, lidx);
4747 
4748  if (owner_rank == exec_nodes->comm_rank){ // I am the owner.
4749  dst_addr = _XMP_get_array_addr(array, lidx);
4750  memcpy(dst_addr, scalar, type_size);
4751  }
4752  }
4753  else if (mode == _XMP_N_GMOVE_OUT){
4754 #ifdef _XMP_MPI3_ONESIDED
4755  _XMP_gmove_inout_scalar(scalar, gmv_desc_leftp, _XMP_N_COARRAY_PUT);
4756 #else
4757  _XMP_fatal("Not supported gmove in/out on non-MPI3 environments");
4758 #endif
4759  }
4760  else {
4761  _XMP_fatal("_XMPF_gmove_garray_scalar: wrong gmove mode");
4762  }
4763 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_gmove_INOUT_SCALAR()

void _XMP_gmove_INOUT_SCALAR ( _XMP_array_t dst_array,
void *  scalar,
  ... 
)
4642  {
4643 
4644  _XMP_gmv_desc_t gmv_desc_leftp;
4645 
4646  va_list args;
4647  va_start(args, scalar);
4648 
4649  // get dst info
4650  //unsigned long long dst_total_elmts = 1;
4651  int dst_dim = dst_array->dim;
4652  int dst_l[dst_dim], dst_u[dst_dim], dst_s[dst_dim];
4653  //unsigned long long dst_d[dst_dim];
4654  for (int i = 0; i < dst_dim; i++) {
4655  dst_l[i] = va_arg(args, int);
4656  int size = va_arg(args, int);
4657  dst_s[i] = va_arg(args, int);
4658  dst_u[i] = dst_l[i] + (size - 1) * dst_s[i];
4659  //dst_d[i] = va_arg(args, unsigned long long);
4660  va_arg(args, unsigned long long); // skip the argument
4661  _XMP_normalize_array_section(&gmv_desc_leftp, i, &(dst_l[i]), &(dst_u[i]), &(dst_s[i]));
4662  //if (dst_s[i] != 0) dst_total_elmts *= _XMP_M_COUNT_TRIPLETi(dst_l[i], dst_u[i], dst_s[i]);
4663  }
4664 
4665  //int mode = va_arg(args, int);
4666 
4667  va_end(args);
4668 
4670 
4671 #ifdef _XMP_MPI3_ONESIDED
4672  _XMP_gmv_desc_t gmv_desc;
4673  int kind[_XMP_N_MAX_DIM];
4674  for (int i = 0; i < dst_dim; i++){
4675  kind[i] = (dst_s[i] == 0) ? XMP_N_GMOVE_INDEX : XMP_N_GMOVE_RANGE;
4676  }
4677 
4678  gmv_desc.is_global = true;
4679  gmv_desc.ndims = dst_dim;
4680 
4681  gmv_desc.a_desc = dst_array;
4682 
4683  gmv_desc.local_data = NULL;
4684  gmv_desc.a_lb = NULL;
4685  gmv_desc.a_ub = NULL;
4686 
4687  gmv_desc.kind = kind;
4688  gmv_desc.lb = dst_l;
4689  gmv_desc.ub = dst_u;
4690  gmv_desc.st = dst_s;
4691 
4692  _XMP_gmove_inout_scalar(scalar, &gmv_desc, _XMP_N_COARRAY_PUT);
4693 #else
4694  _XMP_fatal("Not supported gmove in/out on non-MPI3 environments");
4695 #endif
4696 
4697 }
Here is the call graph for this function:

◆ _XMP_gmove_larray_garray()

void _XMP_gmove_larray_garray ( _XMP_gmv_desc_t gmv_desc_leftp,
_XMP_gmv_desc_t gmv_desc_rightp,
int  mode 
)
5006 {
5007  _XMP_array_t *src_array = gmv_desc_rightp->a_desc;
5008  size_t type_size = src_array->type_size;
5009 
5010  // get dst info
5011  unsigned long long dst_total_elmts = 1;
5012  int dst_dim = gmv_desc_leftp->ndims;
5013  int dst_l[dst_dim], dst_u[dst_dim], dst_s[dst_dim];
5014  unsigned long long dst_d[dst_dim];
5015  int dst_scalar_flag = 1;
5016  for(int i=0;i<dst_dim;i++){
5017  dst_l[i] = gmv_desc_leftp->lb[i];
5018  dst_u[i] = gmv_desc_leftp->ub[i];
5019  dst_s[i] = gmv_desc_leftp->st[i];
5020  dst_d[i] = (i == 0)? 1 : dst_d[i-1]*(gmv_desc_leftp->a_ub[i] - gmv_desc_leftp->a_lb[i]+1);
5021  _XMP_normalize_array_section(gmv_desc_leftp, i, &(dst_l[i]), &(dst_u[i]), &(dst_s[i]));
5022  if (dst_s[i] != 0)
5023  dst_total_elmts *= _XMP_M_COUNT_TRIPLETi(dst_l[i], dst_u[i], dst_s[i]);
5024  dst_scalar_flag &= (dst_s[i] == 0);
5025  }
5026 
5027  // get src info
5028  unsigned long long src_total_elmts = 1;
5029  int src_dim = src_array->dim;
5030  int src_l[src_dim], src_u[src_dim], src_s[src_dim];
5031  unsigned long long src_d[src_dim];
5032  int src_scalar_flag = 1;
5033  for(int i=0;i<src_dim;i++) {
5034  src_l[i] = gmv_desc_rightp->lb[i];
5035  src_u[i] = gmv_desc_rightp->ub[i];
5036  src_s[i] = gmv_desc_rightp->st[i];
5037  src_d[i] = src_array->info[i].dim_acc;
5038  _XMP_normalize_array_section(gmv_desc_rightp, i, &(src_l[i]), &(src_u[i]), &(src_s[i]));
5039  if (src_s[i] != 0) src_total_elmts *= _XMP_M_COUNT_TRIPLETi(src_l[i], src_u[i], src_s[i]);
5040  src_scalar_flag &= (src_s[i] == 0);
5041  }
5042 
5043  if(dst_total_elmts != src_total_elmts && !src_scalar_flag){
5044  _XMP_fatal("wrong assign statement for gmove");
5045  }
5046 
5047  if (mode == _XMP_N_GMOVE_NORMAL){
5048  if (dst_scalar_flag && src_scalar_flag){
5049  char *dst_addr = (char *)gmv_desc_leftp->local_data;
5050  for (int i = 0; i < dst_dim; i++)
5051  dst_addr += ((dst_l[i] - gmv_desc_leftp->a_lb[i]) * dst_d[i]) * type_size;
5052  _XMP_gmove_BCAST_GSCALAR(dst_addr, src_array, src_l);
5053  return;
5054  }
5055  else if (!dst_scalar_flag && src_scalar_flag){
5056  char *tmp = _XMP_alloc(src_array->type_size);
5057  _XMP_gmove_BCAST_GSCALAR(tmp, src_array, src_l);
5058  char *dst_addr = (char *)gmv_desc_leftp->local_data;
5059 
5060  // to 0-based
5061  for (int i = 0; i < dst_dim; i++) {
5062  dst_l[i] -= gmv_desc_leftp->a_lb[i];
5063  dst_u[i] -= gmv_desc_leftp->a_lb[i];
5064  }
5065 
5066  _XMP_gmove_lsection_scalar(dst_addr, dst_dim, dst_l, dst_u, dst_s, dst_d, tmp, type_size);
5067  _XMP_free(tmp);
5068  return;
5069  }
5070  else {
5071  //gmove_total_elmts = dst_total_elmts;
5072  }
5073  }
5074 
5075  _XMP_gmove_array_array_common(gmv_desc_leftp, gmv_desc_rightp,
5076  dst_l, dst_u, dst_s, dst_d,
5077  src_l, src_u, src_s, src_d,
5078  mode);
5079 }
Here is the call graph for this function:

◆ _XMP_gmove_localcopy_ARRAY()

void _XMP_gmove_localcopy_ARRAY ( int  type,
int  type_size,
void *  dst_addr,
int  dst_dim,
int *  dst_l,
int *  dst_u,
int *  dst_s,
unsigned long long *  dst_d,
void *  src_addr,
int  src_dim,
int *  src_l,
int *  src_u,
int *  src_s,
unsigned long long *  src_d 
)
326  {
327  unsigned long long dst_buffer_elmts = 1;
328  for (int i = 0; i < dst_dim; i++) {
329  dst_buffer_elmts *= _XMP_M_COUNT_TRIPLETi(dst_l[i], dst_u[i], dst_s[i]);
330  }
331 
332  unsigned long long src_buffer_elmts = 1;
333  for (int i = 0; i < src_dim; i++) {
334  src_buffer_elmts *= _XMP_M_COUNT_TRIPLETi(src_l[i], src_u[i], src_s[i]);
335  }
336 
337  if (dst_buffer_elmts != src_buffer_elmts) {
338  _XMP_fatal("wrong assign statement for gmove");
339  }
340 
341  void *buffer = _XMP_alloc(dst_buffer_elmts * type_size);
342  (*_xmp_pack_array)(buffer, src_addr, type, type_size, src_dim, src_l, src_u, src_s, src_d);
343  (*_xmp_unpack_array)(dst_addr, buffer, type, type_size, dst_dim, dst_l, dst_u, dst_s, dst_d);
344  _XMP_free(buffer);
345 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_gmove_scalar_garray()

void _XMP_gmove_scalar_garray ( void *  scalar,
_XMP_gmv_desc_t gmv_desc_rightp,
int  mode 
)
4705 {
4706  if (mode == _XMP_N_GMOVE_NORMAL){
4707  _XMP_array_t *array = gmv_desc_rightp->a_desc;
4708  int ndims = gmv_desc_rightp->ndims;
4709  int ridx[ndims];
4710 
4711  for (int i=0;i<ndims;i++)
4712  ridx[i] = gmv_desc_rightp->lb[i];
4713 
4714  _XMP_gmove_BCAST_GSCALAR(scalar, array, ridx);
4715  }
4716  else if (mode == _XMP_N_GMOVE_IN){
4717 #ifdef _XMP_MPI3_ONESIDED
4718  _XMP_gmove_inout_scalar(scalar, gmv_desc_rightp, _XMP_N_COARRAY_GET);
4719 #else
4720  _XMP_fatal("Not supported gmove in/out on non-MPI3 environments");
4721 #endif
4722  }
4723  else {
4724  _XMP_fatal("_XMPF_gmove_scalar_garray: wrong gmove mode");
4725  }
4726 }
Here is the call graph for this function:

◆ _XMP_gmove_SENDRECV_GSCALAR()

void _XMP_gmove_SENDRECV_GSCALAR ( void *  dst_addr,
void *  src_addr,
_XMP_array_t dst_array,
_XMP_array_t src_array,
int  dst_ref_index[],
int  src_ref_index[] 
)
891  {
892  _XMP_ASSERT(dst_array->type_size == src_array->type_size);
893  size_t type_size = dst_array->type_size;
894 
895  if(_XMP_IS_SINGLE) {
896  memcpy(dst_addr, src_addr, type_size);
897  return;
898  }
899 
900  _XMP_nodes_ref_t *dst_ref = _XMP_create_gmove_nodes_ref_SCALAR(dst_array, dst_ref_index);
901  _XMP_nodes_ref_t *src_ref = _XMP_create_gmove_nodes_ref_SCALAR(src_array, src_ref_index);
902 
903  _XMP_nodes_t *exec_nodes = _XMP_get_execution_nodes();
904  _XMP_ASSERT(exec_nodes->is_member);
905 
906  int exec_rank = exec_nodes->comm_rank;
907  MPI_Comm *exec_comm = exec_nodes->comm;
908 
909  // calc dst_ranks
910  int dst_shrink_nodes_size = dst_ref->shrink_nodes_size;
911  int dst_ranks[dst_shrink_nodes_size];
912  if (dst_shrink_nodes_size == 1) {
913  dst_ranks[0] = _XMP_calc_linear_rank(dst_ref->nodes, dst_ref->ref);
914  } else {
915  _XMP_translate_nodes_rank_array_to_ranks(dst_ref->nodes, dst_ranks, dst_ref->ref, dst_shrink_nodes_size);
916  }
917 
918  // calc src_ranks
919  int src_shrink_nodes_size = src_ref->shrink_nodes_size;
920  int src_ranks[src_shrink_nodes_size];
921  if (src_shrink_nodes_size == 1) {
922  src_ranks[0] = _XMP_calc_linear_rank(src_ref->nodes, src_ref->ref);
923  } else {
924  _XMP_translate_nodes_rank_array_to_ranks(src_ref->nodes, src_ranks, src_ref->ref, src_shrink_nodes_size);
925  }
926 
927  int wait_recv = _XMP_N_INT_FALSE;
928  MPI_Request gmove_request;
929  for (int i = 0; i < dst_shrink_nodes_size; i++) {
930  if (dst_ranks[i] == exec_rank) {
931  wait_recv = _XMP_N_INT_TRUE;
932 
933  int src_rank;
934  if ((dst_shrink_nodes_size == src_shrink_nodes_size) ||
935  (dst_shrink_nodes_size < src_shrink_nodes_size)) {
936  src_rank = src_ranks[i];
937  } else {
938  src_rank = src_ranks[i % src_shrink_nodes_size];
939  }
940 
941  MPI_Irecv(dst_addr, type_size, MPI_BYTE, src_rank, _XMP_N_MPI_TAG_GMOVE, *exec_comm, &gmove_request);
942  }
943  }
944 
945  for (int i = 0; i < src_shrink_nodes_size; i++) {
946  if (src_ranks[i] == exec_rank) {
947  if ((dst_shrink_nodes_size == src_shrink_nodes_size) ||
948  (dst_shrink_nodes_size < src_shrink_nodes_size)) {
949  if (i < dst_shrink_nodes_size) {
950  MPI_Send(src_addr, type_size, MPI_BYTE, dst_ranks[i], _XMP_N_MPI_TAG_GMOVE, *exec_comm);
951  }
952  } else {
953  int request_size = _XMP_M_COUNT_TRIPLETi(i, dst_shrink_nodes_size, src_shrink_nodes_size);
954  MPI_Request *requests = _XMP_alloc(sizeof(MPI_Request) * request_size);
955 
956  int request_count = 0;
957  for (int j = i; j < dst_shrink_nodes_size; j += src_shrink_nodes_size) {
958  MPI_Isend(src_addr, type_size, MPI_BYTE, dst_ranks[j], _XMP_N_MPI_TAG_GMOVE, *exec_comm, requests + request_count);
959  request_count++;
960  }
961 
962  MPI_Waitall(request_size, requests, MPI_STATUSES_IGNORE);
963  _XMP_free(requests);
964  }
965  }
966  }
967 
968  if (wait_recv) {
969  MPI_Wait(&gmove_request, MPI_STATUS_IGNORE);
970  }
971 
972  _XMP_free(dst_ref);
973  _XMP_free(src_ref);
974 }
Here is the caller graph for this function:

◆ _XMP_gtol_array_ref_triplet()

void _XMP_gtol_array_ref_triplet ( _XMP_array_t array,
int  dim_index,
int *  lower,
int *  upper,
int *  stride 
)
115  {
116  _XMP_array_info_t *array_info = &(array->info[dim_index]);
117  if (array_info->shadow_type == _XMP_N_SHADOW_FULL) {
118  return;
119  }
120 
121  _XMP_template_t *align_template = array->align_template;
122 
123  int align_template_index = array_info->align_template_index;
124  if (align_template_index != _XMP_N_NO_ALIGN_TEMPLATE) {
125  int align_subscript = array_info->align_subscript;
126 
127  int t_lower = *lower - align_subscript,
128  t_upper = *upper - align_subscript,
129  t_stride = *stride;
130 
131  _XMP_template_info_t *ti = &(align_template->info[align_template_index]);
132  int template_ser_lower = ti->ser_lower;
133 
134  _XMP_template_chunk_t *tc = &(align_template->chunk[align_template_index]);
135  int template_par_width = tc->par_width;
136  int template_par_nodes_size = (tc->onto_nodes_info)->size;
137  int template_par_chunk_width = tc->par_chunk_width;
138 
139  switch (tc->dist_manner) {
141  // do nothing
142  break;
143  case _XMP_N_DIST_BLOCK:
144  {
145  t_lower = _XMP_SM_GTOL_BLOCK(t_lower, template_ser_lower, template_par_chunk_width);
146  t_upper = _XMP_SM_GTOL_BLOCK(t_upper, template_ser_lower, template_par_chunk_width);
147  // t_stride does not change
148  } break;
149  case _XMP_N_DIST_CYCLIC:
150  {
151  t_lower = _XMP_SM_GTOL_CYCLIC(t_lower, template_ser_lower, template_par_nodes_size);
152  t_upper = _XMP_SM_GTOL_CYCLIC(t_upper, template_ser_lower, template_par_nodes_size);
153  t_stride = t_stride / template_par_nodes_size;
154  } break;
156  {
157  t_lower = _XMP_SM_GTOL_BLOCK_CYCLIC(template_par_width, t_lower, template_ser_lower, template_par_nodes_size);
158  t_upper = _XMP_SM_GTOL_BLOCK_CYCLIC(template_par_width, t_upper, template_ser_lower, template_par_nodes_size);
159  // t_stride does not change (in block)
160  } break;
161  default:
162  _XMP_fatal("wrong distribute manner for normal shadow");
163  }
164 
165  *lower = t_lower + align_subscript;
166  *upper = t_upper + align_subscript;
167  *stride = t_stride;
168  }
169 
171  switch (array_info->shadow_type) {
172  case _XMP_N_SHADOW_NONE:
173  // do nothing
174  break;
176  switch (array_info->align_manner) {
179  case _XMP_N_ALIGN_BLOCK:
180  {
181  int shadow_size_lo = array_info->shadow_size_lo;
182  *lower += shadow_size_lo;
183  *upper += shadow_size_lo;
184  } break;
185  case _XMP_N_ALIGN_CYCLIC:
186  // FIXME not supported now
188  // FIXME not supported now
189  default:
190  _XMP_fatal("gmove does not support shadow region for cyclic or block-cyclic distribution");
191  } break;
192  default:
193  _XMP_fatal("unknown shadow type");
194  }
195 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_gtol_calc_offset()

unsigned long long _XMP_gtol_calc_offset ( _XMP_array_t a,
int  g_idx[] 
)
2992 {
2993  int ndims = a->dim;
2994  int l_idx[ndims];
2995  for(int i=0;i<ndims;i++)
2996  xmp_array_gtol(a, i+1, g_idx[i], &l_idx[i]);
2997 
2998  //xmp_dbg_printf("g0 = %d, g1 = %d, l0 = %d, l1 = %d\n", g_idx[0], g_idx[1], l_idx[0], l_idx[1]);
2999 
3000  unsigned long long offset = 0;
3001 
3002  for (int i = 0; i < a->dim; i++){
3003  offset += (l_idx[i] * a->info[i].dim_acc * a->type_size);
3004  //xmp_dbg_printf("(%d) acc = %llu, type_size = %d\n", i, a->info[i].dim_acc, a->type_size);
3005  }
3006 
3007  return offset;
3008 
3009 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_sendrecv_ARRAY()

void _XMP_sendrecv_ARRAY ( int  type,
int  type_size,
MPI_Datatype *  mpi_datatype,
_XMP_array_t dst_array,
int *  dst_array_nodes_ref,
int *  dst_lower,
int *  dst_upper,
int *  dst_stride,
unsigned long long *  dst_dim_acc,
_XMP_array_t src_array,
int *  src_array_nodes_ref,
int *  src_lower,
int *  src_upper,
int *  src_stride,
unsigned long long *  src_dim_acc 
)
482  {
483  _XMP_nodes_t *dst_array_nodes = dst_array->array_nodes;
484  _XMP_nodes_t *src_array_nodes = src_array->array_nodes;
485  void *dst_addr = dst_array->array_addr_p;
486  void *src_addr = src_array->array_addr_p;
487  int dst_dim = dst_array->dim;
488  int src_dim = src_array->dim;
489 
490  _XMP_nodes_t *exec_nodes = _XMP_get_execution_nodes();
491  _XMP_ASSERT(exec_nodes->is_member);
492 
493  int exec_rank = exec_nodes->comm_rank;
494  MPI_Comm *exec_comm = exec_nodes->comm;
495 
496  // calc dst_ranks
497  _XMP_nodes_ref_t *dst_ref = _XMP_create_nodes_ref_for_target_nodes(dst_array_nodes, dst_array_nodes_ref, exec_nodes);
498  int dst_shrink_nodes_size = dst_ref->shrink_nodes_size;
499  int dst_ranks[dst_shrink_nodes_size];
500  if (dst_shrink_nodes_size == 1) {
501  dst_ranks[0] = _XMP_calc_linear_rank(dst_ref->nodes, dst_ref->ref);
502  } else {
503  _XMP_translate_nodes_rank_array_to_ranks(dst_ref->nodes, dst_ranks, dst_ref->ref, dst_shrink_nodes_size);
504  }
505 
506  unsigned long long total_elmts = 1;
507  for (int i = 0; i < dst_dim; i++) {
508  total_elmts *= _XMP_M_COUNT_TRIPLETi(dst_lower[i], dst_upper[i], dst_stride[i]);
509  }
510 
511  // calc src_ranks
512  _XMP_nodes_ref_t *src_ref = _XMP_create_nodes_ref_for_target_nodes(src_array_nodes, src_array_nodes_ref, exec_nodes);
513  int src_shrink_nodes_size = src_ref->shrink_nodes_size;
514  int src_ranks[src_shrink_nodes_size];
515  if (src_shrink_nodes_size == 1) {
516  src_ranks[0] = _XMP_calc_linear_rank(src_ref->nodes, src_ref->ref);
517  } else {
518  _XMP_translate_nodes_rank_array_to_ranks(src_ref->nodes, src_ranks, src_ref->ref, src_shrink_nodes_size);
519  }
520 
521  unsigned long long src_total_elmts = 1;
522  for (int i = 0; i < src_dim; i++) {
523  src_total_elmts *= _XMP_M_COUNT_TRIPLETi(src_lower[i], src_upper[i], src_stride[i]);
524  }
525 
526  if (total_elmts != src_total_elmts) {
527  _XMP_fatal("wrong assign statement for gmove");
528  }
529 
530  // recv phase
531  void *recv_buffer = NULL;
532  void *recv_alloc = NULL;
533  int wait_recv = _XMP_N_INT_FALSE;
534  MPI_Request gmove_request;
535  for (int i = 0; i < dst_shrink_nodes_size; i++) {
536  if (dst_ranks[i] == exec_rank) {
537  wait_recv = _XMP_N_INT_TRUE;
538 
539  int src_rank;
540  if ((dst_shrink_nodes_size == src_shrink_nodes_size) ||
541  (dst_shrink_nodes_size < src_shrink_nodes_size)) {
542  src_rank = src_ranks[i];
543  } else {
544  src_rank = src_ranks[i % src_shrink_nodes_size];
545  }
546 
547  for (int i = 0; i < dst_dim; i++) {
548  _XMP_gtol_array_ref_triplet(dst_array, i, &(dst_lower[i]), &(dst_upper[i]), &(dst_stride[i]));
549  }
550  if(dst_dim == 1 && dst_stride[0] == 1){
551  recv_buffer = (char*)dst_addr + type_size * dst_lower[0];
552  }else{
553  recv_alloc = _XMP_alloc(total_elmts * type_size);
554  recv_buffer = recv_alloc;
555  }
556  MPI_Irecv(recv_buffer, total_elmts, *mpi_datatype, src_rank, _XMP_N_MPI_TAG_GMOVE, *exec_comm, &gmove_request);
557  // fprintf(stderr, "DEBUG: Proc(%d), Irecv(src=%d, total_elmnt=%llu)\n", exec_rank, src_rank, total_elmts);
558  }
559  }
560 
561  // send phase
562  for (int i = 0; i < src_shrink_nodes_size; i++) {
563  if (src_ranks[i] == exec_rank) {
564  void *send_buffer = NULL;
565  void *send_alloc = NULL;
566  for (int j = 0; j < src_dim; j++) {
567  _XMP_gtol_array_ref_triplet(src_array, j, &(src_lower[j]), &(src_upper[j]), &(src_stride[j]));
568  }
569  if(src_dim == 1 && src_stride[0] == 1){
570  send_buffer = (char*)src_addr + type_size * src_lower[0];
571  }else{
572  send_alloc = _XMP_alloc(total_elmts * type_size);
573  send_buffer = send_alloc;
574  (*_xmp_pack_array)(send_buffer, src_addr, type, type_size, src_dim, src_lower, src_upper, src_stride, src_dim_acc);
575  }
576  if ((dst_shrink_nodes_size == src_shrink_nodes_size) ||
577  (dst_shrink_nodes_size < src_shrink_nodes_size)) {
578  if (i < dst_shrink_nodes_size) {
579  MPI_Send(send_buffer, total_elmts, *mpi_datatype, dst_ranks[i], _XMP_N_MPI_TAG_GMOVE, *exec_comm);
580  // fprintf(stderr, "DEBUG: Proc(%d), Send(dst=%d, total_elmnt=%llu)\n", exec_rank, dst_ranks[i], total_elmts);
581  }
582  } else {
583  int request_size = _XMP_M_COUNT_TRIPLETi(i, dst_shrink_nodes_size - 1, src_shrink_nodes_size);
584  MPI_Request *requests = _XMP_alloc(sizeof(MPI_Request) * request_size);
585 
586  int request_count = 0;
587  for (int j = i; j < dst_shrink_nodes_size; j += src_shrink_nodes_size) {
588  MPI_Isend(send_buffer, total_elmts, *mpi_datatype, dst_ranks[j], _XMP_N_MPI_TAG_GMOVE, *exec_comm,
589  requests + request_count);
590  // fprintf(stderr, "DEBUG: Proc(%d), Isend(dst=%d, total_elmnt=%llu)\n", exec_rank, dst_ranks[j], total_elmts);
591  request_count++;
592  }
593 
594  MPI_Waitall(request_size, requests, MPI_STATUSES_IGNORE);
595  _XMP_free(requests);
596  }
597 
598  _XMP_free(send_alloc);
599  }
600  }
601 
602  // wait recv phase
603  if (wait_recv) {
604  MPI_Wait(&gmove_request, MPI_STATUS_IGNORE);
605  if(! (dst_dim == 1 && dst_stride[0] == 1)){
606  (*_xmp_unpack_array)(dst_addr, recv_buffer, type, type_size, dst_dim, dst_lower, dst_upper, dst_stride, dst_dim_acc);
607  }
608  _XMP_free(recv_alloc);
609  }
610 
611  _XMP_finalize_nodes_ref(dst_ref);
612  _XMP_finalize_nodes_ref(src_ref);
613 }
Here is the call graph for this function:

◆ xmp_calc_gmove_array_owner_linear_rank_scalar_()

int xmp_calc_gmove_array_owner_linear_rank_scalar_ ( _XMP_array_t **  a,
int *  ref_index 
)
226  {
227  _XMP_array_t *array = *a;
228  _XMP_nodes_t *array_nodes = array->array_nodes;
229  int array_nodes_dim = array_nodes->dim;
230  int rank_array[array_nodes_dim];
231 
232  _XMP_calc_gmove_rank_array_SCALAR(array, ref_index, rank_array);
233 
234  return _XMP_calc_linear_rank_on_target_nodes(array_nodes, rank_array, _XMP_get_execution_nodes());
235 }
Here is the call graph for this function:

Variable Documentation

◆ _alloc_size

int(* _alloc_size)[_XMP_N_MAX_DIM]

◆ _dim_alloc_size

int _dim_alloc_size

◆ _XMP_pack_comm_set

void(* _XMP_pack_comm_set) (void *sendbuf, int sendbuf_size, _XMP_array_t *a, _XMP_comm_set_t *comm_set[][_XMP_N_MAX_DIM])

◆ _XMP_unpack_comm_set

void(* _XMP_unpack_comm_set) (void *recvbuf, int recvbuf_size, _XMP_array_t *a, _XMP_comm_set_t *comm_set[][_XMP_N_MAX_DIM])

◆ gmv_nodes

_XMP_nodes_t* gmv_nodes

◆ n_gmv_nodes

int n_gmv_nodes
_XMPC_running
int _XMPC_running
Definition: xmp_runtime.c:15
_XMP_array_info_type::align_subscript
long long align_subscript
Definition: xmp_data_struct.h:246
_XMP_calc_template_par_triplet
int _XMP_calc_template_par_triplet(_XMP_template_t *template, int template_index, int nodes_rank, int *template_lower, int *template_upper, int *template_stride)
Definition: xmp_template.c:667
_XMP_template_chunk_type::par_upper
long long par_upper
Definition: xmp_data_struct.h:81
_XMP_array_type::array_nodes
_XMP_nodes_t * array_nodes
Definition: xmp_data_struct.h:306
_XMP_nodes_info_type::size
int size
Definition: xmp_data_struct.h:32
_XMP_nodes_info_type
Definition: xmp_data_struct.h:31
_XMP_template_type::info
_XMP_template_info_t info[1]
Definition: xmp_data_struct.h:115
_XMP_array_info_type::align_template_index
int align_template_index
Definition: xmp_data_struct.h:260
_XMP_N_DIST_BLOCK
#define _XMP_N_DIST_BLOCK
Definition: xmp_constant.h:29
_XMP_gmove_garray_scalar
void _XMP_gmove_garray_scalar(_XMP_gmv_desc_t *gmv_desc_leftp, void *scalar, int mode)
Definition: xmp_gmove.c:4732
_XMP_template_info_type::ser_lower
long long ser_lower
Definition: xmp_data_struct.h:72
xmp_is_async
_Bool xmp_is_async()
Definition: xmp_async.c:20
_XMP_nodes_type::is_member
int is_member
Definition: xmp_data_struct.h:46
_XMP_alloc
void * _XMP_alloc(size_t size)
Definition: xmp_util.c:21
_XMP_gmv_desc_type::lb
int * lb
Definition: xmp_data_struct.h:398
_XMP_N_NO_ONTO_NODES
#define _XMP_N_NO_ONTO_NODES
Definition: xmp_constant.h:24
_XMP_array_info_type
Definition: xmp_data_struct.h:194
_XMP_gmv_desc_type::st
int * st
Definition: xmp_data_struct.h:400
_XMP_template_type::chunk
_XMP_template_chunk_t * chunk
Definition: xmp_data_struct.h:112
_XMP_calc_nodes_index_from_inherit_nodes_index
int _XMP_calc_nodes_index_from_inherit_nodes_index(_XMP_nodes_t *nodes, int inherit_nodes_index)
Definition: xmp_nodes.c:1309
_XMP_N_GMOVE_NORMAL
#define _XMP_N_GMOVE_NORMAL
Definition: xmp_constant.h:69
_XMP_calc_linear_rank
int _XMP_calc_linear_rank(_XMP_nodes_t *n, int *rank_array)
Definition: xmp_nodes.c:1035
_XMP_gmv_desc_type::ub
int * ub
Definition: xmp_data_struct.h:399
_XMP_calc_gmove_array_owner_linear_rank_SCALAR
int _XMP_calc_gmove_array_owner_linear_rank_SCALAR(_XMP_array_t *array, int *ref_index)
Definition: xmp_gmove.c:216
XMP_N_GMOVE_INDEX
#define XMP_N_GMOVE_INDEX
Definition: xmp_constant.h:128
_XMP_template_chunk_type::par_lower
long long par_lower
Definition: xmp_data_struct.h:80
_XMP_translate_nodes_rank_array_to_ranks
void _XMP_translate_nodes_rank_array_to_ranks(_XMP_nodes_t *nodes, int *ranks, int *rank_array, int shrink_nodes_size)
Definition: xmp_nodes.c:1270
_XMP_get_array_addr
void * _XMP_get_array_addr(_XMP_array_t *a, int *gidx)
Definition: xmp_gmove.c:47
_XMP_gmv_desc_type
Definition: xmp_data_struct.h:386
_XMP_array_info_type::shadow_type
int shadow_type
Definition: xmp_data_struct.h:248
_XMP_array_info_type::ser_size
int ser_size
Definition: xmp_data_struct.h:201
_XMP_gmv_desc_type::a_ub
int * a_ub
Definition: xmp_data_struct.h:395
_XMP_N_GMOVE_OUT
#define _XMP_N_GMOVE_OUT
Definition: xmp_constant.h:71
_XMP_template_chunk_type::onto_nodes_info
_XMP_nodes_info_t * onto_nodes_info
Definition: xmp_data_struct.h:94
_XMP_check_template_ref_inclusion
int _XMP_check_template_ref_inclusion(int ref_lower, int ref_upper, int ref_stride, _XMP_template_t *t, int index)
Definition: xmp_template.c:243
_XMP_nodes_ref_type::shrink_nodes_size
int shrink_nodes_size
Definition: xmp_data_struct.h:66
_XMP_calc_global_index_HOMECOPY
int _XMP_calc_global_index_HOMECOPY(_XMP_array_t *dst_array, int dst_dim_index, int *dst_l, int *dst_u, int *dst_s, int *src_l, int *src_u, int *src_s)
Definition: xmp_gmove.c:404
_XMP_calc_linear_rank_on_target_nodes
int _XMP_calc_linear_rank_on_target_nodes(_XMP_nodes_t *n, int *rank_array, _XMP_nodes_t *target_nodes)
Definition: xmp_nodes.c:1049
_XMP_SM_GTOL_CYCLIC
#define _XMP_SM_GTOL_CYCLIC(_i, _m, _P)
Definition: xmp_gmove.c:16
_XMP_gmove_array_array_common
void _XMP_gmove_array_array_common(_XMP_gmv_desc_t *gmv_desc_leftp, _XMP_gmv_desc_t *gmv_desc_rightp, int *dst_l, int *dst_u, int *dst_s, unsigned long long *dst_d, int *src_l, int *src_u, int *src_s, unsigned long long *src_d, int mode)
Definition: xmp_gmove.c:2036
_XMP_template_info_type
Definition: xmp_data_struct.h:70
_XMP_nodes_type::comm_rank
int comm_rank
Definition: xmp_data_struct.h:52
_XMP_template_chunk_type::par_chunk_width
unsigned long long par_chunk_width
Definition: xmp_data_struct.h:86
_XMP_SM_GTOL_BLOCK_CYCLIC
#define _XMP_SM_GTOL_BLOCK_CYCLIC(_b, _i, _m, _P)
Definition: xmp_gmove.c:19
xmp_array_gtol
int xmp_array_gtol(xmp_desc_t d, int dim, int g_idx, int *lidx)
Definition: xmp_lib.c:234
_XMP_N_MPI_TAG_GMOVE
#define _XMP_N_MPI_TAG_GMOVE
Definition: xmp_constant.h:10
_XMP_calc_gmove_rank_array_SCALAR
void _XMP_calc_gmove_rank_array_SCALAR(_XMP_array_t *array, int *ref_index, int *rank_array)
Definition: xmp_gmove.c:197
_XMP_N_GMOVE_IN
#define _XMP_N_GMOVE_IN
Definition: xmp_constant.h:70
_XMP_array_info_type::shadow_size_lo
int shadow_size_lo
Definition: xmp_data_struct.h:249
_XMP_gmove_localcopy_ARRAY
void _XMP_gmove_localcopy_ARRAY(int type, int type_size, void *dst_addr, int dst_dim, int *dst_l, int *dst_u, int *dst_s, unsigned long long *dst_d, void *src_addr, int src_dim, int *src_l, int *src_u, int *src_s, unsigned long long *src_d)
Definition: xmp_gmove.c:322
_XMP_array_type::align_template
_XMP_template_t * align_template
Definition: xmp_data_struct.h:312
_XMP_nodes_ref_type::ref
int * ref
Definition: xmp_data_struct.h:65
_XMP_gmv_desc_type::kind
int * kind
Definition: xmp_data_struct.h:397
_XMP_array_info_type::align_manner
int align_manner
Definition: xmp_data_struct.h:197
_XMP_IS_SINGLE
#define _XMP_IS_SINGLE
Definition: xmp_internal.h:50
_XMP_N_INT_FALSE
#define _XMP_N_INT_FALSE
Definition: xmp_constant.h:5
_XMP_N_DIST_CYCLIC
#define _XMP_N_DIST_CYCLIC
Definition: xmp_constant.h:30
_XMP_template_type
Definition: xmp_data_struct.h:98
_XMP_N_SHADOW_FULL
#define _XMP_N_SHADOW_FULL
Definition: xmp_constant.h:66
_XMP_nodes_ref_type
Definition: xmp_data_struct.h:63
_XMP_SM_GTOL_BLOCK
#define _XMP_SM_GTOL_BLOCK(_i, _m, _w)
Definition: xmp_gmove.c:13
_XMP_array_info_type::ser_lower
int ser_lower
Definition: xmp_data_struct.h:199
_XMP_create_nodes_ref_for_target_nodes
_XMP_nodes_ref_t * _XMP_create_nodes_ref_for_target_nodes(_XMP_nodes_t *n, int *rank_array, _XMP_nodes_t *target_nodes)
Definition: xmp_nodes.c:1234
_XMP_template_chunk_type
Definition: xmp_data_struct.h:78
_XMP_template_chunk_type::onto_nodes_index
int onto_nodes_index
Definition: xmp_data_struct.h:92
_XMP_N_DIST_BLOCK_CYCLIC
#define _XMP_N_DIST_BLOCK_CYCLIC
Definition: xmp_constant.h:31
_XMP_N_DIST_DUPLICATION
#define _XMP_N_DIST_DUPLICATION
Definition: xmp_constant.h:28
_XMP_gmv_desc_type::local_data
void * local_data
Definition: xmp_data_struct.h:393
_XMP_array_type
Definition: xmp_data_struct.h:266
_XMP_N_NO_ALIGN_TEMPLATE
#define _XMP_N_NO_ALIGN_TEMPLATE
Definition: xmp_constant.h:23
_XMP_template_chunk_type::par_stride
int par_stride
Definition: xmp_data_struct.h:85
_XMP_gmv_desc_type::ndims
int ndims
Definition: xmp_data_struct.h:389
_XMP_template_chunk_type::dist_manner
int dist_manner
Definition: xmp_data_struct.h:87
_XMP_N_COARRAY_PUT
#define _XMP_N_COARRAY_PUT
Definition: xmp_constant.h:112
_XMP_array_type::type_size
size_t type_size
Definition: xmp_data_struct.h:274
_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_array_type::info
_XMP_array_info_t info[1]
Definition: xmp_data_struct.h:313
_XMP_free
void _XMP_free(void *p)
Definition: xmp_util.c:37
_XMP_ASSERT
#define _XMP_ASSERT(_flag)
Definition: xmp_internal.h:34
_XMP_gmv_desc_type::is_global
_Bool is_global
Definition: xmp_data_struct.h:388
_XMP_finalize_nodes_ref
void _XMP_finalize_nodes_ref(_XMP_nodes_ref_t *nodes_ref)
Definition: xmp_nodes.c:1228
_XMP_N_ALIGN_NOT_ALIGNED
#define _XMP_N_ALIGN_NOT_ALIGNED
Definition: xmp_constant.h:35
_XMP_check_gmove_array_ref_inclusion_SCALAR
int _XMP_check_gmove_array_ref_inclusion_SCALAR(_XMP_array_t *array, int array_index, int ref_index)
Definition: xmp_gmove.c:309
_XMP_array_type::array_addr_p
void * array_addr_p
Definition: xmp_data_struct.h:279
_XMP_N_ALIGN_BLOCK_CYCLIC
#define _XMP_N_ALIGN_BLOCK_CYCLIC
Definition: xmp_constant.h:39
_XMP_nodes_type::comm
_XMP_comm_t * comm
Definition: xmp_data_struct.h:53
_XMP_array_info_type::par_stride
int par_stride
Definition: xmp_data_struct.h:206
_XMP_array_type::dim
int dim
Definition: xmp_data_struct.h:272
_XMP_fatal
void _XMP_fatal(char *msg)
Definition: xmp_util.c:42
XMP_N_GMOVE_RANGE
#define XMP_N_GMOVE_RANGE
Definition: xmp_constant.h:129
_XMP_gmv_desc_type::a_desc
_XMP_array_t * a_desc
Definition: xmp_data_struct.h:391
_XMP_array_info_type::dim_acc
unsigned long long dim_acc
Definition: xmp_data_struct.h:242
_XMP_array_info_type::par_size
int par_size
Definition: xmp_data_struct.h:207
_XMP_template_chunk_type::par_width
unsigned long long par_width
Definition: xmp_data_struct.h:82
_XMP_N_COARRAY_GET
#define _XMP_N_COARRAY_GET
Definition: xmp_constant.h:111
_XMP_nodes_ref_type::nodes
_XMP_nodes_t * nodes
Definition: xmp_data_struct.h:64
_XMP_gmove_SENDRECV_GSCALAR
void _XMP_gmove_SENDRECV_GSCALAR(void *dst_addr, void *src_addr, _XMP_array_t *dst_array, _XMP_array_t *src_array, int dst_ref_index[], int src_ref_index[])
Definition: xmp_gmove.c:889
_XMP_nodes_type::dim
int dim
Definition: xmp_data_struct.h:47
_XMP_nodes_type
Definition: xmp_data_struct.h:40
_XMP_M_COUNT_TRIPLETi
#define _XMP_M_COUNT_TRIPLETi(l_, u_, s_)
Definition: xmp_gpu_func.hpp:25
_XMP_N_MAX_DIM
#define _XMP_N_MAX_DIM
Definition: xmp_constant.h:6
_XMP_calc_template_owner_SCALAR
int _XMP_calc_template_owner_SCALAR(_XMP_template_t *ref_template, int dim_index, long long ref_index)
Definition: xmp_template.c:632
_XMP_gtol_calc_offset
unsigned long long _XMP_gtol_calc_offset(_XMP_array_t *a, int g_idx[])
Definition: xmp_gmove.c:2991
_XMP_gmove_BCAST_GSCALAR
void _XMP_gmove_BCAST_GSCALAR(void *dst_addr, _XMP_array_t *array, int ref_index[])
Definition: xmp_gmove.c:687
_XMP_N_SHADOW_NORMAL
#define _XMP_N_SHADOW_NORMAL
Definition: xmp_constant.h:65
_XMP_gmove_inout_scalar
void _XMP_gmove_inout_scalar(void *scalar, _XMP_gmv_desc_t *gmv_desc, int rdma_type)
_XMP_N_SHADOW_NONE
#define _XMP_N_SHADOW_NONE
Definition: xmp_constant.h:64
_XMPF_running
int _XMPF_running
Definition: xmp_runtime.c:16
_XMP_array_type::type
int type
Definition: xmp_data_struct.h:273
_XMP_N_INT_TRUE
#define _XMP_N_INT_TRUE
Definition: xmp_constant.h:4
_XMP_N_ALIGN_GBLOCK
#define _XMP_N_ALIGN_GBLOCK
Definition: xmp_constant.h:40
_XMP_normalize_array_section
void _XMP_normalize_array_section(int *lower, int *upper, int *stride)
_XMP_template_chunk_type::mapping_array
long long * mapping_array
Definition: xmp_data_struct.h:88
_XMP_array_type::is_allocated
_Bool is_allocated
Definition: xmp_data_struct.h:270
_XMP_gtol_array_ref_triplet
void _XMP_gtol_array_ref_triplet(_XMP_array_t *array, int dim_index, int *lower, int *upper, int *stride)
Definition: xmp_gmove.c:114
_XMP_get_execution_nodes
void * _XMP_get_execution_nodes(void)
Definition: xmp_nodes_stack.c:46
_XMP_comm_set_type
Definition: xmp_data_struct.h:439
_XMP_array_info_type::par_lower
int par_lower
Definition: xmp_data_struct.h:204
_XMP_gmv_desc_type::a_lb
int * a_lb
Definition: xmp_data_struct.h:394