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

Macros

#define MPI_PORTABLE_PLATFORM_H
 

Functions

void xmp_dbg_printf (char *fmt,...)
 
_XMP_nodes_t_XMP_create_temporary_nodes (_XMP_nodes_t *n)
 
_XMP_nodes_t_XMP_init_nodes_struct_GLOBAL (int dim, int *dim_size, int is_static)
 
_XMP_nodes_t_XMP_init_nodes_struct_EXEC (int dim, int *dim_size, int is_static)
 
_XMP_nodes_t_XMP_init_nodes_struct_NODES_NUMBER (int dim, int ref_lower, int ref_upper, int ref_stride, int *dim_size, int is_static)
 
_XMP_nodes_t_XMP_init_nodes_struct_NODES_NAMED (int dim, _XMP_nodes_t *ref_nodes, int *shrink, int *ref_lower, int *ref_upper, int *ref_stride, int *dim_size, int is_static)
 
void _XMP_init_nodes_STATIC_GLOBAL (_XMP_nodes_t **nodes, int dim,...)
 
void _XMP_init_nodes_DYNAMIC_GLOBAL (_XMP_nodes_t **nodes, int dim,...)
 
void _XMP_init_nodes_STATIC_EXEC (_XMP_nodes_t **nodes, int dim,...)
 
void _XMP_init_nodes_DYNAMIC_EXEC (_XMP_nodes_t **nodes, int dim,...)
 
void _XMP_init_nodes_STATIC_NODES_NUMBER (_XMP_nodes_t **nodes, int dim, int ref_lower, int ref_upper, int ref_stride,...)
 
void _XMP_init_nodes_DYNAMIC_NODES_NUMBER (_XMP_nodes_t **nodes, int dim, int ref_lower, int ref_upper, int ref_stride,...)
 
void _XMP_init_nodes_STATIC_NODES_NAMED (_XMP_nodes_t **nodes, int dim, _XMP_nodes_t *ref_nodes,...)
 
void _XMP_init_nodes_DYNAMIC_NODES_NAMED (_XMP_nodes_t **nodes, int dim, _XMP_nodes_t *ref_nodes,...)
 
void _XMP_finalize_nodes (_XMP_nodes_t *nodes)
 
int _XMP_exec_task_GLOBAL_PART (_XMP_task_desc_t **task_desc, int ref_lower, int ref_upper, int ref_stride)
 
int _XMP_exec_task_NODES_ENTIRE (_XMP_task_desc_t **task_desc, _XMP_nodes_t *ref_nodes)
 
int _XMP_exec_task_NODES_ENTIRE_nocomm (_XMP_nodes_t *ref_nodes)
 
void _XMP_exec_task_NODES_FINALIZE (_XMP_task_desc_t *task_desc)
 
int _XMP_exec_task_NODES_PART (_XMP_task_desc_t **task_desc, _XMP_nodes_t *ref_nodes,...)
 
int _XMP_exec_task_NODES_PART_nocomm (_XMP_nodes_t *ref_nodes,...)
 
_XMP_nodes_t_XMP_create_nodes_by_comm (int is_member, _XMP_comm_t *comm)
 
void _XMP_calc_rank_array (_XMP_nodes_t *n, int *rank_array, int linear_rank)
 
int _XMP_calc_linear_rank (_XMP_nodes_t *n, int *rank_array)
 
int _XMP_calc_linear_rank_on_target_nodes (_XMP_nodes_t *n, int *rank_array, _XMP_nodes_t *target_nodes)
 
_Bool _XMP_calc_coord_on_target_nodes2 (_XMP_nodes_t *n, int *ncoord, _XMP_nodes_t *target_n, int *target_ncoord)
 
_Bool _XMP_calc_coord_on_target_nodes (_XMP_nodes_t *n, int *ncoord, _XMP_nodes_t *target_n, int *target_ncoord)
 
_XMP_nodes_ref_t_XMP_init_nodes_ref (_XMP_nodes_t *n, int *rank_array)
 
void _XMP_finalize_nodes_ref (_XMP_nodes_ref_t *nodes_ref)
 
_XMP_nodes_ref_t_XMP_create_nodes_ref_for_target_nodes (_XMP_nodes_t *n, int *rank_array, _XMP_nodes_t *target_nodes)
 
void _XMP_translate_nodes_rank_array_to_ranks (_XMP_nodes_t *nodes, int *ranks, int *rank_array, int shrink_nodes_size)
 
int _XMP_get_next_rank (_XMP_nodes_t *nodes, int *rank_array)
 
int _XMP_calc_nodes_index_from_inherit_nodes_index (_XMP_nodes_t *nodes, int inherit_nodes_index)
 

Macro Definition Documentation

◆ MPI_PORTABLE_PLATFORM_H

#define MPI_PORTABLE_PLATFORM_H

Function Documentation

◆ _XMP_calc_coord_on_target_nodes()

_Bool _XMP_calc_coord_on_target_nodes ( _XMP_nodes_t n,
int *  ncoord,
_XMP_nodes_t target_n,
int *  target_ncoord 
)
1152 {
1153  if(_XMP_calc_coord_on_target_nodes2(n, ncoord, target_n, target_ncoord))
1154  return true;
1155 
1156  _XMP_nodes_t *p = n->inherit_nodes;
1157  if (p){
1158  int pcoord[_XMP_N_MAX_DIM];
1159  int rank = _XMP_calc_linear_rank(n, ncoord);
1160  //_XMP_calc_rank_array(p, pcoord, rank);
1161 
1162  _XMP_nodes_inherit_info_t *inherit_info = n->inherit_info;
1163  int multiplier[_XMP_N_MAX_DIM];
1164 
1165  int dim = p->dim;
1166  multiplier[0] = 1;
1167  for(int i=1;i<dim;i++)
1168  multiplier[i] = multiplier[i-1] * inherit_info[i].size;
1169 
1170  int j = rank;
1171  for (int i=dim;i>= 0;i--){
1172  if(inherit_info[i].shrink){
1173  pcoord[i] = p->info[i].rank;
1174  }
1175  else{
1176  int rank_dim = j / multiplier[i];
1177  pcoord[i] = inherit_info[i].stride * rank_dim + inherit_info[i].lower;
1178  j = j % multiplier[i];
1179  }
1180  }
1181 
1182  if(_XMP_calc_coord_on_target_nodes(p, pcoord, target_n, target_ncoord))
1183  return true;
1184  }
1185 
1186  return false;
1187 }
Here is the call graph for this function:

◆ _XMP_calc_coord_on_target_nodes2()

_Bool _XMP_calc_coord_on_target_nodes2 ( _XMP_nodes_t n,
int *  ncoord,
_XMP_nodes_t target_n,
int *  target_ncoord 
)
1092 {
1093  if (n == target_n){
1094  //printf("%d, %d\n", n->dim, target_n->dim);
1095  memcpy(target_ncoord, ncoord, sizeof(int) * n->dim);
1096  return true;
1097  }
1098  else if (n->attr == _XMP_ENTIRE_NODES && target_n->attr == _XMP_ENTIRE_NODES){
1099  int rank = _XMP_calc_linear_rank(n, ncoord);
1100  _XMP_calc_rank_array(target_n, target_ncoord, rank);
1101  //xmp_dbg_printf("ncoord = %d, target_ncoord = %d\n", ncoord[0], target_ncoord[0]);
1102  return true;
1103  }
1104 
1105  _XMP_nodes_t *target_p = target_n->inherit_nodes;
1106  if (target_p){
1107  int target_pcoord[_XMP_N_MAX_DIM];
1108  if (_XMP_calc_coord_on_target_nodes2(n, ncoord, target_p, target_pcoord)){
1109  //int target_prank = _XMP_calc_linear_rank(target_p, target_pcoord);
1110  /* printf("dim = %d, m0 = %d, m1 = %d\n", */
1111  /* target_n->dim, target_n->info[0].multiplier, target_n->info[1].multiplier); */
1112  //_XMP_calc_rank_array(target_n, target_ncoord, target_prank);
1113 
1114  _XMP_nodes_inherit_info_t *inherit_info = target_n->inherit_info;
1115  int target_rank = 0;
1116  int multiplier = 1;
1117 
1118  for (int i = 0; i < target_p->dim; i++){
1119 
1120  if (inherit_info[i].shrink){
1121  ;
1122  }
1123  else if (target_pcoord[i] < inherit_info[i].lower || target_pcoord[i] > inherit_info[i].upper){
1124  for (int i = 0; i < target_n->dim; i++){
1125  target_ncoord[i] = -1;
1126  }
1127  return true;
1128  }
1129  else {
1130  int target_rank_dim = (target_pcoord[i] - inherit_info[i].lower) / inherit_info[i].stride;
1131  target_rank += multiplier * target_rank_dim;
1132  //multiplier *= inherit_info[i].size;
1133  multiplier *= _XMP_M_COUNT_TRIPLETi(inherit_info[i].lower, inherit_info[i].upper, inherit_info[i].stride);
1134  }
1135  }
1136 
1137  _XMP_calc_rank_array(target_n, target_ncoord, target_rank);
1138 
1139  return true;
1140  }
1141  }
1142 
1143  return false;
1144 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_calc_linear_rank()

int _XMP_calc_linear_rank ( _XMP_nodes_t n,
int *  rank_array 
)
1036 {
1037  int acc_rank = 0;
1038  int acc_nodes_size = 1;
1039  int nodes_dim = n->dim;
1040 
1041  for(int i=0;i<nodes_dim;i++){
1042  acc_rank += (rank_array[i]) * acc_nodes_size;
1043  acc_nodes_size *= n->info[i].size;
1044  }
1045 
1046  return acc_rank;
1047 }
Here is the caller graph for this function:

◆ _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 
)
1050 {
1051  if(_XMP_compare_nodes(n, target_nodes)){
1052  return _XMP_calc_linear_rank(n, rank_array);
1053  }
1054  else{
1055  _XMP_nodes_t *inherit_nodes = n->inherit_nodes;
1056  if(inherit_nodes == NULL){
1057  // FIXME implement
1058  _XMP_fatal("unsupported case: gmove");
1059  return _XMP_N_INVALID_RANK; // XXX dummy;
1060  }
1061  else{
1062  int inherit_nodes_dim = inherit_nodes->dim;
1063  int *new_rank_array = _XMP_alloc(sizeof(int) * inherit_nodes_dim);
1064  _XMP_nodes_inherit_info_t *inherit_info = n->inherit_info;
1065 
1066  int j = 0;
1067  for(int i=0;i<inherit_nodes_dim;i++){
1068  if(inherit_info[i].shrink){
1069  // FIXME how implement ???
1070  new_rank_array[i] = 0;
1071  }
1072  else{
1073  new_rank_array[i] = ((inherit_info[i].stride) * rank_array[j]) + (inherit_info[i].lower);
1074  j++;
1075  }
1076  }
1077 
1078  int ret = _XMP_calc_linear_rank_on_target_nodes(inherit_nodes, new_rank_array, target_nodes);
1079  _XMP_free(new_rank_array);
1080  return ret;
1081  }
1082  }
1083 }
Here is the caller graph for this function:

◆ _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 
)
1310 {
1311  _XMP_nodes_t *inherit_nodes = nodes->inherit_nodes;
1312  if(inherit_nodes == NULL)
1313  _XMP_fatal("inherit nodes is NULL");
1314 
1315  int nodes_index = 0;
1316  int inherit_nodes_index_count = 0;
1317  int inherit_nodes_dim = inherit_nodes->dim;
1318  for(int i=0;i<inherit_nodes_dim;i++,inherit_nodes_index_count++){
1319  if(inherit_nodes_index_count == inherit_nodes_index){
1320  return nodes_index;
1321  }
1322  else{
1323  if(!nodes->inherit_info[i].shrink)
1324  nodes_index++;
1325  }
1326  }
1327 
1328  _XMP_fatal("the function does not reach here");
1329  return 0;
1330 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_calc_rank_array()

void _XMP_calc_rank_array ( _XMP_nodes_t n,
int *  rank_array,
int  linear_rank 
)
1027 {
1028  int j = linear_rank;
1029  for(int i=n->dim-1;i>=0;i--){
1030  rank_array[i] = j / n->info[i].multiplier;
1031  j = j % n->info[i].multiplier;
1032  }
1033 }
Here is the caller graph for this function:

◆ _XMP_create_nodes_by_comm()

_XMP_nodes_t* _XMP_create_nodes_by_comm ( int  is_member,
_XMP_comm_t comm 
)
1008 {
1009  int size;
1010  MPI_Comm_size(*((MPI_Comm *)comm), &size);
1011 
1012  _XMP_nodes_t *n = _XMP_create_new_nodes(is_member, 1, size, comm);
1013 
1014  n->info[0].size = size;
1015  if(is_member)
1016  MPI_Comm_rank(*((MPI_Comm *)comm), &(n->info[0].rank));
1017 
1018  // calc inherit info
1019  n->inherit_nodes = NULL;
1020  n->inherit_info = NULL;
1021  n->info[0].multiplier = 1;
1022 
1023  return n;
1024 }

◆ _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 
)
1235 {
1236  if(_XMP_compare_nodes(n, target_nodes)){
1237  return _XMP_init_nodes_ref(n, rank_array);
1238  }
1239  else{
1240  _XMP_nodes_t *inherit_nodes = n->inherit_nodes;
1241  if(inherit_nodes == NULL){
1242  // FIXME implement
1243  _XMP_fatal("unsupported case: gmove");
1244  return NULL; // XXX dummy;
1245  // return _XMP_init_nodes_ref(n, rank_array);
1246  }
1247  else{
1248  int inherit_nodes_dim = inherit_nodes->dim;
1249  int *new_rank_array = _XMP_alloc(sizeof(int) * inherit_nodes_dim);
1250  _XMP_nodes_inherit_info_t *inherit_info = n->inherit_info;
1251 
1252  int j = 0;
1253  for(int i = 0; i < inherit_nodes_dim; i++){
1254  if(inherit_info[i].shrink){
1255  new_rank_array[i] = _XMP_N_UNSPECIFIED_RANK;
1256  }
1257  else{
1258  new_rank_array[i] = ((inherit_info[i].stride) * rank_array[j]) + (inherit_info[i].lower);
1259  j++;
1260  }
1261  }
1262 
1263  _XMP_nodes_ref_t *ret = _XMP_create_nodes_ref_for_target_nodes(inherit_nodes, new_rank_array, target_nodes);
1264  _XMP_free(new_rank_array);
1265  return ret;
1266  }
1267  }
1268 }
Here is the caller graph for this function:

◆ _XMP_create_temporary_nodes()

_XMP_nodes_t* _XMP_create_temporary_nodes ( _XMP_nodes_t n)
258 {
259  int is_member = n->is_member;
260  int onto_nodes_dim = n->dim;
261  int onto_nodes_size = n->comm_size;
262  int dim_size[onto_nodes_dim];
263  _XMP_nodes_t *new_node = _XMP_create_new_nodes(is_member, onto_nodes_dim, onto_nodes_size, n->comm);
264  _XMP_nodes_inherit_info_t *inherit_info = _XMP_alloc(sizeof(_XMP_nodes_inherit_info_t) * onto_nodes_dim);
265 
266  for(int i=0;i<onto_nodes_dim;i++){
267  inherit_info[i].shrink = _XMP_N_INT_FALSE;
268  inherit_info[i].lower = 1;
269  inherit_info[i].upper = n->info[i].size;
270  inherit_info[i].stride = 1;
271  dim_size[i] = n->info[i].size;
272  }
273 
274  new_node->inherit_nodes = new_node;
275  new_node->inherit_info = inherit_info;
276 
277  _XMP_init_nodes_info(new_node, dim_size, _XMP_N_INT_TRUE);
278 
279  new_node->info[0].multiplier = 1;
280  for(int i=1;i<onto_nodes_dim;i++)
281  new_node->info[i].multiplier = new_node->info[i-1].multiplier * dim_size[i-1];
282 
283  return new_node;
284 }
Here is the caller graph for this function:

◆ _XMP_exec_task_GLOBAL_PART()

int _XMP_exec_task_GLOBAL_PART ( _XMP_task_desc_t **  task_desc,
int  ref_lower,
int  ref_upper,
int  ref_stride 
)
836 {
837  int lower[1], upper[1], stride[1];
838 
839  lower[0] = ref_lower;
840  upper[0] = ref_upper;
841  stride[0] = ref_stride;
842 
843  _XMP_task_desc_t *desc = NULL;
844  if(*task_desc == NULL){
845  desc = (_XMP_task_desc_t *)_XMP_alloc(sizeof(_XMP_task_desc_t));
846  *task_desc = desc;
847  }
848  else{
849  desc = *task_desc;
850  if(_XMP_compare_task_exec_cond(desc, _XMP_world_nodes, lower, upper, stride)){
851  if(desc->execute){
852  _XMP_push_nodes(desc->nodes);
853  return _XMP_N_INT_TRUE;
854  }
855  else {
856  return _XMP_N_INT_FALSE;
857  }
858  }
859  else{
860  if(desc->nodes != NULL){
861  _XMP_finalize_nodes(desc->nodes);
862  }
863  }
864  }
865 
866  _XMP_nodes_t *n = NULL;
867  _XMP_init_nodes_STATIC_NODES_NUMBER(&n, 1, ref_lower, ref_upper, ref_stride, _XMP_M_COUNT_TRIPLETi(ref_lower, ref_upper, ref_stride));
868  _XMP_set_task_desc(desc, n, n->is_member, _XMP_world_nodes, lower, upper, stride);
869 
870  if(n->is_member){
871  _XMP_push_nodes(n);
872  return _XMP_N_INT_TRUE;
873  }
874  else{
875  return _XMP_N_INT_FALSE;
876  }
877 }
Here is the call graph for this function:

◆ _XMP_exec_task_NODES_ENTIRE()

int _XMP_exec_task_NODES_ENTIRE ( _XMP_task_desc_t **  task_desc,
_XMP_nodes_t ref_nodes 
)
880 {
881  if(ref_nodes->is_member){
882  _XMP_push_nodes(ref_nodes);
883  return true;
884  }
885  else{
886  return false;
887  }
888 }
Here is the call graph for this function:

◆ _XMP_exec_task_NODES_ENTIRE_nocomm()

int _XMP_exec_task_NODES_ENTIRE_nocomm ( _XMP_nodes_t ref_nodes)
891 {
892  return ref_nodes->is_member;
893 }

◆ _XMP_exec_task_NODES_FINALIZE()

void _XMP_exec_task_NODES_FINALIZE ( _XMP_task_desc_t task_desc)
896 {
897  if(task_desc == NULL) return;
898 
899  if(xmp_is_async()){
900  // Keep a node descriptor. After wait_async directive,
901  // the node descriptor is freed.
903  _XMP_free(task_desc);
904  return;
905  }
906 
907  _XMP_finalize_nodes(task_desc->nodes);
908  _XMP_free(task_desc);
909 }
Here is the call graph for this function:

◆ _XMP_exec_task_NODES_PART()

int _XMP_exec_task_NODES_PART ( _XMP_task_desc_t **  task_desc,
_XMP_nodes_t ref_nodes,
  ... 
)
912 {
913  va_list args;
914  va_start(args, ref_nodes);
915  int ref_dim = ref_nodes->dim;
916  int shrink[ref_dim], ref_lower[ref_dim], ref_upper[ref_dim], ref_stride[ref_dim];
917 
918  int acc_dim_size = 1;
919  for(int i=0;i<ref_dim;i++){
920  shrink[i] = va_arg(args, int);
921  if(!shrink[i]){
922  ref_lower[i] = va_arg(args, int);
923  ref_upper[i] = va_arg(args, int);
924  ref_stride[i] = va_arg(args, int);
925  acc_dim_size *= _XMP_M_COUNT_TRIPLETi(ref_lower[i], ref_upper[i], ref_stride[i]);
926  }
927  else{
928  ref_lower[i] = 0;
929  ref_upper[i] = 0;
930  ref_stride[i] = -1;
931  }
932  }
933  va_end(args);
934 
935  _XMP_task_desc_t *desc = NULL;
936  if(*task_desc == NULL){
937  desc = (_XMP_task_desc_t *)_XMP_alloc(sizeof(_XMP_task_desc_t));
938  *task_desc = desc;
939  // printf("communicator created.\n");
940  }
941  else{
942  desc = *task_desc;
943  if(_XMP_compare_task_exec_cond(desc, ref_nodes, ref_lower, ref_upper, ref_stride)){
944  if(desc->execute){
945  _XMP_push_nodes(desc->nodes);
946  // printf("communicator reused.\n");
947  return _XMP_N_INT_TRUE;
948  }
949  else{
950  return _XMP_N_INT_FALSE;
951  }
952  }
953  else{
954  if(desc->nodes != NULL){
955  _XMP_finalize_nodes(desc->nodes);
956  // printf("communicator freed.\n");
957  }
958  }
959  }
960 
961  _XMP_nodes_t *n = _XMP_init_nodes_struct_NODES_NAMED(1, ref_nodes, shrink, ref_lower, ref_upper, ref_stride,
962  &acc_dim_size, _XMP_N_INT_TRUE);
963  _XMP_set_task_desc(desc, n, n->is_member, ref_nodes, ref_lower, ref_upper, ref_stride);
964 
965  if(n->is_member){
966  _XMP_push_nodes(n);
967  return _XMP_N_INT_TRUE;
968  }
969  else{
970  return _XMP_N_INT_FALSE;
971  }
972 }
Here is the call graph for this function:

◆ _XMP_exec_task_NODES_PART_nocomm()

int _XMP_exec_task_NODES_PART_nocomm ( _XMP_nodes_t ref_nodes,
  ... 
)
975 {
976  if (!ref_nodes->is_member) return _XMP_N_INT_FALSE;
977 
978  int ref_dim = ref_nodes->dim;
979  int shrink[ref_dim], ref_lower[ref_dim], ref_upper[ref_dim], ref_stride[ref_dim];
980 
981  va_list args;
982  va_start(args, ref_nodes);
983 
984  for(int i=0;i<ref_dim;i++){
985  shrink[i] = va_arg(args, int);
986  if(!shrink[i]){
987  ref_lower[i] = va_arg(args, int);
988  ref_upper[i] = va_arg(args, int);
989  ref_stride[i] = va_arg(args, int);
990  }
991  }
992  va_end(args);
993 
994  for(int i=0;i<ref_dim;i++){
995  if(shrink[i]) continue;
996 
997  int me = ref_nodes->info[i].rank + 1;
998 
999  if (me < ref_lower[i] || ref_upper[i] < me) return _XMP_N_INT_FALSE;
1000  if ((me-ref_lower[i]) % ref_stride[i] != 0) return _XMP_N_INT_FALSE;
1001  }
1002 
1003  return _XMP_N_INT_TRUE;
1004 }

◆ _XMP_finalize_nodes()

void _XMP_finalize_nodes ( _XMP_nodes_t nodes)
817 {
818  if(!nodes->use_subcomm)
819  _XMP_finalize_comm(nodes->comm);
820 
821  if(nodes->subcomm != NULL){
822  int dim = nodes->dim;
823  int num_comms = 1<<dim;
824  for(int i=0;i<num_comms;i++){
825  MPI_Comm *comm = (MPI_Comm *)nodes->subcomm;
826  MPI_Comm_free(&comm[i]);
827  }
828  free(nodes->subcomm);
829  }
830 
831  _XMP_free(nodes->inherit_info);
832  _XMP_free(nodes);
833 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_finalize_nodes_ref()

void _XMP_finalize_nodes_ref ( _XMP_nodes_ref_t nodes_ref)
1229 {
1230  _XMP_free(nodes_ref->ref);
1231  _XMP_free(nodes_ref);
1232 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_get_next_rank()

int _XMP_get_next_rank ( _XMP_nodes_t nodes,
int *  rank_array 
)
1292 {
1293  int i, dim = nodes->dim;
1294  for(i=0;i<dim;i++){
1295  int size = nodes->info[i].size;
1296  (rank_array[i])++;
1297  if(rank_array[i] == size)
1298  rank_array[i] = 0;
1299  else
1300  break;
1301  }
1302 
1303  if(i == dim)
1304  return _XMP_N_INT_FALSE;
1305  else
1306  return _XMP_N_INT_TRUE;
1307 }

◆ _XMP_init_nodes_DYNAMIC_EXEC()

void _XMP_init_nodes_DYNAMIC_EXEC ( _XMP_nodes_t **  nodes,
int  dim,
  ... 
)
657  {
658  int dim_size[dim - 1];
659 
660  va_list args;
661  va_start(args, dim);
662  for(int i=0;i<dim-1;i++){
663  int dim_size_temp = va_arg(args, int);
664  if(dim_size_temp <= 0)
665  _XMP_fatal("<nodes-size> should be less or equal to zero");
666 
667  dim_size[i] = dim_size_temp;
668  }
669 
671  *nodes = n;
672 
673  int *last_dim_size_p = va_arg(args, int *);
674  *last_dim_size_p = n->info[dim - 1].size;
675  for(int i=0;i<dim;i++){
676  int *rank_p = va_arg(args, int *);
677  *rank_p = n->info[i].rank;
678  }
679  va_end(args);
680 }
Here is the call graph for this function:

◆ _XMP_init_nodes_DYNAMIC_GLOBAL()

void _XMP_init_nodes_DYNAMIC_GLOBAL ( _XMP_nodes_t **  nodes,
int  dim,
  ... 
)
598 {
599  int dim_size[dim];
600  int *dim_size_p[dim];
601 
602  va_list args;
603  va_start(args, dim);
604  for(int i=0;i<dim;i++){
605  dim_size_p[i] = NULL;
606  int dim_size_temp = va_arg(args, int);
607  if(dim_size_temp == -1){
608  dim_size_p[i] = va_arg(args, int *);
609  }
610  else if (dim_size_temp <= 0){
611  _XMP_fatal("<nodes-size> should be more than zero");
612  }
613 
614  dim_size[i] = dim_size_temp;
615  }
616 
618  *nodes = n;
619 
620  // int *last_dim_size_p = va_arg(args, int *);
621  // *last_dim_size_p = n->info[dim - 1].size;
622  for(int i=0;i<dim;i++){
623  if(dim_size_p[i])
624  *dim_size_p[i] = n->info[i].size;
625  }
626 
627  for(int i=0;i<dim;i++){
628  int *rank_p = va_arg(args, int *);
629  *rank_p = n->info[i].rank;
630  }
631  va_end(args);
632 }
Here is the call graph for this function:

◆ _XMP_init_nodes_DYNAMIC_NODES_NAMED()

void _XMP_init_nodes_DYNAMIC_NODES_NAMED ( _XMP_nodes_t **  nodes,
int  dim,
_XMP_nodes_t ref_nodes,
  ... 
)
777 {
778  int ref_dim = ref_nodes->dim;
779  int shrink[ref_dim];
780  int ref_lower[ref_dim], ref_upper[ref_dim], ref_stride[ref_dim];
781  int dim_size[dim-1];
782 
783  va_list args;
784  va_start(args, ref_nodes);
785  for(int i=0;i<ref_dim;i++){
786  shrink[i] = va_arg(args, int);
787  if(!shrink[i]){
788  ref_lower[i] = va_arg(args, int);
789  ref_upper[i] = va_arg(args, int);
790  ref_stride[i] = va_arg(args, int);
791  }
792  }
793 
794  for(int i=0;i<dim-1;i++){
795  int dim_size_temp = va_arg(args, int);
796  if(dim_size_temp <= 0)
797  _XMP_fatal("<nodes-size> should be less or equal to zero");
798  else
799  dim_size[i] = dim_size_temp;
800  }
801 
802  _XMP_nodes_t *n = _XMP_init_nodes_struct_NODES_NAMED(dim, ref_nodes, shrink, ref_lower, ref_upper,
803  ref_stride, dim_size, _XMP_N_INT_FALSE);
804  *nodes = n;
805 
806  int *last_dim_size_p = va_arg(args, int *);
807  *last_dim_size_p = n->info[dim - 1].size;
808 
809  for(int i=0;i<dim;i++){
810  int *rank_p = va_arg(args, int *);
811  *rank_p = n->info[i].rank;
812  }
813  va_end(args);
814 }
Here is the call graph for this function:

◆ _XMP_init_nodes_DYNAMIC_NODES_NUMBER()

void _XMP_init_nodes_DYNAMIC_NODES_NUMBER ( _XMP_nodes_t **  nodes,
int  dim,
int  ref_lower,
int  ref_upper,
int  ref_stride,
  ... 
)
710 {
711  int dim_size[dim - 1];
712 
713  va_list args;
714  va_start(args, ref_stride);
715  for(int i=0;i<dim-1;i++){
716  int dim_size_temp = va_arg(args, int);
717  if(dim_size_temp <= 0)
718  _XMP_fatal("<nodes-size> should be less or equal to zero");
719 
720  dim_size[i] = dim_size_temp;
721  }
722 
723  _XMP_nodes_t *n = _XMP_init_nodes_struct_NODES_NUMBER(dim, ref_lower, ref_upper,
724  ref_stride, dim_size, _XMP_N_INT_FALSE);
725  *nodes = n;
726 
727  int *last_dim_size_p = va_arg(args, int *);
728  *last_dim_size_p = n->info[dim - 1].size;
729 
730  for(int i=0;i<dim;i++){
731  int *rank_p = va_arg(args, int *);
732  *rank_p = n->info[i].rank;
733  }
734  va_end(args);
735 }
Here is the call graph for this function:

◆ _XMP_init_nodes_ref()

_XMP_nodes_ref_t* _XMP_init_nodes_ref ( _XMP_nodes_t n,
int *  rank_array 
)
1208 {
1209  _XMP_nodes_ref_t *nodes_ref = _XMP_alloc(sizeof(_XMP_nodes_ref_t));
1210  int dim = n->dim;
1211  int *new_rank_array = _XMP_alloc(sizeof(int) * dim);
1212 
1213  int shrink_nodes_size = 1;
1214  for(int i=0;i<dim;i++){
1215  new_rank_array[i] = rank_array[i];
1216  if(new_rank_array[i] == _XMP_N_UNSPECIFIED_RANK){
1217  shrink_nodes_size *= (n->info[i].size);
1218  }
1219  }
1220 
1221  nodes_ref->nodes = n;
1222  nodes_ref->ref = new_rank_array;
1223  nodes_ref->shrink_nodes_size = shrink_nodes_size;
1224 
1225  return nodes_ref;
1226 }
Here is the call graph for this function:

◆ _XMP_init_nodes_STATIC_EXEC()

void _XMP_init_nodes_STATIC_EXEC ( _XMP_nodes_t **  nodes,
int  dim,
  ... 
)
634  {
635  int dim_size[dim];
636 
637  va_list args;
638  va_start(args, dim);
639  for (int i=0;i<dim;i++){
640  int dim_size_temp = va_arg(args, int);
641  if(dim_size_temp <= 0)
642  _XMP_fatal("<nodes-size> should be less or equal to zero");
643 
644  dim_size[i] = dim_size_temp;
645  }
646 
648  *nodes = n;
649 
650  for (int i=0;i<dim;i++) {
651  int *rank_p = va_arg(args, int *);
652  *rank_p = n->info[i].rank;
653  }
654  va_end(args);
655 }
Here is the call graph for this function:

◆ _XMP_init_nodes_STATIC_GLOBAL()

void _XMP_init_nodes_STATIC_GLOBAL ( _XMP_nodes_t **  nodes,
int  dim,
  ... 
)
574 {
575  int dim_size[dim];
576 
577  va_list args;
578  va_start(args, dim);
579  for(int i=0;i<dim;i++){
580  int dim_size_temp = va_arg(args, int);
581  if(dim_size_temp <= 0){
582  _XMP_fatal("<nodes-size> should be more than zero");
583  }
584  dim_size[i] = dim_size_temp;
585  }
586 
588  *nodes = n;
589 
590  for(int i=0;i<dim;i++){
591  int *rank_p = va_arg(args, int *);
592  *rank_p = n->info[i].rank;
593  }
594  va_end(args);
595 }
Here is the call graph for this function:

◆ _XMP_init_nodes_STATIC_NODES_NAMED()

void _XMP_init_nodes_STATIC_NODES_NAMED ( _XMP_nodes_t **  nodes,
int  dim,
_XMP_nodes_t ref_nodes,
  ... 
)
738 {
739  int ref_dim = ref_nodes->dim;
740  int shrink[ref_dim];
741  int ref_lower[ref_dim], ref_upper[ref_dim], ref_stride[ref_dim];
742  int dim_size[dim];
743 
744  va_list args;
745  va_start(args, ref_nodes);
746  for(int i=0;i<ref_dim;i++){
747  shrink[i] = va_arg(args, int);
748  if(!shrink[i]){
749  ref_lower[i] = va_arg(args, int);
750  ref_upper[i] = va_arg(args, int);
751  ref_stride[i] = va_arg(args, int);
752  }
753  }
754 
755  for(int i=0;i<dim;i++){
756  int dim_size_temp = va_arg(args, int);
757  if(dim_size_temp <= 0){
758  _XMP_fatal("<nodes-size> should be less or equal to zero");
759  }
760  else{
761  dim_size[i] = dim_size_temp;
762  }
763  }
764 
765  _XMP_nodes_t *n = _XMP_init_nodes_struct_NODES_NAMED(dim, ref_nodes, shrink, ref_lower, ref_upper,
766  ref_stride, dim_size, _XMP_N_INT_TRUE);
767  *nodes = n;
768 
769  for(int i=0;i<dim;i++){
770  int *rank_p = va_arg(args, int *);
771  *rank_p = n->info[i].rank;
772  }
773  va_end(args);
774 }
Here is the call graph for this function:

◆ _XMP_init_nodes_STATIC_NODES_NUMBER()

void _XMP_init_nodes_STATIC_NODES_NUMBER ( _XMP_nodes_t **  nodes,
int  dim,
int  ref_lower,
int  ref_upper,
int  ref_stride,
  ... 
)
684 {
685  int dim_size[dim];
686 
687  va_list args;
688  va_start(args, ref_stride);
689  for(int i=0;i<dim;i++){
690  int dim_size_temp = va_arg(args, int);
691  if(dim_size_temp <= 0)
692  _XMP_fatal("<nodes-size> should be less or equal to zero");
693 
694  dim_size[i] = dim_size_temp;
695  }
696 
697  _XMP_nodes_t *n = _XMP_init_nodes_struct_NODES_NUMBER(dim, ref_lower, ref_upper,
698  ref_stride, dim_size, _XMP_N_INT_TRUE);
699  *nodes = n;
700 
701  for(int i=0;i<dim;i++){
702  int *rank_p = va_arg(args, int *);
703  *rank_p = n->info[i].rank;
704  }
705  va_end(args);
706 }
Here is the call graph for this function:

◆ _XMP_init_nodes_struct_EXEC()

_XMP_nodes_t* _XMP_init_nodes_struct_EXEC ( int  dim,
int *  dim_size,
int  is_static 
)
393 {
394  _XMP_nodes_t *exec_nodes = _XMP_get_execution_nodes();
395  int size = exec_nodes->comm_size;
396 
397  MPI_Comm *comm = _XMP_alloc(sizeof(MPI_Comm));
398  MPI_Comm_dup(*((MPI_Comm *)exec_nodes->comm), comm);
399 
400  _XMP_nodes_t *n = _XMP_create_new_nodes(_XMP_N_INT_TRUE, dim, size, (_XMP_comm_t *)comm);
401 
402  // calc inherit info
403  _XMP_nodes_t *inherit_nodes = _XMP_get_execution_nodes();
404  n->inherit_nodes = inherit_nodes;
405  n->inherit_info = _XMP_calc_inherit_info(inherit_nodes);
407 
408  // calc info
409  _XMP_init_nodes_info(n, dim_size, is_static);
410 
411  n->info[0].multiplier = 1;
412  for(int i=1;i<dim;i++)
413  n->info[i].multiplier = n->info[i-1].multiplier * dim_size[i-1];
414 
415  return n;
416 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _XMP_init_nodes_struct_GLOBAL()

_XMP_nodes_t* _XMP_init_nodes_struct_GLOBAL ( int  dim,
int *  dim_size,
int  is_static 
)
345 {
346  MPI_Comm *comm = _XMP_alloc(sizeof(MPI_Comm));
347  MPI_Comm_dup(MPI_COMM_WORLD, comm);
348 
349  _XMP_nodes_t *n = _XMP_create_new_nodes(_XMP_N_INT_TRUE, dim, _XMP_world_size, (_XMP_comm_t *)comm);
350 
351  // calc inherit info
352  n->inherit_nodes = NULL;
353  n->inherit_info = NULL;
354  n->attr = _XMP_ENTIRE_NODES;
355 
356  // set dim_size if XMP_NODE_SIZEn is set.
357  if(!is_static){
358  is_static = _XMP_N_INT_TRUE;
359  for(int i=0;i<dim;i++){
360  if(dim_size[i] == -1){
361  char name[20];
362  sprintf(name, "XMP_NODE_SIZE%d", i);
363  char *size = getenv(name);
364 
365  if(!size){
366  if(i == dim - 1){
367  is_static = _XMP_N_INT_FALSE;
368  break;
369  }
370  else _XMP_fatal("XMP_NODE_SIZE not specified although '*' is in the dimension of a node array\n");
371  }
372  else{
373  dim_size[i] = atoi(size);
374  if(dim_size[i] <= 0 || dim_size[i] > _XMP_world_size)
375  _XMP_fatal("Wrong value in XMP_NODE_SIZE\n");
376  }
377  }
378  }
379  }
380 
381  _XMP_init_nodes_info(n, dim_size, is_static);
382 
383  n->info[0].multiplier = 1;
384  for(int i=1;i<dim;i++)
385  n->info[i].multiplier = n->info[i-1].multiplier * dim_size[i-1];
386 
387  n->subcomm = create_subcomm(n);
388 
389  return n;
390 }
Here is the caller graph for this function:

◆ _XMP_init_nodes_struct_NODES_NAMED()

_XMP_nodes_t* _XMP_init_nodes_struct_NODES_NAMED ( int  dim,
_XMP_nodes_t ref_nodes,
int *  shrink,
int *  ref_lower,
int *  ref_upper,
int *  ref_stride,
int *  dim_size,
int  is_static 
)
500 {
501  int ref_dim = ref_nodes->dim;
502  int is_ref_member = ref_nodes->is_member;
503  int color = 1;
504  int is_member = _XMP_N_INT_TRUE;
505 
506  if(is_ref_member){
507  int acc_nodes_size = 1;
508  for(int i=0;i<ref_dim;i++){
509  int size = ref_nodes->info[i].size;
510  int rank = ref_nodes->info[i].rank;
511 
512  if(shrink[i]){
513  color += (acc_nodes_size * rank);
514  }
515  else{
516  _XMP_validate_nodes_ref(&ref_lower[i], &ref_upper[i], &ref_stride[i], size);
517  is_member = is_member && _XMP_check_nodes_ref_inclusion(ref_lower[i], ref_upper[i], ref_stride[i], size, rank);
518  }
519 
520  acc_nodes_size *= size;
521  }
522 
523  if(!is_member){
524  color = 0;
525  }
526  }
527  else{
528  is_member = _XMP_N_INT_FALSE;
529  color = 0;
530  }
531 
532  int comm_size = 1;
533  for(int i=0;i<ref_dim;i++)
534  if(!shrink[i])
535  comm_size *= _XMP_M_COUNT_TRIPLETi(ref_lower[i], ref_upper[i], ref_stride[i]);
536 
537  MPI_Comm *comm;
538  int use_subcomm;
539  if(check_subcomm(ref_nodes, ref_lower, ref_upper, ref_stride, shrink)){
540  use_subcomm = _XMP_N_INT_TRUE;
541  comm = (MPI_Comm *)get_subcomm(ref_nodes, ref_lower, ref_upper, ref_stride, shrink);
542  // MPI_Comm_dup(*(MPI_Comm *)get_subcomm(ref_nodes, ref_lower, ref_upper, ref_stride, shrink), comm);
543  }
544  else if(comm_size == 1){
545  use_subcomm = _XMP_N_INT_FALSE;
546  comm = _XMP_alloc(sizeof(MPI_Comm));
547  MPI_Comm_dup(MPI_COMM_SELF, comm);
548  }
549  else{
550  use_subcomm = _XMP_N_INT_FALSE;
551  comm = _XMP_alloc(sizeof(MPI_Comm));
552  MPI_Comm_split(*((MPI_Comm *)ref_nodes->comm), color, _XMP_world_rank, comm);
553  }
554 
555  _XMP_nodes_t *n = _XMP_create_new_nodes(is_member, dim, comm_size, (_XMP_comm_t *)comm);
556 
557  // calc inherit info
558  n->inherit_nodes = ref_nodes;
559  n->inherit_info = _XMP_calc_inherit_info_by_ref(ref_nodes, shrink, ref_lower, ref_upper, ref_stride);
561  n->use_subcomm = use_subcomm;
562 
563  // calc info
564  _XMP_init_nodes_info(n, dim_size, is_static);
565 
566  n->info[0].multiplier = 1;
567  for(int i=1;i<dim;i++)
568  n->info[i].multiplier = n->info[i-1].multiplier * dim_size[i-1];
569 
570  return n;
571 }
Here is the caller graph for this function:

◆ _XMP_init_nodes_struct_NODES_NUMBER()

_XMP_nodes_t* _XMP_init_nodes_struct_NODES_NUMBER ( int  dim,
int  ref_lower,
int  ref_upper,
int  ref_stride,
int *  dim_size,
int  is_static 
)
420 {
421  _XMP_validate_nodes_ref(&ref_lower, &ref_upper, &ref_stride, _XMP_world_size);
422  int is_member = _XMP_check_nodes_ref_inclusion(ref_lower, ref_upper, ref_stride, _XMP_world_size, _XMP_world_rank);
423 
424  MPI_Comm *comm = _XMP_alloc(sizeof(MPI_Comm));
425  MPI_Comm_split(*((MPI_Comm *)(_XMP_get_execution_nodes())->comm), is_member, _XMP_world_rank, comm);
426 
427  _XMP_nodes_t *n = _XMP_create_new_nodes(is_member, dim, _XMP_M_COUNT_TRIPLETi(ref_lower, ref_upper, ref_stride),
428  (_XMP_comm_t *)comm);
429 
430  // calc inherit info
431  int shrink[1] = {_XMP_N_INT_FALSE};
432  int l[1] = {ref_lower};
433  int u[1] = {ref_upper};
434  int s[1] = {ref_stride};
436  n->inherit_info = _XMP_calc_inherit_info_by_ref(_XMP_world_nodes, shrink, l, u, s);
437  n->attr = _XMP_ENTIRE_NODES;
438 
439  // calc info
440  _XMP_init_nodes_info(n, dim_size, is_static);
441 
442  n->info[0].multiplier = 1;
443  for(int i=1;i<dim;i++)
444  n->info[i].multiplier = n->info[i-1].multiplier * dim_size[i-1];
445 
446  n->subcomm = create_subcomm(n);
447 
448  return n;
449 }
Here is the caller graph for this function:

◆ _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 
)
1271 {
1272  int calc_flag = _XMP_N_INT_TRUE;
1273  int nodes_dim = nodes->dim;
1274 
1275  for(int i=0;i<nodes_dim;i++){
1276  if(rank_array[i] == _XMP_N_UNSPECIFIED_RANK){
1277  calc_flag = _XMP_N_INT_FALSE;
1278  int nodes_size = nodes->info[i].size;
1279  int new_shrink_nodes_size = shrink_nodes_size / nodes_size;
1280  for(int j=0;j<nodes_size;j++){
1281  rank_array[i] = j;
1282  _XMP_translate_nodes_rank_array_to_ranks(nodes, ranks + (j * new_shrink_nodes_size), rank_array, new_shrink_nodes_size);
1283  }
1284  }
1285  }
1286 
1287  if(calc_flag)
1288  *ranks = _XMP_calc_linear_rank(nodes, rank_array);
1289 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ xmp_dbg_printf()

void xmp_dbg_printf ( char *  fmt,
  ... 
)
38 {
39  char buf[512];
40  va_list args;
41 
42  va_start(args,fmt);
43  vsprintf(buf,fmt,args);
44  va_end(args);
45 
46  printf("[%d] %s",_XMP_world_rank, buf);
47  fflush(stdout);
48 }
_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_nodes_info_type::size
int size
Definition: xmp_data_struct.h:32
_XMP_nodes_type::info
_XMP_nodes_info_t info[1]
Definition: xmp_data_struct.h:60
_XMP_nodes_inherit_info_type::shrink
int shrink
Definition: xmp_data_struct.h:21
_XMP_nodes_type::comm_size
int comm_size
Definition: xmp_data_struct.h:48
_XMP_init_nodes_struct_NODES_NAMED
_XMP_nodes_t * _XMP_init_nodes_struct_NODES_NAMED(int dim, _XMP_nodes_t *ref_nodes, int *shrink, int *ref_lower, int *ref_upper, int *ref_stride, int *dim_size, int is_static)
Definition: xmp_nodes.c:498
_XMP_nodes_type::use_subcomm
int use_subcomm
Definition: xmp_data_struct.h:55
_XMP_task_desc_type
Definition: xmp_data_struct.h:316
_XMP_nodes_type::subcomm
_XMP_comm_t * subcomm
Definition: xmp_data_struct.h:54
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_init_nodes_struct_EXEC
_XMP_nodes_t * _XMP_init_nodes_struct_EXEC(int dim, int *dim_size, int is_static)
Definition: xmp_nodes.c:392
_XMP_nodes_type::inherit_nodes
struct _XMP_nodes_type * inherit_nodes
Definition: xmp_data_struct.h:57
_XMP_EQUIVALENCE_NODES
#define _XMP_EQUIVALENCE_NODES
Definition: xmp_constant.h:124
_XMP_init_nodes_ref
_XMP_nodes_ref_t * _XMP_init_nodes_ref(_XMP_nodes_t *n, int *rank_array)
Definition: xmp_nodes.c:1207
_XMP_world_size
int _XMP_world_size
Definition: xmp_world.c:8
_XMP_world_nodes
void * _XMP_world_nodes
Definition: xmp_world.c:10
_XMP_calc_coord_on_target_nodes2
_Bool _XMP_calc_coord_on_target_nodes2(_XMP_nodes_t *n, int *ncoord, _XMP_nodes_t *target_n, int *target_ncoord)
Definition: xmp_nodes.c:1090
_XMP_nodes_ref_type::shrink_nodes_size
int shrink_nodes_size
Definition: xmp_data_struct.h:66
_XMP_calc_linear_rank
int _XMP_calc_linear_rank(_XMP_nodes_t *n, int *rank_array)
Definition: xmp_nodes.c:1035
_XMP_finalize_comm
void _XMP_finalize_comm(void *comm)
_XMP_world_rank
int _XMP_world_rank
Definition: xmp_world.c:9
_XMP_nodes_dealloc_after_wait_async
void _XMP_nodes_dealloc_after_wait_async(_XMP_nodes_t *n)
Definition: xmp_async.c:272
_XMP_nodes_inherit_info_type
Definition: xmp_data_struct.h:20
_XMP_nodes_ref_type::ref
int * ref
Definition: xmp_data_struct.h:65
_XMP_ENTIRE_NODES
#define _XMP_ENTIRE_NODES
Definition: xmp_constant.h:121
_XMP_N_INT_FALSE
#define _XMP_N_INT_FALSE
Definition: xmp_constant.h:5
_XMP_N_UNSPECIFIED_RANK
#define _XMP_N_UNSPECIFIED_RANK
Definition: xmp_constant.h:22
_XMP_nodes_inherit_info_type::stride
int stride
Definition: xmp_data_struct.h:25
_XMP_init_nodes_struct_NODES_NUMBER
_XMP_nodes_t * _XMP_init_nodes_struct_NODES_NUMBER(int dim, int ref_lower, int ref_upper, int ref_stride, int *dim_size, int is_static)
Definition: xmp_nodes.c:418
_XMP_nodes_ref_type
Definition: xmp_data_struct.h:63
_XMP_nodes_type::inherit_info
_XMP_nodes_inherit_info_t * inherit_info
Definition: xmp_data_struct.h:59
_XMP_comm_t
#define _XMP_comm_t
Definition: xmp_data_struct.h:17
_XMP_task_desc_type::nodes
_XMP_nodes_t * nodes
Definition: xmp_data_struct.h:317
_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_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_nodes_info_type::rank
int rank
Definition: xmp_data_struct.h:35
_XMP_init_nodes_struct_GLOBAL
_XMP_nodes_t * _XMP_init_nodes_struct_GLOBAL(int dim, int *dim_size, int is_static)
Definition: xmp_nodes.c:344
_XMP_N_INVALID_RANK
#define _XMP_N_INVALID_RANK
Definition: xmp_constant.h:21
_XMP_EXECUTING_NODES
#define _XMP_EXECUTING_NODES
Definition: xmp_constant.h:122
_XMP_free
void _XMP_free(void *p)
Definition: xmp_util.c:37
_XMP_nodes_inherit_info_type::lower
int lower
Definition: xmp_data_struct.h:23
_XMP_nodes_info_type::multiplier
int multiplier
Definition: xmp_data_struct.h:37
_XMP_calc_coord_on_target_nodes
_Bool _XMP_calc_coord_on_target_nodes(_XMP_nodes_t *n, int *ncoord, _XMP_nodes_t *target_n, int *target_ncoord)
Definition: xmp_nodes.c:1150
_XMP_nodes_type::comm
_XMP_comm_t * comm
Definition: xmp_data_struct.h:53
_XMP_fatal
void _XMP_fatal(char *msg)
Definition: xmp_util.c:42
_XMP_calc_rank_array
void _XMP_calc_rank_array(_XMP_nodes_t *n, int *rank_array, int linear_rank)
Definition: xmp_nodes.c:1026
_XMP_nodes_ref_type::nodes
_XMP_nodes_t * nodes
Definition: xmp_data_struct.h:64
_XMP_push_nodes
void _XMP_push_nodes(void *nodes)
_XMP_nodes_type::dim
int dim
Definition: xmp_data_struct.h:47
_XMP_nodes_type::attr
int attr
Definition: xmp_data_struct.h:49
_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_finalize_nodes
void _XMP_finalize_nodes(_XMP_nodes_t *nodes)
Definition: xmp_nodes.c:816
_XMP_N_INT_TRUE
#define _XMP_N_INT_TRUE
Definition: xmp_constant.h:4
_XMP_init_nodes_STATIC_NODES_NUMBER
void _XMP_init_nodes_STATIC_NODES_NUMBER(_XMP_nodes_t **nodes, int dim, int ref_lower, int ref_upper, int ref_stride,...)
Definition: xmp_nodes.c:682
_XMP_task_desc_type::execute
int execute
Definition: xmp_data_struct.h:318
_XMP_get_execution_nodes
void * _XMP_get_execution_nodes(void)
Definition: xmp_nodes_stack.c:46
_XMP_nodes_inherit_info_type::upper
int upper
Definition: xmp_data_struct.h:24