2 /* Copyright (c) 2007-2014. The SimGrid Team.
3 * All rights reserved. */
5 /* This program is free software; you can redistribute it and/or modify it
6 * under the terms of the license (GNU LGPL) which comes with this package. */
9 #include "smpi_mpi_dt_private.h"
11 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi,
12 "Logging specific to SMPI (pmpi)");
15 //this function need to be here because of the calls to smpi_bench
16 void TRACE_smpi_set_category(const char *category)
18 //need to end bench otherwise categories for execution tasks are wrong
20 TRACE_internal_smpi_set_category (category);
21 //begin bench after changing process's category
26 /* PMPI User level calls */
28 int PMPI_Init(int *argc, char ***argv)
30 smpi_process_init(argc, argv);
31 smpi_process_mark_as_initialized();
33 int rank = smpi_process_index();
34 TRACE_smpi_init(rank);
35 TRACE_smpi_computing_init(rank);
36 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
37 extra->type = TRACING_INIT;
38 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
39 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
45 int PMPI_Finalize(void)
49 int rank = smpi_process_index();
50 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
51 extra->type = TRACING_FINALIZE;
52 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
54 smpi_process_finalize();
56 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
57 TRACE_smpi_finalize(smpi_process_index());
59 smpi_process_destroy();
63 int PMPI_Finalized(int* flag)
65 *flag=smpi_process_finalized();
69 int PMPI_Get_version (int *version,int *subversion){
70 *version = MPI_VERSION;
71 *subversion= MPI_SUBVERSION;
75 int PMPI_Get_library_version (char *version,int *len){
76 int retval = MPI_SUCCESS;
78 snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2014",SIMGRID_VERSION_MAJOR,
79 SIMGRID_VERSION_MINOR);
80 *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
85 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
87 if (provided != NULL) {
88 *provided = MPI_THREAD_MULTIPLE;
90 return MPI_Init(argc, argv);
93 int PMPI_Query_thread(int *provided)
97 if (provided == NULL) {
100 *provided = MPI_THREAD_MULTIPLE;
101 retval = MPI_SUCCESS;
106 int PMPI_Is_thread_main(int *flag)
111 retval = MPI_ERR_ARG;
113 *flag = smpi_process_index() == 0;
114 retval = MPI_SUCCESS;
119 int PMPI_Abort(MPI_Comm comm, int errorcode)
122 smpi_process_destroy();
123 // FIXME: should kill all processes in comm instead
124 simcall_process_kill(SIMIX_process_self());
128 double PMPI_Wtime(void)
131 if (smpi_process_initialized() && !smpi_process_finalized() && !smpi_process_get_sampling()) {
133 time = SIMIX_get_clock();
136 time = SIMIX_get_clock();
141 extern double sg_maxmin_precision;
142 double PMPI_Wtick(void)
144 return sg_maxmin_precision;
147 int PMPI_Address(void *location, MPI_Aint * address)
152 retval = MPI_ERR_ARG;
154 *address = (MPI_Aint) location;
155 retval = MPI_SUCCESS;
160 int PMPI_Get_address(void *location, MPI_Aint * address)
162 return PMPI_Address(location, address);
165 int PMPI_Type_free(MPI_Datatype * datatype)
168 /* Free a predefined datatype is an error according to the standard, and
169 should be checked for */
170 if (*datatype == MPI_DATATYPE_NULL) {
171 retval = MPI_ERR_ARG;
173 smpi_datatype_free(datatype);
174 retval = MPI_SUCCESS;
179 int PMPI_Type_size(MPI_Datatype datatype, int *size)
183 if (datatype == MPI_DATATYPE_NULL) {
184 retval = MPI_ERR_TYPE;
185 } else if (size == NULL) {
186 retval = MPI_ERR_ARG;
188 *size = (int) smpi_datatype_size(datatype);
189 retval = MPI_SUCCESS;
194 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
198 if (datatype == MPI_DATATYPE_NULL) {
199 retval = MPI_ERR_TYPE;
200 } else if (lb == NULL || extent == NULL) {
201 retval = MPI_ERR_ARG;
203 retval = smpi_datatype_extent(datatype, lb, extent);
208 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
210 return PMPI_Type_get_extent(datatype, lb, extent);
213 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
217 if (datatype == MPI_DATATYPE_NULL) {
218 retval = MPI_ERR_TYPE;
219 } else if (extent == NULL) {
220 retval = MPI_ERR_ARG;
222 *extent = smpi_datatype_get_extent(datatype);
223 retval = MPI_SUCCESS;
228 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
232 if (datatype == MPI_DATATYPE_NULL) {
233 retval = MPI_ERR_TYPE;
234 } else if (disp == NULL) {
235 retval = MPI_ERR_ARG;
237 *disp = smpi_datatype_lb(datatype);
238 retval = MPI_SUCCESS;
243 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
247 if (datatype == MPI_DATATYPE_NULL) {
248 retval = MPI_ERR_TYPE;
249 } else if (disp == NULL) {
250 retval = MPI_ERR_ARG;
252 *disp = smpi_datatype_ub(datatype);
253 retval = MPI_SUCCESS;
258 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
262 if (function == NULL || op == NULL) {
263 retval = MPI_ERR_ARG;
265 *op = smpi_op_new(function, commute);
266 retval = MPI_SUCCESS;
271 int PMPI_Op_free(MPI_Op * op)
276 retval = MPI_ERR_ARG;
277 } else if (*op == MPI_OP_NULL) {
280 smpi_op_destroy(*op);
282 retval = MPI_SUCCESS;
287 int PMPI_Group_free(MPI_Group * group)
292 retval = MPI_ERR_ARG;
294 smpi_group_destroy(*group);
295 *group = MPI_GROUP_NULL;
296 retval = MPI_SUCCESS;
301 int PMPI_Group_size(MPI_Group group, int *size)
305 if (group == MPI_GROUP_NULL) {
306 retval = MPI_ERR_GROUP;
307 } else if (size == NULL) {
308 retval = MPI_ERR_ARG;
310 *size = smpi_group_size(group);
311 retval = MPI_SUCCESS;
316 int PMPI_Group_rank(MPI_Group group, int *rank)
320 if (group == MPI_GROUP_NULL) {
321 retval = MPI_ERR_GROUP;
322 } else if (rank == NULL) {
323 retval = MPI_ERR_ARG;
325 *rank = smpi_group_rank(group, smpi_process_index());
326 retval = MPI_SUCCESS;
331 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
332 MPI_Group group2, int *ranks2)
334 int retval, i, index;
335 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
336 retval = MPI_ERR_GROUP;
338 for (i = 0; i < n; i++) {
339 if(ranks1[i]==MPI_PROC_NULL){
340 ranks2[i]=MPI_PROC_NULL;
342 index = smpi_group_index(group1, ranks1[i]);
343 ranks2[i] = smpi_group_rank(group2, index);
346 retval = MPI_SUCCESS;
351 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
355 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
356 retval = MPI_ERR_GROUP;
357 } else if (result == NULL) {
358 retval = MPI_ERR_ARG;
360 *result = smpi_group_compare(group1, group2);
361 retval = MPI_SUCCESS;
366 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
367 MPI_Group * newgroup)
369 int retval, i, proc1, proc2, size, size2;
371 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
372 retval = MPI_ERR_GROUP;
373 } else if (newgroup == NULL) {
374 retval = MPI_ERR_ARG;
376 size = smpi_group_size(group1);
377 size2 = smpi_group_size(group2);
378 for (i = 0; i < size2; i++) {
379 proc2 = smpi_group_index(group2, i);
380 proc1 = smpi_group_rank(group1, proc2);
381 if (proc1 == MPI_UNDEFINED) {
386 *newgroup = MPI_GROUP_EMPTY;
388 *newgroup = smpi_group_new(size);
389 size2 = smpi_group_size(group1);
390 for (i = 0; i < size2; i++) {
391 proc1 = smpi_group_index(group1, i);
392 smpi_group_set_mapping(*newgroup, proc1, i);
394 for (i = size2; i < size; i++) {
395 proc2 = smpi_group_index(group2, i - size2);
396 smpi_group_set_mapping(*newgroup, proc2, i);
399 retval = MPI_SUCCESS;
404 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
405 MPI_Group * newgroup)
407 int retval, i, proc1, proc2, size;
409 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
410 retval = MPI_ERR_GROUP;
411 } else if (newgroup == NULL) {
412 retval = MPI_ERR_ARG;
414 size = smpi_group_size(group2);
415 for (i = 0; i < size; i++) {
416 proc2 = smpi_group_index(group2, i);
417 proc1 = smpi_group_rank(group1, proc2);
418 if (proc1 == MPI_UNDEFINED) {
423 *newgroup = MPI_GROUP_EMPTY;
425 *newgroup = smpi_group_new(size);
427 for (i = 0; i < smpi_group_size(group2); i++) {
428 proc2 = smpi_group_index(group2, i);
429 proc1 = smpi_group_rank(group1, proc2);
430 if (proc1 != MPI_UNDEFINED) {
431 smpi_group_set_mapping(*newgroup, proc2, j);
436 retval = MPI_SUCCESS;
441 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
443 int retval, i, proc1, proc2, size, size2;
445 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
446 retval = MPI_ERR_GROUP;
447 } else if (newgroup == NULL) {
448 retval = MPI_ERR_ARG;
450 size = size2 = smpi_group_size(group1);
451 for (i = 0; i < size2; i++) {
452 proc1 = smpi_group_index(group1, i);
453 proc2 = smpi_group_rank(group2, proc1);
454 if (proc2 != MPI_UNDEFINED) {
459 *newgroup = MPI_GROUP_EMPTY;
461 *newgroup = smpi_group_new(size);
462 for (i = 0; i < size2; i++) {
463 proc1 = smpi_group_index(group1, i);
464 proc2 = smpi_group_rank(group2, proc1);
465 if (proc2 == MPI_UNDEFINED) {
466 smpi_group_set_mapping(*newgroup, proc1, i);
470 retval = MPI_SUCCESS;
475 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
477 int retval, i, index;
479 if (group == MPI_GROUP_NULL) {
480 retval = MPI_ERR_GROUP;
481 } else if (newgroup == NULL) {
482 retval = MPI_ERR_ARG;
485 *newgroup = MPI_GROUP_EMPTY;
486 } else if (n == smpi_group_size(group)) {
488 if(group!= smpi_comm_group(MPI_COMM_WORLD)
489 && group != MPI_GROUP_NULL
490 && group != smpi_comm_group(MPI_COMM_SELF)
491 && group != MPI_GROUP_EMPTY)
492 smpi_group_use(group);
494 *newgroup = smpi_group_new(n);
495 for (i = 0; i < n; i++) {
496 index = smpi_group_index(group, ranks[i]);
497 smpi_group_set_mapping(*newgroup, index, i);
500 retval = MPI_SUCCESS;
505 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
507 int retval, i, j, newsize, oldsize, index;
509 if (group == MPI_GROUP_NULL) {
510 retval = MPI_ERR_GROUP;
511 } else if (newgroup == NULL) {
512 retval = MPI_ERR_ARG;
516 if(group!= smpi_comm_group(MPI_COMM_WORLD)
517 && group != MPI_GROUP_NULL
518 && group != smpi_comm_group(MPI_COMM_SELF)
519 && group != MPI_GROUP_EMPTY)
520 smpi_group_use(group);
521 } else if (n == smpi_group_size(group)) {
522 *newgroup = MPI_GROUP_EMPTY;
524 oldsize=smpi_group_size(group);
525 newsize = oldsize - n;
526 *newgroup = smpi_group_new(newsize);
528 int* to_exclude=xbt_new0(int, smpi_group_size(group));
529 for(i=0; i<oldsize; i++)
532 to_exclude[ranks[i]]=1;
535 for(i=0; i<oldsize; i++){
536 if(to_exclude[i]==0){
537 index = smpi_group_index(group, i);
538 smpi_group_set_mapping(*newgroup, index, j);
543 xbt_free(to_exclude);
545 retval = MPI_SUCCESS;
550 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
551 MPI_Group * newgroup)
553 int retval, i, j, rank, size, index;
555 if (group == MPI_GROUP_NULL) {
556 retval = MPI_ERR_GROUP;
557 } else if (newgroup == NULL) {
558 retval = MPI_ERR_ARG;
561 *newgroup = MPI_GROUP_EMPTY;
564 for (i = 0; i < n; i++) {
565 for (rank = ranges[i][0]; /* First */
566 rank >= 0; /* Last */
570 rank += ranges[i][2]; /* Stride */
571 if (ranges[i][0]<ranges[i][1]){
572 if(rank > ranges[i][1])
575 if(rank < ranges[i][1])
581 *newgroup = smpi_group_new(size);
583 for (i = 0; i < n; i++) {
584 for (rank = ranges[i][0]; /* First */
585 rank >= 0; /* Last */
587 index = smpi_group_index(group, rank);
588 smpi_group_set_mapping(*newgroup, index, j);
590 rank += ranges[i][2]; /* Stride */
591 if (ranges[i][0]<ranges[i][1]){
592 if(rank > ranges[i][1])
595 if(rank < ranges[i][1])
601 retval = MPI_SUCCESS;
606 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
607 MPI_Group * newgroup)
609 int retval, i, rank, newrank,oldrank, size, index, add;
611 if (group == MPI_GROUP_NULL) {
612 retval = MPI_ERR_GROUP;
613 } else if (newgroup == NULL) {
614 retval = MPI_ERR_ARG;
618 if(group!= smpi_comm_group(MPI_COMM_WORLD)
619 && group != MPI_GROUP_NULL
620 && group != smpi_comm_group(MPI_COMM_SELF)
621 && group != MPI_GROUP_EMPTY)
622 smpi_group_use(group);
624 size = smpi_group_size(group);
625 for (i = 0; i < n; i++) {
626 for (rank = ranges[i][0]; /* First */
627 rank >= 0; /* Last */
631 rank += ranges[i][2]; /* Stride */
632 if (ranges[i][0]<ranges[i][1]){
633 if(rank > ranges[i][1])
636 if(rank < ranges[i][1])
642 *newgroup = MPI_GROUP_EMPTY;
644 *newgroup = smpi_group_new(size);
647 while (newrank < size) {
649 for (i = 0; i < n; i++) {
650 for (rank = ranges[i][0];rank >= 0;){
656 rank += ranges[i][2]; /* Stride */
658 if (ranges[i][0]<ranges[i][1]){
659 if(rank > ranges[i][1])
662 if(rank < ranges[i][1])
668 index = smpi_group_index(group, oldrank);
669 smpi_group_set_mapping(*newgroup, index, newrank);
677 retval = MPI_SUCCESS;
682 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
685 if (comm == MPI_COMM_NULL) {
686 retval = MPI_ERR_COMM;
687 } else if (rank == NULL) {
688 retval = MPI_ERR_ARG;
690 *rank = smpi_comm_rank(comm);
691 retval = MPI_SUCCESS;
696 int PMPI_Comm_size(MPI_Comm comm, int *size)
699 if (comm == MPI_COMM_NULL) {
700 retval = MPI_ERR_COMM;
701 } else if (size == NULL) {
702 retval = MPI_ERR_ARG;
704 *size = smpi_comm_size(comm);
705 retval = MPI_SUCCESS;
710 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
714 if (comm == MPI_COMM_NULL) {
715 retval = MPI_ERR_COMM;
716 } else if (name == NULL || len == NULL) {
717 retval = MPI_ERR_ARG;
719 smpi_comm_get_name(comm, name, len);
720 retval = MPI_SUCCESS;
725 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
729 if (comm == MPI_COMM_NULL) {
730 retval = MPI_ERR_COMM;
731 } else if (group == NULL) {
732 retval = MPI_ERR_ARG;
734 *group = smpi_comm_group(comm);
735 if(*group!= smpi_comm_group(MPI_COMM_WORLD)
736 && *group != MPI_GROUP_NULL
737 && *group != smpi_comm_group(MPI_COMM_SELF)
738 && *group != MPI_GROUP_EMPTY)
739 smpi_group_use(*group);
740 retval = MPI_SUCCESS;
745 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
749 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
750 retval = MPI_ERR_COMM;
751 } else if (result == NULL) {
752 retval = MPI_ERR_ARG;
754 if (comm1 == comm2) { /* Same communicators means same groups */
758 smpi_group_compare(smpi_comm_group(comm1),
759 smpi_comm_group(comm2));
760 if (*result == MPI_IDENT) {
761 *result = MPI_CONGRUENT;
764 retval = MPI_SUCCESS;
769 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
773 if (comm == MPI_COMM_NULL) {
774 retval = MPI_ERR_COMM;
775 } else if (newcomm == NULL) {
776 retval = MPI_ERR_ARG;
778 *newcomm = smpi_comm_new(smpi_comm_group(comm), smpi_comm_topo(comm));
779 retval = MPI_SUCCESS;
784 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
788 if (comm == MPI_COMM_NULL) {
789 retval = MPI_ERR_COMM;
790 } else if (group == MPI_GROUP_NULL) {
791 retval = MPI_ERR_GROUP;
792 } else if (newcomm == NULL) {
793 retval = MPI_ERR_ARG;
794 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
795 *newcomm= MPI_COMM_NULL;
796 retval = MPI_SUCCESS;
799 *newcomm = smpi_comm_new(group, NULL);
800 retval = MPI_SUCCESS;
805 int PMPI_Comm_free(MPI_Comm * comm)
810 retval = MPI_ERR_ARG;
811 } else if (*comm == MPI_COMM_NULL) {
812 retval = MPI_ERR_COMM;
814 smpi_comm_destroy(*comm);
815 *comm = MPI_COMM_NULL;
816 retval = MPI_SUCCESS;
821 int PMPI_Comm_disconnect(MPI_Comm * comm)
823 /* TODO: wait until all communication in comm are done */
827 retval = MPI_ERR_ARG;
828 } else if (*comm == MPI_COMM_NULL) {
829 retval = MPI_ERR_COMM;
831 smpi_comm_destroy(*comm);
832 *comm = MPI_COMM_NULL;
833 retval = MPI_SUCCESS;
838 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
843 if (comm_out == NULL) {
844 retval = MPI_ERR_ARG;
845 } else if (comm == MPI_COMM_NULL) {
846 retval = MPI_ERR_COMM;
848 *comm_out = smpi_comm_split(comm, color, key);
849 retval = MPI_SUCCESS;
856 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
857 int tag, MPI_Comm comm, MPI_Request * request)
862 if (request == NULL) {
863 retval = MPI_ERR_ARG;
864 } else if (comm == MPI_COMM_NULL) {
865 retval = MPI_ERR_COMM;
866 } else if (!is_datatype_valid(datatype)) {
867 retval = MPI_ERR_TYPE;
868 } else if (dst == MPI_PROC_NULL) {
869 retval = MPI_SUCCESS;
871 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
872 retval = MPI_SUCCESS;
875 if (retval != MPI_SUCCESS && request)
876 *request = MPI_REQUEST_NULL;
880 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
881 int tag, MPI_Comm comm, MPI_Request * request)
886 if (request == NULL) {
887 retval = MPI_ERR_ARG;
888 } else if (comm == MPI_COMM_NULL) {
889 retval = MPI_ERR_COMM;
890 } else if (!is_datatype_valid(datatype)) {
891 retval = MPI_ERR_TYPE;
892 } else if (src == MPI_PROC_NULL) {
893 retval = MPI_SUCCESS;
895 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
896 retval = MPI_SUCCESS;
899 if (retval != MPI_SUCCESS && request)
900 *request = MPI_REQUEST_NULL;
904 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype,
905 int dst, int tag, MPI_Comm comm, MPI_Request* request)
910 if (request == NULL) {
911 retval = MPI_ERR_ARG;
912 } else if (comm == MPI_COMM_NULL) {
913 retval = MPI_ERR_COMM;
914 } else if (!is_datatype_valid(datatype)) {
915 retval = MPI_ERR_TYPE;
916 } else if (dst == MPI_PROC_NULL) {
917 retval = MPI_SUCCESS;
919 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
920 retval = MPI_SUCCESS;
923 if (retval != MPI_SUCCESS && request)
924 *request = MPI_REQUEST_NULL;
928 int PMPI_Start(MPI_Request * request)
933 if (request == NULL || *request == MPI_REQUEST_NULL) {
934 retval = MPI_ERR_REQUEST;
936 smpi_mpi_start(*request);
937 retval = MPI_SUCCESS;
943 int PMPI_Startall(int count, MPI_Request * requests)
948 if (requests == NULL) {
949 retval = MPI_ERR_ARG;
951 retval = MPI_SUCCESS;
952 for (i = 0 ; i < count ; i++) {
953 if(requests[i] == MPI_REQUEST_NULL) {
954 retval = MPI_ERR_REQUEST;
957 if(retval != MPI_ERR_REQUEST) {
958 smpi_mpi_startall(count, requests);
965 int PMPI_Request_free(MPI_Request * request)
970 if (*request == MPI_REQUEST_NULL) {
971 retval = MPI_ERR_ARG;
973 smpi_mpi_request_free(request);
974 retval = MPI_SUCCESS;
980 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
981 int tag, MPI_Comm comm, MPI_Request * request)
987 if (request == NULL) {
988 retval = MPI_ERR_ARG;
989 } else if (comm == MPI_COMM_NULL) {
990 retval = MPI_ERR_COMM;
991 } else if (src == MPI_PROC_NULL) {
992 *request = MPI_REQUEST_NULL;
993 retval = MPI_SUCCESS;
994 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
995 retval = MPI_ERR_COMM;
996 } else if (count < 0) {
997 retval = MPI_ERR_COUNT;
998 } else if (buf==NULL && count > 0) {
999 retval = MPI_ERR_COUNT;
1000 } else if (!is_datatype_valid(datatype)) {
1001 retval = MPI_ERR_TYPE;
1002 } else if(tag<0 && tag != MPI_ANY_TAG){
1003 retval = MPI_ERR_TAG;
1007 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1008 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1010 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1011 extra->type = TRACING_IRECV;
1012 extra->send_size = count;
1013 extra->src = src_traced;
1015 extra->datatype1 = encode_datatype(datatype);
1016 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1019 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1020 retval = MPI_SUCCESS;
1023 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1024 (*request)->recv = 1;
1029 if (retval != MPI_SUCCESS && request)
1030 *request = MPI_REQUEST_NULL;
1035 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1036 int tag, MPI_Comm comm, MPI_Request * request)
1041 if (request == NULL) {
1042 retval = MPI_ERR_ARG;
1043 } else if (comm == MPI_COMM_NULL) {
1044 retval = MPI_ERR_COMM;
1045 } else if (dst == MPI_PROC_NULL) {
1046 *request = MPI_REQUEST_NULL;
1047 retval = MPI_SUCCESS;
1048 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1049 retval = MPI_ERR_RANK;
1050 } else if (count < 0) {
1051 retval = MPI_ERR_COUNT;
1052 } else if (buf==NULL && count > 0) {
1053 retval = MPI_ERR_COUNT;
1054 } else if (!is_datatype_valid(datatype)) {
1055 retval = MPI_ERR_TYPE;
1056 } else if(tag<0 && tag != MPI_ANY_TAG){
1057 retval = MPI_ERR_TAG;
1061 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1062 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1064 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1065 extra->type = TRACING_ISEND;
1066 extra->send_size = count;
1068 extra->dst = dst_traced;
1069 extra->datatype1 = encode_datatype(datatype);
1070 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1071 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1074 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1075 retval = MPI_SUCCESS;
1078 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1079 (*request)->send = 1;
1084 if (retval != MPI_SUCCESS && request)
1085 *request = MPI_REQUEST_NULL;
1089 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype,
1090 int dst, int tag, MPI_Comm comm, MPI_Request* request)
1095 if (request == NULL) {
1096 retval = MPI_ERR_ARG;
1097 } else if (comm == MPI_COMM_NULL) {
1098 retval = MPI_ERR_COMM;
1099 } else if (dst == MPI_PROC_NULL) {
1100 *request = MPI_REQUEST_NULL;
1101 retval = MPI_SUCCESS;
1102 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1103 retval = MPI_ERR_RANK;
1104 } else if (count < 0) {
1105 retval = MPI_ERR_COUNT;
1106 } else if (buf==NULL && count > 0) {
1107 retval = MPI_ERR_COUNT;
1108 } else if (!is_datatype_valid(datatype)) {
1109 retval = MPI_ERR_TYPE;
1110 } else if(tag<0 && tag != MPI_ANY_TAG){
1111 retval = MPI_ERR_TAG;
1115 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1116 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1117 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1118 extra->type = TRACING_ISSEND;
1119 extra->send_size = count;
1121 extra->dst = dst_traced;
1122 extra->datatype1 = encode_datatype(datatype);
1123 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1124 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1127 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1128 retval = MPI_SUCCESS;
1131 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1132 (*request)->send = 1;
1137 if (retval != MPI_SUCCESS && request)
1138 *request = MPI_REQUEST_NULL;
1142 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1143 MPI_Comm comm, MPI_Status * status)
1148 if (comm == MPI_COMM_NULL) {
1149 retval = MPI_ERR_COMM;
1150 } else if (src == MPI_PROC_NULL) {
1151 smpi_empty_status(status);
1152 status->MPI_SOURCE = MPI_PROC_NULL;
1153 retval = MPI_SUCCESS;
1154 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1155 retval = MPI_ERR_RANK;
1156 } else if (count < 0) {
1157 retval = MPI_ERR_COUNT;
1158 } else if (buf==NULL && count > 0) {
1159 retval = MPI_ERR_COUNT;
1160 } else if (!is_datatype_valid(datatype)) {
1161 retval = MPI_ERR_TYPE;
1162 } else if(tag<0 && tag != MPI_ANY_TAG){
1163 retval = MPI_ERR_TAG;
1166 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1167 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1168 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1169 extra->type = TRACING_RECV;
1170 extra->send_size = count;
1171 extra->src = src_traced;
1173 extra->datatype1 = encode_datatype(datatype);
1174 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1177 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1178 retval = MPI_SUCCESS;
1181 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1182 if(status!=MPI_STATUS_IGNORE){
1183 src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1184 TRACE_smpi_recv(rank, src_traced, rank);
1186 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1194 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1201 if (comm == MPI_COMM_NULL) {
1202 retval = MPI_ERR_COMM;
1203 } else if (dst == MPI_PROC_NULL) {
1204 retval = MPI_SUCCESS;
1205 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1206 retval = MPI_ERR_RANK;
1207 } else if (count < 0) {
1208 retval = MPI_ERR_COUNT;
1209 } else if (buf==NULL && count > 0) {
1210 retval = MPI_ERR_COUNT;
1211 } else if (!is_datatype_valid(datatype)) {
1212 retval = MPI_ERR_TYPE;
1213 } else if(tag<0 && tag != MPI_ANY_TAG){
1214 retval = MPI_ERR_TAG;
1218 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1219 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1220 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1221 extra->type = TRACING_SEND;
1222 extra->send_size = count;
1224 extra->dst = dst_traced;
1225 extra->datatype1 = encode_datatype(datatype);
1226 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1227 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1230 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1231 retval = MPI_SUCCESS;
1234 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1244 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1249 if (comm == MPI_COMM_NULL) {
1250 retval = MPI_ERR_COMM;
1251 } else if (dst == MPI_PROC_NULL) {
1252 retval = MPI_SUCCESS;
1253 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1254 retval = MPI_ERR_RANK;
1255 } else if (count < 0) {
1256 retval = MPI_ERR_COUNT;
1257 } else if (buf==NULL && count > 0) {
1258 retval = MPI_ERR_COUNT;
1259 } else if (!is_datatype_valid(datatype)){
1260 retval = MPI_ERR_TYPE;
1261 } else if(tag<0 && tag != MPI_ANY_TAG){
1262 retval = MPI_ERR_TAG;
1266 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1267 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1268 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1269 extra->type = TRACING_SSEND;
1270 extra->send_size = count;
1272 extra->dst = dst_traced;
1273 extra->datatype1 = encode_datatype(datatype);
1274 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra); TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1277 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1278 retval = MPI_SUCCESS;
1281 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1289 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1290 int dst, int sendtag, void *recvbuf, int recvcount,
1291 MPI_Datatype recvtype, int src, int recvtag,
1292 MPI_Comm comm, MPI_Status * status)
1298 if (comm == MPI_COMM_NULL) {
1299 retval = MPI_ERR_COMM;
1300 } else if (!is_datatype_valid(sendtype)
1301 || !is_datatype_valid(recvtype)) {
1302 retval = MPI_ERR_TYPE;
1303 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1304 smpi_empty_status(status);
1305 status->MPI_SOURCE = MPI_PROC_NULL;
1306 retval = MPI_SUCCESS;
1307 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1308 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1309 retval = MPI_ERR_RANK;
1310 } else if (sendcount < 0 || recvcount<0) {
1311 retval = MPI_ERR_COUNT;
1312 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1313 retval = MPI_ERR_COUNT;
1314 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1315 retval = MPI_ERR_TAG;
1319 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1320 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1321 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1322 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1323 extra->type = TRACING_SENDRECV;
1324 extra->send_size = sendcount;
1325 extra->recv_size = recvcount;
1326 extra->src = src_traced;
1327 extra->dst = dst_traced;
1328 extra->datatype1 = encode_datatype(sendtype);
1329 extra->datatype2 = encode_datatype(recvtype);
1331 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1332 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1336 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1337 recvcount, recvtype, src, recvtag, comm, status);
1338 retval = MPI_SUCCESS;
1341 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1342 TRACE_smpi_recv(rank, src_traced, rank);
1351 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1352 int dst, int sendtag, int src, int recvtag,
1353 MPI_Comm comm, MPI_Status * status)
1355 //TODO: suboptimal implementation
1358 if (!is_datatype_valid(datatype)) {
1359 retval = MPI_ERR_TYPE;
1360 } else if (count < 0) {
1361 retval = MPI_ERR_COUNT;
1363 int size = smpi_datatype_get_extent(datatype) * count;
1364 recvbuf = xbt_new0(char, size);
1366 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1367 datatype, src, recvtag, comm, status);
1368 if(retval==MPI_SUCCESS){
1369 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1377 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1381 if (request == NULL || flag == NULL) {
1382 retval = MPI_ERR_ARG;
1383 } else if (*request == MPI_REQUEST_NULL) {
1385 smpi_empty_status(status);
1386 retval = MPI_ERR_REQUEST;
1389 int rank = request && (*request)->comm != MPI_COMM_NULL
1390 ? smpi_process_index()
1393 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1394 extra->type = TRACING_TEST;
1395 TRACE_smpi_testing_in(rank, extra);
1397 *flag = smpi_mpi_test(request, status);
1399 TRACE_smpi_testing_out(rank);
1401 retval = MPI_SUCCESS;
1407 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1408 MPI_Status * status)
1413 if (index == NULL || flag == NULL) {
1414 retval = MPI_ERR_ARG;
1416 *flag = smpi_mpi_testany(count, requests, index, status);
1417 retval = MPI_SUCCESS;
1423 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1429 retval = MPI_ERR_ARG;
1431 *flag = smpi_mpi_testall(count, requests, statuses);
1432 retval = MPI_SUCCESS;
1438 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1442 if (status == NULL) {
1443 retval = MPI_ERR_ARG;
1444 } else if (comm == MPI_COMM_NULL) {
1445 retval = MPI_ERR_COMM;
1446 } else if (source == MPI_PROC_NULL) {
1447 smpi_empty_status(status);
1448 status->MPI_SOURCE = MPI_PROC_NULL;
1449 retval = MPI_SUCCESS;
1451 smpi_mpi_probe(source, tag, comm, status);
1452 retval = MPI_SUCCESS;
1459 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1464 retval = MPI_ERR_ARG;
1465 } else if (status == NULL) {
1466 retval = MPI_ERR_ARG;
1467 } else if (comm == MPI_COMM_NULL) {
1468 retval = MPI_ERR_COMM;
1469 } else if (source == MPI_PROC_NULL) {
1471 smpi_empty_status(status);
1472 status->MPI_SOURCE = MPI_PROC_NULL;
1473 retval = MPI_SUCCESS;
1475 smpi_mpi_iprobe(source, tag, comm, flag, status);
1476 retval = MPI_SUCCESS;
1482 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1488 smpi_empty_status(status);
1490 if (request == NULL) {
1491 retval = MPI_ERR_ARG;
1492 } else if (*request == MPI_REQUEST_NULL) {
1493 retval = MPI_ERR_REQUEST;
1497 int rank = request && (*request)->comm != MPI_COMM_NULL
1498 ? smpi_process_index()
1501 int src_traced = (*request)->src;
1502 int dst_traced = (*request)->dst;
1503 MPI_Comm comm = (*request)->comm;
1504 int is_wait_for_receive = (*request)->recv;
1505 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1506 extra->type = TRACING_WAIT;
1507 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1510 smpi_mpi_wait(request, status);
1511 retval = MPI_SUCCESS;
1514 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1515 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1516 if (is_wait_for_receive) {
1517 if(src_traced==MPI_ANY_SOURCE)
1518 src_traced = (status!=MPI_STATUS_IGNORE) ?
1519 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1521 TRACE_smpi_recv(rank, src_traced, dst_traced);
1531 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1538 //save requests information for tracing
1540 int *srcs = xbt_new0(int, count);
1541 int *dsts = xbt_new0(int, count);
1542 int *recvs = xbt_new0(int, count);
1543 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1545 for (i = 0; i < count; i++) {
1546 MPI_Request req = requests[i]; //already received requests are no longer valid
1550 recvs[i] = req->recv;
1551 comms[i] = req->comm;
1554 int rank_traced = smpi_process_index();
1555 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1556 extra->type = TRACING_WAITANY;
1557 extra->send_size=count;
1558 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1561 *index = smpi_mpi_waitany(count, requests, status);
1563 if(*index!=MPI_UNDEFINED){
1564 int src_traced = srcs[*index];
1565 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1566 int dst_traced = dsts[*index];
1567 int is_wait_for_receive = recvs[*index];
1568 if (is_wait_for_receive) {
1569 if(srcs[*index]==MPI_ANY_SOURCE)
1570 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1571 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) :
1573 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1575 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1587 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1592 //save information from requests
1594 int *srcs = xbt_new0(int, count);
1595 int *dsts = xbt_new0(int, count);
1596 int *recvs = xbt_new0(int, count);
1597 int *valid = xbt_new0(int, count);
1598 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1600 //int valid_count = 0;
1601 for (i = 0; i < count; i++) {
1602 MPI_Request req = requests[i];
1603 if(req!=MPI_REQUEST_NULL){
1606 recvs[i] = req->recv;
1607 comms[i] = req->comm;
1613 int rank_traced = smpi_process_index();
1614 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1615 extra->type = TRACING_WAITALL;
1616 extra->send_size=count;
1617 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1619 int retval = smpi_mpi_waitall(count, requests, status);
1621 for (i = 0; i < count; i++) {
1623 //int src_traced = srcs[*index];
1624 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1625 int src_traced = srcs[i];
1626 int dst_traced = dsts[i];
1627 int is_wait_for_receive = recvs[i];
1628 if (is_wait_for_receive) {
1629 if(src_traced==MPI_ANY_SOURCE)
1630 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1631 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) :
1633 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1637 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1649 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1650 int *indices, MPI_Status status[])
1655 if (outcount == NULL) {
1656 retval = MPI_ERR_ARG;
1658 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1659 retval = MPI_SUCCESS;
1665 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1666 int* indices, MPI_Status status[])
1671 if (outcount == NULL) {
1672 retval = MPI_ERR_ARG;
1674 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1675 retval = MPI_SUCCESS;
1682 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1688 if (comm == MPI_COMM_NULL) {
1689 retval = MPI_ERR_COMM;
1692 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1693 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1695 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1696 extra->type = TRACING_BCAST;
1697 extra->send_size = count;
1698 extra->root = root_traced;
1699 extra->datatype1 = encode_datatype(datatype);
1700 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1703 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1704 retval = MPI_SUCCESS;
1706 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1714 int PMPI_Barrier(MPI_Comm comm)
1720 if (comm == MPI_COMM_NULL) {
1721 retval = MPI_ERR_COMM;
1724 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1725 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1726 extra->type = TRACING_BARRIER;
1727 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1729 mpi_coll_barrier_fun(comm);
1730 retval = MPI_SUCCESS;
1732 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1740 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1741 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1742 int root, MPI_Comm comm)
1748 if (comm == MPI_COMM_NULL) {
1749 retval = MPI_ERR_COMM;
1750 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1751 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1752 retval = MPI_ERR_TYPE;
1753 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1754 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1755 retval = MPI_ERR_COUNT;
1758 char* sendtmpbuf = (char*) sendbuf;
1759 int sendtmpcount = sendcount;
1760 MPI_Datatype sendtmptype = sendtype;
1761 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1763 sendtmptype=recvtype;
1766 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1767 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1768 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1769 extra->type = TRACING_GATHER;
1770 extra->send_size = sendtmpcount;
1771 extra->recv_size = recvcount;
1772 extra->root = root_traced;
1773 extra->datatype1 = encode_datatype(sendtmptype);
1774 extra->datatype2 = encode_datatype(recvtype);
1776 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1778 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1779 recvtype, root, comm);
1782 retval = MPI_SUCCESS;
1784 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1792 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1793 void *recvbuf, int *recvcounts, int *displs,
1794 MPI_Datatype recvtype, int root, MPI_Comm comm)
1800 if (comm == MPI_COMM_NULL) {
1801 retval = MPI_ERR_COMM;
1802 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1803 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1804 retval = MPI_ERR_TYPE;
1805 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1806 retval = MPI_ERR_COUNT;
1807 } else if (recvcounts == NULL || displs == NULL) {
1808 retval = MPI_ERR_ARG;
1810 char* sendtmpbuf = (char*) sendbuf;
1811 int sendtmpcount = sendcount;
1812 MPI_Datatype sendtmptype = sendtype;
1813 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1815 sendtmptype=recvtype;
1819 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1820 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1822 int size = smpi_comm_size(comm);
1823 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1824 extra->type = TRACING_GATHERV;
1825 extra->send_size = sendtmpcount;
1826 extra->recvcounts= xbt_malloc(size*sizeof(int));
1827 for(i=0; i< size; i++)//copy data to avoid bad free
1828 extra->recvcounts[i] = recvcounts[i];
1829 extra->num_processes = size;
1830 extra->root = root_traced;
1831 extra->datatype1 = encode_datatype(sendtmptype);
1832 extra->datatype2 = encode_datatype(recvtype);
1834 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1836 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1837 displs, recvtype, root, comm);
1838 retval = MPI_SUCCESS;
1840 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1848 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1849 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1856 if (comm == MPI_COMM_NULL) {
1857 retval = MPI_ERR_COMM;
1858 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1859 (recvtype == MPI_DATATYPE_NULL)){
1860 retval = MPI_ERR_TYPE;
1861 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1863 retval = MPI_ERR_COUNT;
1865 if(sendbuf == MPI_IN_PLACE) {
1866 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1867 sendcount=recvcount;
1871 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1872 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1873 extra->type = TRACING_ALLGATHER;
1874 extra->send_size = sendcount;
1875 extra->recv_size = recvcount;
1876 extra->datatype1 = encode_datatype(sendtype);
1877 extra->datatype2 = encode_datatype(recvtype);
1879 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1881 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1883 retval = MPI_SUCCESS;
1886 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1893 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1894 void *recvbuf, int *recvcounts, int *displs,
1895 MPI_Datatype recvtype, MPI_Comm comm)
1901 if (comm == MPI_COMM_NULL) {
1902 retval = MPI_ERR_COMM;
1903 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1904 (recvtype == MPI_DATATYPE_NULL)){
1905 retval = MPI_ERR_TYPE;
1906 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1907 retval = MPI_ERR_COUNT;
1908 } else if (recvcounts == NULL || displs == NULL) {
1909 retval = MPI_ERR_ARG;
1912 if(sendbuf == MPI_IN_PLACE) {
1913 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1914 sendcount=recvcounts[smpi_comm_rank(comm)];
1918 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1920 int size = smpi_comm_size(comm);
1921 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1922 extra->type = TRACING_ALLGATHERV;
1923 extra->send_size = sendcount;
1924 extra->recvcounts= xbt_malloc(size*sizeof(int));
1925 for(i=0; i< size; i++)//copy data to avoid bad free
1926 extra->recvcounts[i] = recvcounts[i];
1927 extra->num_processes = size;
1928 extra->datatype1 = encode_datatype(sendtype);
1929 extra->datatype2 = encode_datatype(recvtype);
1931 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
1933 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1934 displs, recvtype, comm);
1935 retval = MPI_SUCCESS;
1937 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1945 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1946 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1947 int root, MPI_Comm comm)
1953 if (comm == MPI_COMM_NULL) {
1954 retval = MPI_ERR_COMM;
1955 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1956 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1957 retval = MPI_ERR_TYPE;
1960 if (recvbuf == MPI_IN_PLACE) {
1962 recvcount=sendcount;
1965 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1966 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1967 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1968 extra->type = TRACING_SCATTER;
1969 extra->send_size = sendcount;
1970 extra->recv_size= recvcount;
1971 extra->root = root_traced;
1972 extra->datatype1 = encode_datatype(sendtype);
1973 extra->datatype2 = encode_datatype(recvtype);
1975 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1977 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1978 recvtype, root, comm);
1979 retval = MPI_SUCCESS;
1981 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1989 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1990 MPI_Datatype sendtype, void *recvbuf, int recvcount,
1991 MPI_Datatype recvtype, int root, MPI_Comm comm)
1997 if (comm == MPI_COMM_NULL) {
1998 retval = MPI_ERR_COMM;
1999 } else if (sendcounts == NULL || displs == NULL) {
2000 retval = MPI_ERR_ARG;
2001 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
2002 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
2003 retval = MPI_ERR_TYPE;
2005 if (recvbuf == MPI_IN_PLACE) {
2007 recvcount=sendcounts[smpi_comm_rank(comm)];
2010 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2011 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2013 int size = smpi_comm_size(comm);
2014 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2015 extra->type = TRACING_SCATTERV;
2016 extra->recv_size = recvcount;
2017 extra->sendcounts= xbt_malloc(size*sizeof(int));
2018 for(i=0; i< size; i++)//copy data to avoid bad free
2019 extra->sendcounts[i] = sendcounts[i];
2020 extra->num_processes = size;
2021 extra->root = root_traced;
2022 extra->datatype1 = encode_datatype(sendtype);
2023 extra->datatype2 = encode_datatype(recvtype);
2025 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2028 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
2029 recvcount, recvtype, root, comm);
2030 retval = MPI_SUCCESS;
2032 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2040 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
2041 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
2047 if (comm == MPI_COMM_NULL) {
2048 retval = MPI_ERR_COMM;
2049 } else if (!is_datatype_valid(datatype) || op == MPI_OP_NULL) {
2050 retval = MPI_ERR_ARG;
2053 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2054 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2055 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2056 extra->type = TRACING_REDUCE;
2057 extra->send_size = count;
2058 extra->datatype1 = encode_datatype(datatype);
2059 extra->root = root_traced;
2061 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2063 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2065 retval = MPI_SUCCESS;
2067 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2075 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
2076 MPI_Datatype datatype, MPI_Op op){
2080 if (!is_datatype_valid(datatype) || op == MPI_OP_NULL) {
2081 retval = MPI_ERR_ARG;
2083 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2090 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
2091 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2097 if (comm == MPI_COMM_NULL) {
2098 retval = MPI_ERR_COMM;
2099 } else if (!is_datatype_valid(datatype)) {
2100 retval = MPI_ERR_TYPE;
2101 } else if (op == MPI_OP_NULL) {
2102 retval = MPI_ERR_OP;
2105 char* sendtmpbuf = (char*) sendbuf;
2106 if( sendbuf == MPI_IN_PLACE ) {
2107 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2108 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2111 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2112 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2113 extra->type = TRACING_ALLREDUCE;
2114 extra->send_size = count;
2115 extra->datatype1 = encode_datatype(datatype);
2117 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2119 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2121 if( sendbuf == MPI_IN_PLACE ) {
2122 xbt_free(sendtmpbuf);
2125 retval = MPI_SUCCESS;
2127 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2135 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2136 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2142 if (comm == MPI_COMM_NULL) {
2143 retval = MPI_ERR_COMM;
2144 } else if (!is_datatype_valid(datatype)) {
2145 retval = MPI_ERR_TYPE;
2146 } else if (op == MPI_OP_NULL) {
2147 retval = MPI_ERR_OP;
2150 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2151 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2152 extra->type = TRACING_SCAN;
2153 extra->send_size = count;
2154 extra->datatype1 = encode_datatype(datatype);
2156 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2158 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2159 retval = MPI_SUCCESS;
2161 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2169 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2170 MPI_Op op, MPI_Comm comm){
2175 if (comm == MPI_COMM_NULL) {
2176 retval = MPI_ERR_COMM;
2177 } else if (!is_datatype_valid(datatype)) {
2178 retval = MPI_ERR_TYPE;
2179 } else if (op == MPI_OP_NULL) {
2180 retval = MPI_ERR_OP;
2183 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2184 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2185 extra->type = TRACING_EXSCAN;
2186 extra->send_size = count;
2187 extra->datatype1 = encode_datatype(datatype);
2189 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2191 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2192 retval = MPI_SUCCESS;
2194 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2202 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2203 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2208 if (comm == MPI_COMM_NULL) {
2209 retval = MPI_ERR_COMM;
2210 } else if (!is_datatype_valid(datatype)) {
2211 retval = MPI_ERR_TYPE;
2212 } else if (op == MPI_OP_NULL) {
2213 retval = MPI_ERR_OP;
2214 } else if (recvcounts == NULL) {
2215 retval = MPI_ERR_ARG;
2218 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2220 int size = smpi_comm_size(comm);
2221 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2222 extra->type = TRACING_REDUCE_SCATTER;
2223 extra->send_size = 0;
2224 extra->recvcounts= xbt_malloc(size*sizeof(int));
2225 for(i=0; i< size; i++)//copy data to avoid bad free
2226 extra->recvcounts[i] = recvcounts[i];
2227 extra->num_processes = size;
2228 extra->datatype1 = encode_datatype(datatype);
2230 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2232 void* sendtmpbuf=sendbuf;
2233 if(sendbuf==MPI_IN_PLACE){
2237 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2238 datatype, op, comm);
2239 retval = MPI_SUCCESS;
2241 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2249 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2250 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2255 if (comm == MPI_COMM_NULL) {
2256 retval = MPI_ERR_COMM;
2257 } else if (!is_datatype_valid(datatype)) {
2258 retval = MPI_ERR_TYPE;
2259 } else if (op == MPI_OP_NULL) {
2260 retval = MPI_ERR_OP;
2261 } else if (recvcount < 0) {
2262 retval = MPI_ERR_ARG;
2264 int count=smpi_comm_size(comm);
2267 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2268 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2269 extra->type = TRACING_REDUCE_SCATTER;
2270 extra->send_size = 0;
2271 extra->recvcounts= xbt_malloc(count*sizeof(int));
2272 for(i=0; i< count; i++)//copy data to avoid bad free
2273 extra->recvcounts[i] = recvcount;
2274 extra->num_processes = count;
2275 extra->datatype1 = encode_datatype(datatype);
2277 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2279 int* recvcounts=(int*)xbt_malloc(count);
2280 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2281 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2282 datatype, op, comm);
2283 xbt_free(recvcounts);
2284 retval = MPI_SUCCESS;
2286 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2294 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2295 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2302 if (comm == MPI_COMM_NULL) {
2303 retval = MPI_ERR_COMM;
2304 } else if (sendtype == MPI_DATATYPE_NULL
2305 || recvtype == MPI_DATATYPE_NULL) {
2306 retval = MPI_ERR_TYPE;
2309 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2310 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2311 extra->type = TRACING_ALLTOALL;
2312 extra->send_size = sendcount;
2313 extra->recv_size = recvcount;
2314 extra->datatype1 = encode_datatype(sendtype);
2315 extra->datatype2 = encode_datatype(recvtype);
2317 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2319 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2321 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2329 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2330 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2331 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2337 if (comm == MPI_COMM_NULL) {
2338 retval = MPI_ERR_COMM;
2339 } else if (sendtype == MPI_DATATYPE_NULL
2340 || recvtype == MPI_DATATYPE_NULL) {
2341 retval = MPI_ERR_TYPE;
2342 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2343 || recvdisps == NULL) {
2344 retval = MPI_ERR_ARG;
2347 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2349 int size = smpi_comm_size(comm);
2350 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2351 extra->type = TRACING_ALLTOALLV;
2352 extra->send_size = 0;
2353 extra->recv_size = 0;
2354 extra->recvcounts= xbt_malloc(size*sizeof(int));
2355 extra->sendcounts= xbt_malloc(size*sizeof(int));
2357 for(i=0; i< size; i++){//copy data to avoid bad free
2358 extra->send_size += sendcounts[i];
2359 extra->recv_size += recvcounts[i];
2361 extra->sendcounts[i] = sendcounts[i];
2362 extra->recvcounts[i] = recvcounts[i];
2364 extra->num_processes = size;
2366 extra->datatype1 = encode_datatype(sendtype);
2367 extra->datatype2 = encode_datatype(recvtype);
2369 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2372 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2373 recvbuf, recvcounts, recvdisps, recvtype,
2376 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2385 int PMPI_Get_processor_name(char *name, int *resultlen)
2387 int retval = MPI_SUCCESS;
2389 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2390 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2391 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2392 MPI_MAX_PROCESSOR_NAME - 1 );
2395 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2400 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2402 int retval = MPI_SUCCESS;
2405 if (status == NULL || count == NULL) {
2406 retval = MPI_ERR_ARG;
2407 } else if (!is_datatype_valid(datatype)) {
2408 retval = MPI_ERR_TYPE;
2410 size = smpi_datatype_size(datatype);
2413 } else if (status->count % size != 0) {
2414 retval = MPI_UNDEFINED;
2416 *count = smpi_mpi_get_count(status, datatype);
2422 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2425 if (old_type == MPI_DATATYPE_NULL) {
2426 retval = MPI_ERR_TYPE;
2427 } else if (count<0){
2428 retval = MPI_ERR_COUNT;
2430 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2435 int PMPI_Type_commit(MPI_Datatype* datatype) {
2438 if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2439 retval = MPI_ERR_TYPE;
2441 smpi_datatype_commit(datatype);
2442 retval = MPI_SUCCESS;
2448 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2451 if (old_type == MPI_DATATYPE_NULL) {
2452 retval = MPI_ERR_TYPE;
2453 } else if (count<0 || blocklen<0){
2454 retval = MPI_ERR_COUNT;
2456 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2461 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2464 if (old_type == MPI_DATATYPE_NULL) {
2465 retval = MPI_ERR_TYPE;
2466 } else if (count<0 || blocklen<0){
2467 retval = MPI_ERR_COUNT;
2469 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2474 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2475 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2478 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2481 if (old_type == MPI_DATATYPE_NULL) {
2482 retval = MPI_ERR_TYPE;
2483 } else if (count<0){
2484 retval = MPI_ERR_COUNT;
2486 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2491 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2494 if (old_type == MPI_DATATYPE_NULL) {
2495 retval = MPI_ERR_TYPE;
2496 } else if (count<0){
2497 retval = MPI_ERR_COUNT;
2499 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2504 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2507 if (old_type == MPI_DATATYPE_NULL) {
2508 retval = MPI_ERR_TYPE;
2509 } else if (count<0){
2510 retval = MPI_ERR_COUNT;
2512 int* blocklens=(int*)xbt_malloc(blocklength*count);
2513 for (i=0; i<count;i++)blocklens[i]=blocklength;
2514 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2515 xbt_free(blocklens);
2521 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2524 if (old_type == MPI_DATATYPE_NULL) {
2525 retval = MPI_ERR_TYPE;
2526 } else if (count<0){
2527 retval = MPI_ERR_COUNT;
2529 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2534 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2535 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2538 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2541 if (old_type == MPI_DATATYPE_NULL) {
2542 retval = MPI_ERR_TYPE;
2543 } else if (count<0){
2544 retval = MPI_ERR_COUNT;
2546 int* blocklens=(int*)xbt_malloc(blocklength*count);
2547 for (i=0; i<count;i++)blocklens[i]=blocklength;
2548 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2549 xbt_free(blocklens);
2555 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2559 retval = MPI_ERR_COUNT;
2561 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2566 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2567 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2571 int PMPI_Error_class(int errorcode, int* errorclass) {
2572 // assume smpi uses only standard mpi error codes
2573 *errorclass=errorcode;
2578 int PMPI_Initialized(int* flag) {
2579 *flag=smpi_process_initialized();
2583 /* The topo part of MPI_COMM_WORLD should always be NULL. When other topologies
2584 * will be implemented, not only should we check if the topology is NULL, but
2585 * we should check if it is the good topology type (so we have to add a
2586 * MPIR_Topo_Type field, and replace the MPI_Topology field by an union)*/
2588 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periodic, int reorder, MPI_Comm* comm_cart) {
2591 if (comm_old == MPI_COMM_NULL){
2592 return MPI_ERR_COMM;
2594 else if (ndims < 0 ||
2595 (ndims > 0 && (dims == NULL ||
2596 periodic == NULL)) ||
2597 comm_cart == NULL) {
2600 retval = smpi_mpi_cart_create(comm_old, ndims, dims, periodic, reorder, comm_cart);
2607 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2608 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2609 return MPI_ERR_TOPOLOGY;
2611 if (coords == NULL) {
2614 return smpi_mpi_cart_rank(comm, coords, rank);
2617 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2618 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2619 return MPI_ERR_TOPOLOGY;
2621 if (source == NULL || dest == NULL || direction < 0 ) {
2624 return smpi_mpi_cart_shift(comm, direction, displ, source, dest);
2627 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2628 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2629 return MPI_ERR_TOPOLOGY;
2631 if (rank < 0 || rank >= smpi_comm_size(comm)) {
2632 return MPI_ERR_RANK;
2637 if(coords == NULL) {
2640 return smpi_mpi_cart_coords(comm, rank, maxdims, coords);
2643 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2644 if(comm == NULL || smpi_comm_topo(comm) == NULL) {
2645 return MPI_ERR_TOPOLOGY;
2647 if(maxdims <= 0 || dims == NULL || periods == NULL || coords == NULL) {
2650 return smpi_mpi_cart_get(comm, maxdims, dims, periods, coords);
2653 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2654 if (comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2655 return MPI_ERR_TOPOLOGY;
2657 if (ndims == NULL) {
2660 return smpi_mpi_cartdim_get(comm, ndims);
2663 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2667 if (ndims < 1 || nnodes < 1) {
2668 return MPI_ERR_DIMS;
2671 return smpi_mpi_dims_create(nnodes, ndims, dims);
2674 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2675 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2676 return MPI_ERR_TOPOLOGY;
2678 if (comm_new == NULL) {
2681 return smpi_mpi_cart_sub(comm, remain_dims, comm_new);
2684 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2685 if(oldtype == MPI_DATATYPE_NULL) {
2686 return MPI_ERR_TYPE;
2688 int blocks[3] = { 1, 1, 1 };
2689 MPI_Aint disps[3] = { lb, 0, lb+extent };
2690 MPI_Datatype types[3] = { MPI_LB, oldtype, MPI_UB };
2692 s_smpi_mpi_struct_t* subtype = smpi_datatype_struct_create( blocks,
2697 smpi_datatype_create(newtype,oldtype->size, lb, lb + extent, 1 , subtype, DT_FLAG_VECTOR);
2699 (*newtype)->flags &= ~DT_FLAG_COMMITED;
2704 /* The following calls are not yet implemented and will fail at runtime. */
2705 /* Once implemented, please move them above this notice. */
2707 #define NOT_YET_IMPLEMENTED { \
2708 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__); \
2709 return MPI_SUCCESS; \
2712 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2716 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2721 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2726 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2731 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2736 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2740 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2744 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2748 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2752 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2756 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2760 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2764 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2768 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2772 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2776 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2780 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2784 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2788 int PMPI_Comm_get_errhandler(MPI_Comm comm, MPI_Errhandler* errhandler) {
2792 int PMPI_Cancel(MPI_Request* request) {
2796 int PMPI_Buffer_attach(void* buffer, int size) {
2800 int PMPI_Buffer_detach(void* buffer, int* size) {
2804 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2808 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2813 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2818 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2823 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2828 int PMPI_Comm_free_keyval(int* keyval) {
2832 int PMPI_Pcontrol(const int level )
2837 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2841 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2846 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2851 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2856 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2861 int PMPI_Type_free_keyval(int* keyval) {
2865 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2869 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2873 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2877 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2881 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2885 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2889 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2893 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2897 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2901 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2905 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2909 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2913 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2917 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2921 int PMPI_Keyval_free(int* keyval) {
2925 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2929 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2933 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2937 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2941 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2945 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2949 int PMPI_Win_fence( int assert, MPI_Win win){
2953 int PMPI_Win_free( MPI_Win* win){
2957 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2961 int PMPI_Info_create( MPI_Info *info){
2965 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2969 int PMPI_Info_free( MPI_Info *info){
2973 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2974 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2978 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2979 int *num_addresses, int *num_datatypes, int *combiner){
2983 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2984 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2985 MPI_Datatype* array_of_datatypes){
2989 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2990 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2991 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2995 int PMPI_Type_create_subarray(int ndims,int *array_of_sizes, int *array_of_subsizes, int *array_of_starts, int order, MPI_Datatype oldtype, MPI_Datatype *newtype){
2999 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
3003 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
3004 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
3009 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
3013 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
3017 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
3021 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
3025 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
3029 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
3033 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
3037 int PMPI_Add_error_class( int *errorclass){
3041 int PMPI_Add_error_code( int errorclass, int *errorcode){
3045 int PMPI_Add_error_string( int errorcode, char *string){
3049 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
3053 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
3057 int PMPI_Info_delete(MPI_Info info, char *key){
3061 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
3065 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
3069 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
3073 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
3077 int PMPI_Grequest_start( MPI_Grequest_query_function *query_fn, MPI_Grequest_free_function *free_fn, MPI_Grequest_cancel_function *cancel_fn, void *extra_state, MPI_Request *request){
3081 int PMPI_Grequest_complete( MPI_Request request){
3085 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
3089 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
3093 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3097 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
3101 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
3105 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
3109 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
3113 int PMPI_Open_port( MPI_Info info, char *port_name){
3117 int PMPI_Close_port(char *port_name){
3121 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3125 int PMPI_Comm_spawn( char *command, char **argv, int maxprocs, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
3129 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
3130 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
3131 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
3135 int PMPI_Comm_get_parent( MPI_Comm *parent){